Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 100363 2016-10-19 09:24:47Z 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 
64 class G4VEmModel;
65 class G4PhysicsVector;
66 class G4IonTable;
69 class G4Pow;
70 
72 {
73 
74 public:
75 
76  explicit G4EmCorrections(G4int verb);
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 = nullptr,
164  G4VEmModel* m2 = nullptr);
165 
166  inline G4int GetNumberOfStoppingVectors() const;
167 
168  inline void SetVerbose(G4int verb);
169 
170 private:
171 
172  void Initialise();
173 
174  void BuildCorrectionVector();
175 
176  void SetupKinematics(const G4ParticleDefinition*,
177  const G4Material*,
178  G4double kineticEnergy);
179 
180  G4double KShell(G4double theta, G4double eta);
181 
182  G4double LShell(G4double theta, G4double eta);
183 
184  G4int Index(G4double x, const G4double* y, G4int n) const;
185 
186  G4double Value(G4double xv, G4double x1, G4double x2,
187  G4double y1, G4double y2) const;
188 
189  G4double Value2(G4double xv, G4double yv, G4double x1, G4double x2,
190  G4double y1, G4double y2, G4double z11, G4double z21,
191  G4double z12, G4double z22) const;
192 
193  G4double NuclearStoppingPower(G4double e, G4double z1, G4double z2,
194  G4double m1, G4double m2);
195 
196  // hide assignment operator
197  G4EmCorrections & operator=(const G4EmCorrections &right) = delete;
198  G4EmCorrections(const G4EmCorrections&) = delete;
199 
200  G4Pow* g4calc;
201 
202  static const G4double inveplus;
203  static const G4double ZD[11];
204  static const G4double UK[20];
205  static const G4double VK[20];
206  static G4double ZK[20];
207  static const G4double Eta[29];
208  static G4double CK[20][29];
209  static G4double CL[26][28];
210  static const G4double UL[26];
211  static G4double VL[26];
212 
213  static G4LPhysicsFreeVector* BarkasCorr;
214  static G4LPhysicsFreeVector* ThetaK;
215  static G4LPhysicsFreeVector* ThetaL;
216 
217  G4double theZieglerFactor;
218  G4double alpha2;
219  G4bool lossFlucFlag;
220 
221  std::vector<const G4Material*> currmat;
222  std::map< G4int, std::vector<G4double> > thcorr;
223  size_t ncouples;
224 
225  const G4ParticleDefinition* particle;
226  const G4ParticleDefinition* curParticle;
227  const G4Material* material;
228  const G4Material* curMaterial;
229  const G4ElementVector* theElementVector;
230  const G4double* atomDensity;
231 
232  G4PhysicsVector* curVector;
233 
234  G4IonTable* ionTable;
235  G4VEmModel* ionLEModel;
236  G4VEmModel* ionHEModel;
237 
238  G4double kinEnergy;
239  G4double mass;
240  G4double massFactor;
241  G4double eth;
242  G4double tau;
243  G4double gamma;
244  G4double bg2;
245  G4double beta2;
246  G4double beta;
247  G4double ba2;
248  G4double tmax;
249  G4double charge;
250  G4double q2;
251  G4double eCorrMin;
252  G4double eCorrMax;
253 
254  G4int verbose;
255 
256  G4int nK;
257  G4int nL;
258  G4int nEtaK;
259  G4int nEtaL;
260 
261  G4int nbinCorr;
262  G4int numberOfElements;
263 
264  // Ion stopping data
265  G4int nIons;
266  G4int idx;
267  G4int currentZ;
268  std::vector<G4int> Zion;
269  std::vector<G4int> Aion;
270  std::vector<G4String> materialName;
271 
272  std::vector<const G4ParticleDefinition*> ionList;
273 
274  std::vector<const G4Material*> materialList;
275  std::vector<G4PhysicsVector*> stopData;
276 
277  G4bool isMaster;
278  G4ionEffectiveCharge effCharge;
279 };
280 
281 inline G4int
282 G4EmCorrections::Index(G4double x, const G4double* y, G4int n) const
283 {
284  G4int iddd = n-1;
285  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
286  do {--iddd;} while (iddd>0 && x<y[iddd]);
287  return iddd;
288 }
289 
290 inline G4double G4EmCorrections::Value(G4double xv, G4double x1, G4double x2,
291  G4double y1, G4double y2) const
292 {
293  return y1 + (y2 - y1)*(xv - x1)/(x2 - x1);
294 }
295 
296 inline G4double G4EmCorrections::Value2(G4double xv, G4double yv,
297  G4double x1, G4double x2,
298  G4double y1, G4double y2,
299  G4double z11, G4double z21,
300  G4double z12, G4double z22) const
301 {
302  return (z11*(x2-xv)*(y2-yv) + z22*(xv-x1)*(yv-y1) +
303  0.5*(z12*((x2-xv)*(yv-y1)+(xv-x1)*(y2-yv))+
304  z21*((xv-x1)*(y2-yv)+(yv-y1)*(x2-xv))))
305  / ((x2-x1)*(y2-y1));
306 }
307 
308 inline void
310 {
311  if(mod1) { ionLEModel = mod1; }
312  if(mod2) { ionHEModel = mod2; }
313 }
314 
316 {
317  return nIons;
318 }
319 
320 inline G4double
322  const G4Material* mat,
323  G4double kineticEnergy)
324 {
325  return effCharge.EffectiveCharge(p,mat,kineticEnergy);
326 }
327 
328 inline G4double
330  const G4Material* mat,
331  G4double kineticEnergy)
332 {
333  return effCharge.EffectiveChargeSquareRatio(p,mat,kineticEnergy);
334 }
335 
336 inline void G4EmCorrections::SetupKinematics(const G4ParticleDefinition* p,
337  const G4Material* mat,
338  G4double kineticEnergy)
339 {
340  if(kineticEnergy != kinEnergy || p != particle) {
341  particle = p;
342  kinEnergy = kineticEnergy;
343  mass = p->GetPDGMass();
344  tau = kineticEnergy / mass;
345  gamma = 1.0 + tau;
346  bg2 = tau * (tau+2.0);
347  beta2 = bg2/(gamma*gamma);
348  beta = std::sqrt(beta2);
349  ba2 = beta2/alpha2;
350  G4double ratio = CLHEP::electron_mass_c2/mass;
351  tmax = 2.0*CLHEP::electron_mass_c2*bg2
352  /(1. + 2.0*gamma*ratio + ratio*ratio);
353  charge = p->GetPDGCharge()*inveplus;
354  if(charge > 1.5) { charge = effCharge.EffectiveCharge(p,mat,kinEnergy); }
355  q2 = charge*charge;
356  }
357  if(mat != material) {
358  material = mat;
359  theElementVector = material->GetElementVector();
360  atomDensity = material->GetAtomicNumDensityVector();
361  numberOfElements = material->GetNumberOfElements();
362  }
363 }
364 
366 {
367  verbose = verb;
368 }
369 
370 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
371 
372 #endif
G4double HighOrderCorrections(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy, G4double cutEnergy)
G4double DensityCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
std::vector< G4Element * > G4ElementVector
void SetVerbose(G4int verb)
Definition: G4Pow.hh:56
G4int GetNumberOfStoppingVectors() const
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)
G4double BlochCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
int G4int
Definition: G4Types.hh:78
G4double SpinCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
void SetIonisationModels(G4VEmModel *m1=nullptr, G4VEmModel *m2=nullptr)
void AddStoppingData(G4int Z, G4int A, const G4String &materialName, G4PhysicsVector *dVector)
static constexpr double electron_mass_c2
G4EmCorrections(G4int verb)
double A(double temperature)
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)
G4double BarkasCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double IonHighOrderCorrections(const G4ParticleDefinition *, const G4MaterialCutsCouple *, G4double kineticEnergy)
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:216
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:186
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
static constexpr double m2
Definition: G4SIunits.hh:130
virtual ~G4EmCorrections()
G4double Bethe(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)