Geant4-11
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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class header file
30//
31//
32// File name: G4EmCalculator
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 27.06.2004
37//
38// Modifications:
39// 17.11.2004 Change signature of methods, add new methods (V.Ivanchenko)
40// 11.01.2006 Add GetCSDARange (V.Ivanchenko)
41// 26.01.2006 Rename GetRange -> GetRangeFromRestricteDEDX (V.Ivanchenko)
42// 22.03.2006 Add ComputeElectronicDEDX and ComputeTotalDEDX (V.Ivanchenko)
43// 29.09.2006 Add member loweModel (V.Ivanchenko)
44// 15.03.2007 Add ComputeEnergyCutFromRangeCut methods (V.Ivanchenko)
45// 02.02.2018 Add parameter to FindLambdaTable to store process type (M. Novak)
46//
47// Class Description:
48//
49// Provide access to dE/dx and cross sections
50
51// -------------------------------------------------------------------
52//
53
54#ifndef G4EmCalculator_h
55#define G4EmCalculator_h 1
56
57#include <vector>
58#include "globals.hh"
60
62class G4NistManager;
63class G4Material;
66class G4PhysicsTable;
67class G4VEmModel;
69class G4VEmProcess;
71class G4VProcess;
73class G4Region;
74class G4Element;
75class G4EmCorrections;
76class G4EmParameters;
77class G4IonTable;
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 = nullptr);
96 inline 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 = nullptr);
105 const G4String& part,
106 const G4String& mat,
107 const G4String& s = "world");
108
110 const G4Material*,
111 const G4Region* r = nullptr);
112 inline G4double GetCSDARange(G4double kinEnergy, const G4String& part,
113 const G4String& mat,
114 const G4String& s = "world");
115
117 const G4Material*,
118 const G4Region* r = nullptr);
119 inline G4double GetRange(G4double kinEnergy, const G4String& part,
120 const G4String& mat,
121 const G4String& s = "world");
122
124 const G4Material*,
125 const G4Region* r = nullptr);
126 inline 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 = nullptr);
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 = nullptr);
146 inline G4double GetMeanFreePath(G4double kinEnergy, const G4String& part,
147 const G4String& proc, const G4String& mat,
148 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 inline G4double ComputeDEDX(G4double kinEnergy, const G4String& part,
165 const G4String& proc,
166 const G4String& mat, G4double cut = DBL_MAX);
167
170 const G4Material* mat, G4double cut = DBL_MAX);
171 inline G4double ComputeElectronicDEDX(G4double kinEnergy, const G4String& part,
172 const G4String& mat, G4double cut = DBL_MAX);
173
176 const G4Material* mat, G4double rangecut = DBL_MAX);
177 inline G4double ComputeDEDXForCutInRange(G4double kinEnergy, const G4String& part,
178 const G4String& mat,
179 G4double rangecut = DBL_MAX);
180
182 const G4Material*);
183 inline G4double ComputeNuclearDEDX(G4double kinEnergy, const G4String& part,
184 const G4String& mat);
185
187 const G4Material*, G4double cut = DBL_MAX);
188 inline G4double ComputeTotalDEDX(G4double kinEnergy, const G4String& part,
189 const G4String& mat, G4double cut = DBL_MAX);
190
192 G4double kinEnergy, const G4ParticleDefinition*,
193 const G4String& processName, const G4Material*,
194 G4double cut = 0.0);
196 G4double kinEnergy, const G4String& part,
197 const G4String& proc,
198 const G4String& mat, G4double cut = 0.0);
199
201 G4double kinEnergy, const G4ParticleDefinition*,
202 const G4String& processName, G4double Z, G4double A,
203 G4double cut = 0.0);
205 G4double kinEnergy, const G4String& part,
206 const G4String& processName, const G4Element*,
207 G4double cut = 0.0);
208
210 G4double kinEnergy, const G4ParticleDefinition*,
211 const G4String& processName, G4int Z, G4int shellIdx,
212 G4double cut = 0.0);
214 G4double kinEnergy, const G4String& part,
215 const G4String& processName, const G4Element*,
216 G4int shellIdx,
217 G4double cut = 0.0);
218
220 const G4Material*);
221
223 const G4String& part, G4int Z,
225 G4double kinEnergy,
226 const G4Material* mat = nullptr);
227
229 G4double kinEnergy, const G4ParticleDefinition*,
230 const G4String& processName, const G4Material*,
231 G4double cut = 0.0);
233 G4double kinEnergy, const G4String&, const G4String&,
234 const G4String& processName, G4double cut = 0.0);
235
237 G4double range, const G4ParticleDefinition*,
238 const G4Material*);
240 G4double range, const G4String&,
241 const G4String&);
242
243 //===========================================================================
244 // Methods to access particles, materials, regions, processes
245 //===========================================================================
246
248
250
251 const G4Material* FindMaterial(const G4String&);
252
253 const G4Region* FindRegion(const G4String&);
254
256 const G4Region* r = nullptr);
257
259 const G4String& processName);
260
261 void SetupMaterial(const G4Material*);
262
263 void SetupMaterial(const G4String&);
264
265 void SetVerbose(G4int val);
266
267 // hide copy and assign
268 G4EmCalculator & operator=(const G4EmCalculator &right) = delete;
270
271private:
272
274
276
278 const G4String& processName,
279 G4double kinEnergy, G4int& proctype);
280
282 const G4String& processName,
283 G4double kinEnergy);
284
286
288 const G4String& processName);
289
291 const G4String& processName);
292
294 const G4String& processName);
295
297 G4VProcess* proc);
298
299 void CheckMaterial(G4int Z);
300
306
307 // cache
309 const G4Material* currentMaterial = nullptr;
310 const G4Material* cutMaterial = nullptr;
315
321
324
330
334
335 G4bool isIon = false;
337
338 std::vector<const G4Material*> localMaterials;
339 std::vector<const G4MaterialCutsCouple*> localCouples;
340 std::vector<G4double> localCuts;
341
347};
348
349//....oooOO0OOooo.......oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
350//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
351
352inline
354 const G4String& material, const G4String& reg)
355{
356 return GetDEDX(kinEnergy,FindParticle(particle),
358}
359
360//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
361
362inline
364 const G4String& particle,
365 const G4String& material,
366 const G4String& reg)
367{
368 return GetRangeFromRestricteDEDX(kinEnergy,FindParticle(particle),
370}
371
372//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
373
374inline
376 const G4String& particle,
377 const G4String& material,
378 const G4String& reg)
379{
380 return GetCSDARange(kinEnergy,FindParticle(particle),
382}
383
384//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
385
386inline
388 const G4String& particle,
389 const G4String& material,
390 const G4String& reg)
391{
392 return GetRange(kinEnergy,FindParticle(particle),
394}
395
396//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
397
398inline
400 const G4String& material, const G4String& reg)
401{
402 return GetKinEnergy(range,FindParticle(particle),
404}
405
406//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
407
408inline
410 const G4String& particle,
411 const G4String& processName,
412 const G4String& material,
413 const G4String& reg)
414{
415 return GetCrossSectionPerVolume(kinEnergy,FindParticle(particle),processName,
417}
418
419//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
420
421inline
423 const G4String& particle,
424 const G4String& processName,
425 const G4String& material,
426 const G4String& reg)
427{
428 return GetMeanFreePath(kinEnergy,FindParticle(particle),processName,
430}
431
432//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
433
434inline G4double
436 const G4String& mat, G4double cut)
437{
438 return
439 ComputeElectronicDEDX(kinEnergy,FindParticle(part),FindMaterial(mat),cut);
440}
441
442//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
443
444inline G4double
446 const G4String& part,
447 const G4String& mat,
448 G4double rangecut)
449{
450 return ComputeDEDXForCutInRange(kinEnergy,FindParticle(part),
451 FindMaterial(mat), rangecut);
452}
453
454//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
455
456inline
458 const G4String& part,
459 const G4String& mat,
460 G4double cut)
461{
462 return ComputeTotalDEDX(kinEnergy,FindParticle(part),FindMaterial(mat),cut);
463}
464
465//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
466
467inline
469 const G4String& particle,
470 const G4String& processName,
471 const G4String& material,
472 G4double cut)
473{
474 return ComputeDEDX(kinEnergy,FindParticle(particle),processName,
476}
477
478//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
479
480inline
482 const G4String& particle,
483 const G4String& material)
484{
485 return ComputeNuclearDEDX(kinEnergy,FindParticle(particle),
487}
488
489//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
490
491inline
493 G4double kinEnergy,
494 const G4String& particle,
495 const G4String& processName,
496 const G4String& material,
497 G4double cut)
498{
499 return ComputeCrossSectionPerVolume(kinEnergy,FindParticle(particle),
500 processName,
502}
503
504//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
505
506inline
508 const G4String& particle,
509 const G4String& processName,
510 const G4Element* elm,
511 G4double cut)
512{
513 return ComputeCrossSectionPerAtom(kinEnergy,FindParticle(particle),
514 processName,
515 elm->GetZ(),elm->GetN(),cut);
516}
517
518//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
519
521 G4double kinEnergy, const G4String& part,
522 const G4String& processName, const G4Element* elm,
523 G4int shellIdx, G4double cut)
524{
525 return ComputeCrossSectionPerShell(kinEnergy, FindParticle(part),
526 processName, elm->GetZasInt(),
527 shellIdx, cut);
528}
529
530//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
531
532inline
534 G4double range,
535 const G4String& particle,
536 const G4String& material)
537{
538 return ComputeEnergyCutFromRangeCut(range,FindParticle(particle),
540}
541
542//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
543
544inline
546 const G4String& particle,
547 const G4String& processName,
548 const G4String& material,
549 G4double cut)
550{
551 return ComputeMeanFreePath(kinEnergy,FindParticle(particle),processName,
553}
554
555//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
556
557#endif
G4AtomicShellEnumerator
static const G4double reg
static constexpr double s
Definition: G4SIunits.hh:154
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
G4double GetZ() const
Definition: G4Element.hh:131
G4int GetZasInt() const
Definition: G4Element.hh:132
G4double GetN() const
Definition: G4Element.hh:135
const G4ParticleDefinition * FindParticle(const G4String &)
G4double ComputeMeanFreePath(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=0.0)
G4String currentName
G4VEmModel * loweModel
G4double GetShellIonisationCrossSectionPerAtom(const G4String &part, G4int Z, G4AtomicShellEnumerator shell, G4double kinEnergy)
G4double GetDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
std::vector< const G4Material * > localMaterials
G4EmCalculator(const G4EmCalculator &)=delete
G4double ComputeTotalDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, G4double cut=DBL_MAX)
G4double GetCSDARange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
const G4ParticleDefinition * theGenericIon
G4double cutenergy[3]
G4double GetRange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
G4EmCalculator & operator=(const G4EmCalculator &right)=delete
G4VEnergyLossProcess * FindEnergyLossProcess(const G4ParticleDefinition *)
G4double ComputeNuclearDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *)
const G4ParticleDefinition * currentParticle
const G4ParticleDefinition * FindIon(G4int Z, G4int A)
G4LossTableManager * manager
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)
G4VProcess * curProcess
G4EmParameters * theParameters
G4bool UpdateParticle(const G4ParticleDefinition *, G4double kinEnergy)
G4double ComputeEnergyCutFromRangeCut(G4double range, const G4ParticleDefinition *, const G4Material *)
G4double GetKinEnergy(G4double range, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
G4EmCorrections * corr
G4double GetRangeFromRestricteDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
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)
G4VEnergyLossProcess * FindEnLossProcess(const G4ParticleDefinition *, const G4String &processName)
G4String currentProcessName
G4bool FindEmModel(const G4ParticleDefinition *, const G4String &processName, G4double kinEnergy)
G4VEmModel * currentModel
G4VMultipleScattering * FindMscProcess(const G4ParticleDefinition *, const G4String &processName)
void FindLambdaTable(const G4ParticleDefinition *, const G4String &processName, G4double kinEnergy, G4int &proctype)
void PrintDEDXTable(const G4ParticleDefinition *)
const G4Region * FindRegion(const G4String &)
G4VProcess * FindProcess(const G4ParticleDefinition *part, const G4String &processName)
G4String currentMaterialName
G4double GetCrossSectionPerVolume(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, const G4Region *r=nullptr)
G4double ComputeCrossSectionPerShell(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, G4int Z, G4int shellIdx, G4double cut=0.0)
G4bool ActiveForParticle(const G4ParticleDefinition *part, G4VProcess *proc)
G4VEnergyLossProcess * currentProcess
G4DynamicParticle * dynParticle
const G4ParticleDefinition * lambdaParticle
const G4Material * currentMaterial
G4bool UpdateCouple(const G4Material *, G4double cut)
G4double ComputeShellIonisationCrossSectionPerAtom(const G4String &part, G4int Z, G4AtomicShellEnumerator shell, G4double kinEnergy, const G4Material *mat=nullptr)
const G4MaterialCutsCouple * FindCouple(const G4Material *, const G4Region *r=nullptr)
void SetupMaterial(const G4Material *)
void PrintRangeTable(const G4ParticleDefinition *)
G4String currentParticleName
void PrintInverseRangeTable(const G4ParticleDefinition *)
G4IonTable * ionTable
G4VEmProcess * FindDiscreteProcess(const G4ParticleDefinition *, const G4String &processName)
void CheckMaterial(G4int Z)
const G4MaterialCutsCouple * currentCouple
std::vector< const G4MaterialCutsCouple * > localCouples
G4double ComputeDEDXForCutInRange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *mat, G4double rangecut=DBL_MAX)
G4double GetMeanFreePath(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, const G4Region *r=nullptr)
const G4PhysicsTable * currentLambda
G4double chargeSquare
G4ionEffectiveCharge * ionEffCharge
G4double ComputeElectronicDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *mat, G4double cut=DBL_MAX)
const G4Material * FindMaterial(const G4String &)
G4NistManager * nist
const G4ParticleDefinition * baseParticle
std::vector< G4double > localCuts
const G4Material * cutMaterial
string material
Definition: eplot.py:19
#define DBL_MAX
Definition: templates.hh:62