Geant4  10.02.p01
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 
59 #include <CLHEP/Units/SystemOfUnits.h>
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 static const G4double c_highland = 13.6*CLHEP::MeV ;
71 
72 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
73 
75 {
76 
77 public:
78 
79  G4UrbanMscModel(const G4String& nam = "UrbanMsc");
80 
81  virtual ~G4UrbanMscModel();
82 
83  void Initialise(const G4ParticleDefinition*, const G4DataVector&);
84 
85  void StartTracking(G4Track*);
86 
88  G4double KineticEnergy,
89  G4double AtomicNumber,
90  G4double AtomicWeight=0.,
91  G4double cut =0.,
93 
95 
97  G4double& currentMinimalStep);
98 
99  G4double ComputeGeomPathLength(G4double truePathLength);
100 
101  G4double ComputeTrueStepLength(G4double geomStepLength);
102 
103  G4double ComputeTheta0(G4double truePathLength, G4double KineticEnergy);
104 
105  inline void SetNewDisplacementFlag(G4bool);
106 
107 private:
108 
109  G4double SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy);
110 
111  void SampleDisplacement(G4double sinTheta, G4double phi);
112 
113  void SampleDisplacementNew(G4double sinTheta, G4double phi);
114 
115  inline void SetParticle(const G4ParticleDefinition*);
116 
117  inline void UpdateCache();
118 
119  inline G4double Randomizetlimit();
120 
121  inline G4double SimpleScattering(G4double xmeanth, G4double x2meanth);
122 
123  // hide assignment operator
126 
127  CLHEP::HepRandomEngine* rndmEngineMod;
128 
132 
135 
139 
148 
154 
156 
162 
164 
169 
171 
176 
179 
182 
185 
186 };
187 
188 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
189 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
190 
191 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
192 
194 {
195  displacementFlag = val;
196 }
197 
198 inline
200 {
201  if (p != particle) {
202  particle = p;
203  mass = p->GetPDGMass();
206  }
207 }
208 
209 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
210 inline
212 {
213  G4double temptlimit = tlimit;
214  if(tlimit > tlimitmin)
215  {
216  do {
217  temptlimit = G4RandGauss::shoot(rndmEngineMod,tlimit,0.3*tlimit);
218  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
219  } while ((temptlimit < tlimitmin) ||
220  (temptlimit > 2.*tlimit-tlimitmin));
221  }
222  else temptlimit = tlimitmin;
223 
224  return temptlimit;
225 }
226 
227 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
228 inline
230 {
231  lnZ = G4Log(Zeff);
232  // correction in theta0 formula
233  G4double w = G4Exp(lnZ/6.);
234  G4double facz = 0.990395+w*(-0.168386+w*0.093286) ;
235  coeffth1 = facz*(1. - 8.7780e-2/Zeff);
236  coeffth2 = facz*(4.0780e-2 + 1.7315e-4*Zeff);
237 
238  // tail parameters
239  G4double Z13 = w*w;
240  coeffc1 = 2.3785 - Z13*(4.1981e-1 - Z13*6.3100e-2);
241  coeffc2 = 4.7526e-1 + Z13*(1.7694 - Z13*3.3885e-1);
242  coeffc3 = 2.3683e-1 - Z13*(1.8111 - Z13*3.2774e-1);
243  coeffc4 = 1.7888e-2 + Z13*(1.9659e-2 - Z13*2.6664e-3);
244 
245  Z2 = Zeff*Zeff;
246  Z23 = Z13*Z13;
247 
248  Zold = Zeff;
249 }
250 
251 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
252 
253 inline
255 {
256  // 'large angle scattering'
257  // 2 model functions with correct xmean and x2mean
258  G4double a = (2.*xmeanth+9.*x2meanth-3.)/(2.*xmeanth-3.*x2meanth+1.);
259  G4double prob = (a+2.)*xmeanth/a;
260 
261  // sampling
262  G4double cth = 1.;
263  if(rndmEngineMod->flat() < prob) {
264  cth = -1.+2.*G4Exp(G4Log(rndmEngineMod->flat())/(a+1.));
265  } else {
266  cth = -1.+2.*rndmEngineMod->flat();
267  }
268  return cth;
269 }
270 
271 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
272 
273 
274 #endif
275 
G4ParticleChangeForMSC * fParticleChange
ThreeVector shoot(const G4int Ap, const G4int Af)
CLHEP::HepRandomEngine * rndmEngineMod
static const double MeV
Definition: G4SIunits.hh:211
G4UrbanMscModel(const G4String &nam="UrbanMsc")
CLHEP::Hep3Vector G4ThreeVector
virtual ~G4UrbanMscModel()
G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety)
void SetNewDisplacementFlag(G4bool)
const G4double w[NPOINTSGL]
G4UrbanMscModel & operator=(const G4UrbanMscModel &right)
G4double a
Definition: TRTMaterials.hh:39
const G4ParticleDefinition * particle
G4double currentRadLength
int G4int
Definition: G4Types.hh:78
G4double SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy)
const G4ParticleDefinition * positron
void SampleDisplacementNew(G4double sinTheta, G4double phi)
G4double Randomizetlimit()
bool G4bool
Definition: G4Types.hh:79
G4double ComputeGeomPathLength(G4double truePathLength)
const G4MaterialCutsCouple * couple
G4double ComputeTrueStepLength(G4double geomStepLength)
static const G4double emax
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
void SetParticle(const G4ParticleDefinition *)
void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4LossTableManager * theManager
G4double GetPDGMass() const
void SampleDisplacement(G4double sinTheta, G4double phi)
static const G4double c_highland
G4double currentKinEnergy
double G4double
Definition: G4Types.hh:76
G4double ComputeTheta0(G4double truePathLength, G4double KineticEnergy)
static const double eplus
Definition: G4SIunits.hh:196
G4double GetPDGCharge() const
#define DBL_MAX
Definition: templates.hh:83
G4double SimpleScattering(G4double xmeanth, G4double x2meanth)
void StartTracking(G4Track *)
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *particle, G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight=0., G4double cut=0., G4double emax=DBL_MAX)
G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &currentMinimalStep)