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

Last change on this file since 1317 was 1315, checked in by garnier, 15 years ago

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

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