Geant4  10.01.p03
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 88979 2015-03-17 10:10:21Z 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  }
156  if(cosTetMaxNuc < 1.0) {
157  G4double cost = wokvi->SetupTarget(G4lrint(Z), cutEnergy);
159  /*
160  if(p->GetParticleName() == "e-")
161  G4cout << "G4WentzelVIRelModel::CS: Z= " << G4int(Z)
162  << " e(MeV)= " << kinEnergy
163  << " 1-cosN= " << 1 - cosTetMaxNuc << " cross(bn)= " << cross/barn
164  << " " << particle->GetParticleName() << G4endl;
165  */
166  }
167  return cross;
168 }
169 
170 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
171 
173 {
175  inside = false;
177 }
178 
179 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
180 
182  const G4Track& track,
183  G4double& currentMinimalStep)
184 {
185  G4double tlimit = currentMinimalStep;
186  const G4DynamicParticle* dp = track.GetDynamicParticle();
187  G4StepPoint* sp = track.GetStep()->GetPreStepPoint();
188  G4StepStatus stepStatus = sp->GetStepStatus();
189  singleScatteringMode = false;
190  //G4cout << "G4WentzelVIRelModel::ComputeTruePathLengthLimit stepStatus= "
191  // << stepStatus << G4endl;
192 
193 
194  // initialisation for each step, lambda may be computed from scratch
200 
201  // extra check for abnormal situation
202  // this check needed to run MSC with eIoni and eBrem inactivated
203  if(tlimit > currentRange) { tlimit = currentRange; }
204 
205  // stop here if small range particle
206  if(inside || tlimit < tlimitminfix) {
207  return ConvertTrueToGeom(tlimit, currentMinimalStep);
208  }
209 
210  // pre step
211  G4double presafety = sp->GetSafety();
212  // far from geometry boundary
213  if(currentRange < presafety) {
214  inside = true;
215  return ConvertTrueToGeom(tlimit, currentMinimalStep);
216  }
217 
218  // compute presafety again if presafety <= 0 and no boundary
219  // i.e. when it is needed for optimization purposes
220  if(stepStatus != fGeomBoundary && presafety < tlimitminfix) {
221  presafety = ComputeSafety(sp->GetPosition(), tlimit);
222  if(currentRange < presafety) {
223  inside = true;
224  return ConvertTrueToGeom(tlimit, currentMinimalStep);
225  }
226  }
227  /*
228  G4cout << "e(MeV)= " << preKinEnergy/MeV
229  << " " << particle->GetParticleName()
230  << " CurLimit(mm)= " << tlimit/mm <<" safety(mm)= " << presafety/mm
231  << " R(mm)= " <<currentRange/mm
232  << " L0(mm^-1)= " << lambdaeff*mm
233  <<G4endl;
234  */
235 
236  // natural limit for high energy
238  0.7*(1.0 - cosTetMaxNuc)*lambdaeff);
239 
240  // low-energy e-
241  if(cosThetaMax > cosTetMaxNuc) {
242  rlimit = std::min(rlimit, facsafety*presafety);
243  }
244 
245  // cut correction
247  //G4cout << "rcut= " << rcut << " rlimit= " << rlimit
248  // << " presafety= " << presafety
249  // << " 1-cosThetaMax= " <<1-cosThetaMax
250  // << " 1-cosTetMaxNuc= " << 1-cosTetMaxNuc
251  // << G4endl;
252  if(rcut > rlimit) { rlimit = std::min(rlimit, rcut*sqrt(rlimit/rcut)); }
253 
254  if(rlimit < tlimit) { tlimit = rlimit; }
255 
256  tlimit = std::max(tlimit, tlimitminfix);
257 
258  // step limit in infinite media
259  tlimit = std::min(tlimit, 50*currentMaterial->GetRadlen()/facgeom);
260 
261  //compute geomlimit and force few steps within a volume
263  {
264  G4double geomlimit = ComputeGeomLimit(track, presafety, currentRange);
265  tlimit = std::min(tlimit, geomlimit/facgeom);
266  }
267 
268  /*
269  G4cout << particle->GetParticleName() << " e= " << preKinEnergy
270  << " L0= " << lambdaeff << " R= " << currentRange
271  << "tlimit= " << tlimit
272  << " currentMinimalStep= " << currentMinimalStep << G4endl;
273  */
274  return ConvertTrueToGeom(tlimit, currentMinimalStep);
275 }
276 
277 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
278 
280 {
281  tPathLength = truelength;
283 
284  if(lambdaeff > 0.0 && lambdaeff != DBL_MAX) {
286  //G4cout << "ComputeGeomPathLength: tLength= " << tPathLength
287  // << " Leff= " << lambdaeff << " tau= " << tau << G4endl;
288  // small step
289  if(tau < numlimit) {
290  zPathLength *= (1.0 - 0.5*tau + tau*tau/6.0);
291 
292  // medium step
293  } else {
294  G4double e1 = 0.0;
295  if(currentRange > tPathLength) {
297  }
298  e1 = 0.5*(e1 + preKinEnergy);
302  }
303  } else { lambdaeff = DBL_MAX; }
304  //G4cout<<"Comp.geom: zLength= "<<zPathLength
305  // <<" tLength= "<<tPathLength<<G4endl;
306  return zPathLength;
307 }
308 
309 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
310 
312 {
313  // initialisation of single scattering x-section
314  xtsec = 0.0;
316 
317  //G4cout << "Step= " << geomStepLength << " Lambda= " << lambdaeff
318  // << " 1-cosThetaMaxNuc= " << 1 - cosTetMaxNuc << G4endl;
319  // pathalogical case
320  if(lambdaeff == DBL_MAX) {
321  singleScatteringMode = true;
322  zPathLength = geomStepLength;
323  tPathLength = geomStepLength;
324  cosThetaMin = 1.0;
325 
326  // normal case
327  } else {
328 
329  // small step use only single scattering
330  static const G4double singleScatLimit = 1.0e-7;
331  if(geomStepLength < lambdaeff*singleScatLimit*(1.0 - cosTetMaxNuc)) {
332  singleScatteringMode = true;
333  cosThetaMin = 1.0;
334  zPathLength = geomStepLength;
335  tPathLength = geomStepLength;
336 
337  // step defined by transportation
338  } else if(geomStepLength != zPathLength) {
339 
340  // step defined by transportation
341  zPathLength = geomStepLength;
342  G4double tau = geomStepLength/lambdaeff;
343  tPathLength = zPathLength*(1.0 + 0.5*tau + tau*tau/3.0);
344 
345  // energy correction for a big step
346  if(tau > numlimit) {
347  G4double e1 = 0.0;
348  if(currentRange > tPathLength) {
350  }
351  e1 = 0.5*(e1 + preKinEnergy);
354  tau = zPathLength/lambdaeff;
355 
356  if(tau < 0.999999) { tPathLength = -lambdaeff*G4Log(1.0 - tau); }
357  else { tPathLength = currentRange; }
358  }
359  }
360  }
361 
362  // check of step length
363  // define threshold angle between single and multiple scattering
365 
366  // recompute transport cross section - do not change energy
367  // anymore - cannot be applied for big steps
368  if(cosThetaMin > cosTetMaxNuc) {
369 
370  // new computation
372  //G4cout << "%%%% cross= " << cross << " xtsec= " << xtsec << G4endl;
373  if(cross <= 0.0) {
374  singleScatteringMode = true;
376  lambdaeff = DBL_MAX;
377  } else if(xtsec > 0.0) {
378 
379  lambdaeff = 1./cross;
380  G4double tau = zPathLength*cross;
381  if(tau < numlimit) {
382  tPathLength = zPathLength*(1.0 + 0.5*tau + tau*tau/3.0);
383  }
384  else if(tau < 0.999999) { tPathLength = -lambdaeff*G4Log(1.0 - tau); }
385  else { tPathLength = currentRange; }
386 
388  }
389  }
390 
391  /*
392  G4cout <<"Comp.true: zLength= "<<zPathLength<<" tLength= "<<tPathLength
393  <<" Leff(mm)= "<<lambdaeff/mm<<" sig0(1/mm)= " << xtsec <<G4endl;
394  G4cout << particle->GetParticleName() << " 1-cosThetaMin= " << 1-cosThetaMin
395  << " 1-cosTetMaxNuc= " << 1-cosTetMaxNuc
396  << " e(MeV)= " << preKinEnergy/MeV << " " << singleScatteringMode
397  << G4endl;
398  */
399  return tPathLength;
400 }
401 
402 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
403 
406  G4double safety)
407 {
408  fDisplacement.set(0.0,0.0,0.0);
409  //G4cout << "!##! G4WentzelVIRelModel::SampleScattering for "
410  // << particle->GetParticleName() << G4endl;
411 
412  // ignore scattering for zero step length and energy below the limit
413  if(preKinEnergy < lowEnergyLimit || tPathLength <= 0.0)
414  { return fDisplacement; }
415 
416  G4double invlambda = 0.0;
417  if(lambdaeff < DBL_MAX) { invlambda = 0.5/lambdaeff; }
418 
419  // use average kinetic energy over the step
420  G4double cut = (*currentCuts)[currentMaterialIndex];
421  /*
422  G4cout <<"SampleScat: E0(MeV)= "<< preKinEnergy/MeV
423  << " Leff= " << lambdaeff <<" sig0(1/mm)= " << xtsec
424  << " xmsc= " << tPathLength*invlambda
425  << " safety= " << safety << G4endl;
426  */
427 
428  // step limit due msc
429  G4double x0 = tPathLength;
430  // const G4double thinlimit = 0.5;
431  static const G4double thinlimit = 0.1;
432  G4int nMscSteps = 1;
433  // large scattering angle case - two step approach
434  if(tPathLength*invlambda > thinlimit && safety > tlimitminfix) {
435  x0 *= 0.5;
436  nMscSteps = 2;
437  }
438 
439  // step limit due to single scattering
440  G4double x1 = 2*tPathLength;
441  if(0.0 < xtsec) { x1 = -G4Log(G4UniformRand())/xtsec; }
442 
443  const G4ElementVector* theElementVector =
446 
447  // geometry
448  G4double sint, cost, phi;
449  G4ThreeVector temp(0.0,0.0,1.0);
450 
451  // current position and direction relative to the end point
452  // because of magnetic field geometry is computed relatively to the
453  // end point of the step
454  G4ThreeVector dir(0.0,0.0,1.0);
455  fDisplacement.set(0.0,0.0,-zPathLength);
457 
458  // start a loop
459  G4double x2 = x0;
460  G4double step, z;
461  G4bool singleScat;
462  /*
463  G4cout << "Start of the loop x1(mm)= " << x1 << " x2(mm)= " << x2
464  << " 1-cost1= " << 1 - cosThetaMin << " " << singleScatteringMode
465  << " xtsec= " << xtsec << G4endl;
466  */
467  do {
468 
469  // single scattering case
470  if(x1 < x2) {
471  step = x1;
472  singleScat = true;
473  } else {
474  step = x2;
475  singleScat = false;
476  }
477 
478  // new position
479  fDisplacement += step*mscfac*dir;
480 
481  if(singleScat) {
482 
483  // select element
484  G4int i = 0;
485  if(nelm > 1) {
486  G4double qsec = G4UniformRand()*xtsec;
487  for (; i<nelm; ++i) { if(xsecn[i] >= qsec) { break; } }
488  }
489  G4double cosTetM =
490  wokvi->SetupTarget(G4lrint((*theElementVector)[i]->GetZ()), cut);
491  //G4cout << "!!! " << cosThetaMin << " " << cosTetM
492  // << " " << prob[i] << G4endl;
493  temp = wokvi->SampleSingleScattering(cosThetaMin, cosTetM, prob[i]);
494 
495  // direction is changed
496  temp.rotateUz(dir);
497  dir = temp;
498 
499  // new proposed step length
500  x1 = -G4Log(G4UniformRand())/xtsec;
501  x2 -= step;
502  if(x2 <= 0.0) { --nMscSteps; }
503 
504  // multiple scattering
505  } else {
506  --nMscSteps;
507  x1 -= step;
508  x2 = x0;
509 
510  // for multiple scattering x0 should be used as a step size
511  if(!singleScatteringMode) {
512  G4double z0 = x0*invlambda;
513 
514  // correction to keep first moment
515 
516  // sample z in interval 0 - 1
517  if(z0 > 5.0) { z = G4UniformRand(); }
518  else {
519  G4double zzz = 0.0;
520  if(z0 > 0.01) { zzz = G4Exp(-1.0/z0); }
521  z = -z0*G4Log(1.0 - (1.0 - zzz)*G4UniformRand());
522  // /(1.0 - (1.0/z0 + 1.0)*zzz);
523  }
524 
525  cost = 1.0 - 2.0*z/*factCM*/;
526  if(cost > 1.0) { cost = 1.0; }
527  else if(cost < -1.0) { cost =-1.0; }
528  sint = sqrt((1.0 - cost)*(1.0 + cost));
529  phi = twopi*G4UniformRand();
530  G4double vx1 = sint*cos(phi);
531  G4double vy1 = sint*sin(phi);
532 
533  // lateral displacement
534  if (latDisplasment && x0 > tlimitminfix) {
535  G4double rms = invsqrt12*sqrt(2*z0);
536  G4double r = x0*mscfac;
537  G4double dx = r*(0.5*vx1 + rms*G4RandGauss::shoot(0.0,1.0));
538  G4double dy = r*(0.5*vy1 + rms*G4RandGauss::shoot(0.0,1.0));
539  G4double dz;
540  G4double d = r*r - dx*dx - dy*dy;
541  if(d >= 0.0) { dz = sqrt(d) - r; }
542  else { dx = dy = dz = 0.0; }
543 
544  // change position
545  temp.set(dx,dy,dz);
546  temp.rotateUz(dir);
547  fDisplacement += temp;
548  }
549  // change direction
550  temp.set(vx1,vy1,cost);
551  temp.rotateUz(dir);
552  dir = temp;
553  }
554  }
555  } while (0 < nMscSteps);
556 
557  dir.rotateUz(oldDirection);
558 
559  //G4cout << "G4WentzelVIRelModel sampling of scattering is done" << G4endl;
560  // end of sampling -------------------------------
561 
563 
564  // lateral displacement
565  fDisplacement.rotateUz(oldDirection);
566 
567  /*
568  G4cout << " r(mm)= " << fDisplacement.mag()
569  << " safety= " << safety
570  << " trueStep(mm)= " << tPathLength
571  << " geomStep(mm)= " << zPathLength
572  << " x= " << fDisplacement.x()
573  << " y= " << fDisplacement.y()
574  << " z= " << fDisplacement.z()
575  << G4endl;
576  */
577 
578  //G4cout<< "G4WentzelVIRelModel::SampleScattering end NewDir= "
579  // << dir<< G4endl;
580  return fDisplacement;
581 }
582 
583 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
584 
586 {
587  // prepare recomputation of x-sections
588  const G4ElementVector* theElementVector = currentMaterial->GetElementVector();
589  const G4double* theAtomNumDensityVector =
592  if(nelm > nelments) {
593  nelments = nelm;
594  xsecn.resize(nelm);
595  prob.resize(nelm);
596  }
597  G4double cut = (*currentCuts)[currentMaterialIndex];
598  // cosTetMaxNuc = wokvi->GetCosThetaNuc();
599 
600  // check consistency
601  xtsec = 0.0;
602  if(cosTetMaxNuc > cosThetaMin) { return 0.0; }
603 
604  // loop over elements
605  G4double xs = 0.0;
606  for (G4int i=0; i<nelm; ++i) {
607  G4double costm =
608  wokvi->SetupTarget(G4lrint((*theElementVector)[i]->GetZ()), cut);
609  G4double density = theAtomNumDensityVector[i];
610 
611  G4double esec = 0.0;
612  if(costm < cosThetaMin) {
613 
614  // recompute the transport x-section
615  if(1.0 > cosThetaMin) {
617  }
618  // recompute the total x-section
621  nucsec += esec;
622  if(nucsec > 0.0) { esec /= nucsec; }
623  xtsec += nucsec*density;
624  }
625  xsecn[i] = xtsec;
626  prob[i] = esec;
627  //G4cout << i << " xs= " << xs << " xtsec= " << xtsec
628  // << " 1-cosThetaMin= " << 1-cosThetaMin
629  // << " 1-cosTetMaxNuc2= " <<1-cosTetMaxNuc2<< G4endl;
630  }
631 
632  //G4cout << "ComputeXS result: xsec(1/mm)= " << xs
633  // << " txsec(1/mm)= " << xtsec <<G4endl;
634  return xs;
635 }
636 
637 //....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)
G4ThreeVector SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio=0.0)
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
const G4double pi
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()
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
G4double SetupTarget(G4int Z, G4double cut=DBL_MAX)
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:309
G4double density
Definition: TRTMaterials.hh:39
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:206
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:434
#define G4UniformRand()
Definition: Randomize.hh:93
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
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
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:346
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:194
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:644
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
const G4Material * currentMaterial
G4ProductionCuts * GetProductionCuts() const
#define DBL_MAX
Definition: templates.hh:83
static const double mm
Definition: G4SIunits.hh:102