Geant4  10.02.p01
G4WentzelVIRelModel.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: G4WentzelVIRelModel.cc 91726 2015-08-03 15:41:36Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4WentzelVIRelModel
34 //
35 // Author: V.Ivanchenko
36 //
37 // Creation date: 08.06.2012 from G4WentzelVIRelModel
38 //
39 // Modifications:
40 //
41 // Class Description:
42 //
43 // Implementation of the model of multiple scattering based on
44 // G.Wentzel, Z. Phys. 40 (1927) 590.
45 // H.W.Lewis, Phys Rev 78 (1950) 526.
46 // J.M. Fernandez-Varea et al., NIM B73 (1993) 447.
47 // L.Urban, CERN-OPEN-2006-077.
48 
49 // -------------------------------------------------------------------
50 //
51 
52 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
54 
55 #include "G4WentzelVIRelModel.hh"
56 #include "G4PhysicalConstants.hh"
57 #include "G4SystemOfUnits.hh"
58 #include "Randomize.hh"
60 #include "G4PhysicsTableHelper.hh"
61 #include "G4ElementVector.hh"
62 #include "G4ProductionCutsTable.hh"
63 #include "G4LossTableManager.hh"
64 #include "G4Pow.hh"
65 #include "G4NistManager.hh"
66 #include "G4Log.hh"
67 #include "G4Exp.hh"
68 
69 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
70 
71 using namespace std;
72 
74  G4VMscModel("WentzelVIRel"),
75  numlimit(0.1),
76  currentCouple(0),
77  cosThetaMin(1.0),
78  isCombined(combined),
79  inside(false),
80  singleScatteringMode(false)
81 {
82  invsqrt12 = 1./sqrt(12.);
83  tlimitminfix = 1.e-6*mm;
84  lowEnergyLimit = 1.0*eV;
85  particle = 0;
86  nelments = 5;
87  xsecn.resize(nelments);
88  prob.resize(nelments);
92  wokvi = new G4WentzelVIRelXSection(combined);
93 
95  xtsec = 0;
97  cosThetaMax = cosTetMaxNuc = 1.0;
98 
99  fParticleChange = 0;
100  currentCuts = 0;
101  currentMaterial = 0;
102 }
103 
104 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
105 
107 {
108  delete wokvi;
109 }
110 
111 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
112 
114  const G4DataVector& cuts)
115 {
116  // reset parameters
117  SetupParticle(p);
118  currentRange = 0.0;
119 
120  if(isCombined) {
121  G4double tet = PolarAngleLimit();
122  if(tet >= pi) { cosThetaMax = -1.0; }
123  else if(tet > 0.0) { cosThetaMax = cos(tet); }
124  }
125 
127  /*
128  G4cout << "G4WentzelVIRelModel: " << particle->GetParticleName()
129  << " 1-cos(ThetaLimit)= " << 1 - cosThetaMax
130  << G4endl;
131  */
132  currentCuts = &cuts;
133 
134  // set values of some data members
136 }
137 
138 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
139 
141  const G4ParticleDefinition* p,
142  G4double kinEnergy,
143  G4double Z, G4double,
144  G4double cutEnergy, G4double)
145 {
146  G4double cross = 0.0;
147  if(p != particle) { SetupParticle(p); }
148  if(kinEnergy < lowEnergyLimit) { return cross; }
149  if(!CurrentCouple()) {
150  G4Exception("G4WentzelVIRelModel::ComputeCrossSectionPerAtom", "em0011",
151  FatalException, " G4MaterialCutsCouple is not defined");
152  return 0.0;
153  }
155  G4int iz = G4int(Z);
156  G4double tmass = proton_mass_c2;
157  if(1 < iz) {
158  tmass = fNistManager->GetAtomicMassAmu(iz)*amu_c2;
159  }
160  cosTetMaxNuc = wokvi->SetupKinematic(kinEnergy, currentMaterial, cutEnergy, tmass);
161  if(cosTetMaxNuc < 1.0) {
162  G4double cost = wokvi->SetupTarget(G4lrint(Z), cutEnergy);
164  /*
165  if(p->GetParticleName() == "e-")
166  G4cout << "G4WentzelVIRelModel::CS: Z= " << G4int(Z)
167  << " e(MeV)= " << kinEnergy
168  << " 1-cosN= " << 1 - cosTetMaxNuc << " cross(bn)= " << cross/barn
169  << " " << particle->GetParticleName() << G4endl;
170  */
171  }
172  return cross;
173 }
174 
175 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
176 
178 {
180  inside = false;
182 }
183 
184 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
185 
187  const G4Track& track,
188  G4double& currentMinimalStep)
189 {
190  G4double tlimit = currentMinimalStep;
191  const G4DynamicParticle* dp = track.GetDynamicParticle();
192  G4StepPoint* sp = track.GetStep()->GetPreStepPoint();
193  G4StepStatus stepStatus = sp->GetStepStatus();
194  singleScatteringMode = false;
195  //G4cout << "G4WentzelVIRelModel::ComputeTruePathLengthLimit stepStatus= "
196  // << stepStatus << G4endl;
197 
198 
199  // initialisation for each step, lambda may be computed from scratch
204 
207  rcut, proton_mass_c2);
208 
209  // extra check for abnormal situation
210  // this check needed to run MSC with eIoni and eBrem inactivated
211  if(tlimit > currentRange) { tlimit = currentRange; }
212 
213  // stop here if small range particle
214  if(inside || tlimit < tlimitminfix) {
215  return ConvertTrueToGeom(tlimit, currentMinimalStep);
216  }
217 
218  // pre step
219  G4double presafety = sp->GetSafety();
220  // far from geometry boundary
221  if(currentRange < presafety) {
222  inside = true;
223  return ConvertTrueToGeom(tlimit, currentMinimalStep);
224  }
225 
226  // compute presafety again if presafety <= 0 and no boundary
227  // i.e. when it is needed for optimization purposes
228  if(stepStatus != fGeomBoundary && presafety < tlimitminfix) {
229  presafety = ComputeSafety(sp->GetPosition(), tlimit);
230  if(currentRange < presafety) {
231  inside = true;
232  return ConvertTrueToGeom(tlimit, currentMinimalStep);
233  }
234  }
235  /*
236  G4cout << "e(MeV)= " << preKinEnergy/MeV
237  << " " << particle->GetParticleName()
238  << " CurLimit(mm)= " << tlimit/mm <<" safety(mm)= " << presafety/mm
239  << " R(mm)= " <<currentRange/mm
240  << " L0(mm^-1)= " << lambdaeff*mm
241  <<G4endl;
242  */
243 
244  // natural limit for high energy
246  0.7*(1.0 - cosTetMaxNuc)*lambdaeff);
247 
248  // low-energy e-
249  if(cosThetaMax > cosTetMaxNuc) {
250  rlimit = std::min(rlimit, facsafety*presafety);
251  }
252 
253  // cut correction
254  //G4cout << "rcut= " << rcut << " rlimit= " << rlimit
255  // << " presafety= " << presafety
256  // << " 1-cosThetaMax= " <<1-cosThetaMax
257  // << " 1-cosTetMaxNuc= " << 1-cosTetMaxNuc
258  // << G4endl;
259  if(rcut > rlimit) { rlimit = std::min(rlimit, rcut*sqrt(rlimit/rcut)); }
260 
261  if(rlimit < tlimit) { tlimit = rlimit; }
262 
263  tlimit = std::max(tlimit, tlimitminfix);
264 
265  // step limit in infinite media
266  tlimit = std::min(tlimit, 50*currentMaterial->GetRadlen()/facgeom);
267 
268  //compute geomlimit and force few steps within a volume
270  {
271  G4double geomlimit = ComputeGeomLimit(track, presafety, currentRange);
272  tlimit = std::min(tlimit, geomlimit/facgeom);
273  }
274 
275  /*
276  G4cout << particle->GetParticleName() << " e= " << preKinEnergy
277  << " L0= " << lambdaeff << " R= " << currentRange
278  << "tlimit= " << tlimit
279  << " currentMinimalStep= " << currentMinimalStep << G4endl;
280  */
281  return ConvertTrueToGeom(tlimit, currentMinimalStep);
282 }
283 
284 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
285 
287 {
288  tPathLength = truelength;
290 
291  if(lambdaeff > 0.0 && lambdaeff != DBL_MAX) {
293  //G4cout << "ComputeGeomPathLength: tLength= " << tPathLength
294  // << " Leff= " << lambdaeff << " tau= " << tau << G4endl;
295  // small step
296  if(tau < numlimit) {
297  zPathLength *= (1.0 - 0.5*tau + tau*tau/6.0);
298 
299  // medium step
300  } else {
301  G4double e1 = 0.0;
302  if(currentRange > tPathLength) {
304  }
305  e1 = 0.5*(e1 + preKinEnergy);
306  cosTetMaxNuc = wokvi->SetupKinematic(e1, currentMaterial, 0.0, proton_mass_c2);
309  }
310  } else { lambdaeff = DBL_MAX; }
311  //G4cout<<"Comp.geom: zLength= "<<zPathLength
312  // <<" tLength= "<<tPathLength<<G4endl;
313  return zPathLength;
314 }
315 
316 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
317 
319 {
320  // initialisation of single scattering x-section
321  xtsec = 0.0;
323 
324  //G4cout << "Step= " << geomStepLength << " Lambda= " << lambdaeff
325  // << " 1-cosThetaMaxNuc= " << 1 - cosTetMaxNuc << G4endl;
326  // pathalogical case
327  if(lambdaeff == DBL_MAX) {
328  singleScatteringMode = true;
329  zPathLength = geomStepLength;
330  tPathLength = geomStepLength;
331  cosThetaMin = 1.0;
332 
333  // normal case
334  } else {
335 
336  // small step use only single scattering
337  static const G4double singleScatLimit = 1.0e-7;
338  if(geomStepLength < lambdaeff*singleScatLimit*(1.0 - cosTetMaxNuc)) {
339  singleScatteringMode = true;
340  cosThetaMin = 1.0;
341  zPathLength = geomStepLength;
342  tPathLength = geomStepLength;
343 
344  // step defined by transportation
345  } else if(geomStepLength != zPathLength) {
346 
347  // step defined by transportation
348  zPathLength = geomStepLength;
349  G4double tau = geomStepLength/lambdaeff;
350  tPathLength = zPathLength*(1.0 + 0.5*tau + tau*tau/3.0);
351 
352  // energy correction for a big step
353  if(tau > numlimit) {
354  G4double e1 = 0.0;
355  if(currentRange > tPathLength) {
357  }
358  e1 = 0.5*(e1 + preKinEnergy);
359  cosTetMaxNuc = wokvi->SetupKinematic(e1, currentMaterial, 0.0, proton_mass_c2);
361  tau = zPathLength/lambdaeff;
362 
363  if(tau < 0.999999) { tPathLength = -lambdaeff*G4Log(1.0 - tau); }
364  else { tPathLength = currentRange; }
365  }
366  }
367  }
368 
369  // check of step length
370  // define threshold angle between single and multiple scattering
372 
373  // recompute transport cross section - do not change energy
374  // anymore - cannot be applied for big steps
375  if(cosThetaMin > cosTetMaxNuc) {
376 
377  // new computation
379  //G4cout << "%%%% cross= " << cross << " xtsec= " << xtsec << G4endl;
380  if(cross <= 0.0) {
381  singleScatteringMode = true;
383  lambdaeff = DBL_MAX;
384  } else if(xtsec > 0.0) {
385 
386  lambdaeff = 1./cross;
387  G4double tau = zPathLength*cross;
388  if(tau < numlimit) {
389  tPathLength = zPathLength*(1.0 + 0.5*tau + tau*tau/3.0);
390  }
391  else if(tau < 0.999999) { tPathLength = -lambdaeff*G4Log(1.0 - tau); }
392  else { tPathLength = currentRange; }
393 
395  }
396  }
397 
398  /*
399  G4cout <<"Comp.true: zLength= "<<zPathLength<<" tLength= "<<tPathLength
400  <<" Leff(mm)= "<<lambdaeff/mm<<" sig0(1/mm)= " << xtsec <<G4endl;
401  G4cout << particle->GetParticleName() << " 1-cosThetaMin= " << 1-cosThetaMin
402  << " 1-cosTetMaxNuc= " << 1-cosTetMaxNuc
403  << " e(MeV)= " << preKinEnergy/MeV << " " << singleScatteringMode
404  << G4endl;
405  */
406  return tPathLength;
407 }
408 
409 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
410 
413  G4double safety)
414 {
415  fDisplacement.set(0.0,0.0,0.0);
416  //G4cout << "!##! G4WentzelVIRelModel::SampleScattering for "
417  // << particle->GetParticleName() << G4endl;
418 
419  // ignore scattering for zero step length and energy below the limit
420  if(preKinEnergy < lowEnergyLimit || tPathLength <= 0.0)
421  { return fDisplacement; }
422 
423  G4double invlambda = 0.0;
424  if(lambdaeff < DBL_MAX) { invlambda = 0.5/lambdaeff; }
425 
426  // use average kinetic energy over the step
427  G4double cut = (*currentCuts)[currentMaterialIndex];
428  /*
429  G4cout <<"SampleScat: E0(MeV)= "<< preKinEnergy/MeV
430  << " Leff= " << lambdaeff <<" sig0(1/mm)= " << xtsec
431  << " xmsc= " << tPathLength*invlambda
432  << " safety= " << safety << G4endl;
433  */
434 
435  // step limit due msc
436  G4double x0 = tPathLength;
437  // const G4double thinlimit = 0.5;
438  static const G4double thinlimit = 0.1;
439  G4int nMscSteps = 1;
440  // large scattering angle case - two step approach
441  if(tPathLength*invlambda > thinlimit && safety > tlimitminfix) {
442  x0 *= 0.5;
443  nMscSteps = 2;
444  }
445 
446  // step limit due to single scattering
447  G4double x1 = 2*tPathLength;
448  if(0.0 < xtsec) { x1 = -G4Log(G4UniformRand())/xtsec; }
449 
450  const G4ElementVector* theElementVector =
453 
454  // geometry
455  G4double sint, cost, phi;
456  G4ThreeVector temp(0.0,0.0,1.0);
457 
458  // current position and direction relative to the end point
459  // because of magnetic field geometry is computed relatively to the
460  // end point of the step
461  G4ThreeVector dir(0.0,0.0,1.0);
462  fDisplacement.set(0.0,0.0,-zPathLength);
464 
465  // start a loop
466  G4double x2 = x0;
467  G4double step, z;
468  G4bool singleScat;
469  /*
470  G4cout << "Start of the loop x1(mm)= " << x1 << " x2(mm)= " << x2
471  << " 1-cost1= " << 1 - cosThetaMin << " " << singleScatteringMode
472  << " xtsec= " << xtsec << G4endl;
473  */
474  do {
475 
476  // single scattering case
477  if(x1 < x2) {
478  step = x1;
479  singleScat = true;
480  } else {
481  step = x2;
482  singleScat = false;
483  }
484 
485  // new position
486  fDisplacement += step*mscfac*dir;
487 
488  if(singleScat) {
489 
490  // select element
491  G4int i = 0;
492  if(nelm > 1) {
493  G4double qsec = G4UniformRand()*xtsec;
494  for (; i<nelm; ++i) { if(xsecn[i] >= qsec) { break; } }
495  }
496  G4double cosTetM =
497  wokvi->SetupTarget(G4lrint((*theElementVector)[i]->GetZ()), cut);
498  //G4cout << "!!! " << cosThetaMin << " " << cosTetM
499  // << " " << prob[i] << G4endl;
500  temp = wokvi->SampleSingleScattering(cosThetaMin, cosTetM, prob[i]);
501 
502  // direction is changed
503  temp.rotateUz(dir);
504  dir = temp;
505 
506  // new proposed step length
507  x1 = -G4Log(G4UniformRand())/xtsec;
508  x2 -= step;
509  if(x2 <= 0.0) { --nMscSteps; }
510 
511  // multiple scattering
512  } else {
513  --nMscSteps;
514  x1 -= step;
515  x2 = x0;
516 
517  // for multiple scattering x0 should be used as a step size
518  if(!singleScatteringMode) {
519  G4double z0 = x0*invlambda;
520 
521  // correction to keep first moment
522 
523  // sample z in interval 0 - 1
524  if(z0 > 5.0) { z = G4UniformRand(); }
525  else {
526  G4double zzz = 0.0;
527  if(z0 > 0.01) { zzz = G4Exp(-1.0/z0); }
528  z = -z0*G4Log(1.0 - (1.0 - zzz)*G4UniformRand());
529  // /(1.0 - (1.0/z0 + 1.0)*zzz);
530  }
531 
532  cost = 1.0 - 2.0*z/*factCM*/;
533  if(cost > 1.0) { cost = 1.0; }
534  else if(cost < -1.0) { cost =-1.0; }
535  sint = sqrt((1.0 - cost)*(1.0 + cost));
536  phi = twopi*G4UniformRand();
537  G4double vx1 = sint*cos(phi);
538  G4double vy1 = sint*sin(phi);
539 
540  // lateral displacement
541  if (latDisplasment && x0 > tlimitminfix) {
542  G4double rms = invsqrt12*sqrt(2*z0);
543  G4double r = x0*mscfac;
544  G4double dx = r*(0.5*vx1 + rms*G4RandGauss::shoot(0.0,1.0));
545  G4double dy = r*(0.5*vy1 + rms*G4RandGauss::shoot(0.0,1.0));
546  G4double dz;
547  G4double d = r*r - dx*dx - dy*dy;
548  if(d >= 0.0) { dz = sqrt(d) - r; }
549  else { dx = dy = dz = 0.0; }
550 
551  // change position
552  temp.set(dx,dy,dz);
553  temp.rotateUz(dir);
554  fDisplacement += temp;
555  }
556  // change direction
557  temp.set(vx1,vy1,cost);
558  temp.rotateUz(dir);
559  dir = temp;
560  }
561  }
562  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
563  } while (0 < nMscSteps);
564 
565  dir.rotateUz(oldDirection);
566 
567  //G4cout << "G4WentzelVIRelModel sampling of scattering is done" << G4endl;
568  // end of sampling -------------------------------
569 
571 
572  // lateral displacement
573  fDisplacement.rotateUz(oldDirection);
574 
575  /*
576  G4cout << " r(mm)= " << fDisplacement.mag()
577  << " safety= " << safety
578  << " trueStep(mm)= " << tPathLength
579  << " geomStep(mm)= " << zPathLength
580  << " x= " << fDisplacement.x()
581  << " y= " << fDisplacement.y()
582  << " z= " << fDisplacement.z()
583  << G4endl;
584  */
585 
586  //G4cout<< "G4WentzelVIRelModel::SampleScattering end NewDir= "
587  // << dir<< G4endl;
588  return fDisplacement;
589 }
590 
591 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
592 
594 {
595  // prepare recomputation of x-sections
596  const G4ElementVector* theElementVector = currentMaterial->GetElementVector();
597  const G4double* theAtomNumDensityVector =
600  if(nelm > nelments) {
601  nelments = nelm;
602  xsecn.resize(nelm);
603  prob.resize(nelm);
604  }
605  G4double cut = (*currentCuts)[currentMaterialIndex];
606  // cosTetMaxNuc = wokvi->GetCosThetaNuc();
607 
608  // check consistency
609  xtsec = 0.0;
610  if(cosTetMaxNuc > cosThetaMin) { return 0.0; }
611 
612  // loop over elements
613  G4double xs = 0.0;
614  for (G4int i=0; i<nelm; ++i) {
615  G4double costm =
616  wokvi->SetupTarget(G4lrint((*theElementVector)[i]->GetZ()), cut);
617  G4double density = theAtomNumDensityVector[i];
618 
619  G4double esec = 0.0;
620  if(costm < cosThetaMin) {
621 
622  // recompute the transport x-section
623  if(1.0 > cosThetaMin) {
625  }
626  // recompute the total x-section
629  nucsec += esec;
630  if(nucsec > 0.0) { esec /= nucsec; }
631  xtsec += nucsec*density;
632  }
633  xsecn[i] = xtsec;
634  prob[i] = esec;
635  //G4cout << i << " xs= " << xs << " xtsec= " << xtsec
636  // << " 1-cosThetaMin= " << 1-cosThetaMin
637  // << " 1-cosTetMaxNuc2= " <<1-cosTetMaxNuc2<< G4endl;
638  }
639 
640  //G4cout << "ComputeXS result: xsec(1/mm)= " << xs
641  // << " txsec(1/mm)= " << xtsec <<G4endl;
642  return xs;
643 }
644 
645 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4double facgeom
Definition: G4VMscModel.hh:178
void StartTracking(G4Track *)
ThreeVector shoot(const G4int Ap, const G4int Af)
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &currentMinimalStep)
virtual void StartTracking(G4Track *)
Definition: G4VEmModel.cc:284
static G4LossTableManager * Instance()
virtual G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety)
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
G4double ConvertTrueToGeom(G4double &tLength, G4double &gLength)
Definition: G4VMscModel.hh:247
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4bool latDisplasment
Definition: G4VMscModel.hh:190
virtual G4double ComputeGeomPathLength(G4double truePathLength)
G4ParticleDefinition * GetDefinition() const
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight=0., G4double cut=DBL_MAX, G4double emax=DBL_MAX)
G4WentzelVIRelXSection * wokvi
G4double ComputeTransportCrossSectionPerAtom(G4double CosThetaMax)
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
static G4NistManager * Instance()
G4StepStatus
Definition: G4StepStatus.hh:49
void DefineMaterial(const G4MaterialCutsCouple *)
G4StepPoint * GetPreStepPoint() const
G4double GetEnergy(const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.hh:315
G4double density
Definition: TRTMaterials.hh:39
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:206
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:452
#define G4UniformRand()
Definition: Randomize.hh:97
const G4ThreeVector & GetPosition() const
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.hh:295
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat, G4double cut, G4double tmass)
G4double SetupTarget(G4int Z, G4double cut)
G4ParticleChangeForMSC * GetParticleChangeForMSC(const G4ParticleDefinition *p=0)
Definition: G4VMscModel.cc:90
bool G4bool
Definition: G4Types.hh:79
G4double iz
Definition: TRTMaterials.hh:39
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
static const double twopi
Definition: G4SIunits.hh:75
std::vector< G4double > prob
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:352
G4double GetRadlen() const
Definition: G4Material.hh:220
const G4DataVector * currentCuts
G4NistManager * fNistManager
static const G4double e1
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4WentzelVIRelModel(G4bool combined=true)
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4LossTableManager * theManager
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
static const double eV
Definition: G4SIunits.hh:212
static const double pi
Definition: G4SIunits.hh:74
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
virtual G4double ComputeTrueStepLength(G4double geomStepLength)
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:662
G4double GetAtomicMassAmu(const G4String &symb) const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double GetSafety() const
G4MscStepLimitType steppingAlgorithm
Definition: G4VMscModel.hh:187
std::vector< G4double > xsecn
G4ParticleChangeForMSC * fParticleChange
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
const G4ParticleDefinition * particle
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
const G4MaterialCutsCouple * currentCouple
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)
void SetupParticle(const G4ParticleDefinition *)
double G4double
Definition: G4Types.hh:76
G4ThreeVector & SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio)
const G4Material * currentMaterial
G4ProductionCuts * GetProductionCuts() const
#define DBL_MAX
Definition: templates.hh:83
static const double mm
Definition: G4SIunits.hh:114