Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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$
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"
63
65class G4Material;
68class G4PhysicsTable;
69class G4VEmModel;
71class G4VEmProcess;
73class G4VProcess;
75class G4Region;
76class G4Element;
77class G4EmCorrections;
78
80{
81
82public:
83
85
87
88 //===========================================================================
89 // Methods to access precalculated dE/dx and cross sections
90 // Materials should exist in the list of the G4MaterialCutsCouple
91 //===========================================================================
92
94 const G4Material*,
95 const G4Region* r = 0);
96 G4double GetDEDX(G4double kinEnergy, const G4String& part,
97 const G4String& mat,
98 const G4String& s = "world");
99
101 const G4ParticleDefinition*,
102 const G4Material*,
103 const G4Region* r = 0);
104 G4double GetRangeFromRestricteDEDX(G4double kinEnergy, const G4String& part,
105 const G4String& mat,
106 const G4String& s = "world");
107
109 const G4Material*,
110 const G4Region* r = 0);
111 G4double GetCSDARange(G4double kinEnergy, const G4String& part,
112 const G4String& mat,
113 const G4String& s = "world");
114
116 const G4Material*,
117 const G4Region* r = 0);
118 G4double GetRange(G4double kinEnergy, const G4String& part,
119 const G4String& mat,
120 const G4String& s = "world");
121
123 const G4Material*,
124 const G4Region* r = 0);
125 G4double GetKinEnergy(G4double range, const G4String& part,
126 const G4String& mat,
127 const G4String& s = "world");
128
130 G4double kinEnergy, const G4ParticleDefinition*,
131 const G4String& processName, const G4Material*,
132 const G4Region* r = 0);
134 G4double kinEnergy, const G4String& part, const G4String& proc,
135 const G4String& mat, const G4String& s = "world");
136
138 const G4String& part, G4int Z,
140 G4double kinEnergy);
141
143 const G4String& processName, const G4Material*,
144 const G4Region* r = 0);
145 G4double GetMeanFreePath(G4double kinEnergy, const G4String& part,
146 const G4String& proc,
147 const G4String& mat, const G4String& s = "world");
148
150
152
154
155 //===========================================================================
156 // Methods to calculate dE/dx and cross sections "on fly"
157 // Existing tables and G4MaterialCutsCouples are not used
158 //===========================================================================
159
161 const G4String& processName, const G4Material*,
162 G4double cut = DBL_MAX);
163 G4double ComputeDEDX(G4double kinEnergy, const G4String& part,
164 const G4String& proc,
165 const G4String& mat, G4double cut = DBL_MAX);
166
169 const G4Material* mat, G4double cut = DBL_MAX);
170 G4double ComputeElectronicDEDX(G4double kinEnergy, const G4String& part,
171 const G4String& mat, G4double cut = DBL_MAX);
172
174 const G4Material*);
175 G4double ComputeNuclearDEDX(G4double kinEnergy, const G4String& part,
176 const G4String& mat);
177
179 const G4Material*, G4double cut = DBL_MAX);
180 G4double ComputeTotalDEDX(G4double kinEnergy, const G4String& part,
181 const G4String& mat, G4double cut = DBL_MAX);
182
184 G4double kinEnergy, const G4ParticleDefinition*,
185 const G4String& processName, const G4Material*,
186 G4double cut = 0.0);
188 G4double kinEnergy, const G4String& part,
189 const G4String& proc,
190 const G4String& mat, G4double cut = 0.0);
191
193 G4double kinEnergy, const G4ParticleDefinition*,
194 const G4String& processName, G4double Z, G4double A,
195 G4double cut = 0.0);
197 G4double kinEnergy, const G4String& part,
198 const G4String& processName, const G4Element*,
199 G4double cut = 0.0);
200
202 const G4Material*);
203
205 const G4String& part, G4int Z,
207 G4double kinEnergy,
208 const G4Material* mat = 0);
209
211 G4double kinEnergy, const G4ParticleDefinition*,
212 const G4String& processName, const G4Material*,
213 G4double cut = 0.0);
215 G4double kinEnergy, const G4String&, const G4String&,
216 const G4String& processName, G4double cut = 0.0);
217
219 G4double range, const G4ParticleDefinition*,
220 const G4Material*);
222 G4double range, const G4String&,
223 const G4String&);
224
225 //===========================================================================
226 // Methods to access particles, materials, regions
227 //===========================================================================
228
230
232
233 const G4Material* FindMaterial(const G4String&);
234
235 const G4Region* FindRegion(const G4String&);
236
238 const G4Region* r = 0);
239
240 void SetVerbose(G4int val);
241
242 //===========================================================================
243 // Private methods
244 //===========================================================================
245
246private:
247
248 G4bool UpdateParticle(const G4ParticleDefinition*, G4double kinEnergy);
249
250 G4bool UpdateCouple(const G4Material*, G4double cut);
251
252 void FindLambdaTable(const G4ParticleDefinition*,
253 const G4String& processName,
254 G4double kinEnergy);
255
256 G4bool FindEmModel(const G4ParticleDefinition*,
257 const G4String& processName,
258 G4double kinEnergy);
259
260 G4VEnergyLossProcess* FindEnergyLossProcess(const G4ParticleDefinition*);
261
262 G4VEnergyLossProcess* FindEnLossProcess(const G4ParticleDefinition*,
263 const G4String& processName);
264
265 G4VEmProcess* FindDiscreteProcess(const G4ParticleDefinition*,
266 const G4String& processName);
267
268 G4VMultipleScattering* FindMscProcess(const G4ParticleDefinition*,
269 const G4String& processName);
270
271 G4bool ActiveForParticle(const G4ParticleDefinition* part,
272 G4VProcess* proc);
273
274 G4EmCalculator & operator=(const G4EmCalculator &right);
276
277 std::vector<const G4Material*> localMaterials;
278 std::vector<const G4MaterialCutsCouple*> localCouples;
279
280 G4LossTableManager* manager;
281 G4EmCorrections* corr;
282 G4DataVector localCuts;
283 G4int nLocalMaterials;
284
285 G4int verbose;
286
287 // cash
288 G4int currentCoupleIndex;
289 const G4MaterialCutsCouple* currentCouple;
290 const G4Material* currentMaterial;
291 const G4ParticleDefinition* currentParticle;
292 const G4ParticleDefinition* lambdaParticle;
293 const G4ParticleDefinition* baseParticle;
294 const G4PhysicsTable* currentLambda;
295 G4VEmModel* currentModel;
296 G4VEmModel* loweModel;
297 G4VEnergyLossProcess* currentProcess;
298
299 const G4ParticleDefinition* theGenericIon;
300 G4ionEffectiveCharge* ionEffCharge;
301 G4DynamicParticle dynParticle;
302
303 G4String currentName;
304 G4String lambdaName;
305 G4double currentCut;
306 G4double chargeSquare;
307 G4double massRatio;
308 G4double mass;
309 G4bool isIon;
310 G4bool isApplicable;
311
312 G4String currentParticleName;
313 G4String currentMaterialName;
314 G4String currentProcessName;
315};
316
317//....oooOO0OOooo.......oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
318
319#endif
G4AtomicShellEnumerator
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
const G4ParticleDefinition * FindParticle(const G4String &)
G4double ComputeMeanFreePath(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=0.0)
G4double GetShellIonisationCrossSectionPerAtom(const G4String &part, G4int Z, G4AtomicShellEnumerator shell, G4double kinEnergy)
G4double GetRangeFromRestricteDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
G4double ComputeTotalDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, G4double cut=DBL_MAX)
G4double ComputeNuclearDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *)
const G4ParticleDefinition * FindIon(G4int Z, G4int A)
G4double ComputeGammaAttenuationLength(G4double kinEnergy, const G4Material *)
G4double ComputeDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=DBL_MAX)
void SetVerbose(G4int val)
G4double ComputeEnergyCutFromRangeCut(G4double range, const G4ParticleDefinition *, const G4Material *)
const G4MaterialCutsCouple * FindCouple(const G4Material *, const G4Region *r=0)
G4double ComputeCrossSectionPerVolume(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=0.0)
G4double ComputeCrossSectionPerAtom(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, G4double Z, G4double A, G4double cut=0.0)
G4double GetCrossSectionPerVolume(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, const G4Region *r=0)
void PrintDEDXTable(const G4ParticleDefinition *)
const G4Region * FindRegion(const G4String &)
G4double GetCSDARange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
G4double GetKinEnergy(G4double range, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
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)
void PrintRangeTable(const G4ParticleDefinition *)
void PrintInverseRangeTable(const G4ParticleDefinition *)
G4double GetMeanFreePath(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, const G4Region *r=0)
G4double GetRange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
G4double ComputeElectronicDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *mat, G4double cut=DBL_MAX)
const G4Material * FindMaterial(const G4String &)
#define DBL_MAX
Definition: templates.hh:83