Geant4  10.03
G4UrbanAdjointMscModel.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: G4UrbanAdjointMscModel
36 //
37 // Author: Laszlo Urban
38 //
39 // Creation date: 19.02.2013
40 //
41 // Created from G4UrbanAdjointMscModel96
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 G4UrbanAdjointMscModel_h
55 #define G4UrbanAdjointMscModel_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 #include "G4Electron.hh"
66 
67 
69 class G4SafetyHelper;
70 class G4LossTableManager;
71 
72 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
73 
75 {
76 
77 public:
78 
79  explicit G4UrbanAdjointMscModel(const G4String& nam = "UrbanMsc");
80 
81  virtual ~G4UrbanAdjointMscModel();
82 
83  virtual void Initialise(const G4ParticleDefinition*,
84  const G4DataVector&) override;
85 
86  virtual void StartTracking(G4Track*) override;
87 
88  virtual G4double
90  G4double KineticEnergy,
91  G4double AtomicNumber,
92  G4double AtomicWeight=0.,
93  G4double cut =0.,
94  G4double emax=DBL_MAX) override;
95 
97  G4double safety) override;
98 
99  virtual G4double
100  ComputeTruePathLengthLimit(const G4Track& track,
101  G4double& currentMinimalStep) override;
102 
103  virtual G4double ComputeGeomPathLength(G4double truePathLength) override;
104 
105  virtual G4double ComputeTrueStepLength(G4double geomStepLength) override;
106 
107  G4double ComputeTheta0(G4double truePathLength, G4double KineticEnergy);
108 
109  inline void SetNewDisplacementFlag(G4bool);
110 
111 private:
112 
113  G4double SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy);
114 
115  void SampleDisplacement(G4double sinTheta, G4double phi);
116 
117  void SampleDisplacementNew(G4double sinTheta, G4double phi);
118 
119  inline void SetParticle(const G4ParticleDefinition*);
120 
121  inline void UpdateCache();
122 
123  inline G4double Randomizetlimit();
124 
125  inline G4double SimpleScattering(G4double xmeanth, G4double x2meanth);
126 
127  // hide assignment operator
130 
131  CLHEP::HepRandomEngine* rndmEngineMod;
132 
136 
139 
143 
152 
158 
160 
166 
168 
173 
175 
180 
183 
186 
189 
190 };
191 
192 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
193 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
194 
195 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
196 
198 {
199  displacementFlag = val;
200 }
201 
202 inline
204 { const G4ParticleDefinition* p1 =p;
205 
206  if (p->GetParticleName() =="adj_e-") p1= G4Electron::Electron();
207 
208  if (p1 != particle) {
209  particle = p1;
210  mass = p1->GetPDGMass();
213  }
214 }
215 
216 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
217 
219 {
220  G4double temptlimit = tlimit;
221  if(tlimit > tlimitmin)
222  {
223  G4double delta = tlimit-tlimitmin;
224  do {
225  temptlimit = G4RandGauss::shoot(rndmEngineMod,tlimit,0.1*delta);
226  // Loop checking, 10-Apr-2016, Laszlo Urban
227  } while ((temptlimit < tlimit-delta) ||
228  (temptlimit > tlimit+delta));
229  }
230  else { temptlimit = tlimitmin; }
231 
232  return temptlimit;
233 }
234 
235 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
236 
238 {
239  lnZ = G4Log(Zeff);
240  // correction in theta0 formula
241  G4double w = G4Exp(lnZ/6.);
242  G4double facz = 0.990395+w*(-0.168386+w*0.093286) ;
243  coeffth1 = facz*(1. - 8.7780e-2/Zeff);
244  coeffth2 = facz*(4.0780e-2 + 1.7315e-4*Zeff);
245 
246  // tail parameters
247  G4double Z13 = w*w;
248  coeffc1 = 2.3785 - Z13*(4.1981e-1 - Z13*6.3100e-2);
249  coeffc2 = 4.7526e-1 + Z13*(1.7694 - Z13*3.3885e-1);
250  coeffc3 = 2.3683e-1 - Z13*(1.8111 - Z13*3.2774e-1);
251  coeffc4 = 1.7888e-2 + Z13*(1.9659e-2 - Z13*2.6664e-3);
252 
253  Z2 = Zeff*Zeff;
254  Z23 = Z13*Z13;
255 
256  Zold = Zeff;
257 }
258 
259 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
260 
261 inline
263 {
264  // 'large angle scattering'
265  // 2 model functions with correct xmean and x2mean
266  G4double a = (2.*xmeanth+9.*x2meanth-3.)/(2.*xmeanth-3.*x2meanth+1.);
267  G4double prob = (a+2.)*xmeanth/a;
268 
269  // sampling
270  G4double cth = 1.;
271  if(rndmEngineMod->flat() < prob) {
272  cth = -1.+2.*G4Exp(G4Log(rndmEngineMod->flat())/(a+1.));
273  } else {
274  cth = -1.+2.*rndmEngineMod->flat();
275  }
276  return cth;
277 }
278 
279 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
280 
281 
282 #endif
283 
void SampleDisplacementNew(G4double sinTheta, G4double phi)
ThreeVector shoot(const G4int Ap, const G4int Af)
G4double SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy)
CLHEP::Hep3Vector G4ThreeVector
CLHEP::HepRandomEngine * rndmEngineMod
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &currentMinimalStep) override
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
G4UrbanAdjointMscModel(const G4String &nam="UrbanMsc")
int G4int
Definition: G4Types.hh:78
G4ParticleChangeForMSC * fParticleChange
const G4String & GetParticleName() const
G4LossTableManager * theManager
virtual G4double ComputeTrueStepLength(G4double geomStepLength) override
const G4MaterialCutsCouple * couple
virtual void StartTracking(G4Track *) override
bool G4bool
Definition: G4Types.hh:79
const G4ParticleDefinition * particle
static constexpr double eplus
Definition: G4SIunits.hh:199
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *particle, G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight=0., G4double cut=0., G4double emax=DBL_MAX) override
void SetParticle(const G4ParticleDefinition *)
G4double SimpleScattering(G4double xmeanth, G4double x2meanth)
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
const G4ParticleDefinition * positron
G4double GetPDGMass() const
G4double ComputeTheta0(G4double truePathLength, G4double KineticEnergy)
void SampleDisplacement(G4double sinTheta, G4double phi)
virtual G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety) override
static G4Electron * Electron()
Definition: G4Electron.cc:94
virtual G4double ComputeGeomPathLength(G4double truePathLength) override
double G4double
Definition: G4Types.hh:76
G4UrbanAdjointMscModel & operator=(const G4UrbanAdjointMscModel &right)=delete
G4double GetPDGCharge() const
#define DBL_MAX
Definition: templates.hh:83
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override