Geant4  10.01.p03
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 83007 2014-07-24 14:46: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 G4Material;
68 class G4PhysicsTable;
69 class G4VEmModel;
71 class G4VEmProcess;
73 class G4VProcess;
75 class G4Region;
76 class G4Element;
77 class G4EmCorrections;
78 class G4IonTable;
79 
81 {
82 
83 public:
84 
86 
88 
89  //===========================================================================
90  // Methods to access precalculated dE/dx and cross sections
91  // Materials should exist in the list of the G4MaterialCutsCouple
92  //===========================================================================
93 
94  G4double GetDEDX(G4double kinEnergy, const G4ParticleDefinition*,
95  const G4Material*,
96  const G4Region* r = 0);
97  G4double GetDEDX(G4double kinEnergy, const G4String& part,
98  const G4String& mat,
99  const G4String& s = "world");
100 
102  const G4ParticleDefinition*,
103  const G4Material*,
104  const G4Region* r = 0);
105  G4double GetRangeFromRestricteDEDX(G4double kinEnergy, const G4String& part,
106  const G4String& mat,
107  const G4String& s = "world");
108 
110  const G4Material*,
111  const G4Region* r = 0);
112  G4double GetCSDARange(G4double kinEnergy, const G4String& part,
113  const G4String& mat,
114  const G4String& s = "world");
115 
116  G4double GetRange(G4double kinEnergy, const G4ParticleDefinition*,
117  const G4Material*,
118  const G4Region* r = 0);
119  G4double GetRange(G4double kinEnergy, const G4String& part,
120  const G4String& mat,
121  const G4String& s = "world");
122 
124  const G4Material*,
125  const G4Region* r = 0);
126  G4double GetKinEnergy(G4double range, const G4String& part,
127  const G4String& mat,
128  const G4String& s = "world");
129 
131  G4double kinEnergy, const G4ParticleDefinition*,
132  const G4String& processName, const G4Material*,
133  const G4Region* r = 0);
135  G4double kinEnergy, const G4String& part, const G4String& proc,
136  const G4String& mat, const G4String& s = "world");
137 
139  const G4String& part, G4int Z,
141  G4double kinEnergy);
142 
144  const G4String& processName, const G4Material*,
145  const G4Region* r = 0);
146  G4double GetMeanFreePath(G4double kinEnergy, const G4String& part,
147  const G4String& proc,
148  const G4String& mat, const G4String& s = "world");
149 
151 
153 
155 
156  //===========================================================================
157  // Methods to calculate dE/dx and cross sections "on fly"
158  // Existing tables and G4MaterialCutsCouples are not used
159  //===========================================================================
160 
162  const G4String& processName, const G4Material*,
163  G4double cut = DBL_MAX);
164  G4double ComputeDEDX(G4double kinEnergy, const G4String& part,
165  const G4String& proc,
166  const G4String& mat, G4double cut = DBL_MAX);
167 
169  const G4ParticleDefinition*,
170  const G4Material* mat, G4double cut = DBL_MAX);
171  G4double ComputeElectronicDEDX(G4double kinEnergy, const G4String& part,
172  const G4String& mat, G4double cut = DBL_MAX);
173 
175  const G4Material*);
176  G4double ComputeNuclearDEDX(G4double kinEnergy, const G4String& part,
177  const G4String& mat);
178 
180  const G4Material*, G4double cut = DBL_MAX);
181  G4double ComputeTotalDEDX(G4double kinEnergy, const G4String& part,
182  const G4String& mat, G4double cut = DBL_MAX);
183 
185  G4double kinEnergy, const G4ParticleDefinition*,
186  const G4String& processName, const G4Material*,
187  G4double cut = 0.0);
189  G4double kinEnergy, const G4String& part,
190  const G4String& proc,
191  const G4String& mat, G4double cut = 0.0);
192 
194  G4double kinEnergy, const G4ParticleDefinition*,
195  const G4String& processName, G4double Z, G4double A,
196  G4double cut = 0.0);
198  G4double kinEnergy, const G4String& part,
199  const G4String& processName, const G4Element*,
200  G4double cut = 0.0);
201 
203  const G4Material*);
204 
206  const G4String& part, G4int Z,
208  G4double kinEnergy,
209  const G4Material* mat = 0);
210 
212  G4double kinEnergy, const G4ParticleDefinition*,
213  const G4String& processName, const G4Material*,
214  G4double cut = 0.0);
216  G4double kinEnergy, const G4String&, const G4String&,
217  const G4String& processName, G4double cut = 0.0);
218 
220  G4double range, const G4ParticleDefinition*,
221  const G4Material*);
223  G4double range, const G4String&,
224  const G4String&);
225 
226  //===========================================================================
227  // Methods to access particles, materials, regions, processes
228  //===========================================================================
229 
231 
233 
234  const G4Material* FindMaterial(const G4String&);
235 
236  const G4Region* FindRegion(const G4String&);
237 
239  const G4Region* r = 0);
240 
242  const G4String& processName);
243 
244  void SetVerbose(G4int val);
245 
246  //===========================================================================
247  // Private methods
248  //===========================================================================
249 
250 private:
251 
253 
254  G4bool UpdateCouple(const G4Material*, G4double cut);
255 
257  const G4String& processName,
258  G4double kinEnergy);
259 
261  const G4String& processName,
262  G4double kinEnergy);
263 
265 
267  const G4String& processName);
268 
270  const G4String& processName);
271 
273  const G4String& processName);
274 
276  G4VProcess* proc);
277 
280 
281  std::vector<const G4Material*> localMaterials;
282  std::vector<const G4MaterialCutsCouple*> localCouples;
283 
289 
291 
292  // cash
303 
307 
316 
320 };
321 
322 //....oooOO0OOooo.......oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
323 
324 #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 &)
std::vector< const G4MaterialCutsCouple * > localCouples
const G4ParticleDefinition * currentParticle
void PrintRangeTable(const G4ParticleDefinition *)
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 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 double s
Definition: G4SIunits.hh:150
G4double GetCSDARange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
G4double GetShellIonisationCrossSectionPerAtom(const G4String &part, G4int Z, G4AtomicShellEnumerator shell, G4double kinEnergy)
G4VEnergyLossProcess * FindEnergyLossProcess(const G4ParticleDefinition *)
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)
static const G4double A[nN]
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
G4double ComputeMeanFreePath(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, 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
G4AtomicShellEnumerator