Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4GoudsmitSaundersonMscModel.hh
이 파일의 문서화 페이지로 가기
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.hh 108305 2018-02-02 13:08:43Z gcosmo $
27 //
28 // ----------------------------------------------------------------------------
29 //
30 // GEANT4 Class header 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 PLERIMINARYY version of updated class.
42 // All algorithms of the class were 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 fgIsUsePWATotalXsecData=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 fgIsUsePWATotalXsecData=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 // 02.02.2018 M. Novak: implemented CrossSectionPerVolume interface method (used only for testing)
102 //
103 // Class description:
104 // Kawrakow-Bielajew Goudsmit-Saunderson MSC model based on the screened Rutherford DCS
105 // for elastic scattering of e-/e+. Option, to include (Mott) correction (see above), is
106 // also available now (SetOptionMottCorrection(true)). An EGSnrc like error-free stepping
107 // algorithm (UseSafety) is available beyond the usual Geant4 step limitation algorithms
108 // and true to geomerty and geometry to true step length computations that were adopted
109 // from the Urban model[5]. The most accurate setting: error-free stepping (UseSafety)
110 // with Mott-correction (SetOptionMottCorrection(true)).
111 //
112 // References:
113 // [1] A.F.Bielajew, NIMB 111 (1996) 195-208
114 // [2] I.Kawrakow, A.F.Bielajew, NIMB 134(1998) 325-336
115 // [3] I.Kawrakow, E.Mainegra-Hing, D.W.O.Rogers, F.Tessier,B.R.B.Walters, NRCC
116 // Report PIRS-701 (2013)
117 // [4] F.Salvat, A.Jablonski, C.J. Powell, CPC 165(2005) 157-190
118 // [5] L.Urban, Preprint CERN-OPEN-2006-077 (2006)
119 //
120 // -----------------------------------------------------------------------------
121 
122 #ifndef G4GoudsmitSaundersonMscModel_h
123 #define G4GoudsmitSaundersonMscModel_h 1
124 
126 
127 #include "G4VMscModel.hh"
128 #include "G4PhysicsTable.hh"
129 #include "G4MaterialCutsCouple.hh"
130 #include "globals.hh"
131 
132 
133 class G4DataVector;
135 class G4LossTableManager;
137 class G4GSPWACorrections;
138 
140 {
141 public:
142 
143  G4GoudsmitSaundersonMscModel(const G4String& nam = "GoudsmitSaunderson");
144 
146 
147  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
148 
149  virtual void InitialiseLocal(const G4ParticleDefinition* p, G4VEmModel* masterModel);
150 
151 
152  virtual G4ThreeVector& SampleScattering(const G4ThreeVector&, G4double safety);
153 
154  virtual G4double ComputeTruePathLengthLimit(const G4Track& track, G4double& currentMinimalStep);
155 
156  virtual G4double ComputeGeomPathLength(G4double truePathLength);
157 
158  virtual G4double ComputeTrueStepLength(G4double geomStepLength);
159 
160  // method to compute first transport cross section per Volume (i.e. macroscropic first transport cross section; this
161  // method is used only for testing and not during a normal simulation)
163  G4double cutEnergy = 0.0, G4double maxEnergy = DBL_MAX);
164 
165  void StartTracking(G4Track*);
166 
167  void SampleMSC();
168 
170 
172 
174 
176 
178 
180 
182 
183 private:
184  inline void SetParticle(const G4ParticleDefinition* p);
185 
186  inline G4double GetLambda(G4double);
187 
188  // hide assignment operator
192 
193  inline G4double Randomizetlimit();
194 
195 private:
197  //
200  //
208  //
216  //
217  //
220  //
223  //
225  //
230 
233 
236  //
237  G4double fLambda0; // elastic mean free path
238  G4double fLambda1; // first transport mean free path
239  G4double fScrA; // screening parameter
240  G4double fG1; // first transport coef.
241  // in case of Mott-correction
245  //
249  //
252  //
253  G4bool fIsEndedUpOnBoundary; // step ended up on boundary i.e. transportation is the winer
262  //
265 };
266 
268 inline
270 {
271  if (p != particle) {
272  particle = p;
274  mass = p->GetPDGMass();
275  }
276 }
277 
278 
279 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
280 inline
282 {
283  G4double temptlimit = tlimit;
284  do {
285  temptlimit = G4RandGauss::shoot(rndmEngineMod,tlimit,0.1*tlimit);
286  } while ( (temptlimit<0.) || (temptlimit>2.*tlimit));
287 
288  return temptlimit;
289 }
290 
291 
292 
293 #endif
const G4ParticleDefinition * particle
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4double GetTransportMeanFreePath(const G4ParticleDefinition *, G4double)
virtual void InitialiseLocal(const G4ParticleDefinition *p, G4VEmModel *masterModel)
const char * p
Definition: xmltok.h:285
G4double GetPDGCharge() const
virtual G4double ComputeTrueStepLength(G4double geomStepLength)
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &currentMinimalStep)
const G4MaterialCutsCouple * currentCouple
G4double GetPDGMass() const
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4double GetTransportMeanFreePathOnly(const G4ParticleDefinition *, G4double)
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
ThreeVector shoot(const G4int Ap, const G4int Af)
virtual G4double ComputeGeomPathLength(G4double truePathLength)
G4GoudsmitSaundersonTable * GetGSTable()
G4GoudsmitSaundersonMscModel(const G4String &nam="GoudsmitSaunderson")
int G4int
Definition: G4Types.hh:78
void SetParticle(const G4ParticleDefinition *p)
virtual G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety)
G4double GetLambda(G4double)
G4GoudsmitSaundersonMscModel & operator=(const G4GoudsmitSaundersonMscModel &right)
#define DBL_MAX
Definition: templates.hh:83
static constexpr double eplus