Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4GoudsmitSaundersonMscModel.cc
이 파일의 문서화 페이지로 가기
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: G4GoudsmitSaundersonMscModel.cc 108305 2018-02-02 13:08:43Z gcosmo $
27 //
28 // ----------------------------------------------------------------------------
29 //
30 // GEANT4 Class implementation file
31 //
32 // File name: G4GoudsmitSaundersonMscModel
33 //
34 // Author: Mihaly Novak / (Omrane Kadri)
35 //
36 // Creation date: 20.02.2009
37 //
38 // Modifications:
39 // 04.03.2009 V.Ivanchenko cleanup and format according to Geant4 EM style
40 // 12.05.2010 O.Kadri: adding Qn1 and Qn12 as private doubles
41 // 18.05.2015 M. Novak provide PRELIMINARY verison of the revised class
42 // This class has been revised and updated, new methods added.
43 // A new version of Kawrakow-Bielajew Goudsmit-Saunderson MSC model
44 // based on the screened Rutherford DCS for elastic scattering of
45 // electrons/positrons has been introduced[1,2]. The corresponding MSC
46 // angular distributions over a 2D parameter grid have been recomputed
47 // and the CDFs are now stored in a variable transformed (smooth) form[2,3]
48 // together with the corresponding rational interpolation parameters.
49 // These angular distributions are handled by the new
50 // G4GoudsmitSaundersonTable class that is responsible to sample if
51 // it was no, single, few or multiple scattering case and delivers the
52 // angular deflection (i.e. cos(theta) and sin(theta)).
53 // Two screening options are provided:
54 // - if fIsUsePWATotalXsecData=TRUE i.e. SetOptionPWAScreening(TRUE)
55 // was called before initialisation: screening parameter value A is
56 // determined such that the first transport coefficient G1(A)
57 // computed according to the screened Rutherford DCS for elastic
58 // scattering will reproduce the one computed from the PWA elastic
59 // and first transport mean free paths[4].
60 // - if fIsUsePWATotalXsecData=FALSE i.e. default value or
61 // SetOptionPWAScreening(FALSE) was called before initialisation:
62 // screening parameter value A is computed according to Moliere's
63 // formula (by using material dependent parameters \chi_cc2 and b_c
64 // precomputed for each material used at initialization in
65 // G4GoudsmitSaundersonTable) [3]
66 // Elastic and first trasport mean free paths are used consistently.
67 // The new version is self-consistent, several times faster, more
68 // robust and accurate compared to the earlier version.
69 // Spin effects as well as a more accurate energy loss correction and
70 // computations of Lewis moments will be implemented later on.
71 // 02.09.2015 M. Novak: first version of new step limit is provided.
72 // fUseSafetyPlus corresponds to Urban fUseSafety (default)
73 // fUseDistanceToBoundary corresponds to Urban fUseDistanceToBoundary
74 // fUseSafety corresponds to EGSnrc error-free stepping algorithm
75 // Range factor can be significantly higher at each case than in Urban.
76 // 23.08.2017 M. Novak: added corrections to account spin effects (Mott-correction).
77 // It can be activated by setting the fIsMottCorrection flag to be true
78 // before initialization using the SetOptionMottCorrection() public method.
79 // The fMottCorrection member is responsible to handle pre-computed Mott
80 // correction (rejection) functions obtained by numerically computing
81 // Goudsmit-Saunderson agnular distributions based on a DCS accounting spin
82 // effects and screening corrections. The DCS used to compute the accurate
83 // GS angular distributions is: DCS_{cor} = DCS_{SR}x[ DCS_{R}/DCS_{Mott}] where :
84 // # DCS_{SR} is the relativistic Screened-Rutherford DCS (first Born approximate
85 // solution of the Klein-Gordon i.e. relativistic Schrodinger equation =>
86 // scattering of spinless e- on exponentially screened Coulomb potential)
87 // note: the default (without using Mott-correction) GS angular distributions
88 // are based on this DCS_{SR} with Moliere's screening parameter!
89 // # DCS_{R} is the Rutherford DCS which is the same as above but without
90 // screening
91 // # DCS_{Mott} is the Mott DCS i.e. solution of the Dirac equation with a bare
92 // Coulomb potential i.e. scattering of particles with spin (e- or e+) on a
93 // point-like unscreened Coulomb potential
94 // # moreover, the screening parameter of the DCS_{cor} was determined such that
95 // the DCS_{cor} with this corrected screening parameter reproduce the first
96 // transport cross sections obtained from the corresponding most accurate DCS
97 // (i.e. from elsepa [4])
98 // Unlike the default GS, the Mott-corrected angular distributions are particle type
99 // (different for e- and e+ <= the DCS_{Mott} and the screening correction) and target
100 // (Z and material) dependent.
101 // 27.10.2017 M. Novak:
102 // - Mott-correction flag is set now through the G4EmParameters
103 // - new form of PWA correction to integrated quantities and screening (default)
104 // - changed step limit flag conventions:
105 // # fUseSafety corresponds to Urban's fUseSafety
106 // # fUseDistanceToBoundary corresponds to Urban's fUseDistanceToBoundary
107 // # fUseSafetyPlus corresponds to the error-free stepping algorithm
108 // 02.02.2018 M. Novak: implemented CrossSectionPerVolume interface method (used only for testing)
109 //
110 // Class description:
111 // Kawrakow-Bielajew Goudsmit-Saunderson MSC model based on the screened Rutherford DCS
112 // for elastic scattering of e-/e+. Option, to include (Mott) correction (see above), is
113 // also available now (SetOptionMottCorrection(true)). An EGSnrc like error-free stepping
114 // algorithm (UseSafety) is available beyond the usual Geant4 step limitation algorithms
115 // and true to geomerty and geometry to true step length computations that were adopted
116 // from the Urban model[5]. The most accurate setting: error-free stepping i.e. the
117 // UseSafetyPlus MSC step limit with Mott-correction (SetOptionMottCorrection(true)). Both
118 // are expected to be set through the G4EmParameters singleton before initialisation:
119 // # G4EmParameters::Instance()->SetMscStepLimitType(fUseSafetyPlus);
120 // # G4EmParameters::Instance()->SetUseMottCorrection(true);
121 //
122 //
123 // References:
124 // [1] A.F.Bielajew, NIMB 111 (1996) 195-208
125 // [2] I.Kawrakow, A.F.Bielajew, NIMB 134(1998) 325-336
126 // [3] I.Kawrakow, E.Mainegra-Hing, D.W.O.Rogers, F.Tessier,B.R.B.Walters, NRCC
127 // Report PIRS-701 (2013)
128 // [4] F.Salvat, A.Jablonski, C.J. Powell, CPC 165(2005) 157-190
129 // [5] L.Urban, Preprint CERN-OPEN-2006-077 (2006)
130 //
131 // -----------------------------------------------------------------------------
132 
133 
135 
137 #include "G4GSPWACorrections.hh"
138 
139 #include "G4PhysicalConstants.hh"
140 #include "G4SystemOfUnits.hh"
141 
142 #include "G4ParticleChangeForMSC.hh"
143 #include "G4DynamicParticle.hh"
144 #include "G4Electron.hh"
145 #include "G4Positron.hh"
146 
147 #include "G4LossTableManager.hh"
148 #include "G4EmParameters.hh"
149 #include "G4Track.hh"
150 #include "G4PhysicsTable.hh"
151 #include "Randomize.hh"
152 #include "G4Log.hh"
153 #include "G4Exp.hh"
154 #include "G4Pow.hh"
155 #include <fstream>
156 
157 
158 // set accurate energy loss and dispalcement sampling to be always on now
160 // set the usual optimization to be always active now
162 
163 
165  : G4VMscModel(nam) {
166  charge = 0;
168  //
169  lambdalimit = 1.*mm;
170  fr = 0.1;
171  rangeinit = 1.e+21;
172  geombig = 1.e+50*mm;
173  geomlimit = geombig;
174  tgeom = geombig;
175  tlimit = 1.e+10*mm;
176  presafety = 0.*mm;
177  //
178  particle = 0;
180  firstStep = true;
181  currentKinEnergy = 0.0;
182  currentRange = 0.0;
183  //
184  tlimitminfix2 = 1.*nm;
185  tausmall = 1.e-16;
187  taulim = 1.e-6;
188  //
189  facsafety = 0.6;
190 
191  currentCouple = nullptr;
192  fParticleChange = nullptr;
193  //
194  fZeff = 1.;
195  //
196  par1 = 0.;
197  par2 = 0.;
198  par3 = 0.;
199  //
200  // Moliere screeing parameter will be used and (by default) corrections are
201  // appalied to the integrated quantities (screeing parameter, elastic mfp, first
202  // and second moments) derived from the corresponding PWA quantities
203  // this PWA correction is ignored if Mott-correction is set to true because
204  // Mott-correction contains all these corrections as well
205  fIsUsePWACorrection = true;
206  //
207  fIsUseMottCorrection = false;
208  //
209  fLambda0 = 0.0; // elastic mean free path
210  fLambda1 = 0.0; // first transport mean free path
211  fScrA = 0.0; // screening parameter
212  fG1 = 0.0; // first transport coef.
213  //
214  fMCtoScrA = 1.0;
215  fMCtoQ1 = 1.0;
216  fMCtoG2PerG1 = 1.0;
217  //
218  fTheTrueStepLenght = 0.;
220  fTheZPathLenght = 0.;
221  //
222  fTheDisplacementVector.set(0.,0.,0.);
223  fTheNewDirection.set(0.,0.,1.);
224  //
225  fIsEverythingWasDone = false;
226  fIsMultipleSacettring = false;
227  fIsSingleScattering = false;
228  fIsEndedUpOnBoundary = false;
229  fIsNoScatteringInMSC = false;
230  fIsNoDisplace = false;
231  fIsInsideSkin = false;
232  fIsWasOnBoundary = false;
233  fIsFirstRealStep = false;
234  rndmEngineMod = G4Random::getTheEngine();
235  //
236  fGSTable = nullptr;
237  fPWACorrection = nullptr;
238 }
239 
240 
242  if (IsMaster()) {
243  if (fGSTable) {
244  delete fGSTable;
245  fGSTable = nullptr;
246  }
247  if (fPWACorrection) {
248  delete fPWACorrection;
249  fPWACorrection = nullptr;
250  }
251  }
252 }
253 
254 
256  SetParticle(p);
257  // -create GoudsmitSaundersonTable and init its Mott-correction member if
258  // Mott-correction was required
259  if (IsMaster()) {
260  // get the Mott-correction flag from EmParameters
261  if (G4EmParameters::Instance()->UseMottCorrection()) {
262  fIsUseMottCorrection = true;
263  }
264  // Mott-correction includes other way of PWA x-section corrections so deactivate it even if it was true
265  // when Mott-correction is activated by the user
266  if (fIsUseMottCorrection) {
267  fIsUsePWACorrection = false;
268  }
269  // clear GS-table
270  if (fGSTable) {
271  delete fGSTable;
272  fGSTable = nullptr;
273  }
274  // clear PWA corrections table if any
275  if (fPWACorrection) {
276  delete fPWACorrection;
277  fPWACorrection = nullptr;
278  }
279  // create GS-table
280  G4bool isElectron = true;
281  if (p->GetPDGCharge()>0.) {
282  isElectron = false;
283  }
284  fGSTable = new G4GoudsmitSaundersonTable(isElectron);
285  // G4GSTable will be initialised:
286  // - Screened-Rutherford DCS based GS angular distributions will be loaded only if they are not there yet
287  // - Mott-correction will be initialised if Mott-correction was requested to be used
289  // - set PWA correction (correction to integrated quantites from Dirac-PWA)
291  // init
293  // create PWA corrections table if it was requested (and not disactivated because active Mott-correction)
294  if (fIsUsePWACorrection) {
295  fPWACorrection = new G4GSPWACorrections(isElectron);
297  }
298  }
300 }
301 
302 
304  fGSTable = static_cast<G4GoudsmitSaundersonMscModel*>(masterModel)->GetGSTable();
307  fPWACorrection = static_cast<G4GoudsmitSaundersonMscModel*>(masterModel)->GetPWACorrection();
308 }
309 
310 
311 // computes macroscopic first transport cross section: used only in testing not during mc transport
313  const G4ParticleDefinition*,
315  G4double,
316  G4double) {
317  G4double xsecTr1 = 0.; // cross section per volume i.e. macroscopic 1st transport cross section
318  G4double efEnergy = kineticEnergy;
319  //
320  fLambda0 = 0.0; // elastic mean free path
321  fLambda1 = 0.0; // first transport mean free path
322  fScrA = 0.0; // screening parameter
323  fG1 = 0.0; // first transport coef.
324  // use Moliere's screening (with Mott-corretion if it was requested)
325  if (efEnergy<10.*CLHEP::eV) efEnergy = 10.*CLHEP::eV;
326  // total mometum square
327  G4double pt2 = efEnergy*(efEnergy+2.0*electron_mass_c2);
328  // beta square
329  G4double beta2 = pt2/(pt2+electron_mass_c2*electron_mass_c2);
330  // current material index
331  G4int matindx = mat->GetIndex();
332  // Moliere's b_c
333  G4double bc = fGSTable->GetMoliereBc(matindx);
334  // get the Mott-correcton factors if Mott-correcton was requested by the user
335  fMCtoScrA = 1.0;
336  fMCtoQ1 = 1.0;
337  fMCtoG2PerG1 = 1.0;
338  G4double scpCor = 1.0;
339  if (fIsUseMottCorrection) {
340  fGSTable->GetMottCorrectionFactors(G4Log(efEnergy), beta2, matindx, fMCtoScrA, fMCtoQ1, fMCtoG2PerG1);
341  // ! no scattering power correction since the current couple is not set before this interface method is called
342  // scpCor = fGSTable->ComputeScatteringPowerCorrection(currentCouple, efEnergy);
343  } else if (fIsUsePWACorrection) {
345  // scpCor = fGSTable->ComputeScatteringPowerCorrection(currentCouple, efEnergy);
346  }
347  // screening parameter:
348  // - if Mott-corretioncorrection: the Screened-Rutherford times Mott-corretion DCS with this
349  // screening parameter gives back the (elsepa) PWA first transport cross section
350  // - if PWA correction: he Screened-Rutherford DCS with this screening parameter
351  // gives back the (elsepa) PWA first transport cross section
352  fScrA = fGSTable->GetMoliereXc2(matindx)/(4.0*pt2*bc)*fMCtoScrA;
353  // elastic mean free path in Geant4 internal lenght units: the neglected (1+screening parameter) term is corrected
354  // (if Mott-corretion: the corrected screening parameter is used for this (1+A) correction + Moliere b_c is also
355  // corrected with the screening parameter correction)
356  fLambda0 = beta2*(1.+fScrA)*fMCtoScrA/bc/scpCor;
357  // first transport coefficient (if Mott-corretion: the corrected screening parameter is used (it will be fully
358  // consistent with the one used during the pre-computation of the Mott-correted GS angular distributions))
359  fG1 = 2.0*fScrA*((1.0+fScrA)*G4Log(1.0/fScrA+1.0)-1.0);
360  // first transport mean free path
362  xsecTr1 = 1./fLambda1;
363  return xsecTr1;
364 }
365 
366 
367 // gives back the first transport mean free path in internal G4 units
368 G4double
371  // kinetic energy is assumed to be in Geant4 internal energy unit which is MeV
372  G4double efEnergy = kineticEnergy;
373  //
375  //
376  fLambda0 = 0.0; // elastic mean free path
377  fLambda1 = 0.0; // first transport mean free path
378  fScrA = 0.0; // screening parameter
379  fG1 = 0.0; // first transport coef.
380 
381  // use Moliere's screening (with Mott-corretion if it was requested)
382  if (efEnergy<10.*CLHEP::eV) efEnergy = 10.*CLHEP::eV;
383  // total mometum square
384  G4double pt2 = efEnergy*(efEnergy+2.0*electron_mass_c2);
385  // beta square
386  G4double beta2 = pt2/(pt2+electron_mass_c2*electron_mass_c2);
387  // current material index
388  G4int matindx = mat->GetIndex();
389  // Moliere's b_c
390  G4double bc = fGSTable->GetMoliereBc(matindx);
391  // get the Mott-correcton factors if Mott-correcton was requested by the user
392  fMCtoScrA = 1.0;
393  fMCtoQ1 = 1.0;
394  fMCtoG2PerG1 = 1.0;
395  G4double scpCor = 1.0;
396  if (fIsUseMottCorrection) {
397  fGSTable->GetMottCorrectionFactors(G4Log(efEnergy), beta2, matindx, fMCtoScrA, fMCtoQ1, fMCtoG2PerG1);
399  } else if (fIsUsePWACorrection) {
401  // scpCor = fGSTable->ComputeScatteringPowerCorrection(currentCouple, efEnergy);
402  }
403  // screening parameter:
404  // - if Mott-corretioncorrection: the Screened-Rutherford times Mott-corretion DCS with this
405  // screening parameter gives back the (elsepa) PWA first transport cross section
406  // - if PWA correction: he Screened-Rutherford DCS with this screening parameter
407  // gives back the (elsepa) PWA first transport cross section
408  fScrA = fGSTable->GetMoliereXc2(matindx)/(4.0*pt2*bc)*fMCtoScrA;
409  // elastic mean free path in Geant4 internal lenght units: the neglected (1+screening parameter) term is corrected
410  // (if Mott-corretion: the corrected screening parameter is used for this (1+A) correction + Moliere b_c is also
411  // corrected with the screening parameter correction)
412  fLambda0 = beta2*(1.+fScrA)*fMCtoScrA/bc/scpCor;
413  // first transport coefficient (if Mott-corretion: the corrected screening parameter is used (it will be fully
414  // consistent with the one used during the pre-computation of the Mott-correted GS angular distributions))
415  fG1 = 2.0*fScrA*((1.0+fScrA)*G4Log(1.0/fScrA+1.0)-1.0);
416  // first transport mean free path
418 
419  return fLambda1;
420 }
421 
422 
423 G4double
426  // kinetic energy is assumed to be in Geant4 internal energy unit which is MeV
427  G4double efEnergy = kineticEnergy;
428  //
430  //
431  G4double lambda0 = 0.0; // elastc mean free path
432  G4double lambda1 = 0.0; // first transport mean free path
433  G4double scrA = 0.0; // screening parametr
434  G4double g1 = 0.0; // first transport mean free path
435 
436  // use Moliere's screening (with Mott-corretion if it was requested)
437  if (efEnergy<10.*CLHEP::eV) efEnergy = 10.*CLHEP::eV;
438  // total mometum square in Geant4 internal energy2 units which is MeV2
439  G4double pt2 = efEnergy*(efEnergy+2.0*electron_mass_c2);
440  G4double beta2 = pt2/(pt2+electron_mass_c2*electron_mass_c2);
441  G4int matindx = mat->GetIndex();
442  G4double bc = fGSTable->GetMoliereBc(matindx);
443  // get the Mott-correcton factors if Mott-correcton was requested by the user
444  G4double mctoScrA = 1.0;
445  G4double mctoQ1 = 1.0;
446  G4double mctoG2PerG1 = 1.0;
447  G4double scpCor = 1.0;
448  if (fIsUseMottCorrection) {
449  fGSTable->GetMottCorrectionFactors(G4Log(efEnergy), beta2, matindx, mctoScrA, mctoQ1, mctoG2PerG1);
451  } else if (fIsUsePWACorrection) {
452  fPWACorrection->GetPWACorrectionFactors(G4Log(efEnergy), beta2, matindx, mctoScrA, mctoQ1, mctoG2PerG1);
453  // scpCor = fGSTable->ComputeScatteringPowerCorrection(currentCouple, efEnergy);
454  }
455  scrA = fGSTable->GetMoliereXc2(matindx)/(4.0*pt2*bc)*mctoScrA;
456  // total elastic mean free path in Geant4 internal lenght units
457  lambda0 = beta2*(1.+scrA)*mctoScrA/bc/scpCor;
458  g1 = 2.0*scrA*((1.0+scrA)*G4Log(1.0/scrA+1.0)-1.0);
459  lambda1 = lambda0/g1;
460 
461  return lambda1;
462 }
463 
464 
467  firstStep = true;
469  rangeinit = 1.e+21;
470 }
471 
472 
474  G4double& currentMinimalStep) {
475  G4double skindepth = 0.;
476  //
477  const G4DynamicParticle* dp = track.GetDynamicParticle();
478  G4StepPoint* sp = track.GetStep()->GetPreStepPoint();
479  G4StepStatus stepStatus = sp->GetStepStatus();
485  // elastic and first transport mfp, screening parameter and G1 are also set
486  // (Mott-correction will be used if it was requested by the user)
488  // Set initial values:
489  // : lengths are initialised to currentMinimalStep which is the true, minimum
490  // step length from all other physics
491  fTheTrueStepLenght = currentMinimalStep;
492  fTheTransportDistance = currentMinimalStep;
493  fTheZPathLenght = currentMinimalStep; // will need to be converted
494  fTheDisplacementVector.set(0.,0.,0.);
495  fTheNewDirection.set(0.,0.,1.);
496 
497  // Can everything be done in the step limit phase ?
498  fIsEverythingWasDone = false;
499  // Multiple scattering needs to be sample ?
500  fIsMultipleSacettring = false;
501  // Single scattering needs to be sample ?
502  fIsSingleScattering = false;
503  // Was zero deflection in multiple scattering sampling ?
504  fIsNoScatteringInMSC = false;
505  // Do not care about displacement in MSC sampling
506  // ( used only in the case of gIsOptimizationOn = true)
507  fIsNoDisplace = false;
508  // get pre-step point safety
509  presafety = sp->GetSafety();
510  //
512  // distance will take into account max-fluct.
513  G4double distance = currentRange;
514  distance *= (1.20-fZeff*(1.62e-2-9.22e-5*fZeff));
515  //
516  // Possible optimization : if the distance is samller than the safety -> the
517  // particle will never leave this volume -> dispalcement
518  // as the effect of multiple elastic scattering can be skipped
519  // Important : this optimization can cause problems if one does scoring
520  // in a bigger volume since MSC won't be done deep inside the volume when
521  // distance < safety so don't use optimized-mode in such case.
522  if (gIsOptimizationOn && (distance<presafety)) {
523  // Indicate that we need to do MSC after transportation and no dispalcement.
524  fIsMultipleSacettring = true;
525  fIsNoDisplace = true;
527  //Compute geomlimit (and presafety) :
528  // - geomlimit will be:
529  // == the straight line distance to the boundary if currentRange is
530  // longer than that
531  // == a big value [geombig = 1.e50*mm] if currentRange is shorter than
532  // the straight line distance to the boundary
533  // - presafety will be updated as well
534  // So the particle can travell 'gemlimit' distance (along a straight
535  // line!) in its current direction:
536  // (1) before reaching a boundary (geomlimit < geombig) OR
537  // (2) before reaching its current range (geomlimit == geombig)
539  // Record that the particle is on a boundary
540  if ( (stepStatus==fGeomBoundary) || (stepStatus==fUndefined && presafety==0.0)) {
541  fIsWasOnBoundary = true;
542  }
543  // Set skin depth = skin x elastic_mean_free_path
544  skindepth = skin*fLambda0;
545  // Init the flag that indicates that the particle are within a skindepth
546  // distance from a boundary
547  fIsInsideSkin = false;
548  // Check if we can try Single Scattering because we are within skindepth
549  // distance from/to a boundary OR the current minimum true-step-length is
550  // shorter than skindepth. NOTICE: the latest has only efficieny reasons
551  // because the MSC angular sampling is fine for any short steps but much
552  // faster to try single scattering in case of short steps.
553  if ((stepStatus==fGeomBoundary) || (presafety<skindepth) || (fTheTrueStepLenght<skindepth)) {
554  // check if we are within skindepth distance from a boundary
555  if ((stepStatus == fGeomBoundary) || (presafety < skindepth)) {
556  fIsInsideSkin = true;
557  fIsWasOnBoundary = true;
558  }
559  //Try single scattering:
560  // - sample distance to next single scattering interaction (sslimit)
561  // - compare to current minimum length
562  // == if sslimit is the shorter:
563  // - set the step length to sslimit
564  // - indicate that single scattering needs to be done
565  // == else : nothing to do
566  //- in both cases, the step length was very short so geometrical and
567  // true path length are the same
568  G4double sslimit = -1.*fLambda0*G4Log(G4UniformRand());
569  // compare to current minimum step length
570  if (sslimit<fTheTrueStepLenght) {
571  fTheTrueStepLenght = sslimit;
572  fIsSingleScattering = true;
573  }
574  // short step -> true step length equal to geometrical path length
576  // Set taht everything is done in step-limit phase so no MSC call
577  // We will check if we need to perform the single-scattering angular
578  // sampling i.e. if single elastic scattering was the winer!
579  fIsEverythingWasDone = true;
580  } else {
581  // After checking we know that we cannot try single scattering so we will
582  // need to make an MSC step
583  // Indicate that we need to make and MSC step. We do not check if we can
584  // do it now i.e. if presafety>final_true_step_length so we let the
585  // fIsEverythingWasDone = false which indicates that we will perform
586  // MSC after transportation.
587  fIsMultipleSacettring = true;
588  // Init the first-real-step falg: it will indicate if we do the first
589  // non-single scattering step in this volume with this particle
590  fIsFirstRealStep = false;
591  // If previously the partcile was on boundary it was within skin as
592  // well. When it is not within skin anymore it has just left the skin
593  // so we make the first real MSC step with the particle.
595  // reset the 'was on boundary' indicator flag
596  fIsWasOnBoundary = false;
597  fIsFirstRealStep = true;
598  }
599  // If this is the first-real msc step (the partcile has just left the
600  // skin) or this is the first step with the particle (was born or
601  // primary):
602  // - set the initial range that will be used later to limit its step
603  // (only in this volume, because after boundary crossing at the
604  // first-real MSC step we will reset)
605  // - don't let the partcile to cross the volume just in one step
606  if (firstStep || fIsFirstRealStep || rangeinit>1.e+20) {
608  // If geomlimit < geombig than the particle might reach the boundary
609  // along its initial direction before losing its energy (in this step)
610  // Otherwise we can be sure that the particle will lose it energy
611  // before reaching the boundary along a starigth line so there is no
612  // geometrical limit appalied. [However, tgeom is set only in the
613  // first or the first-real MSC step. After the first or first real
614  // MSC step the direction will change tgeom won't guaranty anything!
615  // But we will try to end up within skindepth from the boundary using
616  // the actual value of geomlimit(See later at step reduction close to
617  // boundary).]
618  if (geomlimit<geombig) {
619  // transfrom straight line distance to the boundary to real step
620  // length based on the mean values (using the prestep point
621  // first-transport mean free path i.e. no energy loss correction)
622  if ((1.-geomlimit/fLambda1)> 0.) {
624  }
625  // the 2-different case that could lead us here
626  if (firstStep) {
627  tgeom = 2.*geomlimit/facgeom;
628  } else {
630  }
631  } else {
632  tgeom = geombig;
633  }
634  }
635  // True step length limit from range factor. Noteice, that the initial
636  // range is used that was set at the first step or first-real MSC step
637  // in this volume with this particle.
639  // Take the minimum of the true step length limits coming from
640  // geometrical constraint or range-factor limitation
642  // Step reduction close to boundary: we try to end up within skindepth
643  // from the boundary ( Notice: in case of mag. field it might not work
644  // because geomlimit is the straigth line distance to the boundary in
645  // the currect direction (if geomlimit<geombig) and mag. field can
646  // change the initial direction. So te particle might hit some boundary
647  // before in a different direction. However, here we restrict the true
648  // path length to this (straight line) lenght so the corresponding
649  // transport distance (straight line) will be even shorter than
650  // geomlimit-0.999*skindepth after the change of true->geom.
651  if (geomlimit<geombig) {
652  tlimit = std::min(tlimit, geomlimit-0.999*skindepth);
653  }
654  // randomize 1st step or 1st 'normal' step in volume
655  if (firstStep || fIsFirstRealStep) {
657  } else {
659  }
660  }
661  } else if (steppingAlgorithm==fUseSafetyPlus) { // THE ERROR_FREE stepping alg.
664  // Set skin depth = skin x elastic_mean_free_path
665  skindepth = skin*fLambda0;
666  // Check if we can try Single Scattering because we are within skindepth
667  // distance from/to a boundary OR the current minimum true-step-length is
668  // shorter than skindepth. NOTICE: the latest has only efficieny reasons
669  // because the MSC angular sampling is fine for any short steps but much
670  // faster to try single scattering in case of short steps.
671  if ((stepStatus==fGeomBoundary) || (presafety<skindepth) || (fTheTrueStepLenght<skindepth)) {
672  //Try single scattering:
673  // - sample distance to next single scattering interaction (sslimit)
674  // - compare to current minimum length
675  // == if sslimit is the shorter:
676  // - set the step length to sslimit
677  // - indicate that single scattering needs to be done
678  // == else : nothing to do
679  //- in both cases, the step length was very short so geometrical and
680  // true path length are the same
681  G4double sslimit = -1.*fLambda0*G4Log(G4UniformRand());
682  // compare to current minimum step length
683  if (sslimit<fTheTrueStepLenght) {
684  fTheTrueStepLenght = sslimit;
685  fIsSingleScattering = true;
686  }
687  // short step -> true step length equal to geometrical path length
689  // Set taht everything is done in step-limit phase so no MSC call
690  // We will check if we need to perform the single-scattering angular
691  // sampling i.e. if single elastic scattering was the winer!
692  fIsEverythingWasDone = true;
693  } else {
694  // After checking we know that we cannot try single scattering so we will
695  // need to make an MSC step
696  // Indicate that we need to make and MSC step.
697  fIsMultipleSacettring = true;
698  fIsEverythingWasDone = true;
699  // limit from range factor
701  // never let the particle go further than the safety if we are out of the skin
702  // if we are here we are out of the skin, presafety > 0.
705  }
706  // make sure that we are still within the aplicability of condensed histry model
707  // i.e. true step length is not longer than first transport mean free path.
708  // We schould take into account energy loss along 0.5x lambda_transport1
709  // step length as well. So let it 0.5 x lambda_transport1
711  }
712  } else {
713  // This is the default stepping algorithm: the fastest but the least
714  // accurate that corresponds to fUseSafety in Urban model. Note, that GS
715  // model can handle any short steps so we do not need the minimum limits
716  //
717  // NO single scattering in case of skin or short steps (by defult the MSC
718  // model will be single or even no scattering in case of short steps
719  // compared to the elastic mean free path.)
720  //
721  // indicate that MSC needs to be done (always and always after transportation)
722  fIsMultipleSacettring = true;
723  if (stepStatus!=fGeomBoundary) {
725  }
726  // Far from boundary-> in optimized mode do not sample dispalcement.
727  if ((distance<presafety) && (gIsOptimizationOn)) {
728  fIsNoDisplace = true;
729  } else {
730  // Urban like
731  if (firstStep || (stepStatus==fGeomBoundary) || rangeinit>1.e+20) {
733  fr = facrange;
734 // We don't use this: we won't converge to the single scattering results with
735 // decreasing range-factor.
736 // rangeinit = std::max(rangeinit, fLambda1);
737 // if(fLambda1 > lambdalimit) {
738 // fr *= (0.75+0.25*fLambda1/lambdalimit);
739 // }
740 
741  }
742  //step limit
744  // first step randomization
745  if (firstStep || stepStatus==fGeomBoundary) {
747  } else {
749  }
750  }
751  }
752  //
753  // unset first-step
754  firstStep =false;
755  // performe single scattering, multiple scattering if this later can be done safely here
756  if (fIsEverythingWasDone) {
757  if (fIsSingleScattering) {
758  // sample single scattering
759  //G4double ekin = 0.5*(currentKinEnergy + GetEnergy(particle,currentRange-fTheTrueStepLenght,currentCouple));
763  G4double cost = fGSTable->SingleScattering(1., fScrA, lekin, beta2, currentMaterialIndex);
764  // protection
765  if (cost<-1.) cost = -1.;
766  if (cost> 1.) cost = 1.;
767  // compute sint
768  G4double dum = 1.-cost;
769  G4double sint = std::sqrt(dum*(2.-dum));
771  G4double sinPhi = std::sin(phi);
772  G4double cosPhi = std::cos(phi);
773  fTheNewDirection.set(sint*cosPhi,sint*sinPhi,cost);
774  } else if (fIsMultipleSacettring) {
775  // sample multiple scattering
776  SampleMSC(); // fTheZPathLenght, fTheDisplacementVector and fTheNewDirection will be set
777  } // and if single scattering but it was longer => nothing to do
778  } //else { do nothing here but after transportation
779  //
780  return ConvertTrueToGeom(fTheTrueStepLenght,currentMinimalStep);
781 }
782 
783 
785  // convert true ->geom
786  // It is called from the step limitation ComputeTruePathLengthLimit if
787  // !fIsEverythingWasDone but protect:
788  par1 = -1.;
789  par2 = par3 = 0.;
790  // if fIsEverythingWasDone = TRUE => fTheZPathLenght is already set
791  // so return with the already known value
792  // Otherwise:
793  if (!fIsEverythingWasDone) {
794  // this correction needed to run MSC with eIoni and eBrem inactivated
795  // and makes no harm for a normal run
797  // do the true -> geom transformation
799  // z = t for very small true-path-length
801  return fTheZPathLenght;
802  }
804  if (tau<=tausmall) {
806  } else if (fTheTrueStepLenght<currentRange*dtrl) {
807  if (tau<taulim) fTheZPathLenght = fTheTrueStepLenght*(1.-0.5*tau) ;
808  else fTheZPathLenght = fLambda1*(1.-G4Exp(-tau));
810  par1 = 1./currentRange ; // alpha =1/range_init for Ekin<mass
811  par2 = 1./(par1*fLambda1) ; // 1/(alphaxlambda01)
812  par3 = 1.+par2 ; // 1+1/
814  fTheZPathLenght = 1./(par1*par3) * (1.-std::pow(1.-par1*fTheTrueStepLenght,par3));
815  } else {
816  fTheZPathLenght = 1./(par1*par3);
817  }
818  } else {
822  //
823  par1 = (fLambda1-lambda1)/(fLambda1*fTheTrueStepLenght); // alpha
824  par2 = 1./(par1*fLambda1);
825  par3 = 1.+par2 ;
826  G4Pow *g4calc = G4Pow::GetInstance();
827  fTheZPathLenght = 1./(par1*par3) * (1.-g4calc->powA(1.-par1*fTheTrueStepLenght,par3));
828  }
829  }
831  //
832  return fTheZPathLenght;
833 }
834 
835 
837  // init
838  fIsEndedUpOnBoundary = false;
839  // step defined other than transportation
840  if (geomStepLength==fTheZPathLenght) {
841  return fTheTrueStepLenght;
842  }
843  // else ::
844  // - set the flag that transportation was the winer so DoNothin in DOIT !!
845  // - convert geom -> true by using the mean value
846  fIsEndedUpOnBoundary = true; // OR LAST STEP
847  fTheZPathLenght = geomStepLength;
848  // was a short single scattering step
850  fTheTrueStepLenght = geomStepLength;
851  return fTheTrueStepLenght;
852  }
853  // t = z for very small step
854  if (geomStepLength<tlimitminfix2) {
855  fTheTrueStepLenght = geomStepLength;
856  // recalculation
857  } else {
858  G4double tlength = geomStepLength;
859  if (geomStepLength>fLambda1*tausmall) {
860  if (par1< 0.) {
861  tlength = -fLambda1*G4Log(1.-geomStepLength/fLambda1) ;
862  } else {
863  if (par1*par3*geomStepLength<1.) {
864  G4Pow *g4calc = G4Pow::GetInstance();
865  tlength = (1.-g4calc->powA( 1.-par1*par3*geomStepLength,1./par3))/par1;
866  } else {
867  tlength = currentRange;
868  }
869  }
870  if (tlength<geomStepLength || tlength>fTheTrueStepLenght) {
871  tlength = geomStepLength;
872  }
873  }
874  fTheTrueStepLenght = tlength;
875  }
876  //
877  return fTheTrueStepLenght;
878 }
879 
883  // single scattering was and scattering happend
884  fTheNewDirection.rotateUz(oldDirection);
886  return fTheDisplacementVector;
887  } else if (steppingAlgorithm==fUseSafetyPlus) { // error-free stepping
888  if (fIsEndedUpOnBoundary) { // do nothing on the boundary
889  return fTheDisplacementVector;
890  } else if (fIsEverythingWasDone) { // evrything is done if not optimizations case !!!
891  // check single scattering and see if it happened
892  if (fIsSingleScattering) {
893  fTheNewDirection.rotateUz(oldDirection);
895  return fTheDisplacementVector;
896  }
897  // check if multiple scattering happened and do things only if scattering was really happening
899  fTheNewDirection.rotateUz(oldDirection);
900  fTheDisplacementVector.rotateUz(oldDirection);
902  }
903  // The only thing that could happen if we are here (fUseSafety and fIsEverythingWasDone)
904  // is that single scattering was tried but did not win so scattering did not happen.
905  // So no displacement and no scattering
906  return fTheDisplacementVector;
907  }
908  //
909  // The only thing that could still happen with fUseSafetyPlus is that we are in the
910  // optimization branch: so sample MSC angle here (no displacement)
911  }
912  //else MSC needs to be done here
913  SampleMSC();
914  if (!fIsNoScatteringInMSC) {
915  fTheNewDirection.rotateUz(oldDirection);
917  if (!fIsNoDisplace) {
918  fTheDisplacementVector.rotateUz(oldDirection);
919  }
920  }
921  //
922  return fTheDisplacementVector;
923 }
924 
925 
927  fIsNoScatteringInMSC = false;
928  // kinetic energy is assumed to be in Geant4 internal energy unit which is MeV
930  //
931  // Energy loss correction: 2 version
932  G4double eloss = 0.0;
933 // if (fTheTrueStepLenght > currentRange*dtrl) {
935 // } else {
936 // eloss = fTheTrueStepLenght*GetDEDX(particle,kineticEnergy,currentCouple);
937 // }
938 
939  G4double tau = 0.;// = kineticEnergy/electron_mass_c2; // where kinEnergy is the mean kinetic energy
940  G4double tau2 = 0.;// = tau*tau;
941  G4double eps0 = 0.;// = eloss/kineticEnergy0; // energy loss fraction to the begin step energy
942  G4double epsm = 0.;// = eloss/kineticEnergy; // energy loss fraction to the mean step energy
943 
944  // - init.
945  G4double efEnergy = kineticEnergy;
946  G4double efStep = fTheTrueStepLenght;
947 
948  G4double kineticEnergy0 = kineticEnergy;
949  if (gIsUseAccurate) { // - use accurate energy loss correction
950  kineticEnergy -= 0.5*eloss; // mean energy along the full step
951  // other parameters for energy loss corrections
952  tau = kineticEnergy/electron_mass_c2; // where kinEnergy is the mean kinetic energy
953  tau2 = tau*tau;
954  eps0 = eloss/kineticEnergy0; // energy loss fraction to the begin step energy
955  epsm = eloss/kineticEnergy; // energy loss fraction to the mean step energy
956 
957  efEnergy = kineticEnergy * (1.-epsm*epsm*(6.+10.*tau+5.*tau2)/(24.*tau2+48.*tau+72.));
958  G4double dum = 0.166666*(4.+tau*(6.+tau*(7.+tau*(4.+tau))))*(epsm/((tau+1.)*(tau+2.)))*(epsm/((tau+1.)*(tau+2.)));
959  efStep = fTheTrueStepLenght*(1.-dum);
960  } else { // - take only mean energy
961  kineticEnergy -= 0.5*eloss; // mean energy along the full step
962  efEnergy = kineticEnergy;
963  G4double factor = 1./(1.+0.9784671*kineticEnergy); //0.9784671 = 1/(2*m_e)
964  eps0 = eloss/kineticEnergy0;
965  epsm = eps0/(1.-0.5*eps0);
966  G4double temp = 0.3*(1 -factor*(1.-0.333333*factor))*eps0*eps0;
967  efStep = fTheTrueStepLenght*(1.+temp);
968  }
969  //
970  // compute elastic mfp, first transport mfp, screening parameter, and G1 (with Mott-correction
971  // if it was requested by the user)
973  // s/lambda_el
974  G4double lambdan=0.;
975  if (fLambda0>0.0) {
976  lambdan=efStep/fLambda0;
977  }
978  if (lambdan<=1.0e-12) {
979  if (fIsEverythingWasDone) {
981  }
982  fIsNoScatteringInMSC = true;
983  return;
984  }
985  // first moment: 2.* lambdan *scrA*((1.+scrA)*log(1.+1./scrA)-1.);
986  G4double Qn1 = lambdan *fG1;
987  // sample scattering angles
988  // new direction, relative to the orriginal one is in {uss,vss,wss}
989  G4double cosTheta1 = 1.0, sinTheta1 = 0.0, cosTheta2 = 1.0, sinTheta2 = 0.0;
990  G4double cosPhi1 = 1.0, sinPhi1 = 0.0, cosPhi2 = 1.0, sinPhi2 = 0.0;
991  G4double uss = 0.0, vss = 0.0, wss = 1.0;
992  G4double x_coord = 0.0, y_coord = 0.0, z_coord = 1.0;
993  G4double u2 = 0.0, v2 = 0.0;
994  // if we are above the upper grid limit with lambdaxG1=true-length/first-trans-mfp
995  // => izotropic distribution: lambG1_max =7.992 but set it to 7
996  if (0.5*Qn1 > 7.0){
997  cosTheta1 = 1.-2.*G4UniformRand();
998  sinTheta1 = std::sqrt((1.-cosTheta1)*(1.+cosTheta1));
999  cosTheta2 = 1.-2.*G4UniformRand();
1000  sinTheta2 = std::sqrt((1.-cosTheta2)*(1.+cosTheta2));
1001  } else {
1002  // sample 2 scattering cost1, sint1, cost2 and sint2 for half path
1003  G4double lekin = G4Log(efEnergy);
1004  G4double pt2 = efEnergy*(efEnergy+2.0*CLHEP::electron_mass_c2);
1006  // backup GS angular dtr pointer (kinetic energy and delta index in case of Mott-correction)
1007  // if the first was an msc sampling (the same will be used if the second is also an msc step)
1009  G4int mcEkinIdx = -1;
1010  G4int mcDeltIdx = -1;
1011  G4double transfPar = 0.;
1012  G4bool isMsc = fGSTable->Sampling(0.5*lambdan, 0.5*Qn1, fScrA, cosTheta1, sinTheta1, lekin, beta2,
1013  currentMaterialIndex, &gsDtr, mcEkinIdx, mcDeltIdx, transfPar,
1014  true);
1015  fGSTable->Sampling(0.5*lambdan, 0.5*Qn1, fScrA, cosTheta2, sinTheta2, lekin, beta2,
1016  currentMaterialIndex, &gsDtr, mcEkinIdx, mcDeltIdx, transfPar, !isMsc);
1017  if (cosTheta1+cosTheta2==2.) { // no scattering happened
1020  fIsNoScatteringInMSC = true;
1021  return;
1022  }
1023  }
1024  // sample 2 azimuthal angles
1026  sinPhi1 = std::sin(phi1);
1027  cosPhi1 = std::cos(phi1);
1029  sinPhi2 = std::sin(phi2);
1030  cosPhi2 = std::cos(phi2);
1031 
1032  // compute final direction realtive to z-dir
1033  u2 = sinTheta2*cosPhi2;
1034  v2 = sinTheta2*sinPhi2;
1035  G4double u2p = cosTheta1*u2 + sinTheta1*cosTheta2;
1036  uss = u2p*cosPhi1 - v2*sinPhi1;
1037  vss = u2p*sinPhi1 + v2*cosPhi1;
1038  wss = cosTheta1*cosTheta2 - sinTheta1*u2;
1039 
1040  // set new direction (is scattering frame)
1041  fTheNewDirection.set(uss,vss,wss);
1042 
1043  // set the fTheZPathLenght if we don't sample displacement and
1044  // we should do everything at the step-limit-phase before we return
1047 
1048  // in optimized-mode if the current-safety > current-range we do not use dispalcement
1049  if(fIsNoDisplace)
1050  return;
1051 
1053  // Compute final position
1054  Qn1 *= fMCtoQ1;
1055  if (gIsUseAccurate) {
1056  // correction parameter
1057  G4double par =1.;
1058  if(Qn1<0.7) par = 1.;
1059  else if (Qn1<7.0) par = -0.031376*Qn1+1.01356;
1060  else par = 0.79;
1061 
1062  // Moments with energy loss correction
1063  // --first the uncorrected (for energy loss) values of gamma, eta, a1=a2=0.5*(1-eta), delta
1064  // gamma = G_2/G_1 based on G2 computed from A by using the Wentzel DCS form of G2
1065  G4double loga = G4Log(1.0+1.0/fScrA);
1066  G4double gamma = 6.0*fScrA*(1.0 + fScrA)*(loga*(1.0 + 2.0*fScrA) - 2.0)/fG1;
1067  gamma *= fMCtoG2PerG1;
1068  // sample eta from p(eta)=2*eta i.e. P(eta) = eta_square ;-> P(eta) = rand --> eta = sqrt(rand)
1069  G4double eta = std::sqrt(G4UniformRand());
1070  G4double eta1 = 0.5*(1 - eta); // used more than once
1071  // 0.5 +sqrt(6)/6 = 0.9082483;
1072  // 1/(4*sqrt(6)) = 0.1020621;
1073  // (4-sqrt(6)/(24*sqrt(6))) = 0.026374715
1074  // delta = 0.9082483-(0.1020621-0.0263747*gamma)*Qn1 without energy loss cor.
1075  G4double delta = 0.9082483-(0.1020621-0.0263747*gamma)*Qn1;
1076 
1077  // compute alpha1 and alpha2 for energy loss correction
1078  G4double temp1 = 2.0 + tau;
1079  G4double temp = (2.0+tau*temp1)/((tau+1.0)*temp1);
1080  //Take logarithmic dependence
1081  temp = temp - (tau+1.0)/((tau+2.0)*(loga*(1.0+fScrA)-1.0));
1082  temp = temp * epsm;
1083  temp1 = 1.0 - temp;
1084  delta = delta + 0.40824829*(eps0*(tau+1.0)/((tau+2.0)*
1085  (loga*(1.0+fScrA)-1.0)*(loga*(1.0+2.0*fScrA)-2.0)) - 0.25*temp*temp);
1086  G4double b = eta*delta;
1087  G4double c = eta*(1.0-delta);
1088 
1089  //Calculate transport direction cosines:
1090  // ut,vt,wt is the final position divided by the true step length
1091  G4double w1v2 = cosTheta1*v2;
1092  G4double ut = b*sinTheta1*cosPhi1 + c*(cosPhi1*u2 - sinPhi1*w1v2) + eta1*uss*temp1;
1093  G4double vt = b*sinTheta1*sinPhi1 + c*(sinPhi1*u2 + cosPhi1*w1v2) + eta1*vss*temp1;
1094  G4double wt = eta1*(1+temp) + b*cosTheta1 + c*cosTheta2 + eta1*wss*temp1;
1095 
1096  // long step correction
1097  ut *=par;
1098  vt *=par;
1099  wt *=par;
1100 
1101  // final position relative to the pre-step point in the scattering frame
1102  // ut = x_f/s so needs to multiply by s
1103  x_coord = ut*fTheTrueStepLenght;
1104  y_coord = vt*fTheTrueStepLenght;
1105  z_coord = wt*fTheTrueStepLenght;
1106 
1108  // We sample in the step limit so set fTheZPathLenght = transportDistance
1109  // and lateral displacement (x_coord,y_coord,z_coord-transportDistance)
1110  //Calculate transport distance
1111  G4double transportDistance = std::sqrt(x_coord*x_coord+y_coord*y_coord+z_coord*z_coord);
1112  // protection
1113  if(transportDistance>fTheTrueStepLenght)
1114  transportDistance = fTheTrueStepLenght;
1115  fTheZPathLenght = transportDistance;
1116  }
1117  // else:: we sample in the DoIt so
1118  // the fTheZPathLenght was already set and was taken as transport along zet
1119  fTheDisplacementVector.set(x_coord,y_coord,z_coord-fTheZPathLenght);
1120  } else {
1121  // compute zz = <z>/tPathLength
1122  // s -> true-path-length
1123  // z -> geom-path-length:: when PRESTA is used z =(def.) <z>
1124  // r -> lateral displacement = s/2 sin(theta) => x_f = r cos(phi); y_f = r sin(phi)
1125  G4double zz = 0.0;
1127  // We sample in the step limit so set fTheZPathLenght = transportDistance
1128  // and lateral displacement (x_coord,y_coord,z_coord-transportDistance)
1129  if(Qn1<0.1) { // use 3-order Taylor approximation of (1-exp(-x))/x around x=0
1130  zz = 1.0 - Qn1*(0.5 - Qn1*(0.166666667 - 0.041666667*Qn1)); // 1/6 =0.166..7 ; 1/24=0.041..
1131  } else {
1132  zz = (1.-G4Exp(-Qn1))/Qn1;
1133  }
1134  } else {
1135  // we sample in the DoIt so
1136  // the fTheZPathLenght was already set and was taken as transport along zet
1138  }
1139 
1140  G4double rr = (1.-zz*zz)/(1.-wss*wss); // s^2 >= <z>^2+r^2 :: where r^2 = s^2/4 sin^2(theta)
1141  if(rr >= 0.25) rr = 0.25; // (1-<z>^2/s^2)/sin^2(theta) >= r^2/(s^2 sin^2(theta)) = 1/4 must hold
1142  G4double rperp = fTheTrueStepLenght*std::sqrt(rr); // this is r/sint
1143  x_coord = rperp*uss;
1144  y_coord = rperp*vss;
1145  z_coord = zz*fTheTrueStepLenght;
1146 
1148  G4double transportDistance = std::sqrt(x_coord*x_coord + y_coord*y_coord + z_coord*z_coord);
1149  fTheZPathLenght = transportDistance;
1150  }
1151 
1152  fTheDisplacementVector.set(x_coord,y_coord,z_coord- fTheZPathLenght);
1153  }
1154 }
void set(double x, double y, double z)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
T max(const T t1, const T t2)
brief Return the largest of the two arguments
const G4ParticleDefinition * particle
void GetMottCorrectionFactors(G4double logekin, G4double beta2, G4int matindx, G4double &mcToScr, G4double &mcToQ1, G4double &mcToG2PerG1)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4double facsafety
Definition: G4VMscModel.hh:177
G4double GetTransportMeanFreePath(const G4ParticleDefinition *, G4double)
G4double ConvertTrueToGeom(G4double &tLength, G4double &gLength)
Definition: G4VMscModel.hh:245
G4StepPoint * GetPreStepPoint() const
G4double dtrl
Definition: G4VMscModel.hh:179
static constexpr double nm
Definition: G4SIunits.hh:112
static constexpr double mm
Definition: G4SIunits.hh:115
G4double SingleScattering(G4double lambdaval, G4double scra, G4double lekin, G4double beta2, G4int matindx)
G4double ComputeGeomLimit(const G4Track &, G4double &presafety, G4double limit)
Definition: G4VMscModel.hh:255
G4double ComputeSafety(const G4ThreeVector &position, G4double limit=DBL_MAX)
Definition: G4VMscModel.hh:237
virtual void InitialiseLocal(const G4ParticleDefinition *p, G4VEmModel *masterModel)
size_t GetIndex() const
Definition: G4Material.hh:262
const char * p
Definition: xmltok.h:285
G4bool Sampling(G4double lambdaval, G4double qval, G4double scra, G4double &cost, G4double &sint, G4double lekin, G4double beta2, G4int matindx, GSMSCAngularDtr **gsDtr, G4int &mcekini, G4int &mcdelti, G4double &transfPar, G4bool isfirst)
G4double facgeom
Definition: G4VMscModel.hh:176
G4double GetPDGCharge() const
virtual G4double ComputeTrueStepLength(G4double geomStepLength)
Double_t zz
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &currentMinimalStep)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
const G4MaterialCutsCouple * currentCouple
void Initialise(G4double lownergylimit, G4double highenergylimit)
G4double GetMoliereBc(G4int matindx)
void GetPWACorrectionFactors(G4double logekin, G4double beta2, G4int matindx, G4double &corToScr, G4double &corToQ1, G4double &corToG2PerG1)
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:455
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:616
G4double facrange
Definition: G4VMscModel.hh:175
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4StepStatus GetStepStatus() const
static G4Pow * GetInstance()
Definition: G4Pow.cc:57
G4ParticleDefinition * GetDefinition() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:242
G4double GetMoliereXc2(G4int matindx)
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.hh:283
static constexpr double electron_mass_c2
const G4ThreeVector & GetPosition() const
#define G4UniformRand()
Definition: Randomize.hh:53
G4double GetTransportMeanFreePathOnly(const G4ParticleDefinition *, G4double)
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
G4StepStatus
Definition: G4StepStatus.hh:49
G4double skin
Definition: G4VMscModel.hh:178
static constexpr double eV
Float_t mat
const G4Step * GetStep() const
virtual G4double ComputeGeomPathLength(G4double truePathLength)
G4double ComputeScatteringPowerCorrection(const G4MaterialCutsCouple *matcut, G4double ekin)
static constexpr double twopi
Definition: SystemOfUnits.h:55
G4GoudsmitSaundersonTable * GetGSTable()
G4GoudsmitSaundersonMscModel(const G4String &nam="GoudsmitSaunderson")
int G4int
Definition: G4Types.hh:78
G4ParticleChangeForMSC * GetParticleChangeForMSC(const G4ParticleDefinition *p=nullptr)
Definition: G4VMscModel.cc:91
void SetParticle(const G4ParticleDefinition *p)
virtual G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety)
G4MscStepLimitType steppingAlgorithm
Definition: G4VMscModel.hh:185
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:227
G4double GetSafety() const
Definition: G4Pow.hh:56
G4double GetZeffective() const
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
G4double GetKineticEnergy() const
static G4LossTableManager * Instance()
G4bool IsMaster() const
Definition: G4VEmModel.hh:700
const G4Material * GetMaterial() const
G4double GetEnergy(const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.hh:303
G4bool isElectron(G4int ityp)
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:609
const G4DynamicParticle * GetDynamicParticle() const
static G4EmParameters * Instance()
T min(const T t1, const T t2)
brief Return the smallest of the two arguments