Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4UrbanMscModel.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: $
27 // GEANT4 tag $Name: $
28 //
29 // -------------------------------------------------------------------
30 //
31 //
32 // GEANT4 Class header file
33 //
34 //
35 // File name: G4UrbanMscModel
36 //
37 // Author: Laszlo Urban
38 //
39 // Creation date: 19.02.2013
40 //
41 // Created from G4UrbanMscModel96
42 //
43 // New parametrization for theta0
44 // Correction for very small step length
45 //
46 // Class Description:
47 //
48 // Implementation of the model of multiple scattering based on
49 // H.W.Lewis Phys Rev 78 (1950) 526 and L.Urban model
50 
51 // -------------------------------------------------------------------
52 //
53 
54 #ifndef G4UrbanMscModel_h
55 #define G4UrbanMscModel_h 1
56 
57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
58 
60 
61 #include "G4VMscModel.hh"
62 #include "G4MscStepLimitType.hh"
63 #include "G4Log.hh"
64 #include "G4Exp.hh"
65 
67 class G4SafetyHelper;
68 class G4LossTableManager;
69 
70 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71 
73 {
74 
75 public:
76 
77  explicit G4UrbanMscModel(const G4String& nam = "UrbanMsc");
78 
79  virtual ~G4UrbanMscModel();
80 
81  virtual void Initialise(const G4ParticleDefinition*,
82  const G4DataVector&) override;
83 
84  virtual void StartTracking(G4Track*) override;
85 
86  virtual G4double
88  G4double KineticEnergy,
89  G4double AtomicNumber,
90  G4double AtomicWeight=0.,
91  G4double cut =0.,
92  G4double emax=DBL_MAX) override;
93 
95  G4double safety) override;
96 
97  virtual G4double
99  G4double& currentMinimalStep) override;
100 
101  virtual G4double ComputeGeomPathLength(G4double truePathLength) override;
102 
103  virtual G4double ComputeTrueStepLength(G4double geomStepLength) override;
104 
105  G4double ComputeTheta0(G4double truePathLength, G4double KineticEnergy);
106 
107  inline void SetNewDisplacementFlag(G4bool);
108 
109 private:
110 
111  G4double SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy);
112 
113  void SampleDisplacement(G4double sinTheta, G4double phi);
114 
115  void SampleDisplacementNew(G4double sinTheta, G4double phi);
116 
117  inline void SetParticle(const G4ParticleDefinition*);
118 
119  inline void UpdateCache();
120 
121  inline G4double Randomizetlimit();
122 
123  inline G4double SimpleScattering(G4double xmeanth, G4double x2meanth);
124 
125  // hide assignment operator
126  G4UrbanMscModel & operator=(const G4UrbanMscModel &right) = delete;
127  G4UrbanMscModel(const G4UrbanMscModel&) = delete;
128 
129  CLHEP::HepRandomEngine* rndmEngineMod;
130 
131  const G4ParticleDefinition* particle;
132  const G4ParticleDefinition* positron;
133  G4ParticleChangeForMSC* fParticleChange;
134 
135  const G4MaterialCutsCouple* couple;
136  G4LossTableManager* theManager;
137 
138  G4double mass;
139  G4double charge,ChargeSquare;
140  G4double masslimite,lambdalimit,fr;
141 
142  G4double taubig;
143  G4double tausmall;
144  G4double taulim;
145  G4double currentTau;
146  G4double tlimit;
147  G4double tlimitmin;
148  G4double tlimitminfix,tlimitminfix2;
149  G4double tgeom;
150 
151  G4double geombig;
152  G4double geommin;
153  G4double geomlimit;
154  G4double skindepth;
155  G4double smallstep;
156 
157  G4double presafety;
158 
159  G4double lambda0;
160  G4double lambdaeff;
161  G4double tPathLength;
162  G4double zPathLength;
163  G4double par1,par2,par3;
164 
165  G4double stepmin;
166 
167  G4double currentKinEnergy;
168  G4double currentRange;
169  G4double rangeinit;
170  G4double currentRadLength;
171 
172  G4int currentMaterialIndex;
173 
174  G4double Zold;
175  G4double Zeff,Z2,Z23,lnZ;
176  G4double coeffth1,coeffth2;
177  G4double coeffc1,coeffc2,coeffc3,coeffc4;
178 
179  G4bool firstStep;
180  G4bool insideskin;
181 
182  G4bool latDisplasmentbackup ;
183  G4bool displacementFlag;
184 
185  G4double rangecut;
186  G4double drr,finalr;
187 
188 };
189 
190 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
191 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
192 
193 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
194 
196 {
197  displacementFlag = val;
198 }
199 
200 inline
201 void G4UrbanMscModel::SetParticle(const G4ParticleDefinition* p)
202 {
203  if (p != particle) {
204  particle = p;
205  mass = p->GetPDGMass();
206  charge = p->GetPDGCharge()/CLHEP::eplus;
207  ChargeSquare = charge*charge;
208  }
209 }
210 
211 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
212 
213 inline G4double G4UrbanMscModel::Randomizetlimit()
214 {
215  G4double temptlimit = tlimit;
216  if(tlimit > tlimitmin)
217  {
218  G4double delta = tlimit-tlimitmin;
219  do {
220  temptlimit = G4RandGauss::shoot(rndmEngineMod,tlimit,0.1*delta);
221  // Loop checking, 10-Apr-2016, Laszlo Urban
222  } while ((temptlimit < tlimit-delta) ||
223  (temptlimit > tlimit+delta));
224  }
225  else { temptlimit = tlimitmin; }
226 
227  return temptlimit;
228 }
229 
230 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
231 
232 inline void G4UrbanMscModel::UpdateCache()
233 {
234  lnZ = G4Log(Zeff);
235  // correction in theta0 formula
236  G4double w = G4Exp(lnZ/6.);
237  G4double facz = 0.990395+w*(-0.168386+w*0.093286) ;
238  coeffth1 = facz*(1. - 8.7780e-2/Zeff);
239  coeffth2 = facz*(4.0780e-2 + 1.7315e-4*Zeff);
240 
241  // tail parameters
242  G4double Z13 = w*w;
243  coeffc1 = 2.3785 - Z13*(4.1981e-1 - Z13*6.3100e-2);
244  coeffc2 = 4.7526e-1 + Z13*(1.7694 - Z13*3.3885e-1);
245  coeffc3 = 2.3683e-1 - Z13*(1.8111 - Z13*3.2774e-1);
246  coeffc4 = 1.7888e-2 + Z13*(1.9659e-2 - Z13*2.6664e-3);
247 
248  Z2 = Zeff*Zeff;
249  Z23 = Z13*Z13;
250 
251  Zold = Zeff;
252 }
253 
254 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
255 
256 inline
257 G4double G4UrbanMscModel::SimpleScattering(G4double xmeanth, G4double x2meanth)
258 {
259  // 'large angle scattering'
260  // 2 model functions with correct xmean and x2mean
261  G4double a = (2.*xmeanth+9.*x2meanth-3.)/(2.*xmeanth-3.*x2meanth+1.);
262  G4double prob = (a+2.)*xmeanth/a;
263 
264  // sampling
265  G4double cth = 1.;
266  if(rndmEngineMod->flat() < prob) {
267  cth = -1.+2.*G4Exp(G4Log(rndmEngineMod->flat())/(a+1.));
268  } else {
269  cth = -1.+2.*rndmEngineMod->flat();
270  }
271  return cth;
272 }
273 
274 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
275 
276 
277 #endif
278 
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
ThreeVector shoot(const G4int Ap, const G4int Af)
virtual G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety) override
G4UrbanMscModel(const G4String &nam="UrbanMsc")
virtual G4double ComputeGeomPathLength(G4double truePathLength) override
virtual ~G4UrbanMscModel()
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
void SetNewDisplacementFlag(G4bool)
const char * p
Definition: xmltok.h:285
virtual double flat()=0
int G4int
Definition: G4Types.hh:78
bool G4bool
Definition: G4Types.hh:79
static constexpr double eplus
static const G4double emax
G4double G4Log(G4double x)
Definition: G4Log.hh:230
virtual G4double ComputeTrueStepLength(G4double geomStepLength) override
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &currentMinimalStep) override
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *particle, G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight=0., G4double cut=0., G4double emax=DBL_MAX) override
G4double GetPDGMass() const
double G4double
Definition: G4Types.hh:76
G4double ComputeTheta0(G4double truePathLength, G4double KineticEnergy)
G4double GetPDGCharge() const
virtual void StartTracking(G4Track *) override
#define DBL_MAX
Definition: templates.hh:83