Geant4_10
G4EmCorrections.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: G4EmCorrections.hh 70371 2013-05-29 15:18:07Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4EmCorrections
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 13.01.2005
38 //
39 // Modifications:
40 // 28.04.2006 General cleanup, add finite size corrections (V.Ivanchenko)
41 // 13.05.2006 Add corrections for ion stopping (V.Ivanhcenko)
42 // 20.05.2008 Removed Finite Size correction (V.Ivanchenko)
43 // 12.09.2008 Added inlined interfaces to effective charge (V.Ivanchenko)
44 // 19.04.2012 Fix reproducibility problem (A.Ribon)
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 
58 
59 #include "globals.hh"
60 #include "G4ionEffectiveCharge.hh"
61 #include "G4Material.hh"
62 #include "G4ParticleDefinition.hh"
63 #include "G4NistManager.hh"
64 
65 class G4VEmModel;
66 class G4PhysicsVector;
67 class G4IonTable;
70 
72 {
73 
74 public:
75 
77 
78  virtual ~G4EmCorrections();
79 
81  const G4Material*,
82  G4double kineticEnergy,
83  G4double cutEnergy);
84 
86  const G4MaterialCutsCouple*,
87  G4double kineticEnergy);
88 
90  const G4Material*,
91  G4double kineticEnergy);
92 
94  const G4Material*,
95  G4double kineticEnergy);
96 
98  const G4Material*,
99  G4double kineticEnergy);
100 
102  const G4Material*,
103  G4double kineticEnergy);
104 
106  const G4Material*,
107  G4double kineticEnergy);
108 
110  const G4Material*,
111  G4double kineticEnergy);
112 
114  const G4Material*,
115  G4double kineticEnergy);
116 
118  const G4Material*,
119  G4double kineticEnergy);
120 
122  const G4Material*,
123  G4double kineticEnergy);
124 
126  const G4Material*,
127  G4double kineticEnergy);
128 
130  const G4Material*,
131  G4double kineticEnergy);
132 
134  const G4Material*,
135  G4double kineticEnergy);
136 
138  const G4Material*,
139  G4double kineticEnergy,
140  G4bool fluct = true);
141 
142  void AddStoppingData(G4int Z, G4int A, const G4String& materialName,
143  G4PhysicsVector* dVector);
144 
145  void InitialiseForNewRun();
146 
147  // effective charge correction using stopping power data
149  const G4Material*,
150  G4double kineticEnergy);
151 
152  // effective charge of an ion
154  const G4Material*,
155  G4double kineticEnergy);
156 
157  inline
159  const G4Material*,
160  G4double kineticEnergy);
161 
162  // ionisation models for ions
163  inline void SetIonisationModels(G4VEmModel* m1 = 0, G4VEmModel* m2 = 0);
164 
166 
167 private:
168 
169  void Initialise();
170 
171  void BuildCorrectionVector();
172 
173  void SetupKinematics(const G4ParticleDefinition*,
174  const G4Material*,
175  G4double kineticEnergy);
176 
177  G4double KShell(G4double theta, G4double eta);
178 
179  G4double LShell(G4double theta, G4double eta);
180 
181  G4int Index(G4double x, G4double* y, G4int n);
182 
184 
185  G4double Value2(G4double xv, G4double yv, G4double x1, G4double x2,
186  G4double y1, G4double y2,
187  G4double z11, G4double z21, G4double z12, G4double z22);
188 
189  G4double NuclearStoppingPower(G4double e, G4double z1, G4double z2,
190  G4double m1, G4double m2);
191 
192  // hide assignment operator
193  G4EmCorrections & operator=(const G4EmCorrections &right);
195 
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 Z23[100];
232 
233  G4LPhysicsFreeVector* BarkasCorr;
234  G4LPhysicsFreeVector* ThetaK;
235  G4LPhysicsFreeVector* ThetaL;
236 
237  std::vector<const G4Material*> currmat;
238  std::map< G4int, std::vector<G4double> > thcorr;
239  size_t ncouples;
240 
241  const G4ParticleDefinition* particle;
242  const G4ParticleDefinition* curParticle;
243  const G4Material* material;
244  const G4Material* curMaterial;
245  const G4ElementVector* theElementVector;
246  const G4double* atomDensity;
247 
248  G4int numberOfElements;
249  G4double kinEnergy;
250  G4double mass;
251  G4double massFactor;
252  G4double formfact;
253  G4double eth;
254  G4double tau;
255  G4double gamma;
256  G4double bg2;
257  G4double beta2;
258  G4double beta;
259  G4double ba2;
260  G4double tmax;
261  G4double charge;
262  G4double q2;
263  G4double eCorrMin;
264  G4double eCorrMax;
265  G4int nbinCorr;
266 
267  G4ionEffectiveCharge effCharge;
268 
269  G4NistManager* nist;
270  G4IonTable* ionTable;
271  G4VEmModel* ionLEModel;
272  G4VEmModel* ionHEModel;
273 
274  // Ion stopping data
275  G4int nIons;
276  G4int idx;
277  G4int currentZ;
278  std::vector<G4int> Zion;
279  std::vector<G4int> Aion;
280  std::vector<G4String> materialName;
281 
282  std::vector<const G4ParticleDefinition*> ionList;
283 
284  std::vector<const G4Material*> materialList;
285  std::vector<G4PhysicsVector*> stopData;
286  G4PhysicsVector* curVector;
287 };
288 
289 inline G4int G4EmCorrections::Index(G4double x, G4double* y, G4int n)
290 {
291  G4int iddd = n-1;
292  do {iddd--;} while (iddd>0 && x<y[iddd]);
293  return iddd;
294 }
295 
296 inline G4double G4EmCorrections::Value(G4double xv, G4double x1, G4double x2,
298 {
299  return y1 + (y2 - y1)*(xv - x1)/(x2 - x1);
300 }
301 
302 inline G4double G4EmCorrections::Value2(G4double xv, G4double yv,
305  G4double z11, G4double z21,
306  G4double z12, G4double z22)
307 {
308  return (z11*(x2-xv)*(y2-yv) + z22*(xv-x1)*(yv-y1) +
309  0.5*(z12*((x2-xv)*(yv-y1)+(xv-x1)*(y2-yv))+
310  z21*((xv-x1)*(y2-yv)+(yv-y1)*(x2-xv))))
311  / ((x2-x1)*(y2-y1));
312 }
313 
314 inline void
316 {
317  if(mod1) { ionLEModel = mod1; }
318  if(mod2) { ionHEModel = mod2; }
319 }
320 
322 {
323  return nIons;
324 }
325 
326 inline G4double
328  const G4Material* mat,
329  G4double kineticEnergy)
330 {
331  return effCharge.EffectiveCharge(p,mat,kineticEnergy);
332 }
333 
334 inline G4double
336  const G4Material* mat,
337  G4double kineticEnergy)
338 {
339  return effCharge.EffectiveChargeSquareRatio(p,mat,kineticEnergy);
340 }
341 
342 inline void G4EmCorrections::SetupKinematics(const G4ParticleDefinition* p,
343  const G4Material* mat,
344  G4double kineticEnergy)
345 {
346  if(kineticEnergy != kinEnergy || p != particle) {
347  particle = p;
348  kinEnergy = kineticEnergy;
349  mass = p->GetPDGMass();
350  tau = kineticEnergy / mass;
351  gamma = 1.0 + tau;
352  bg2 = tau * (tau+2.0);
353  beta2 = bg2/(gamma*gamma);
354  beta = std::sqrt(beta2);
355  ba2 = beta2/alpha2;
356  G4double ratio = CLHEP::electron_mass_c2/mass;
357  tmax = 2.0*CLHEP::electron_mass_c2*bg2 /(1. + 2.0*gamma*ratio + ratio*ratio);
358  charge = p->GetPDGCharge()/CLHEP::eplus;
359  //if(charge < 1.5) {q2 = charge*charge;}
360  //else {
361  // q2 = effCharge.EffectiveChargeSquareRatio(p,mat,kinEnergy);
362  // charge = std::sqrt(q2);
363  //}
364  if(charge > 1.5) { charge = effCharge.EffectiveCharge(p,mat,kinEnergy); }
365  q2 = charge*charge;
366  }
367  if(mat != material) {
368  material = mat;
369  theElementVector = material->GetElementVector();
370  atomDensity = material->GetAtomicNumDensityVector();
371  numberOfElements = material->GetNumberOfElements();
372  }
373 }
374 
375 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
376 
377 #endif
Double_t y2[nxs]
Definition: Style.C:21
G4double HighOrderCorrections(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy, G4double cutEnergy)
G4double DensityCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Double_t y1[nxs]
Definition: Style.C:20
std::vector< G4Element * > G4ElementVector
Double_t x2[nxs]
Definition: Style.C:19
void SetIonisationModels(G4VEmModel *m1=0, G4VEmModel *m2=0)
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
const char * p
Definition: xmltok.h:285
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *p, const G4Material *material, G4double kineticEnergy)
G4double MottCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double EffectiveChargeCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
tuple x
Definition: test.py:50
G4double BlochCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
int G4int
Definition: G4Types.hh:78
G4double SpinCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
void AddStoppingData(G4int Z, G4int A, const G4String &materialName, G4PhysicsVector *dVector)
Float_t mat
Definition: plot.C:40
Double_t y
Definition: plot.C:279
Char_t n[5]
Float_t Z
Definition: plot.C:39
G4double NuclearDEDX(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy, G4bool fluct=true)
bool G4bool
Definition: G4Types.hh:79
G4double IonBarkasCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Double_t x1[nxs]
Definition: Style.C:18
G4double BarkasCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double IonHighOrderCorrections(const G4ParticleDefinition *, const G4MaterialCutsCouple *, G4double kineticEnergy)
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:214
G4double EffectiveCharge(const G4ParticleDefinition *p, const G4Material *material, G4double kineticEnergy)
G4double ShellCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double GetPDGMass() const
G4double KShellCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double LShellCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
G4double ComputeIonCorrections(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
double G4double
Definition: G4Types.hh:76
G4double ShellCorrectionSTD(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double GetPDGCharge() const
G4int GetNumberOfStoppingVectors()
virtual ~G4EmCorrections()
G4double Bethe(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)