Geant4 10.7.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//
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class header file
31//
32//
33// File name: G4EmCalculator
34//
35// Author: Vladimir Ivanchenko
36//
37// Creation date: 27.06.2004
38//
39// Modifications:
40// 17.11.2004 Change signature of methods, add new methods (V.Ivanchenko)
41// 11.01.2006 Add GetCSDARange (V.Ivanchenko)
42// 26.01.2006 Rename GetRange -> GetRangeFromRestricteDEDX (V.Ivanchenko)
43// 22.03.2006 Add ComputeElectronicDEDX and ComputeTotalDEDX (V.Ivanchenko)
44// 29.09.2006 Add member loweModel (V.Ivanchenko)
45// 15.03.2007 Add ComputeEnergyCutFromRangeCut methods (V.Ivanchenko)
46// 02.02.2018 Add parameter to FindLambdaTable to store process type (M. Novak)
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 G4NistManager;
66class G4Material;
69class G4PhysicsTable;
70class G4VEmModel;
72class G4VEmProcess;
74class G4VProcess;
76class G4Region;
77class G4Element;
78class G4EmCorrections;
79class G4EmParameters;
80class G4IonTable;
81
83{
84
85public:
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
97 const G4Material*,
98 const G4Region* r = nullptr);
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 = nullptr);
108 const G4String& part,
109 const G4String& mat,
110 const G4String& s = "world");
111
113 const G4Material*,
114 const G4Region* r = nullptr);
115 inline G4double GetCSDARange(G4double kinEnergy, const G4String& part,
116 const G4String& mat,
117 const G4String& s = "world");
118
120 const G4Material*,
121 const G4Region* r = nullptr);
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 = nullptr);
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 = nullptr);
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 = nullptr);
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
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
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 = nullptr);
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 = nullptr);
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
274private:
275
276 G4bool UpdateParticle(const G4ParticleDefinition*, G4double kinEnergy);
277
278 G4bool UpdateCouple(const G4Material*, G4double cut);
279
280 void FindLambdaTable(const G4ParticleDefinition*,
281 const G4String& processName,
282 G4double kinEnergy, G4int& proctype);
283
284 G4bool FindEmModel(const G4ParticleDefinition*,
285 const G4String& processName,
286 G4double kinEnergy);
287
288 G4VEnergyLossProcess* FindEnergyLossProcess(const G4ParticleDefinition*);
289
290 G4VEnergyLossProcess* FindEnLossProcess(const G4ParticleDefinition*,
291 const G4String& processName);
292
293 G4VEmProcess* FindDiscreteProcess(const G4ParticleDefinition*,
294 const G4String& processName);
295
296 G4VMultipleScattering* FindMscProcess(const G4ParticleDefinition*,
297 const G4String& processName);
298
299 G4bool ActiveForParticle(const G4ParticleDefinition* part,
300 G4VProcess* proc);
301
302 void CheckMaterial(G4int Z);
303
304 // hide copy and assign
305 G4EmCalculator & operator=(const G4EmCalculator &right) = delete;
306 G4EmCalculator(const G4EmCalculator&) = delete;
307
308 std::vector<const G4Material*> localMaterials;
309 std::vector<const G4MaterialCutsCouple*> localCouples;
310
311 G4EmParameters* theParameters;
312 G4LossTableManager* manager;
313 G4NistManager* nist;
314 G4IonTable* ionTable;
315 G4EmCorrections* corr;
316 G4DataVector localCuts;
317 G4int nLocalMaterials;
318
319 G4int verbose;
320
321 // cache
322 G4int currentCoupleIndex;
323 const G4MaterialCutsCouple* currentCouple;
324 const G4Material* currentMaterial;
325 const G4Material* cutMaterial;
326 const G4ParticleDefinition* currentParticle;
327 const G4ParticleDefinition* lambdaParticle;
328 const G4ParticleDefinition* baseParticle;
329 const G4PhysicsTable* currentLambda;
330
331 G4VEmModel* currentModel;
332 G4VEmModel* loweModel;
333 G4VEnergyLossProcess* currentProcess;
334 G4VProcess* curProcess;
335
336 const G4ParticleDefinition* theGenericIon;
337 G4ionEffectiveCharge* ionEffCharge;
338 G4DynamicParticle dynParticle;
339
340 G4String currentName;
341 G4String lambdaName;
342 G4double currentCut;
343 G4double chargeSquare;
344 G4double massRatio;
345 G4double mass;
346 G4double cutenergy[3];
347 G4bool isIon;
348 G4bool isApplicable;
349
350 G4String currentParticleName;
351 G4String currentMaterialName;
352 G4String currentProcessName;
353};
354
355//....oooOO0OOooo.......oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
356//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
357
358inline
360 const G4String& material, const G4String& reg)
361{
362 return GetDEDX(kinEnergy,FindParticle(particle),
363 FindMaterial(material),FindRegion(reg));
364}
365
366//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
367
368inline
370 const G4String& particle,
371 const G4String& material,
372 const G4String& reg)
373{
374 return GetRangeFromRestricteDEDX(kinEnergy,FindParticle(particle),
375 FindMaterial(material),FindRegion(reg));
376}
377
378//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
379
380inline
382 const G4String& particle,
383 const G4String& material,
384 const G4String& reg)
385{
386 return GetCSDARange(kinEnergy,FindParticle(particle),
387 FindMaterial(material),FindRegion(reg));
388}
389
390//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
391
392inline
394 const G4String& particle,
395 const G4String& material,
396 const G4String& reg)
397{
398 return GetRange(kinEnergy,FindParticle(particle),
399 FindMaterial(material),FindRegion(reg));
400}
401
402//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
403
404inline
406 const G4String& material, const G4String& reg)
407{
408 return GetKinEnergy(range,FindParticle(particle),
409 FindMaterial(material),FindRegion(reg));
410}
411
412//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
413
414inline
416 const G4String& particle,
417 const G4String& processName,
418 const G4String& material,
419 const G4String& reg)
420{
421 return GetCrossSectionPerVolume(kinEnergy,FindParticle(particle),processName,
422 FindMaterial(material),FindRegion(reg));
423}
424
425//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
426
427inline
429 const G4String& particle,
430 const G4String& processName,
431 const G4String& material,
432 const G4String& reg)
433{
434 return GetMeanFreePath(kinEnergy,FindParticle(particle),processName,
435 FindMaterial(material),FindRegion(reg));
436}
437
438//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
439
440inline G4double
442 const G4String& mat, G4double cut)
443{
444 return
445 ComputeElectronicDEDX(kinEnergy,FindParticle(part),FindMaterial(mat),cut);
446}
447
448//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
449
450inline G4double
452 const G4String& part,
453 const G4String& mat,
454 G4double rangecut)
455{
456 return ComputeDEDXForCutInRange(kinEnergy,FindParticle(part),
457 FindMaterial(mat), rangecut);
458}
459
460//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
461
462inline
464 const G4String& part,
465 const G4String& mat,
466 G4double cut)
467{
468 return ComputeTotalDEDX(kinEnergy,FindParticle(part),FindMaterial(mat),cut);
469}
470
471//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
472
473inline
475 const G4String& particle,
476 const G4String& processName,
477 const G4String& material,
478 G4double cut)
479{
480 return ComputeDEDX(kinEnergy,FindParticle(particle),processName,
481 FindMaterial(material),cut);
482}
483
484//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
485
486inline
488 const G4String& particle,
489 const G4String& material)
490{
491 return ComputeNuclearDEDX(kinEnergy,FindParticle(particle),
492 FindMaterial(material));
493}
494
495//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
496
497inline
499 G4double kinEnergy,
500 const G4String& particle,
501 const G4String& processName,
502 const G4String& material,
503 G4double cut)
504{
505 return ComputeCrossSectionPerVolume(kinEnergy,FindParticle(particle),
506 processName,
507 FindMaterial(material),cut);
508}
509
510//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
511
512inline
514 const G4String& particle,
515 const G4String& processName,
516 const G4Element* elm,
517 G4double cut)
518{
519 return ComputeCrossSectionPerAtom(kinEnergy,FindParticle(particle),
520 processName,
521 elm->GetZ(),elm->GetN(),cut);
522}
523
524//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
525
527 G4double kinEnergy, const G4String& part,
528 const G4String& processName, const G4Element* elm,
529 G4int shellIdx, G4double cut)
530{
531 return ComputeCrossSectionPerShell(kinEnergy, FindParticle(part),
532 processName, elm->GetZasInt(),
533 shellIdx, cut);
534}
535
536//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
537
538inline
540 G4double range,
541 const G4String& particle,
542 const G4String& material)
543{
544 return ComputeEnergyCutFromRangeCut(range,FindParticle(particle),
545 FindMaterial(material));
546}
547
548//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
549
550inline
552 const G4String& particle,
553 const G4String& processName,
554 const G4String& material,
555 G4double cut)
556{
557 return ComputeMeanFreePath(kinEnergy,FindParticle(particle),processName,
558 FindMaterial(material),cut);
559}
560
561//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
562
563#endif
G4AtomicShellEnumerator
double A(double temperature)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4double GetZ() const
Definition: G4Element.hh:130
G4int GetZasInt() const
Definition: G4Element.hh:131
G4double GetN() const
Definition: G4Element.hh:134
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 GetDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
G4double ComputeTotalDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, G4double cut=DBL_MAX)
G4double GetCSDARange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
G4double GetRange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
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 *)
G4double GetKinEnergy(G4double range, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
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)
void PrintDEDXTable(const G4ParticleDefinition *)
const G4Region * FindRegion(const G4String &)
G4VProcess * FindProcess(const G4ParticleDefinition *part, const G4String &processName)
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)
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 *)
void PrintInverseRangeTable(const G4ParticleDefinition *)
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)
G4double ComputeElectronicDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *mat, G4double cut=DBL_MAX)
const G4Material * FindMaterial(const G4String &)
#define DBL_MAX
Definition: templates.hh:62