Geant4  10.01
G4WentzelVIModel.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: G4WentzelVIModel.cc 85306 2014-10-27 14:17:47Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4WentzelVIModel
34 //
35 // Author: V.Ivanchenko
36 //
37 // Creation date: 09.04.2008 from G4MuMscModel
38 //
39 // Modifications:
40 // 27-05-2010 V.Ivanchenko added G4WentzelOKandVIxSection class to
41 // compute cross sections and sample scattering angle
42 //
43 //
44 // Class Description:
45 //
46 // Implementation of the model of multiple scattering based on
47 // G.Wentzel, Z. Phys. 40 (1927) 590.
48 // H.W.Lewis, Phys Rev 78 (1950) 526.
49 // J.M. Fernandez-Varea et al., NIM B73 (1993) 447.
50 // L.Urban, CERN-OPEN-2006-077.
51 
52 // -------------------------------------------------------------------
53 //
54 
55 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57 
58 #include "G4WentzelVIModel.hh"
59 #include "G4PhysicalConstants.hh"
60 #include "G4SystemOfUnits.hh"
61 #include "Randomize.hh"
63 #include "G4PhysicsTableHelper.hh"
64 #include "G4ElementVector.hh"
65 #include "G4ProductionCutsTable.hh"
66 #include "G4Log.hh"
67 #include "G4Exp.hh"
68 
69 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
70 
71 using namespace std;
72 
73 const G4int minNCollisions = 10;
74 
76  G4VMscModel(nam),
77  ssFactor(1.05),
78  invssFactor(1.0),
79  currentCouple(0),
80  inside(false),
81  singleScatteringMode(false),
82  cosThetaMin(1.0),
83  cosThetaMax(-1.0),
84  fSecondMoments(0),
85  idx2(0),
86  numlimit(0.1),
87  isCombined(combined),
88  useSecondMoment(false)
89 {
91  invsqrt12 = 1./sqrt(12.);
92  tlimitminfix = 1.e-6*mm;
93  lowEnergyLimit = 1.0*eV;
94  particle = 0;
95  nelments = 5;
96  xsecn.resize(nelments);
97  prob.resize(nelments);
98  wokvi = new G4WentzelOKandVIxSection(combined);
99  fixedCut = -1.0;
100 
102  = currentRange = xtsec = cosTetMaxNuc = 0.0;
104 
105  fParticleChange = 0;
106  currentCuts = 0;
107  currentMaterial = 0;
108 }
109 
110 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
111 
113 {
114  delete wokvi;
115  if(fSecondMoments && IsMaster()) {
116  delete fSecondMoments;
117  fSecondMoments = 0;
118  }
119 }
120 
121 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
122 
124  const G4DataVector& cuts)
125 {
126  // reset parameters
127  SetupParticle(p);
128  currentRange = 0.0;
129 
130  if(isCombined) {
131  G4double tet = PolarAngleLimit();
132  if(tet <= 0.0) { cosThetaMax = 1.0; }
133  else if(tet < CLHEP::pi) { cosThetaMax = cos(tet); }
134  }
135 
136  //G4cout << "G4WentzelVIModel::Initialise " << p->GetParticleName() << G4endl;
138  /*
139  G4cout << "G4WentzelVIModel: " << particle->GetParticleName()
140  << " 1-cos(ThetaLimit)= " << 1 - cosThetaMax
141  << " SingScatFactor= " << ssFactor
142  << G4endl;
143  */
144  currentCuts = &cuts;
145 
146  // set values of some data members
148 
149  // build second moment table only if transport table is build
151  if(useSecondMoment && IsMaster() && table) {
152  fSecondMoments =
154  // Access to materials
155  const G4ProductionCutsTable* theCoupleTable =
157  size_t numOfCouples = theCoupleTable->GetTableSize();
158 
159  G4bool splineFlag = true;
160  G4PhysicsVector* aVector = 0;
161  for(size_t i=0; i<numOfCouples; ++i) {
162 
163  //G4cout<< "i= " << i << " Flag= " << GetFlag(i) << G4endl;
164  if(table->GetFlag(i)) {
165  DefineMaterial(theCoupleTable->GetMaterialCutsCouple(i));
166 
167  delete (*fSecondMoments)[i];
168  aVector = new G4PhysicsVector((*table)[i]);
169  size_t n = aVector->GetVectorLength();
170 
171  for(size_t j=0; j<n; ++j) {
172  G4double e = aVector->Energy(j);
173  aVector->PutValue(j, ComputeSecondMoment(p, e)*e*e);
174  }
175  if(splineFlag) { aVector->FillSecondDerivatives(); }
177  }
178  }
179  }
180 }
181 
182 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
183 
185  G4VEmModel* masterModel)
186 {
187  fSecondMoments = static_cast<G4WentzelVIModel*>(masterModel)
189 }
190 
191 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
192 
194  const G4ParticleDefinition* p,
195  G4double kinEnergy,
196  G4double Z, G4double,
197  G4double cutEnergy, G4double)
198 {
199  G4double cross = 0.0;
200  if(p != particle) { SetupParticle(p); }
201  if(kinEnergy < lowEnergyLimit) { return cross; }
202  if(!CurrentCouple()) {
203  G4Exception("G4WentzelVIModel::ComputeCrossSectionPerAtom", "em0011",
204  FatalException, " G4MaterialCutsCouple is not defined");
205  return 0.0;
206  }
209  if(cosTetMaxNuc < 1.0) {
210  G4double cut = cutEnergy;
211  if(fixedCut > 0.0) { cut = fixedCut; }
212  G4double cost = wokvi->SetupTarget(G4lrint(Z), cut);
214  /*
215  if(p->GetParticleName() == "e-")
216  G4cout << "G4WentzelVIModel::CS: Z= " << G4int(Z) << " e(MeV)= "<<kinEnergy
217  << " 1-cosN= " << 1 - cosTetMaxNuc << " cross(bn)= " << cross/barn
218  << " " << particle->GetParticleName() << G4endl;
219  */
220  }
221  return cross;
222 }
223 
224 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
225 
227 {
229  inside = false;
230 }
231 
232 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
233 
235  const G4Track& track,
236  G4double& currentMinimalStep)
237 {
238  G4double tlimit = currentMinimalStep;
239  const G4DynamicParticle* dp = track.GetDynamicParticle();
240  G4StepPoint* sp = track.GetStep()->GetPreStepPoint();
241  G4StepStatus stepStatus = sp->GetStepStatus();
242  singleScatteringMode = false;
243 
244  //G4cout << "G4WentzelVIModel::ComputeTruePathLengthLimit stepStatus= "
245  // << stepStatus << " " << track.GetDefinition()->GetParticleName()
246  // << G4endl;
247 
248  // initialisation for each step, lambda may be computed from scratch
255 
256  //G4cout << "lambdaeff= " << lambdaeff << " Range= " << currentRange
257  // << " tlimit= " << tlimit << " 1-cost= " << 1 - cosTetMaxNuc << G4endl;
258 
259  // extra check for abnormal situation
260  // this check needed to run MSC with eIoni and eBrem inactivated
261  if(tlimit > currentRange) { tlimit = currentRange; }
262 
263  // stop here if small range particle
264  if(inside || tlimit < tlimitminfix) {
265  return ConvertTrueToGeom(tlimit, currentMinimalStep);
266  }
267 
268  // pre step
269  G4double presafety = sp->GetSafety();
270  // far from geometry boundary
271  if(currentRange < presafety) {
272  inside = true;
273  return ConvertTrueToGeom(tlimit, currentMinimalStep);
274  }
275 
276  // compute presafety again if presafety <= 0 and no boundary
277  // i.e. when it is needed for optimization purposes
278  if(stepStatus != fGeomBoundary && presafety < tlimitminfix) {
279  presafety = ComputeSafety(sp->GetPosition(), tlimit);
280  if(currentRange < presafety) {
281  inside = true;
282  return ConvertTrueToGeom(tlimit, currentMinimalStep);
283  }
284  }
285  /*
286  G4cout << "e(MeV)= " << preKinEnergy/MeV
287  << " " << particle->GetParticleName()
288  << " CurLimit(mm)= " << tlimit/mm <<" safety(mm)= " << presafety/mm
289  << " R(mm)= " <<currentRange/mm
290  << " L0(mm^-1)= " << lambdaeff*mm
291  << G4endl;
292  */
293  // natural limit for high energy
295  (1.0 - cosTetMaxNuc)*lambdaeff*invssFactor);
296 
297  // low-energy e-
298  if(cosThetaMax > cosTetMaxNuc) {
299  rlimit = std::min(rlimit, facsafety*presafety);
300  }
301 
302  // cut correction
304  //G4cout << "rcut= " << rcut << " rlimit= " << rlimit << " presafety= "
305  // << presafety << " 1-cosThetaMax= " <<1-cosThetaMax
306  //<< " 1-cosTetMaxNuc= " << 1-cosTetMaxNuc << G4endl;
307  if(rcut > rlimit) { rlimit = std::min(rlimit, rcut*sqrt(rlimit/rcut)); }
308 
309  tlimit = std::min(tlimit, rlimit);
310  tlimit = std::max(tlimit, tlimitminfix);
311 
312  // step limit in infinite media
313  tlimit = std::min(tlimit, 50*currentMaterial->GetRadlen()/facgeom);
314 
315  //compute geomlimit and force few steps within a volume
317  && stepStatus == fGeomBoundary) {
318 
319  G4double geomlimit = ComputeGeomLimit(track, presafety, currentRange);
320  tlimit = std::min(tlimit, geomlimit/facgeom);
321  }
322  /*
323  G4cout << particle->GetParticleName() << " e= " << preKinEnergy
324  << " L0= " << lambdaeff << " R= " << currentRange
325  << " tlimit= " << tlimit
326  << " currentMinimalStep= " << currentMinimalStep << G4endl;
327  */
328  return ConvertTrueToGeom(tlimit, currentMinimalStep);
329 }
330 
331 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
332 
334 {
335  zPathLength = tPathLength = truelength;
336 
337  // small step use only single scattering
338  cosThetaMin = 1.0;
340  //G4cout << "xtsec= " << xtsec << " Nav= "
341  // << zPathLength*xtsec << G4endl;
342  if(0.0 >= lambdaeff || G4int(zPathLength*xtsec) < minNCollisions) {
343  singleScatteringMode = true;
344  lambdaeff = DBL_MAX;
345 
346  } else {
347  //G4cout << "ComputeGeomPathLength: tLength= " << tPathLength
348  // << " Leff= " << lambdaeff << G4endl;
349  // small step
352  zPathLength *= (1.0 - 0.5*tau + tau*tau/6.0);
353 
354  // medium step
355  } else {
356  G4double e1 = 0.0;
357  if(currentRange > tPathLength) {
359  }
360  effKinEnergy = 0.5*(e1 + preKinEnergy);
363  //G4cout << " tLength= "<< tPathLength<< " Leff= " << lambdaeff << G4endl;
365  if(tPathLength*numlimit < lambdaeff) {
366  zPathLength *= (1.0 - G4Exp(-tPathLength/lambdaeff));
367  }
368  }
369  }
370  //G4cout << "Comp.geom: zLength= "<<zPathLength<<" tLength= "
371  // << tPathLength<< " Leff= " << lambdaeff << G4endl;
372  return zPathLength;
373 }
374 
375 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
376 
378 {
379  // initialisation of single scattering x-section
380  /*
381  G4cout << "ComputeTrueStepLength: Step= " << geomStepLength
382  << " geomL= " << zPathLength
383  << " Lambda= " << lambdaeff
384  << " 1-cosThetaMaxNuc= " << 1 - cosTetMaxNuc << G4endl;
385  */
387  zPathLength = tPathLength = geomStepLength;
388 
389  } else {
390 
391  // step defined by transportation
392  // change both geom and true step lengths
393  if(geomStepLength < zPathLength) {
394 
395  // single scattering
396  if(G4int(geomStepLength*xtsec) < minNCollisions) {
397  zPathLength = tPathLength = geomStepLength;
398  lambdaeff = DBL_MAX;
399  singleScatteringMode = true;
400 
401  // multiple scattering
402  } else {
403  // small step
404  if(geomStepLength < numlimit*lambdaeff) {
405  G4double tau = geomStepLength/lambdaeff;
406  tPathLength = geomStepLength*(1.0 + 0.5*tau + tau*tau/3.0);
407 
408  // energy correction for a big step
409  } else {
410  tPathLength *= geomStepLength/zPathLength;
411  G4double e1 = 0.0;
412  if(currentRange > tPathLength) {
414  }
415  effKinEnergy = 0.5*(e1 + preKinEnergy);
418  G4double tau = geomStepLength/lambdaeff;
419 
420  if(tau < 0.999999) { tPathLength = -lambdaeff*G4Log(1.0 - tau); }
421  else { tPathLength = currentRange; }
422  }
423  zPathLength = geomStepLength;
424  }
425  }
426  }
427  // check of step length
428  // define threshold angle between single and multiple scattering
429  if(!singleScatteringMode) {
431  xtsec = 0.0;
432 
433  // recompute transport cross section - do not change energy
434  // anymore - cannot be applied for big steps
435  if(cosThetaMin > cosTetMaxNuc) {
436  // new computation
438  //G4cout << "%%%% cross= " << cross << " xtsec= " << xtsec
439  // << " 1-cosTMin= " << 1.0 - cosThetaMin << G4endl;
440  if(cross <= 0.0) {
441  singleScatteringMode = true;
443  lambdaeff = DBL_MAX;
444  cosThetaMin = 1.0;
445  } else if(xtsec > 0.0) {
446 
447  lambdaeff = 1./cross;
448  G4double tau = zPathLength*cross;
449  if(tau < numlimit) {
450  tPathLength = zPathLength*(1.0 + 0.5*tau + tau*tau/3.0);
451  } else if(tau < 0.999999) {
452  tPathLength = -lambdaeff*G4Log(1.0 - tau);
453  } else {
455  }
456  }
457  }
458  }
460  /*
461  G4cout <<"Comp.true: zLength= "<<zPathLength<<" tLength= "<<tPathLength
462  <<" Leff(mm)= "<<lambdaeff/mm<<" sig0(1/mm)= " << xtsec <<G4endl;
463  G4cout << particle->GetParticleName() << " 1-cosThetaMin= " << 1-cosThetaMin
464  << " 1-cosTetMaxNuc= " << 1-cosTetMaxNuc
465  << " e(MeV)= " << preKinEnergy/MeV << " "
466  << " SSmode= " << singleScatteringMode << G4endl;
467  */
468  return tPathLength;
469 }
470 
471 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
472 
475  G4double /*safety*/)
476 {
477  fDisplacement.set(0.0,0.0,0.0);
478  //G4cout << "!##! G4WentzelVIModel::SampleScattering for "
479  // << particle->GetParticleName() << G4endl;
480 
481  // ignore scattering for zero step length and energy below the limit
482  if(preKinEnergy < lowEnergyLimit || tPathLength <= 0.0)
483  { return fDisplacement; }
484 
485  G4double invlambda = 0.0;
486  if(lambdaeff < DBL_MAX) { invlambda = 0.5/lambdaeff; }
487 
488  // use average kinetic energy over the step
489  G4double cut = (*currentCuts)[currentMaterialIndex];
490  if(fixedCut > 0.0) { cut = fixedCut; }
491  /*
492  G4cout <<"SampleScat: E0(MeV)= "<< preKinEnergy/MeV
493  << " Leff= " << lambdaeff <<" sig0(1/mm)= " << xtsec
494  << " xmsc= " << tPathLength*invlambda
495  << " safety= " << safety << G4endl;
496  */
497  // step limit due msc
498  G4int nMscSteps = 1;
499  G4double x0 = tPathLength;
500  G4double z0 = x0*invlambda;
501  //G4double zzz = 0.0;
502  G4double prob2 = 0.0;
503 
504  // large scattering angle case - two step approach
505 
506  if(!singleScatteringMode) {
507  static const G4double zzmin = 0.05;
508  if(useSecondMoment) {
509  G4double z1 = invlambda*invlambda;
511  prob2 = (z2 - z1)/(1.5*z1 - z2);
512  }
513  // if(z0 > zzmin && safety > tlimitminfix) {
514  if(z0 > zzmin) {
515  x0 *= 0.5;
516  z0 *= 0.5;
517  nMscSteps = 2;
518  }
519  //if(z0 > zzmin) { zzz = G4Exp(-1.0/z0); }
520  G4double zzz = 0.0;
521  if(z0 > zzmin) {
522  zzz = G4Exp(-1.0/z0);
523  z0 += zzz;
524  prob2 *= (1 + zzz);
525  }
526  prob2 /= (1 + prob2);
527  }
528 
529  // step limit due to single scattering
530  G4double x1 = 2*tPathLength;
531  if(0.0 < xtsec) { x1 = -G4Log(G4UniformRand())/xtsec; }
532 
533  // no scattering case
534  if(singleScatteringMode && x1 > tPathLength)
535  { return fDisplacement; }
536 
537  const G4ElementVector* theElementVector =
540 
541  // geometry
542  G4double sint, cost, phi;
543  G4ThreeVector temp(0.0,0.0,1.0);
544 
545  // current position and direction relative to the end point
546  // because of magnetic field geometry is computed relatively to the
547  // end point of the step
548  G4ThreeVector dir(0.0,0.0,1.0);
549  fDisplacement.set(0.0,0.0,-zPathLength);
550 
552 
553  // start a loop
554  G4double x2 = x0;
555  G4double step, z;
556  G4bool singleScat;
557  /*
558  G4cout << "Start of the loop x1(mm)= " << x1 << " x2(mm)= " << x2
559  << " 1-cost1= " << 1 - cosThetaMin << " SSmode= " << singleScatteringMode
560  << " xtsec= " << xtsec << " Nst= " << nMscSteps << G4endl;
561  */
562  do {
563 
564  //G4cout << "# x1(mm)= "<< x1<< " x2(mm)= "<< x2 << G4endl;
565  // single scattering case
566  if(singleScatteringMode && x1 > x2) {
567  fDisplacement += x2*mscfac*dir;
568  break;
569  }
570 
571  // what is next single of multiple?
572  if(x1 <= x2) {
573  step = x1;
574  singleScat = true;
575  } else {
576  step = x2;
577  singleScat = false;
578  }
579 
580  //G4cout << "# step(mm)= "<< step<< " singlScat= "<< singleScat << G4endl;
581 
582  // new position
583  fDisplacement += step*mscfac*dir;
584 
585  if(singleScat) {
586 
587  // select element
588  G4int i = 0;
589  if(nelm > 1) {
590  G4double qsec = G4UniformRand()*xtsec;
591  for (; i<nelm; ++i) { if(xsecn[i] >= qsec) { break; } }
592  }
593  G4double cosTetM =
594  wokvi->SetupTarget(G4lrint((*theElementVector)[i]->GetZ()), cut);
595  //G4cout << "!!! " << cosThetaMin << " " << cosTetM << " "
596  // << prob[i] << G4endl;
597  temp = wokvi->SampleSingleScattering(cosThetaMin, cosTetM, prob[i]);
598 
599  // direction is changed
600  temp.rotateUz(dir);
601  dir = temp;
602  //G4cout << dir << G4endl;
603 
604  // new proposed step length
605  x2 -= step;
606  x1 = -G4Log(G4UniformRand())/xtsec;
607 
608  // multiple scattering
609  } else {
610  --nMscSteps;
611  x1 -= step;
612  x2 = x0;
613 
614  // sample z in interval 0 - 1
615  G4bool isFirst = true;
616  if(prob2 > 0.0 && G4UniformRand() < prob2) { isFirst = false; }
617  do {
618  //z = -z0*G4Log(1.0 - (1.0 - zzz)*G4UniformRand());
619  if(isFirst) { z = -G4Log(G4UniformRand()); }
620  else { z = G4RandGamma::shoot(2.0,2.0); }
621  z *= z0;
622  } while(z > 1.0);
623 
624  cost = 1.0 - 2.0*z/*factCM*/;
625  if(cost > 1.0) { cost = 1.0; }
626  else if(cost < -1.0) { cost =-1.0; }
627  sint = sqrt((1.0 - cost)*(1.0 + cost));
628  phi = twopi*G4UniformRand();
629  G4double vx1 = sint*cos(phi);
630  G4double vy1 = sint*sin(phi);
631 
632  // lateral displacement
633  if (latDisplasment) {
634  G4double rms = invsqrt12*sqrt(2*z0);
635  G4double r = x0*mscfac;
636  G4double dx = r*(0.5*vx1 + rms*G4RandGauss::shoot(0.0,1.0));
637  G4double dy = r*(0.5*vy1 + rms*G4RandGauss::shoot(0.0,1.0));
638  G4double dz;
639  G4double d = r*r - dx*dx - dy*dy;
640 
641  // change position
642  if(d >= 0.0) {
643  dz = sqrt(d) - r;
644  temp.set(dx,dy,dz);
645  temp.rotateUz(dir);
646  fDisplacement += temp;
647  }
648  }
649  // change direction
650  temp.set(vx1,vy1,cost);
651  temp.rotateUz(dir);
652  dir = temp;
653  }
654  } while (0 < nMscSteps);
655 
656  dir.rotateUz(oldDirection);
657 
658  //G4cout<<"G4WentzelVIModel sampling is done 1-cost= "<< 1.-dir.z()<<G4endl;
659  // end of sampling -------------------------------
660 
662 
663  // lateral displacement
664  fDisplacement.rotateUz(oldDirection);
665 
666  /*
667  G4cout << " r(mm)= " << fDisplacement.mag()
668  << " safety= " << safety
669  << " trueStep(mm)= " << tPathLength
670  << " geomStep(mm)= " << zPathLength
671  << " x= " << fDisplacement.x()
672  << " y= " << fDisplacement.y()
673  << " z= " << fDisplacement.z()
674  << G4endl;
675  */
676 
677  //G4cout<< "G4WentzelVIModel::SampleScattering end NewDir= " << dir<< G4endl;
678  return fDisplacement;
679 }
680 
681 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
682 
684 {
685  // prepare recomputation of x-sections
686  const G4ElementVector* theElementVector = currentMaterial->GetElementVector();
687  const G4double* theAtomNumDensityVector =
690  if(nelm > nelments) {
691  nelments = nelm;
692  xsecn.resize(nelm);
693  prob.resize(nelm);
694  }
695 
696  // check consistency
697  xtsec = 0.0;
698  if(cosTetMaxNuc >= cosTheta) { return 0.0; }
699 
700  G4double cut = (*currentCuts)[currentMaterialIndex];
701  if(fixedCut > 0.0) { cut = fixedCut; }
702 
703  // loop over elements
704  G4double xs = 0.0;
705  for (G4int i=0; i<nelm; ++i) {
706  G4double costm =
707  wokvi->SetupTarget(G4lrint((*theElementVector)[i]->GetZ()), cut);
708  G4double density = theAtomNumDensityVector[i];
709 
710  G4double esec = 0.0;
711  if(costm < cosTheta) {
712 
713  // recompute the transport x-section
714  if(1.0 > cosTheta) {
715  xs += density*wokvi->ComputeTransportCrossSectionPerAtom(cosTheta);
716  }
717  // recompute the total x-section
718  G4double nucsec = wokvi->ComputeNuclearCrossSection(cosTheta, costm);
719  esec = wokvi->ComputeElectronCrossSection(cosTheta, costm);
720  nucsec += esec;
721  if(nucsec > 0.0) { esec /= nucsec; }
722  xtsec += nucsec*density;
723  }
724  xsecn[i] = xtsec;
725  prob[i] = esec;
726  //G4cout << i << " xs= " << xs << " xtsec= " << xtsec
727  // << " 1-cosTheta= " << 1-cosTheta
728  // << " 1-cosTetMaxNuc2= " <<1-cosTetMaxNuc2<< G4endl;
729  }
730 
731  //G4cout << "ComputeXS result: xsec(1/mm)= " << xs
732  // << " txsec(1/mm)= " << xtsec <<G4endl;
733  return xs;
734 }
735 
736 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
737 
739  G4double kinEnergy)
740 {
741  G4double xs = 0.0;
742 
743  SetupParticle(p);
745 
746  if(cosTetMaxNuc >= 1.0) { return xs; }
747 
748  const G4ElementVector* theElementVector = currentMaterial->GetElementVector();
749  const G4double* theAtomNumDensityVector =
752 
753  G4double cut = (*currentCuts)[currentMaterialIndex];
754  if(fixedCut > 0.0) { cut = fixedCut; }
755 
756  // loop over elements
757  for (G4int i=0; i<nelm; ++i) {
758  G4double costm =
759  wokvi->SetupTarget(G4lrint((*theElementVector)[i]->GetZ()), cut);
760  xs += theAtomNumDensityVector[i]
762  }
763  return xs;
764 }
765 
766 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
767 
769 {
770  if(val > 0.05) {
771  ssFactor = val;
772  invssFactor = 1.0/(val - 0.05);
773  }
774 }
775 
776 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
virtual ~G4WentzelVIModel()
void SetSingleScatteringFactor(G4double)
G4double facgeom
Definition: G4VMscModel.hh:178
ThreeVector shoot(const G4int Ap, const G4int Af)
G4WentzelOKandVIxSection * wokvi
std::vector< G4double > prob
static void SetPhysicsVector(G4PhysicsTable *physTable, size_t idx, G4PhysicsVector *vec)
std::vector< G4Element * > G4ElementVector
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
const G4DynamicParticle * GetDynamicParticle() const
G4double facrange
Definition: G4VMscModel.hh:177
G4double z
Definition: TRTMaterials.hh:39
G4double GetProductionCut(G4int index) const
G4StepStatus GetStepStatus() const
virtual G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety)
G4double ComputeSecondTransportMoment(G4double CosThetaMax)
const G4double pi
void DefineMaterial(const G4MaterialCutsCouple *)
G4double ConvertTrueToGeom(G4double &tLength, G4double &gLength)
Definition: G4VMscModel.hh:247
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4bool latDisplasment
Definition: G4VMscModel.hh:190
const G4MaterialCutsCouple * currentCouple
G4bool IsMaster() const
Definition: G4VEmModel.hh:699
G4PhysicsTable * GetCrossSectionTable()
Definition: G4VEmModel.hh:807
G4ParticleChangeForMSC * fParticleChange
G4ParticleDefinition * GetDefinition() const
size_t GetVectorLength() const
const G4Material * currentMaterial
G4WentzelVIModel(G4bool comb=true, const G4String &nam="WentzelVIUni")
const G4Step * GetStep() const
G4double ComputeSafety(const G4ThreeVector &position, G4double limit=DBL_MAX)
Definition: G4VMscModel.hh:239
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
int G4int
Definition: G4Types.hh:78
void FillSecondDerivatives()
G4double SetupTarget(G4int Z, G4double cut=DBL_MAX)
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &currentMinimalStep)
static G4PhysicsTable * PreparePhysicsTable(G4PhysicsTable *physTable)
G4StepStatus
Definition: G4StepStatus.hh:49
G4StepPoint * GetPreStepPoint() const
G4double GetEnergy(const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.hh:309
G4double density
Definition: TRTMaterials.hh:39
G4double ComputeTransportCrossSectionPerAtom(G4double CosThetaMax)
std::vector< G4double > xsecn
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:433
#define G4UniformRand()
Definition: Randomize.hh:95
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight=0., G4double cut=DBL_MAX, G4double emax=DBL_MAX)
const G4ThreeVector & GetPosition() const
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.hh:289
G4ParticleChangeForMSC * GetParticleChangeForMSC(const G4ParticleDefinition *p=0)
Definition: G4VMscModel.cc:89
bool G4bool
Definition: G4Types.hh:79
G4double SecondMoment(const G4ParticleDefinition *, const G4MaterialCutsCouple *, G4double kineticEnergy)
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
void PutValue(size_t index, G4double theValue)
G4double ComputeTransportXSectionPerVolume(G4double cosTheta)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
const G4int n
G4double Energy(size_t index) const
G4double ComputeGeomLimit(const G4Track &, G4double &presafety, G4double limit)
Definition: G4VMscModel.hh:257
G4ThreeVector fDisplacement
Definition: G4VMscModel.hh:186
G4double GetTransportMeanFreePath(const G4ParticleDefinition *part, G4double kinEnergy)
Definition: G4VMscModel.hh:346
G4double GetRadlen() const
Definition: G4Material.hh:218
static const G4double e1
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
static G4ProductionCutsTable * GetProductionCutsTable()
static const double eV
Definition: G4SIunits.hh:194
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
const G4int minNCollisions
int G4lrint(double ad)
Definition: templates.hh:163
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double facsafety
Definition: G4VMscModel.hh:179
void SetupParticle(const G4ParticleDefinition *)
const G4DataVector * currentCuts
void StartTracking(G4Track *)
G4PhysicsTable * GetSecondMomentTable()
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:643
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
const G4ParticleDefinition * particle
G4double GetSafety() const
G4double ComputeSecondMoment(const G4ParticleDefinition *, G4double kineticEnergy)
G4MscStepLimitType steppingAlgorithm
Definition: G4VMscModel.hh:187
G4PhysicsTable * fSecondMoments
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
G4bool GetFlag(size_t i) const
double G4double
Definition: G4Types.hh:76
virtual G4double ComputeTrueStepLength(G4double geomStepLength)
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
G4ProductionCuts * GetProductionCuts() const
#define DBL_MAX
Definition: templates.hh:83
static const double mm
Definition: G4SIunits.hh:102
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)
virtual G4double ComputeGeomPathLength(G4double truePathLength)
G4ThreeVector & SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio=0.0)