Geant4  10.01.p03
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 90583 2015-06-04 11:16:41Z 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 "G4EmParameters.hh"
67 #include "G4Log.hh"
68 #include "G4Exp.hh"
69 
70 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71 
72 using namespace std;
73 
74 const G4int minNCollisions = 10;
75 
77  G4VMscModel(nam),
78  ssFactor(1.05),
79  invssFactor(1.0),
80  currentCouple(0),
81  inside(false),
82  singleScatteringMode(false),
83  cosThetaMin(1.0),
84  cosThetaMax(-1.0),
85  fSecondMoments(0),
86  idx2(0),
87  numlimit(0.1),
88  isCombined(combined),
89  useSecondMoment(false)
90 {
92  invsqrt12 = 1./sqrt(12.);
93  tlimitminfix = 1.e-6*mm;
94  lowEnergyLimit = 1.0*eV;
95  particle = 0;
96  nelments = 5;
97  xsecn.resize(nelments);
98  prob.resize(nelments);
99  wokvi = new G4WentzelOKandVIxSection(combined);
100  fixedCut = -1.0;
101 
103  = currentRange = xtsec = cosTetMaxNuc = 0.0;
105 
106  fParticleChange = 0;
107  currentCuts = 0;
108  currentMaterial = 0;
109 }
110 
111 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
112 
114 {
115  delete wokvi;
116  if(fSecondMoments && IsMaster()) {
117  delete fSecondMoments;
118  fSecondMoments = 0;
119  }
120 }
121 
122 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
123 
125  const G4DataVector& cuts)
126 {
127  // reset parameters
128  SetupParticle(p);
129  currentRange = 0.0;
130 
131  if(isCombined) {
132  G4double tet = PolarAngleLimit();
133  if(tet <= 0.0) { cosThetaMax = 1.0; }
134  else if(tet < CLHEP::pi) { cosThetaMax = cos(tet); }
135  }
136 
137  //G4cout << "G4WentzelVIModel::Initialise " << p->GetParticleName() << G4endl;
139  /*
140  G4cout << "G4WentzelVIModel: " << particle->GetParticleName()
141  << " 1-cos(ThetaLimit)= " << 1 - cosThetaMax
142  << " SingScatFactor= " << ssFactor
143  << G4endl;
144  */
145  currentCuts = &cuts;
146 
147  // set values of some data members
149 
150  // build second moment table only if transport table is build
152  if(useSecondMoment && IsMaster() && table) {
153 
154  //G4cout << "### G4WentzelVIModel::Initialise: build 2nd moment table "
155  // << table << G4endl;
156  fSecondMoments =
158  // Access to materials
159  const G4ProductionCutsTable* theCoupleTable =
161  size_t numOfCouples = theCoupleTable->GetTableSize();
162 
163  G4bool splineFlag = true;
164  G4PhysicsVector* aVector = 0;
165  G4PhysicsVector* bVector = 0;
168  if(emin < emax) {
170  *G4lrint(std::log10(emax/emin));
171  if(n < 3) { n = 3; }
172 
173  for(size_t i=0; i<numOfCouples; ++i) {
174 
175  //G4cout<< "i= " << i << " Flag= " << fSecondMoments->GetFlag(i)
176  // << G4endl;
177  if(fSecondMoments->GetFlag(i)) {
178  DefineMaterial(theCoupleTable->GetMaterialCutsCouple(i));
179 
180  delete (*fSecondMoments)[i];
181  if(!aVector) {
182  aVector = new G4PhysicsLogVector(emin, emax, n);
183  bVector = aVector;
184  } else {
185  bVector = new G4PhysicsVector(*aVector);
186  }
187  for(size_t j=0; j<n; ++j) {
188  G4double e = bVector->Energy(j);
189  bVector->PutValue(j, ComputeSecondMoment(p, e)*e*e);
190  }
191  if(splineFlag) { bVector->FillSecondDerivatives(); }
192  (*fSecondMoments)[i] = bVector;
193  }
194  }
195  }
196  //G4cout << *fSecondMoments << G4endl;
197  }
198 }
199 
200 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
201 
203  G4VEmModel* masterModel)
204 {
205  fSecondMoments = static_cast<G4WentzelVIModel*>(masterModel)
207 }
208 
209 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
210 
212  const G4ParticleDefinition* p,
213  G4double kinEnergy,
214  G4double Z, G4double,
215  G4double cutEnergy, G4double)
216 {
217  G4double cross = 0.0;
218  if(p != particle) { SetupParticle(p); }
219  if(kinEnergy < lowEnergyLimit) { return cross; }
220  if(!CurrentCouple()) {
221  G4Exception("G4WentzelVIModel::ComputeCrossSectionPerAtom", "em0011",
222  FatalException, " G4MaterialCutsCouple is not defined");
223  return 0.0;
224  }
227  if(cosTetMaxNuc < 1.0) {
228  G4double cut = cutEnergy;
229  if(fixedCut > 0.0) { cut = fixedCut; }
230  G4double cost = wokvi->SetupTarget(G4lrint(Z), cut);
232  /*
233  if(p->GetParticleName() == "e-")
234  G4cout << "G4WentzelVIModel::CS: Z= " << G4int(Z) << " e(MeV)= "<<kinEnergy
235  << " 1-cosN= " << 1 - cosTetMaxNuc << " cross(bn)= " << cross/barn
236  << " " << particle->GetParticleName() << G4endl;
237  */
238  }
239  return cross;
240 }
241 
242 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
243 
245 {
248 }
249 
250 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
251 
253  const G4Track& track,
254  G4double& currentMinimalStep)
255 {
256  G4double tlimit = currentMinimalStep;
257  const G4DynamicParticle* dp = track.GetDynamicParticle();
258  G4StepPoint* sp = track.GetStep()->GetPreStepPoint();
259  G4StepStatus stepStatus = sp->GetStepStatus();
260  singleScatteringMode = false;
261 
262  //G4cout << "G4WentzelVIModel::ComputeTruePathLengthLimit stepStatus= "
263  // << stepStatus << " " << track.GetDefinition()->GetParticleName()
264  // << G4endl;
265 
266  // initialisation for each step, lambda may be computed from scratch
273 
274  //G4cout << "lambdaeff= " << lambdaeff << " Range= " << currentRange
275  // << " tlimit= " << tlimit << " 1-cost= " << 1 - cosTetMaxNuc << G4endl;
276 
277  // extra check for abnormal situation
278  // this check needed to run MSC with eIoni and eBrem inactivated
279  if(tlimit > currentRange) { tlimit = currentRange; }
280 
281  // stop here if small range particle
282  if(tlimit < tlimitminfix) {
283  return ConvertTrueToGeom(tlimit, currentMinimalStep);
284  }
285 
286  // pre step
287  G4double presafety = sp->GetSafety();
288  // far from geometry boundary
289  if(currentRange < presafety) {
290  return ConvertTrueToGeom(tlimit, currentMinimalStep);
291  }
292 
293  // compute presafety again if presafety <= 0 and no boundary
294  // i.e. when it is needed for optimization purposes
295  if(stepStatus != fGeomBoundary && presafety < tlimitminfix) {
296  presafety = ComputeSafety(sp->GetPosition(), tlimit);
297  if(currentRange < presafety) {
298  return ConvertTrueToGeom(tlimit, currentMinimalStep);
299  }
300  }
301  /*
302  G4cout << "e(MeV)= " << preKinEnergy/MeV
303  << " " << particle->GetParticleName()
304  << " CurLimit(mm)= " << tlimit/mm <<" safety(mm)= " << presafety/mm
305  << " R(mm)= " <<currentRange/mm
306  << " L0(mm^-1)= " << lambdaeff*mm
307  << G4endl;
308  */
309  // natural limit for high energy
311  (1.0 - cosTetMaxNuc)*lambdaeff*invssFactor);
312 
313  // low-energy e-
314  if(cosThetaMax > cosTetMaxNuc) {
315  rlimit = std::min(rlimit, facsafety*presafety);
316  }
317 
318  // cut correction
320  //G4cout << "rcut= " << rcut << " rlimit= " << rlimit << " presafety= "
321  // << presafety << " 1-cosThetaMax= " <<1-cosThetaMax
322  //<< " 1-cosTetMaxNuc= " << 1-cosTetMaxNuc << G4endl;
323  if(rcut > rlimit) { rlimit = std::min(rlimit, rcut*sqrt(rlimit/rcut)); }
324 
325  tlimit = std::min(tlimit, rlimit);
326  tlimit = std::max(tlimit, tlimitminfix);
327 
328  // step limit in infinite media
329  tlimit = std::min(tlimit, 50*currentMaterial->GetRadlen()/facgeom);
330 
331  //compute geomlimit and force few steps within a volume
333  && stepStatus == fGeomBoundary) {
334 
335  G4double geomlimit = ComputeGeomLimit(track, presafety, currentRange);
336  tlimit = std::min(tlimit, geomlimit/facgeom);
337  }
338  /*
339  G4cout << particle->GetParticleName() << " e= " << preKinEnergy
340  << " L0= " << lambdaeff << " R= " << currentRange
341  << " tlimit= " << tlimit
342  << " currentMinimalStep= " << currentMinimalStep << G4endl;
343  */
344  return ConvertTrueToGeom(tlimit, currentMinimalStep);
345 }
346 
347 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
348 
350 {
351  zPathLength = tPathLength = truelength;
352 
353  // small step use only single scattering
354  cosThetaMin = 1.0;
356  //G4cout << "xtsec= " << xtsec << " Nav= "
357  // << zPathLength*xtsec << G4endl;
358  if(0.0 >= lambdaeff || G4int(zPathLength*xtsec) < minNCollisions) {
359  singleScatteringMode = true;
360  lambdaeff = DBL_MAX;
361 
362  } else {
363  //G4cout << "ComputeGeomPathLength: tLength= " << tPathLength
364  // << " Leff= " << lambdaeff << G4endl;
365  // small step
368  zPathLength *= (1.0 - 0.5*tau + tau*tau/6.0);
369 
370  // medium step
371  } else {
372  G4double e1 = 0.0;
373  if(currentRange > tPathLength) {
375  }
376  effKinEnergy = 0.5*(e1 + preKinEnergy);
379  //G4cout << " tLength= "<< tPathLength<< " Leff= " << lambdaeff << G4endl;
381  if(tPathLength*numlimit < lambdaeff) {
382  zPathLength *= (1.0 - G4Exp(-tPathLength/lambdaeff));
383  }
384  }
385  }
386  //G4cout << "Comp.geom: zLength= "<<zPathLength<<" tLength= "
387  // << tPathLength<< " Leff= " << lambdaeff << G4endl;
388  return zPathLength;
389 }
390 
391 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
392 
394 {
395  // initialisation of single scattering x-section
396  /*
397  G4cout << "ComputeTrueStepLength: Step= " << geomStepLength
398  << " geomL= " << zPathLength
399  << " Lambda= " << lambdaeff
400  << " 1-cosThetaMaxNuc= " << 1 - cosTetMaxNuc << G4endl;
401  */
403  zPathLength = tPathLength = geomStepLength;
404 
405  } else {
406 
407  // step defined by transportation
408  // change both geom and true step lengths
409  if(geomStepLength < zPathLength) {
410 
411  // single scattering
412  if(G4int(geomStepLength*xtsec) < minNCollisions) {
413  zPathLength = tPathLength = geomStepLength;
414  lambdaeff = DBL_MAX;
415  singleScatteringMode = true;
416 
417  // multiple scattering
418  } else {
419  // small step
420  if(geomStepLength < numlimit*lambdaeff) {
421  G4double tau = geomStepLength/lambdaeff;
422  tPathLength = geomStepLength*(1.0 + 0.5*tau + tau*tau/3.0);
423 
424  // energy correction for a big step
425  } else {
426  tPathLength *= geomStepLength/zPathLength;
427  G4double e1 = 0.0;
428  if(currentRange > tPathLength) {
430  }
431  effKinEnergy = 0.5*(e1 + preKinEnergy);
434  G4double tau = geomStepLength/lambdaeff;
435 
436  if(tau < 0.999999) { tPathLength = -lambdaeff*G4Log(1.0 - tau); }
437  else { tPathLength = currentRange; }
438  }
439  zPathLength = geomStepLength;
440  }
441  }
442  }
443  // check of step length
444  // define threshold angle between single and multiple scattering
445  if(!singleScatteringMode) {
447  xtsec = 0.0;
448 
449  // recompute transport cross section - do not change energy
450  // anymore - cannot be applied for big steps
451  if(cosThetaMin > cosTetMaxNuc) {
452  // new computation
454  //G4cout << "%%%% cross= " << cross << " xtsec= " << xtsec
455  // << " 1-cosTMin= " << 1.0 - cosThetaMin << G4endl;
456  if(cross <= 0.0) {
457  singleScatteringMode = true;
459  lambdaeff = DBL_MAX;
460  cosThetaMin = 1.0;
461  } else if(xtsec > 0.0) {
462 
463  lambdaeff = 1./cross;
464  G4double tau = zPathLength*cross;
465  if(tau < numlimit) {
466  tPathLength = zPathLength*(1.0 + 0.5*tau + tau*tau/3.0);
467  } else if(tau < 0.999999) {
468  tPathLength = -lambdaeff*G4Log(1.0 - tau);
469  } else {
471  }
472  }
473  }
474  }
476  /*
477  G4cout <<"Comp.true: zLength= "<<zPathLength<<" tLength= "<<tPathLength
478  <<" Leff(mm)= "<<lambdaeff/mm<<" sig0(1/mm)= " << xtsec <<G4endl;
479  G4cout << particle->GetParticleName() << " 1-cosThetaMin= " << 1-cosThetaMin
480  << " 1-cosTetMaxNuc= " << 1-cosTetMaxNuc
481  << " e(MeV)= " << preKinEnergy/MeV << " "
482  << " SSmode= " << singleScatteringMode << G4endl;
483  */
484  return tPathLength;
485 }
486 
487 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
488 
491  G4double /*safety*/)
492 {
493  fDisplacement.set(0.0,0.0,0.0);
494  //G4cout << "!##! G4WentzelVIModel::SampleScattering for "
495  // << particle->GetParticleName() << G4endl;
496 
497  // ignore scattering for zero step length and energy below the limit
498  if(preKinEnergy < lowEnergyLimit || tPathLength <= 0.0)
499  { return fDisplacement; }
500 
501  G4double invlambda = 0.0;
502  if(lambdaeff < DBL_MAX) { invlambda = 0.5/lambdaeff; }
503 
504  // use average kinetic energy over the step
505  G4double cut = (*currentCuts)[currentMaterialIndex];
506  if(fixedCut > 0.0) { cut = fixedCut; }
507  /*
508  G4cout <<"SampleScat: E0(MeV)= "<< preKinEnergy/MeV
509  << " Leff= " << lambdaeff <<" sig0(1/mm)= " << xtsec
510  << " xmsc= " << tPathLength*invlambda
511  << " safety= " << safety << G4endl;
512  */
513  // step limit due msc
514  G4int nMscSteps = 1;
515  G4double x0 = tPathLength;
516  G4double z0 = x0*invlambda;
517  //G4double zzz = 0.0;
518  G4double prob2 = 0.0;
519 
520  // large scattering angle case - two step approach
521  if(!singleScatteringMode) {
522  static const G4double zzmin = 0.05;
523  if(useSecondMoment) {
524  G4double z1 = invlambda*invlambda;
526  prob2 = (z2 - z1)/(1.5*z1 - z2);
527  }
528  // if(z0 > zzmin && safety > tlimitminfix) {
529  if(z0 > zzmin) {
530  x0 *= 0.5;
531  z0 *= 0.5;
532  nMscSteps = 2;
533  }
534  //if(z0 > zzmin) { zzz = G4Exp(-1.0/z0); }
535  G4double zzz = 0.0;
536  if(z0 > zzmin) {
537  zzz = G4Exp(-1.0/z0);
538  z0 += zzz;
539  prob2 *= (1 + zzz);
540  }
541  prob2 /= (1 + prob2);
542  }
543 
544  // step limit due to single scattering
545  G4double x1 = 2*tPathLength;
546  if(0.0 < xtsec) { x1 = -G4Log(rndmEngineMod->flat())/xtsec; }
547 
548  // no scattering case
549  if(singleScatteringMode && x1 > tPathLength)
550  { return fDisplacement; }
551 
552  const G4ElementVector* theElementVector =
555 
556  // geometry
557  G4double sint, cost, phi;
558  G4ThreeVector temp(0.0,0.0,1.0);
559 
560  // current position and direction relative to the end point
561  // because of magnetic field geometry is computed relatively to the
562  // end point of the step
563  G4ThreeVector dir(0.0,0.0,1.0);
564  fDisplacement.set(0.0,0.0,-zPathLength);
565 
567 
568  // start a loop
569  G4double x2 = x0;
570  G4double step, z;
571  G4bool singleScat;
572  /*
573  G4cout << "Start of the loop x1(mm)= " << x1 << " x2(mm)= " << x2
574  << " 1-cost1= " << 1 - cosThetaMin << " SSmode= " << singleScatteringMode
575  << " xtsec= " << xtsec << " Nst= " << nMscSteps << G4endl;
576  */
577  do {
578 
579  //G4cout << "# x1(mm)= "<< x1<< " x2(mm)= "<< x2 << G4endl;
580  // single scattering case
581  if(singleScatteringMode && x1 > x2) {
582  fDisplacement += x2*mscfac*dir;
583  break;
584  }
585 
586  // what is next single of multiple?
587  if(x1 <= x2) {
588  step = x1;
589  singleScat = true;
590  } else {
591  step = x2;
592  singleScat = false;
593  }
594 
595  //G4cout << "# step(mm)= "<< step<< " singlScat= "<< singleScat << G4endl;
596 
597  // new position
598  fDisplacement += step*mscfac*dir;
599 
600  if(singleScat) {
601 
602  // select element
603  G4int i = 0;
604  if(nelm > 1) {
605  G4double qsec = rndmEngineMod->flat()*xtsec;
606  for (; i<nelm; ++i) { if(xsecn[i] >= qsec) { break; } }
607  }
608  G4double cosTetM =
609  wokvi->SetupTarget(G4lrint((*theElementVector)[i]->GetZ()), cut);
610  //G4cout << "!!! " << cosThetaMin << " " << cosTetM << " "
611  // << prob[i] << G4endl;
612  temp = wokvi->SampleSingleScattering(cosThetaMin, cosTetM, prob[i]);
613 
614  // direction is changed
615  temp.rotateUz(dir);
616  dir = temp;
617  //G4cout << dir << G4endl;
618 
619  // new proposed step length
620  x2 -= step;
621  x1 = -G4Log(rndmEngineMod->flat())/xtsec;
622 
623  // multiple scattering
624  } else {
625  --nMscSteps;
626  x1 -= step;
627  x2 = x0;
628 
629  // sample z in interval 0 - 1
630  G4bool isFirst = true;
631  if(prob2 > 0.0 && rndmEngineMod->flat() < prob2) { isFirst = false; }
632  do {
633  //z = -z0*G4Log(1.0 - (1.0 - zzz)*rndmEngineMod->flat());
634  if(isFirst) { z = -G4Log(rndmEngineMod->flat()); }
635  else { z = G4RandGamma::shoot(rndmEngineMod, 2.0, 2.0); }
636  z *= z0;
637  } while(z > 1.0);
638 
639  cost = 1.0 - 2.0*z/*factCM*/;
640  if(cost > 1.0) { cost = 1.0; }
641  else if(cost < -1.0) { cost =-1.0; }
642  sint = sqrt((1.0 - cost)*(1.0 + cost));
643  phi = twopi*rndmEngineMod->flat();
644  G4double vx1 = sint*cos(phi);
645  G4double vy1 = sint*sin(phi);
646 
647  // lateral displacement
648  if (latDisplasment) {
649  G4double rms = invsqrt12*sqrt(2*z0);
650  G4double r = x0*mscfac;
651  G4double dx = r*(0.5*vx1 + rms*G4RandGauss::shoot(rndmEngineMod, 0.0, 1.0));
652  G4double dy = r*(0.5*vy1 + rms*G4RandGauss::shoot(rndmEngineMod, 0.0, 1.0));
653  G4double dz;
654  G4double d = r*r - dx*dx - dy*dy;
655 
656  // change position
657  if(d >= 0.0) {
658  dz = sqrt(d) - r;
659  temp.set(dx,dy,dz);
660  temp.rotateUz(dir);
661  fDisplacement += temp;
662  }
663  }
664  // change direction
665  temp.set(vx1,vy1,cost);
666  temp.rotateUz(dir);
667  dir = temp;
668  }
669  } while (0 < nMscSteps);
670 
671  dir.rotateUz(oldDirection);
672 
673  //G4cout<<"G4WentzelVIModel sampling is done 1-cost= "<< 1.-dir.z()<<G4endl;
674  // end of sampling -------------------------------
675 
677 
678  // lateral displacement
679  fDisplacement.rotateUz(oldDirection);
680 
681  /*
682  G4cout << " r(mm)= " << fDisplacement.mag()
683  << " safety= " << safety
684  << " trueStep(mm)= " << tPathLength
685  << " geomStep(mm)= " << zPathLength
686  << " x= " << fDisplacement.x()
687  << " y= " << fDisplacement.y()
688  << " z= " << fDisplacement.z()
689  << G4endl;
690  */
691 
692  //G4cout<< "G4WentzelVIModel::SampleScattering end NewDir= " << dir<< G4endl;
693  return fDisplacement;
694 }
695 
696 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
697 
699 {
700  // prepare recomputation of x-sections
701  const G4ElementVector* theElementVector = currentMaterial->GetElementVector();
702  const G4double* theAtomNumDensityVector =
705  if(nelm > nelments) {
706  nelments = nelm;
707  xsecn.resize(nelm);
708  prob.resize(nelm);
709  }
710 
711  // check consistency
712  xtsec = 0.0;
713  if(cosTetMaxNuc >= cosTheta) { return 0.0; }
714 
715  G4double cut = (*currentCuts)[currentMaterialIndex];
716  if(fixedCut > 0.0) { cut = fixedCut; }
717 
718  // loop over elements
719  G4double xs = 0.0;
720  for (G4int i=0; i<nelm; ++i) {
721  G4double costm =
722  wokvi->SetupTarget(G4lrint((*theElementVector)[i]->GetZ()), cut);
723  G4double density = theAtomNumDensityVector[i];
724 
725  G4double esec = 0.0;
726  if(costm < cosTheta) {
727 
728  // recompute the transport x-section
729  if(1.0 > cosTheta) {
730  xs += density*wokvi->ComputeTransportCrossSectionPerAtom(cosTheta);
731  }
732  // recompute the total x-section
733  G4double nucsec = wokvi->ComputeNuclearCrossSection(cosTheta, costm);
734  esec = wokvi->ComputeElectronCrossSection(cosTheta, costm);
735  nucsec += esec;
736  if(nucsec > 0.0) { esec /= nucsec; }
737  xtsec += nucsec*density;
738  }
739  xsecn[i] = xtsec;
740  prob[i] = esec;
741  //G4cout << i << " xs= " << xs << " xtsec= " << xtsec
742  // << " 1-cosTheta= " << 1-cosTheta
743  // << " 1-cosTetMaxNuc2= " <<1-cosTetMaxNuc2<< G4endl;
744  }
745 
746  //G4cout << "ComputeXS result: xsec(1/mm)= " << xs
747  // << " txsec(1/mm)= " << xtsec <<G4endl;
748  return xs;
749 }
750 
751 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
752 
754  G4double kinEnergy)
755 {
756  G4double xs = 0.0;
757 
758  SetupParticle(p);
760 
761  if(cosTetMaxNuc >= 1.0) { return xs; }
762 
763  const G4ElementVector* theElementVector = currentMaterial->GetElementVector();
764  const G4double* theAtomNumDensityVector =
767 
768  G4double cut = (*currentCuts)[currentMaterialIndex];
769  if(fixedCut > 0.0) { cut = fixedCut; }
770 
771  // loop over elements
772  for (G4int i=0; i<nelm; ++i) {
773  G4double costm =
774  wokvi->SetupTarget(G4lrint((*theElementVector)[i]->GetZ()), cut);
775  xs += theAtomNumDensityVector[i]
777  }
778  return xs;
779 }
780 
781 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
782 
784 {
785  if(val > 0.05) {
786  ssFactor = val;
787  invssFactor = 1.0/(val - 0.05);
788  }
789 }
790 
791 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4double LowEnergyActivationLimit() const
Definition: G4VEmModel.hh:637
G4double HighEnergyActivationLimit() const
Definition: G4VEmModel.hh:630
virtual ~G4WentzelVIModel()
G4int NumberOfBinsPerDecade() const
void SetSingleScatteringFactor(G4double)
G4double facgeom
Definition: G4VMscModel.hh:178
ThreeVector shoot(const G4int Ap, const G4int Af)
G4WentzelOKandVIxSection * wokvi
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:623
std::vector< G4double > prob
virtual void StartTracking(G4Track *)
Definition: G4VEmModel.cc:284
std::vector< G4Element * > G4ElementVector
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:616
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:700
G4PhysicsTable * GetCrossSectionTable()
Definition: G4VEmModel.hh:808
G4ParticleChangeForMSC * fParticleChange
G4ParticleDefinition * GetDefinition() 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:190
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:206
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:434
CLHEP::HepRandomEngine * rndmEngineMod
Definition: G4VEmModel.hh:400
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:90
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:220
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)
static G4EmParameters * Instance()
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:644
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:186
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)