source: trunk/source/processes/electromagnetic/utils/include/G4EmCorrections.hh@ 1057

Last change on this file since 1057 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

File size: 11.0 KB
Line 
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: G4EmCorrections.hh,v 1.24 2008/09/12 14:44:48 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-02 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class header file
32//
33//
34// File name: G4EmCorrections
35//
36// Author: Vladimir Ivanchenko
37//
38// Creation date: 13.01.2005
39//
40// Modifications:
41// 28.04.2006 General cleanup, add finite size corrections (V.Ivanchenko)
42// 13.05.2006 Add corrections for ion stopping (V.Ivanhcenko)
43// 20.05.2008 Removed Finite Size correction (V.Ivanchenko)
44// 12.09.2008 Added inlined interfaces to effective charge (V.Ivanchenko)
45//
46// Class Description:
47//
48// This class provides calculation of EM corrections to ionisation
49//
50
51// -------------------------------------------------------------------
52//
53
54#ifndef G4EmCorrections_h
55#define G4EmCorrections_h 1
56
57#include "globals.hh"
58#include "G4AtomicShells.hh"
59#include "G4ionEffectiveCharge.hh"
60#include "G4Material.hh"
61#include "G4ParticleDefinition.hh"
62#include "G4NistManager.hh"
63
64class G4VEmModel;
65class G4PhysicsVector;
66class G4IonTable;
67class G4MaterialCutsCouple;
68
69class G4EmCorrections
70{
71
72public:
73
74 G4EmCorrections();
75
76 virtual ~G4EmCorrections();
77
78 G4double HighOrderCorrections(const G4ParticleDefinition*,
79 const G4Material*,
80 G4double kineticEnergy,
81 G4double cutEnergy);
82
83 G4double IonHighOrderCorrections(const G4ParticleDefinition*,
84 const G4MaterialCutsCouple*,
85 G4double kineticEnergy);
86
87 G4double ComputeIonCorrections(const G4ParticleDefinition*,
88 const G4Material*,
89 G4double kineticEnergy);
90
91 G4double IonBarkasCorrection(const G4ParticleDefinition*,
92 const G4Material*,
93 G4double kineticEnergy);
94
95 G4double Bethe(const G4ParticleDefinition*,
96 const G4Material*,
97 G4double kineticEnergy);
98
99 G4double SpinCorrection(const G4ParticleDefinition*,
100 const G4Material*,
101 G4double kineticEnergy);
102
103 G4double KShellCorrection(const G4ParticleDefinition*,
104 const G4Material*,
105 G4double kineticEnergy);
106
107 G4double LShellCorrection(const G4ParticleDefinition*,
108 const G4Material*,
109 G4double kineticEnergy);
110
111 G4double ShellCorrection(const G4ParticleDefinition*,
112 const G4Material*,
113 G4double kineticEnergy);
114
115 G4double ShellCorrectionSTD(const G4ParticleDefinition*,
116 const G4Material*,
117 G4double kineticEnergy);
118
119 G4double DensityCorrection(const G4ParticleDefinition*,
120 const G4Material*,
121 G4double kineticEnergy);
122
123 G4double BarkasCorrection(const G4ParticleDefinition*,
124 const G4Material*,
125 G4double kineticEnergy);
126
127 G4double BlochCorrection(const G4ParticleDefinition*,
128 const G4Material*,
129 G4double kineticEnergy);
130
131 G4double MottCorrection(const G4ParticleDefinition*,
132 const G4Material*,
133 G4double kineticEnergy);
134
135 G4double NuclearDEDX(const G4ParticleDefinition*,
136 const G4Material*,
137 G4double kineticEnergy,
138 G4bool fluct = true);
139
140 void AddStoppingData(G4int Z, G4int A, const G4String& materialName,
141 G4PhysicsVector* dVector);
142
143 void InitialiseForNewRun();
144
145 // effective charge correction using stopping power data
146 G4double EffectiveChargeCorrection(const G4ParticleDefinition*,
147 const G4Material*,
148 G4double kineticEnergy);
149
150 // effective charge of an ion
151 inline G4double GetParticleCharge(const G4ParticleDefinition*,
152 const G4Material*,
153 G4double kineticEnergy);
154
155 inline
156 G4double EffectiveChargeSquareRatio(const G4ParticleDefinition*,
157 const G4Material*,
158 G4double kineticEnergy);
159
160 // ionisation models for ions
161 inline void SetIonisationModels(G4VEmModel* m1 = 0, G4VEmModel* m2 = 0);
162
163 inline G4int GetNumberOfStoppingVectors();
164
165private:
166
167 void Initialise();
168
169 void BuildCorrectionVector();
170
171 void SetupKinematics(const G4ParticleDefinition*,
172 const G4Material*,
173 G4double kineticEnergy);
174
175 G4double KShell(G4double theta, G4double eta);
176
177 G4double LShell(G4double theta, G4double eta);
178
179 G4int Index(G4double x, G4double* y, G4int n);
180
181 G4double Value(G4double xv, G4double x1, G4double x2, G4double y1, G4double y2);
182
183 G4double Value2(G4double xv, G4double yv, G4double x1, G4double x2,
184 G4double y1, G4double y2,
185 G4double z11, G4double z21, G4double z12, G4double z22);
186
187 G4double NuclearStoppingPower(G4double e, G4double z1, G4double z2,
188 G4double m1, G4double m2);
189
190 // hide assignment operator
191 G4EmCorrections & operator=(const G4EmCorrections &right);
192 G4EmCorrections(const G4EmCorrections&);
193
194 G4double engBarkas[47];
195 G4double corBarkas[47];
196 G4double ed[104];
197 G4double a[104];
198 G4double theZieglerFactor;
199 G4double alpha2;
200 G4bool lossFlucFlag;
201
202 G4int verbose;
203
204 G4int nK;
205 G4int nL;
206 G4int nEtaK;
207 G4int nEtaL;
208
209 G4double COSEB[14];
210 G4double COSXI[14];
211 G4double ZD[11];
212
213 G4double TheK[20];
214 G4double SK[20];
215 G4double TK[20];
216 G4double UK[20];
217 G4double VK[20];
218 G4double ZK[20];
219
220 G4double TheL[26];
221 G4double SL[26];
222 G4double TL[26];
223 G4double UL[26];
224 G4double VL[26];
225
226 G4double Eta[29];
227 G4double CK[20][29];
228 G4double CL[26][28];
229 G4double HM[53];
230 G4double HN[31];
231 G4double MSH[93];
232 G4double TAU[93];
233 G4double Z23[100];
234
235 std::vector<const G4Material*> currmat;
236 std::vector<G4double> thcorr[100];
237 size_t ncouples;
238
239 const G4ParticleDefinition* particle;
240 const G4ParticleDefinition* curParticle;
241 const G4Material* material;
242 const G4Material* curMaterial;
243 const G4ElementVector* theElementVector;
244 const G4double* atomDensity;
245
246 G4int numberOfElements;
247 G4double kinEnergy;
248 G4double mass;
249 G4double massFactor;
250 G4double formfact;
251 G4double eth;
252 G4double tau;
253 G4double gamma;
254 G4double bg2;
255 G4double beta2;
256 G4double beta;
257 G4double ba2;
258 G4double tmax;
259 G4double charge;
260 G4double q2;
261 G4double eCorrMin;
262 G4double eCorrMax;
263 G4int nbinCorr;
264
265 G4AtomicShells shells;
266 G4ionEffectiveCharge effCharge;
267
268 G4NistManager* nist;
269 const G4IonTable* ionTable;
270 G4VEmModel* ionLEModel;
271 G4VEmModel* ionHEModel;
272
273 // Ion stopping data
274 G4int nIons;
275 G4int idx;
276 G4int currentZ;
277 std::vector<G4int> Zion;
278 std::vector<G4int> Aion;
279 std::vector<G4String> materialName;
280
281 std::vector<const G4ParticleDefinition*> ionList;
282
283 std::vector<const G4Material*> materialList;
284 std::vector<G4PhysicsVector*> stopData;
285 G4PhysicsVector* curVector;
286};
287
288inline G4int G4EmCorrections::Index(G4double x, G4double* y, G4int n)
289{
290 G4int iddd = n-1;
291 do {iddd--;} while (iddd>0 && x<y[iddd]);
292 return iddd;
293}
294
295inline G4double G4EmCorrections::Value(G4double xv, G4double x1, G4double x2,
296 G4double y1, G4double y2)
297{
298 return y1 + (y2 - y1)*(xv - x1)/(x2 - x1);
299}
300
301inline G4double G4EmCorrections::Value2(G4double xv, G4double yv,
302 G4double x1, G4double x2,
303 G4double y1, G4double y2,
304 G4double z11, G4double z21,
305 G4double z12, G4double z22)
306{
307 return (z11*(x2-xv)*(y2-yv) + z22*(xv-x1)*(yv-y1) +
308 0.5*(z12*((x2-xv)*(yv-y1)+(xv-x1)*(y2-yv))+
309 z21*((xv-x1)*(y2-yv)+(yv-y1)*(x2-xv))))
310 / ((x2-x1)*(y2-y1));
311}
312
313inline
314void G4EmCorrections::SetIonisationModels(G4VEmModel* m1, G4VEmModel* m2)
315{
316 if(m1) ionLEModel = m1;
317 if(m2) ionHEModel = m2;
318}
319
320inline G4int G4EmCorrections::GetNumberOfStoppingVectors()
321{
322 return nIons;
323}
324
325inline G4double
326G4EmCorrections::GetParticleCharge(const G4ParticleDefinition* p,
327 const G4Material* mat,
328 G4double kineticEnergy)
329{
330 return effCharge.EffectiveCharge(p,mat,kineticEnergy);
331}
332
333inline G4double
334G4EmCorrections::EffectiveChargeSquareRatio(const G4ParticleDefinition* p,
335 const G4Material* mat,
336 G4double kineticEnergy)
337{
338 return effCharge.EffectiveChargeSquareRatio(p,mat,kineticEnergy);
339}
340
341inline void G4EmCorrections::SetupKinematics(const G4ParticleDefinition* p,
342 const G4Material* mat,
343 G4double kineticEnergy)
344{
345 if(kineticEnergy != kinEnergy || p != particle) {
346 particle = p;
347 kinEnergy = kineticEnergy;
348 mass = p->GetPDGMass();
349 tau = kineticEnergy / mass;
350 gamma = 1.0 + tau;
351 bg2 = tau * (tau+2.0);
352 beta2 = bg2/(gamma*gamma);
353 beta = std::sqrt(beta2);
354 ba2 = beta2/alpha2;
355 G4double ratio = electron_mass_c2/mass;
356 tmax = 2.0*electron_mass_c2*bg2 /(1. + 2.0*gamma*ratio + ratio*ratio);
357 charge = p->GetPDGCharge()/eplus;
358 if(charge < 1.5) {q2 = charge*charge;}
359 else {
360 q2 = effCharge.EffectiveChargeSquareRatio(p,mat,kinEnergy);
361 charge = std::sqrt(q2);
362 }
363 }
364 if(mat != material) {
365 material = mat;
366 theElementVector = material->GetElementVector();
367 atomDensity = material->GetAtomicNumDensityVector();
368 numberOfElements = material->GetNumberOfElements();
369 }
370}
371
372//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
373
374#endif
Note: See TracBrowser for help on using the repository browser.