Geant4  10.02.p02
G4EmCalculator.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: G4EmCalculator.hh 97616 2016-06-06 12:37:57Z gcosmo $
27 //
28 //
29 // -------------------------------------------------------------------
30 //
31 // GEANT4 Class header file
32 //
33 //
34 // File name: G4EmCalculator
35 //
36 // Author: Vladimir Ivanchenko
37 //
38 // Creation date: 27.06.2004
39 //
40 // Modifications:
41 // 17.11.2004 Change signature of methods, add new methods (V.Ivanchenko)
42 // 11.01.2006 Add GetCSDARange (V.Ivanchenko)
43 // 26.01.2006 Rename GetRange -> GetRangeFromRestricteDEDX (V.Ivanchenko)
44 // 22.03.2006 Add ComputeElectronicDEDX and ComputeTotalDEDX (V.Ivanchenko)
45 // 29.09.2006 Add member loweModel (V.Ivanchenko)
46 // 15.03.2007 Add ComputeEnergyCutFromRangeCut methods (V.Ivanchenko)
47 //
48 // Class Description:
49 //
50 // Provide access to dE/dx and cross sections
51 
52 // -------------------------------------------------------------------
53 //
54 
55 #ifndef G4EmCalculator_h
56 #define G4EmCalculator_h 1
57 
58 #include <vector>
59 #include "globals.hh"
60 #include "G4DataVector.hh"
61 #include "G4DynamicParticle.hh"
62 #include "G4VAtomDeexcitation.hh"
63 
64 class G4LossTableManager;
65 class G4NistManager;
66 class G4Material;
69 class G4PhysicsTable;
70 class G4VEmModel;
72 class G4VEmProcess;
74 class G4VProcess;
76 class G4Region;
77 class G4Element;
78 class G4EmCorrections;
79 class G4EmParameters;
80 class G4IonTable;
81 
83 {
84 
85 public:
86 
88 
90 
91  //===========================================================================
92  // Methods to access precalculated dE/dx and cross sections
93  // Materials should exist in the list of the G4MaterialCutsCouple
94  //===========================================================================
95 
96  G4double GetDEDX(G4double kinEnergy, const G4ParticleDefinition*,
97  const G4Material*,
98  const G4Region* r = 0);
99  inline G4double GetDEDX(G4double kinEnergy, const G4String& part,
100  const G4String& mat,
101  const G4String& s = "world");
102 
104  const G4ParticleDefinition*,
105  const G4Material*,
106  const G4Region* r = 0);
107  inline G4double GetRangeFromRestricteDEDX(G4double kinEnergy,
108  const G4String& part,
109  const G4String& mat,
110  const G4String& s = "world");
111 
113  const G4Material*,
114  const G4Region* r = 0);
115  inline G4double GetCSDARange(G4double kinEnergy, const G4String& part,
116  const G4String& mat,
117  const G4String& s = "world");
118 
119  G4double GetRange(G4double kinEnergy, const G4ParticleDefinition*,
120  const G4Material*,
121  const G4Region* r = 0);
122  inline G4double GetRange(G4double kinEnergy, const G4String& part,
123  const G4String& mat,
124  const G4String& s = "world");
125 
127  const G4Material*,
128  const G4Region* r = 0);
129  inline G4double GetKinEnergy(G4double range, const G4String& part,
130  const G4String& mat,
131  const G4String& s = "world");
132 
134  G4double kinEnergy, const G4ParticleDefinition*,
135  const G4String& processName, const G4Material*,
136  const G4Region* r = 0);
138  G4double kinEnergy, const G4String& part, const G4String& proc,
139  const G4String& mat, const G4String& s = "world");
140 
142  const G4String& part, G4int Z,
144  G4double kinEnergy);
145 
147  const G4String& processName, const G4Material*,
148  const G4Region* r = 0);
149  inline G4double GetMeanFreePath(G4double kinEnergy, const G4String& part,
150  const G4String& proc, const G4String& mat,
151  const G4String& s = "world");
152 
154 
156 
158 
159  //===========================================================================
160  // Methods to calculate dE/dx and cross sections "on fly"
161  // Existing tables and G4MaterialCutsCouples are not used
162  //===========================================================================
163 
165  const G4String& processName, const G4Material*,
166  G4double cut = DBL_MAX);
167  inline G4double ComputeDEDX(G4double kinEnergy, const G4String& part,
168  const G4String& proc,
169  const G4String& mat, G4double cut = DBL_MAX);
170 
172  const G4ParticleDefinition*,
173  const G4Material* mat, G4double cut = DBL_MAX);
174  inline G4double ComputeElectronicDEDX(G4double kinEnergy, const G4String& part,
175  const G4String& mat, G4double cut = DBL_MAX);
176 
178  const G4ParticleDefinition*,
179  const G4Material* mat, G4double rangecut = DBL_MAX);
180  inline G4double ComputeDEDXForCutInRange(G4double kinEnergy, const G4String& part,
181  const G4String& mat,
182  G4double rangecut = DBL_MAX);
183 
185  const G4Material*);
186  inline G4double ComputeNuclearDEDX(G4double kinEnergy, const G4String& part,
187  const G4String& mat);
188 
190  const G4Material*, G4double cut = DBL_MAX);
191  inline G4double ComputeTotalDEDX(G4double kinEnergy, const G4String& part,
192  const G4String& mat, G4double cut = DBL_MAX);
193 
195  G4double kinEnergy, const G4ParticleDefinition*,
196  const G4String& processName, const G4Material*,
197  G4double cut = 0.0);
199  G4double kinEnergy, const G4String& part,
200  const G4String& proc,
201  const G4String& mat, G4double cut = 0.0);
202 
204  G4double kinEnergy, const G4ParticleDefinition*,
205  const G4String& processName, G4double Z, G4double A,
206  G4double cut = 0.0);
208  G4double kinEnergy, const G4String& part,
209  const G4String& processName, const G4Element*,
210  G4double cut = 0.0);
211 
213  G4double kinEnergy, const G4ParticleDefinition*,
214  const G4String& processName, G4int Z, G4int shellIdx,
215  G4double cut = 0.0);
217  G4double kinEnergy, const G4String& part,
218  const G4String& processName, const G4Element*,
219  G4int shellIdx,
220  G4double cut = 0.0);
221 
223  const G4Material*);
224 
226  const G4String& part, G4int Z,
228  G4double kinEnergy,
229  const G4Material* mat = 0);
230 
232  G4double kinEnergy, const G4ParticleDefinition*,
233  const G4String& processName, const G4Material*,
234  G4double cut = 0.0);
236  G4double kinEnergy, const G4String&, const G4String&,
237  const G4String& processName, G4double cut = 0.0);
238 
240  G4double range, const G4ParticleDefinition*,
241  const G4Material*);
243  G4double range, const G4String&,
244  const G4String&);
245 
246  //===========================================================================
247  // Methods to access particles, materials, regions, processes
248  //===========================================================================
249 
251 
253 
254  const G4Material* FindMaterial(const G4String&);
255 
256  const G4Region* FindRegion(const G4String&);
257 
259  const G4Region* r = 0);
260 
262  const G4String& processName);
263 
264  void SetupMaterial(const G4Material*);
265 
266  void SetupMaterial(const G4String&);
267 
268  void SetVerbose(G4int val);
269 
270  //===========================================================================
271  // Private methods
272  //===========================================================================
273 
274 private:
275 
277 
278  G4bool UpdateCouple(const G4Material*, G4double cut);
279 
281  const G4String& processName,
282  G4double kinEnergy);
283 
285  const G4String& processName,
286  G4double kinEnergy);
287 
289 
291  const G4String& processName);
292 
294  const G4String& processName);
295 
297  const G4String& processName);
298 
300  G4VProcess* proc);
301 
302  void CheckMaterial(G4int Z);
303 
304  // hide copy and assign
307 
308  std::vector<const G4Material*> localMaterials;
309  std::vector<const G4MaterialCutsCouple*> localCouples;
310 
318 
320 
321  // cache
333 
337 
347 
351 };
352 
353 //....oooOO0OOooo.......oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
354 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
355 
356 inline
358  const G4String& material, const G4String& reg)
359 {
360  return GetDEDX(kinEnergy,FindParticle(particle),
361  FindMaterial(material),FindRegion(reg));
362 }
363 
364 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
365 
366 inline
368  const G4String& particle,
369  const G4String& material,
370  const G4String& reg)
371 {
372  return GetRangeFromRestricteDEDX(kinEnergy,FindParticle(particle),
373  FindMaterial(material),FindRegion(reg));
374 }
375 
376 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
377 
378 inline
380  const G4String& particle,
381  const G4String& material,
382  const G4String& reg)
383 {
384  return GetCSDARange(kinEnergy,FindParticle(particle),
385  FindMaterial(material),FindRegion(reg));
386 }
387 
388 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
389 
390 inline
392  const G4String& particle,
393  const G4String& material,
394  const G4String& reg)
395 {
396  return GetRange(kinEnergy,FindParticle(particle),
397  FindMaterial(material),FindRegion(reg));
398 }
399 
400 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
401 
402 inline
404  const G4String& material, const G4String& reg)
405 {
406  return GetKinEnergy(range,FindParticle(particle),
407  FindMaterial(material),FindRegion(reg));
408 }
409 
410 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
411 
412 inline
414  const G4String& particle,
415  const G4String& processName,
416  const G4String& material,
417  const G4String& reg)
418 {
419  return GetCrossSectionPerVolume(kinEnergy,FindParticle(particle),processName,
420  FindMaterial(material),FindRegion(reg));
421 }
422 
423 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
424 
425 inline
427  const G4String& particle,
428  const G4String& processName,
429  const G4String& material,
430  const G4String& reg)
431 {
432  return GetMeanFreePath(kinEnergy,FindParticle(particle),processName,
433  FindMaterial(material),FindRegion(reg));
434 }
435 
436 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
437 
438 inline G4double
440  const G4String& mat, G4double cut)
441 {
442  return
443  ComputeElectronicDEDX(kinEnergy,FindParticle(part),FindMaterial(mat),cut);
444 }
445 
446 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
447 
448 inline G4double
450  const G4String& part,
451  const G4String& mat,
452  G4double rangecut)
453 {
454  return ComputeDEDXForCutInRange(kinEnergy,FindParticle(part),
455  FindMaterial(mat), rangecut);
456 }
457 
458 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
459 
460 inline
462  const G4String& part,
463  const G4String& mat,
464  G4double cut)
465 {
466  return ComputeTotalDEDX(kinEnergy,FindParticle(part),FindMaterial(mat),cut);
467 }
468 
469 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
470 
471 inline
473  const G4String& particle,
474  const G4String& processName,
475  const G4String& material,
476  G4double cut)
477 {
478  return ComputeDEDX(kinEnergy,FindParticle(particle),processName,
479  FindMaterial(material),cut);
480 }
481 
482 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
483 
484 inline
486  const G4String& particle,
487  const G4String& material)
488 {
489  return ComputeNuclearDEDX(kinEnergy,FindParticle(particle),
490  FindMaterial(material));
491 }
492 
493 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
494 
495 inline
497  G4double kinEnergy,
498  const G4String& particle,
499  const G4String& processName,
500  const G4String& material,
501  G4double cut)
502 {
503  return ComputeCrossSectionPerVolume(kinEnergy,FindParticle(particle),
504  processName,
505  FindMaterial(material),cut);
506 }
507 
508 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
509 
510 inline
512  const G4String& particle,
513  const G4String& processName,
514  const G4Element* elm,
515  G4double cut)
516 {
517  return ComputeCrossSectionPerAtom(kinEnergy,FindParticle(particle),
518  processName,
519  elm->GetZ(),elm->GetN(),cut);
520 }
521 
522 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
523 
525  G4double kinEnergy, const G4String& part,
526  const G4String& processName, const G4Element* elm,
527  G4int shellIdx, G4double cut)
528 {
529  return ComputeCrossSectionPerShell(kinEnergy, FindParticle(part),
530  processName, G4lrint(elm->GetZ()),
531  shellIdx, cut);
532 }
533 
534 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
535 
536 inline
538  G4double range,
539  const G4String& particle,
540  const G4String& material)
541 {
542  return ComputeEnergyCutFromRangeCut(range,FindParticle(particle),
543  FindMaterial(material));
544 }
545 
546 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
547 
548 inline
550  const G4String& particle,
551  const G4String& processName,
552  const G4String& material,
553  G4double cut)
554 {
555  return ComputeMeanFreePath(kinEnergy,FindParticle(particle),processName,
556  FindMaterial(material),cut);
557 }
558 
559 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
560 
561 #endif
G4VEmModel * loweModel
G4VEmProcess * FindDiscreteProcess(const G4ParticleDefinition *, const G4String &processName)
const G4PhysicsTable * currentLambda
void FindLambdaTable(const G4ParticleDefinition *, const G4String &processName, G4double kinEnergy)
const G4ParticleDefinition * theGenericIon
const G4Material * FindMaterial(const G4String &)
const G4Material * cutMaterial
G4double GetN() const
Definition: G4Element.hh:134
std::vector< const G4MaterialCutsCouple * > localCouples
const G4ParticleDefinition * currentParticle
void PrintRangeTable(const G4ParticleDefinition *)
G4double GetZ() const
Definition: G4Element.hh:131
G4bool UpdateCouple(const G4Material *, G4double cut)
G4double ComputeTotalDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, G4double cut=DBL_MAX)
G4double ComputeElectronicDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *mat, G4double cut=DBL_MAX)
G4EmCorrections * corr
G4ionEffectiveCharge * ionEffCharge
const G4MaterialCutsCouple * currentCouple
G4double GetCrossSectionPerVolume(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, const G4Region *r=0)
G4double GetMeanFreePath(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, const G4Region *r=0)
void SetupMaterial(const G4Material *)
void PrintInverseRangeTable(const G4ParticleDefinition *)
G4double GetRange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
int G4int
Definition: G4Types.hh:78
G4DynamicParticle dynParticle
G4double ComputeShellIonisationCrossSectionPerAtom(const G4String &part, G4int Z, G4AtomicShellEnumerator shell, G4double kinEnergy, const G4Material *mat=0)
G4double GetDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
static const G4double reg
static const double s
Definition: G4SIunits.hh:168
G4double GetCSDARange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
G4double cutenergy[3]
G4double GetShellIonisationCrossSectionPerAtom(const G4String &part, G4int Z, G4AtomicShellEnumerator shell, G4double kinEnergy)
G4VEnergyLossProcess * FindEnergyLossProcess(const G4ParticleDefinition *)
double A(double temperature)
G4VEnergyLossProcess * currentProcess
const G4MaterialCutsCouple * FindCouple(const G4Material *, const G4Region *r=0)
G4double ComputeDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=DBL_MAX)
G4String currentName
const G4ParticleDefinition * baseParticle
bool G4bool
Definition: G4Types.hh:79
G4bool FindEmModel(const G4ParticleDefinition *, const G4String &processName, G4double kinEnergy)
G4double ComputeNuclearDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *)
G4double GetRangeFromRestricteDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
G4DataVector localCuts
G4EmCalculator & operator=(const G4EmCalculator &right)
void CheckMaterial(G4int Z)
G4bool ActiveForParticle(const G4ParticleDefinition *part, G4VProcess *proc)
G4VMultipleScattering * FindMscProcess(const G4ParticleDefinition *, const G4String &processName)
std::vector< const G4Material * > localMaterials
G4VEnergyLossProcess * FindEnLossProcess(const G4ParticleDefinition *, const G4String &processName)
const G4ParticleDefinition * FindParticle(const G4String &)
G4double chargeSquare
int G4lrint(double ad)
Definition: templates.hh:163
G4double ComputeMeanFreePath(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=0.0)
G4double ComputeCrossSectionPerShell(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, G4int Z, G4int shellIdx, G4double cut=0.0)
G4double GetKinEnergy(G4double range, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
const G4Region * FindRegion(const G4String &)
G4IonTable * ionTable
G4String currentMaterialName
const G4ParticleDefinition * FindIon(G4int Z, G4int A)
G4String currentProcessName
G4String currentParticleName
G4double ComputeCrossSectionPerVolume(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=0.0)
G4double ComputeGammaAttenuationLength(G4double kinEnergy, const G4Material *)
G4double ComputeCrossSectionPerAtom(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, G4double Z, G4double A, G4double cut=0.0)
G4VProcess * FindProcess(const G4ParticleDefinition *part, const G4String &processName)
void SetVerbose(G4int val)
void PrintDEDXTable(const G4ParticleDefinition *)
const G4ParticleDefinition * lambdaParticle
double G4double
Definition: G4Types.hh:76
G4double ComputeEnergyCutFromRangeCut(G4double range, const G4ParticleDefinition *, const G4Material *)
const G4Material * currentMaterial
G4VEmModel * currentModel
G4bool UpdateParticle(const G4ParticleDefinition *, G4double kinEnergy)
#define DBL_MAX
Definition: templates.hh:83
G4LossTableManager * manager
G4double ComputeDEDXForCutInRange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *mat, G4double rangecut=DBL_MAX)
G4AtomicShellEnumerator
G4EmParameters * theParameters
G4NistManager * nist