Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4GoudsmitSaundersonMscModel.hh
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: G4GoudsmitSaundersonMscModel.hh 94933 2015-12-18 09:22:52Z 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 //
77 // Class description:
78 // Kawrakow-Bielajew Goudsmit-Saunderson MSC model based on the screened
79 // Rutherford DCS for elastic scattering of electrons/positrons. Step limitation
80 // algorithm as well as true to geomerty and geometry to true step length
81 // computations are adopted from Urban model[5].
82 //
83 // References:
84 // [1] A.F.Bielajew, NIMB 111 (1996) 195-208
85 // [2] I.Kawrakow, A.F.Bielajew, NIMB 134(1998) 325-336
86 // [3] I.Kawrakow, E.Mainegra-Hing, D.W.O.Rogers, F.Tessier,B.R.B.Walters, NRCC
87 // Report PIRS-701 (2013)
88 // [4] F.Salvat, A.Jablonski, C.J. Powell, CPC 165(2005) 157-190
89 // [5] L.Urban, Preprint CERN-OPEN-2006-077 (2006)
90 //
91 // -----------------------------------------------------------------------------
92 
93 #ifndef G4GoudsmitSaundersonMscModel_h
94 #define G4GoudsmitSaundersonMscModel_h 1
95 
97 
98 #include "G4VMscModel.hh"
99 #include "G4PhysicsTable.hh"
100 #include "G4MaterialCutsCouple.hh"
101 #include "globals.hh"
102 
103 
104 class G4DataVector;
106 class G4LossTableManager;
108 class G4PWATotalXsecTable;
109 
110 
112 {
113 public:
114 
115  G4GoudsmitSaundersonMscModel(const G4String& nam = "GoudsmitSaunderson");
116 
118 
119  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
120 
121  void StartTracking(G4Track*);
122 
124  void SingleScattering(G4double &cost, G4double &sint);
125  void SampleMSC();
126 
128  G4double safety);
129 
130  virtual G4double ComputeTruePathLengthLimit(const G4Track& track,
131  G4double& currentMinimalStep);
132 
133  virtual G4double ComputeGeomPathLength(G4double truePathLength);
134 
135  virtual G4double ComputeTrueStepLength(G4double geomStepLength);
136 
137  void SetOptionPWAScreening(G4bool opt){fIsUsePWATotalXsecData=opt;}
138 
139 private:
140  inline void SetParticle(const G4ParticleDefinition* p);
141 
142  inline G4double GetLambda(G4double);
143 
144  // hide assignment operator
147  G4double GetTransportMeanFreePathOnly(const G4ParticleDefinition*,G4double);
148 
149  inline G4double Randomizetlimit();
150 
151 private:
152  CLHEP::HepRandomEngine* rndmEngineMod;
153 
154 
155  G4double lowKEnergy;
156  G4double highKEnergy;
157  G4double currentKinEnergy;
158  G4double currentRange;
159 
160  G4double fr,rangeinit,geombig,geomlimit;
161  G4double lambdalimit,tlimit,tgeom;
162  G4int charge,currentMaterialIndex;
163 
164  G4bool firstStep;
165 
166  G4double par1,par2,par3,tlimitminfix2,tausmall,mass,taulim;
167 
168  G4LossTableManager* theManager;
169  const G4ParticleDefinition* particle;
170  G4ParticleChangeForMSC* fParticleChange;
171  const G4MaterialCutsCouple* currentCouple;
172 
173  static G4GoudsmitSaundersonTable* fgGSTable;
174  static G4PWATotalXsecTable* fgPWAXsecTable;
175 
176  G4bool fIsUsePWATotalXsecData;
177 
178  G4double presafety;
179  G4double fZeff;
180 
181  //
182  G4double fLambda0; // elastic mean free path
183  G4double fLambda1; // first transport mean free path
184  G4double fScrA; // screening parameter
185  G4double fG1; // first transport coef.
186  //
187  G4double fTheTrueStepLenght;
188  G4double fTheTransportDistance;
189  G4double fTheZPathLenght;
190  G4ThreeVector fTheDisplacementVector;
191  G4ThreeVector fTheNewDirection;
192  //
193  G4bool fIsEndedUpOnBoundary; // step ended up on boundary i.e. transportation is the winer
194  G4bool fIsMultipleSacettring;
195  G4bool fIsSingleScattering;
196  G4bool fIsEverythingWasDone;
197  G4bool fIsNoScatteringInMSC;
198  G4bool fIsNoDisplace;
199  G4bool fIsInsideSkin;
200  G4bool fIsWasOnBoundary;
201  G4bool fIsFirstRealStep;
202  //
203  static G4bool fgIsUseAccurate;
204  static G4bool fgIsOptimizationOn;
205 };
206 
208 inline
209 void G4GoudsmitSaundersonMscModel::SetParticle(const G4ParticleDefinition* p)
210 {
211  if (p != particle) {
212  particle = p;
213  charge = (G4int)(p->GetPDGCharge()/CLHEP::eplus);
214  mass = p->GetPDGMass();
215  }
216 }
217 
218 
219 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
220 inline
221 G4double G4GoudsmitSaundersonMscModel::Randomizetlimit()
222 {
223  G4double temptlimit = tlimit;
224  do {
225  temptlimit = G4RandGauss::shoot(rndmEngineMod,tlimit,0.3*tlimit);
226  } while ( (temptlimit<0.) || (temptlimit > 2.*tlimit));
227 
228  return temptlimit;
229 }
230 
231 
232 
233 #endif
G4GoudsmitSaundersonMscModel(const G4String &nam="GoudsmitSaunderson")
virtual G4double ComputeTrueStepLength(G4double geomStepLength)
ThreeVector shoot(const G4int Ap, const G4int Af)
const char * p
Definition: xmltok.h:285
int G4int
Definition: G4Types.hh:78
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &currentMinimalStep)
void SingleScattering(G4double &cost, G4double &sint)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
bool G4bool
Definition: G4Types.hh:79
static constexpr double eplus
G4double GetPDGMass() const
G4double GetTransportMeanFreePath(const G4ParticleDefinition *, G4double)
virtual G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety)
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
virtual G4double ComputeGeomPathLength(G4double truePathLength)