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: G4UrbanMscModel90.hh,v 1.4 2008/10/29 14:15:30 vnivanch Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-02 $ |
---|
28 | // |
---|
29 | // ------------------------------------------------------------------- |
---|
30 | // |
---|
31 | // |
---|
32 | // GEANT4 Class header file |
---|
33 | // |
---|
34 | // |
---|
35 | // File name: G4UrbanMscModel90 |
---|
36 | // |
---|
37 | // Author: V.Ivanchenko clone Laszlo Urban model |
---|
38 | // |
---|
39 | // Creation date: 07.12.2007 |
---|
40 | // |
---|
41 | // Modifications: |
---|
42 | // |
---|
43 | // |
---|
44 | // Class Description: |
---|
45 | // |
---|
46 | // Implementation of the model of multiple scattering based on |
---|
47 | // H.W.Lewis Phys Rev 78 (1950) 526 and L.Urban model |
---|
48 | |
---|
49 | // ------------------------------------------------------------------- |
---|
50 | // |
---|
51 | |
---|
52 | #ifndef G4UrbanMscModel90_h |
---|
53 | #define G4UrbanMscModel90_h 1 |
---|
54 | |
---|
55 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
56 | |
---|
57 | #include "G4VMscModel.hh" |
---|
58 | #include "G4PhysicsTable.hh" |
---|
59 | #include "G4MscStepLimitType.hh" |
---|
60 | |
---|
61 | class G4ParticleChangeForMSC; |
---|
62 | class G4SafetyHelper; |
---|
63 | class G4LossTableManager; |
---|
64 | |
---|
65 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
66 | |
---|
67 | class G4UrbanMscModel90 : public G4VMscModel |
---|
68 | { |
---|
69 | |
---|
70 | public: |
---|
71 | |
---|
72 | G4UrbanMscModel90(const G4String& nam = "UrbanMscUni90"); |
---|
73 | |
---|
74 | virtual ~G4UrbanMscModel90(); |
---|
75 | |
---|
76 | void Initialise(const G4ParticleDefinition*, const G4DataVector&); |
---|
77 | |
---|
78 | G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition* particle, |
---|
79 | G4double KineticEnergy, |
---|
80 | G4double AtomicNumber, |
---|
81 | G4double AtomicWeight=0., |
---|
82 | G4double cut =0., |
---|
83 | G4double emax=DBL_MAX); |
---|
84 | |
---|
85 | void SampleSecondaries(std::vector<G4DynamicParticle*>*, |
---|
86 | const G4MaterialCutsCouple*, |
---|
87 | const G4DynamicParticle*, |
---|
88 | G4double, |
---|
89 | G4double); |
---|
90 | |
---|
91 | void SampleScattering(const G4DynamicParticle*, |
---|
92 | G4double safety); |
---|
93 | |
---|
94 | G4double ComputeTruePathLengthLimit(const G4Track& track, |
---|
95 | G4PhysicsTable* theLambdaTable, |
---|
96 | G4double currentMinimalStep); |
---|
97 | |
---|
98 | G4double ComputeGeomPathLength(G4double truePathLength); |
---|
99 | |
---|
100 | G4double ComputeTrueStepLength(G4double geomStepLength); |
---|
101 | |
---|
102 | G4double ComputeTheta0(G4double truePathLength, |
---|
103 | G4double KineticEnergy); |
---|
104 | |
---|
105 | private: |
---|
106 | |
---|
107 | G4double SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy); |
---|
108 | |
---|
109 | G4double SampleDisplacement(); |
---|
110 | |
---|
111 | G4double LatCorrelation(); |
---|
112 | |
---|
113 | void GeomLimit(const G4Track& track); |
---|
114 | |
---|
115 | inline G4double GetLambda(G4double kinEnergy); |
---|
116 | |
---|
117 | inline void SetParticle(const G4ParticleDefinition*); |
---|
118 | |
---|
119 | // hide assignment operator |
---|
120 | G4UrbanMscModel90 & operator=(const G4UrbanMscModel90 &right); |
---|
121 | G4UrbanMscModel90(const G4UrbanMscModel90&); |
---|
122 | |
---|
123 | const G4ParticleDefinition* particle; |
---|
124 | G4ParticleChangeForMSC* fParticleChange; |
---|
125 | |
---|
126 | G4SafetyHelper* safetyHelper; |
---|
127 | G4PhysicsTable* theLambdaTable; |
---|
128 | const G4MaterialCutsCouple* couple; |
---|
129 | G4LossTableManager* theManager; |
---|
130 | |
---|
131 | G4double mass; |
---|
132 | G4double charge; |
---|
133 | |
---|
134 | G4double masslimite,masslimitmu; |
---|
135 | |
---|
136 | G4double taubig; |
---|
137 | G4double tausmall; |
---|
138 | G4double taulim; |
---|
139 | G4double currentTau; |
---|
140 | G4double frscaling1,frscaling2; |
---|
141 | G4double tlimit; |
---|
142 | G4double tlimitmin; |
---|
143 | G4double tlimitminfix; |
---|
144 | |
---|
145 | G4double nstepmax; |
---|
146 | G4double geombig; |
---|
147 | G4double geommin; |
---|
148 | G4double geomlimit; |
---|
149 | G4double skindepth; |
---|
150 | G4double smallstep; |
---|
151 | |
---|
152 | G4double presafety; |
---|
153 | |
---|
154 | G4double lambda0; |
---|
155 | G4double lambdaeff; |
---|
156 | G4double tPathLength; |
---|
157 | G4double zPathLength; |
---|
158 | G4double par1,par2,par3; |
---|
159 | |
---|
160 | G4double stepmin; |
---|
161 | |
---|
162 | G4double currentKinEnergy; |
---|
163 | G4double currentRange; |
---|
164 | G4double currentRadLength; |
---|
165 | |
---|
166 | G4double Zeff; |
---|
167 | |
---|
168 | G4int currentMaterialIndex; |
---|
169 | |
---|
170 | G4bool isInitialized; |
---|
171 | G4bool inside; |
---|
172 | G4bool insideskin; |
---|
173 | |
---|
174 | }; |
---|
175 | |
---|
176 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
177 | |
---|
178 | inline |
---|
179 | G4double G4UrbanMscModel90::GetLambda(G4double e) |
---|
180 | { |
---|
181 | G4double x; |
---|
182 | if(theLambdaTable) { |
---|
183 | G4bool b; |
---|
184 | x = ((*theLambdaTable)[currentMaterialIndex])->GetValue(e, b); |
---|
185 | } else { |
---|
186 | x = CrossSection(couple,particle,e); |
---|
187 | } |
---|
188 | if(x > DBL_MIN) x = 1./x; |
---|
189 | else x = DBL_MAX; |
---|
190 | return x; |
---|
191 | } |
---|
192 | |
---|
193 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
194 | |
---|
195 | inline |
---|
196 | void G4UrbanMscModel90::SetParticle(const G4ParticleDefinition* p) |
---|
197 | { |
---|
198 | if (p != particle) { |
---|
199 | particle = p; |
---|
200 | mass = p->GetPDGMass(); |
---|
201 | charge = p->GetPDGCharge()/eplus; |
---|
202 | } |
---|
203 | } |
---|
204 | |
---|
205 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
206 | |
---|
207 | #endif |
---|
208 | |
---|