Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4UrbanMscModel.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: $
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 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
124  G4UrbanMscModel & operator=(const G4UrbanMscModel &right) = delete;
125  G4UrbanMscModel(const G4UrbanMscModel&) = delete;
126 
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 inline
193 {
194  if (p != particle) {
195  particle = p;
196  mass = p->GetPDGMass();
199  }
200 }
201 
202 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
203 
205 {
206  G4double res = tlimitmin;
207  if(tlimit > tlimitmin)
208  {
210  res = std::max(res, tlimitmin);
211  }
212  return res;
213 }
214 
215 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
216 
218 {
219  lnZ = G4Log(Zeff);
220  // correction in theta0 formula
221  G4double w = G4Exp(lnZ/6.);
222  G4double facz = 0.990395+w*(-0.168386+w*0.093286) ;
223  coeffth1 = facz*(1. - 8.7780e-2/Zeff);
224  coeffth2 = facz*(4.0780e-2 + 1.7315e-4*Zeff);
225 
226  // tail parameters
227  G4double Z13 = w*w;
228  coeffc1 = 2.3785 - Z13*(4.1981e-1 - Z13*6.3100e-2);
229  coeffc2 = 4.7526e-1 + Z13*(1.7694 - Z13*3.3885e-1);
230  coeffc3 = 2.3683e-1 - Z13*(1.8111 - Z13*3.2774e-1);
231  coeffc4 = 1.7888e-2 + Z13*(1.9659e-2 - Z13*2.6664e-3);
232 
233  Z2 = Zeff*Zeff;
234  Z23 = Z13*Z13;
235 
236  Zold = Zeff;
237 }
238 
239 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
240 
241 inline
243 {
244  // 'large angle scattering'
245  // 2 model functions with correct xmean and x2mean
246  G4double a = (2.*xmeanth+9.*x2meanth-3.)/(2.*xmeanth-3.*x2meanth+1.);
247  G4double prob = (a+2.)*xmeanth/a;
248 
249  // sampling
250  G4double rdm = rndmEngineMod->flat();
251  G4double cth = (rndmEngineMod->flat() < prob)
252  ? -1.+2.*G4Exp(G4Log(rdm)/(a+1.)) : -1.+2.*rdm;
253  return cth;
254 }
255 
256 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
257 
258 
259 #endif
260 
G4UrbanMscModel(const G4String &nam="UrbanMsc")
G4double Randomizetlimit()
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
T max(const T t1, const T t2)
brief Return the largest of the two arguments
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
const G4ParticleDefinition * positron
virtual G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety) override
virtual void StartTracking(G4Track *) override
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4LossTableManager * theManager
const char * p
Definition: xmltok.h:285
void SampleDisplacement(G4double sinTheta, G4double phi)
G4double ComputeTheta0(G4double truePathLength, G4double KineticEnergy)
G4double GetPDGCharge() const
static const G4double emax
G4double GetPDGMass() const
G4double G4Log(G4double x)
Definition: G4Log.hh:230
const G4MaterialCutsCouple * couple
G4double SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy)
G4ParticleChangeForMSC * fParticleChange
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
virtual double flat()=0
virtual ~G4UrbanMscModel()
void SetParticle(const G4ParticleDefinition *)
virtual G4double ComputeGeomPathLength(G4double truePathLength) override
CLHEP::HepRandomEngine * rndmEngineMod
ThreeVector shoot(const G4int Ap, const G4int Af)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *particle, G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight=0., G4double cut=0., G4double emax=DBL_MAX) override
int G4int
Definition: G4Types.hh:78
G4double SimpleScattering(G4double xmeanth, G4double x2meanth)
G4double currentRadLength
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &currentMinimalStep) override
G4double currentKinEnergy
void SampleDisplacementNew(G4double sinTheta, G4double phi)
#define DBL_MAX
Definition: templates.hh:83
G4UrbanMscModel & operator=(const G4UrbanMscModel &right)=delete
virtual G4double ComputeTrueStepLength(G4double geomStepLength) override
const G4ParticleDefinition * particle
static constexpr double eplus