Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4EmCalculator.hh
이 파일의 문서화 페이지로 가기
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 108306 2018-02-02 13:10:06Z 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 // 02.02.2018 Add parameter to FindLambdaTable to store process type (M. Novak)
48 //
49 // Class Description:
50 //
51 // Provide access to dE/dx and cross sections
52 
53 // -------------------------------------------------------------------
54 //
55 
56 #ifndef G4EmCalculator_h
57 #define G4EmCalculator_h 1
58 
59 #include <vector>
60 #include "globals.hh"
61 #include "G4DataVector.hh"
62 #include "G4DynamicParticle.hh"
63 #include "G4VAtomDeexcitation.hh"
64 
65 class G4LossTableManager;
66 class G4NistManager;
67 class G4Material;
70 class G4PhysicsTable;
71 class G4VEmModel;
73 class G4VEmProcess;
75 class G4VProcess;
77 class G4Region;
78 class G4Element;
79 class G4EmCorrections;
80 class G4EmParameters;
81 class G4IonTable;
82 
84 {
85 
86 public:
87 
89 
91 
92  //===========================================================================
93  // Methods to access precalculated dE/dx and cross sections
94  // Materials should exist in the list of the G4MaterialCutsCouple
95  //===========================================================================
96 
97  G4double GetDEDX(G4double kinEnergy, const G4ParticleDefinition*,
98  const G4Material*,
99  const G4Region* r = nullptr);
100  inline G4double GetDEDX(G4double kinEnergy, const G4String& part,
101  const G4String& mat,
102  const G4String& s = "world");
103 
105  const G4ParticleDefinition*,
106  const G4Material*,
107  const G4Region* r = nullptr);
108  inline G4double GetRangeFromRestricteDEDX(G4double kinEnergy,
109  const G4String& part,
110  const G4String& mat,
111  const G4String& s = "world");
112 
114  const G4Material*,
115  const G4Region* r = nullptr);
116  inline G4double GetCSDARange(G4double kinEnergy, const G4String& part,
117  const G4String& mat,
118  const G4String& s = "world");
119 
120  G4double GetRange(G4double kinEnergy, const G4ParticleDefinition*,
121  const G4Material*,
122  const G4Region* r = nullptr);
123  inline G4double GetRange(G4double kinEnergy, const G4String& part,
124  const G4String& mat,
125  const G4String& s = "world");
126 
128  const G4Material*,
129  const G4Region* r = nullptr);
130  inline G4double GetKinEnergy(G4double range, const G4String& part,
131  const G4String& mat,
132  const G4String& s = "world");
133 
135  G4double kinEnergy, const G4ParticleDefinition*,
136  const G4String& processName, const G4Material*,
137  const G4Region* r = nullptr);
139  G4double kinEnergy, const G4String& part, const G4String& proc,
140  const G4String& mat, const G4String& s = "world");
141 
143  const G4String& part, G4int Z,
145  G4double kinEnergy);
146 
148  const G4String& processName, const G4Material*,
149  const G4Region* r = nullptr);
150  inline G4double GetMeanFreePath(G4double kinEnergy, const G4String& part,
151  const G4String& proc, const G4String& mat,
152  const G4String& s = "world");
153 
155 
157 
159 
160  //===========================================================================
161  // Methods to calculate dE/dx and cross sections "on fly"
162  // Existing tables and G4MaterialCutsCouples are not used
163  //===========================================================================
164 
166  const G4String& processName, const G4Material*,
167  G4double cut = DBL_MAX);
168  inline G4double ComputeDEDX(G4double kinEnergy, const G4String& part,
169  const G4String& proc,
170  const G4String& mat, G4double cut = DBL_MAX);
171 
173  const G4ParticleDefinition*,
174  const G4Material* mat, G4double cut = DBL_MAX);
175  inline G4double ComputeElectronicDEDX(G4double kinEnergy, const G4String& part,
176  const G4String& mat, G4double cut = DBL_MAX);
177 
179  const G4ParticleDefinition*,
180  const G4Material* mat, G4double rangecut = DBL_MAX);
181  inline G4double ComputeDEDXForCutInRange(G4double kinEnergy, const G4String& part,
182  const G4String& mat,
183  G4double rangecut = DBL_MAX);
184 
186  const G4Material*);
187  inline G4double ComputeNuclearDEDX(G4double kinEnergy, const G4String& part,
188  const G4String& mat);
189 
191  const G4Material*, G4double cut = DBL_MAX);
192  inline G4double ComputeTotalDEDX(G4double kinEnergy, const G4String& part,
193  const G4String& mat, G4double cut = DBL_MAX);
194 
196  G4double kinEnergy, const G4ParticleDefinition*,
197  const G4String& processName, const G4Material*,
198  G4double cut = 0.0);
200  G4double kinEnergy, const G4String& part,
201  const G4String& proc,
202  const G4String& mat, G4double cut = 0.0);
203 
205  G4double kinEnergy, const G4ParticleDefinition*,
206  const G4String& processName, G4double Z, G4double A,
207  G4double cut = 0.0);
209  G4double kinEnergy, const G4String& part,
210  const G4String& processName, const G4Element*,
211  G4double cut = 0.0);
212 
214  G4double kinEnergy, const G4ParticleDefinition*,
215  const G4String& processName, G4int Z, G4int shellIdx,
216  G4double cut = 0.0);
218  G4double kinEnergy, const G4String& part,
219  const G4String& processName, const G4Element*,
220  G4int shellIdx,
221  G4double cut = 0.0);
222 
224  const G4Material*);
225 
227  const G4String& part, G4int Z,
229  G4double kinEnergy,
230  const G4Material* mat = nullptr);
231 
233  G4double kinEnergy, const G4ParticleDefinition*,
234  const G4String& processName, const G4Material*,
235  G4double cut = 0.0);
237  G4double kinEnergy, const G4String&, const G4String&,
238  const G4String& processName, G4double cut = 0.0);
239 
241  G4double range, const G4ParticleDefinition*,
242  const G4Material*);
244  G4double range, const G4String&,
245  const G4String&);
246 
247  //===========================================================================
248  // Methods to access particles, materials, regions, processes
249  //===========================================================================
250 
252 
254 
255  const G4Material* FindMaterial(const G4String&);
256 
257  const G4Region* FindRegion(const G4String&);
258 
260  const G4Region* r = nullptr);
261 
263  const G4String& processName);
264 
265  void SetupMaterial(const G4Material*);
266 
267  void SetupMaterial(const G4String&);
268 
269  void SetVerbose(G4int val);
270 
271  //===========================================================================
272  // Private methods
273  //===========================================================================
274 
275 private:
276 
278 
279  G4bool UpdateCouple(const G4Material*, G4double cut);
280 
282  const G4String& processName,
283  G4double kinEnergy, G4int& proctype);
284 
286  const G4String& processName,
287  G4double kinEnergy);
288 
290 
292  const G4String& processName);
293 
295  const G4String& processName);
296 
298  const G4String& processName);
299 
301  G4VProcess* proc);
302 
303  void CheckMaterial(G4int Z);
304 
305  // hide copy and assign
306  G4EmCalculator & operator=(const G4EmCalculator &right) = delete;
307  G4EmCalculator(const G4EmCalculator&) = delete;
308 
309  std::vector<const G4Material*> localMaterials;
310  std::vector<const G4MaterialCutsCouple*> localCouples;
311 
319 
321 
322  // cache
331 
336 
340 
350 
354 };
355 
356 //....oooOO0OOooo.......oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
357 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
358 
359 inline
361  const G4String& material, const G4String& reg)
362 {
363  return GetDEDX(kinEnergy,FindParticle(particle),
364  FindMaterial(material),FindRegion(reg));
365 }
366 
367 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
368 
369 inline
371  const G4String& particle,
372  const G4String& material,
373  const G4String& reg)
374 {
375  return GetRangeFromRestricteDEDX(kinEnergy,FindParticle(particle),
376  FindMaterial(material),FindRegion(reg));
377 }
378 
379 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
380 
381 inline
383  const G4String& particle,
384  const G4String& material,
385  const G4String& reg)
386 {
387  return GetCSDARange(kinEnergy,FindParticle(particle),
388  FindMaterial(material),FindRegion(reg));
389 }
390 
391 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
392 
393 inline
395  const G4String& particle,
396  const G4String& material,
397  const G4String& reg)
398 {
399  return GetRange(kinEnergy,FindParticle(particle),
400  FindMaterial(material),FindRegion(reg));
401 }
402 
403 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
404 
405 inline
407  const G4String& material, const G4String& reg)
408 {
409  return GetKinEnergy(range,FindParticle(particle),
410  FindMaterial(material),FindRegion(reg));
411 }
412 
413 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
414 
415 inline
417  const G4String& particle,
418  const G4String& processName,
419  const G4String& material,
420  const G4String& reg)
421 {
422  return GetCrossSectionPerVolume(kinEnergy,FindParticle(particle),processName,
423  FindMaterial(material),FindRegion(reg));
424 }
425 
426 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
427 
428 inline
430  const G4String& particle,
431  const G4String& processName,
432  const G4String& material,
433  const G4String& reg)
434 {
435  return GetMeanFreePath(kinEnergy,FindParticle(particle),processName,
436  FindMaterial(material),FindRegion(reg));
437 }
438 
439 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
440 
441 inline G4double
443  const G4String& mat, G4double cut)
444 {
445  return
446  ComputeElectronicDEDX(kinEnergy,FindParticle(part),FindMaterial(mat),cut);
447 }
448 
449 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
450 
451 inline G4double
453  const G4String& part,
454  const G4String& mat,
455  G4double rangecut)
456 {
457  return ComputeDEDXForCutInRange(kinEnergy,FindParticle(part),
458  FindMaterial(mat), rangecut);
459 }
460 
461 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
462 
463 inline
465  const G4String& part,
466  const G4String& mat,
467  G4double cut)
468 {
469  return ComputeTotalDEDX(kinEnergy,FindParticle(part),FindMaterial(mat),cut);
470 }
471 
472 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
473 
474 inline
476  const G4String& particle,
477  const G4String& processName,
478  const G4String& material,
479  G4double cut)
480 {
481  return ComputeDEDX(kinEnergy,FindParticle(particle),processName,
482  FindMaterial(material),cut);
483 }
484 
485 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
486 
487 inline
489  const G4String& particle,
490  const G4String& material)
491 {
492  return ComputeNuclearDEDX(kinEnergy,FindParticle(particle),
493  FindMaterial(material));
494 }
495 
496 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
497 
498 inline
500  G4double kinEnergy,
501  const G4String& particle,
502  const G4String& processName,
503  const G4String& material,
504  G4double cut)
505 {
506  return ComputeCrossSectionPerVolume(kinEnergy,FindParticle(particle),
507  processName,
508  FindMaterial(material),cut);
509 }
510 
511 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
512 
513 inline
515  const G4String& particle,
516  const G4String& processName,
517  const G4Element* elm,
518  G4double cut)
519 {
520  return ComputeCrossSectionPerAtom(kinEnergy,FindParticle(particle),
521  processName,
522  elm->GetZ(),elm->GetN(),cut);
523 }
524 
525 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
526 
528  G4double kinEnergy, const G4String& part,
529  const G4String& processName, const G4Element* elm,
530  G4int shellIdx, G4double cut)
531 {
532  return ComputeCrossSectionPerShell(kinEnergy, FindParticle(part),
533  processName, elm->GetZasInt(),
534  shellIdx, cut);
535 }
536 
537 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
538 
539 inline
541  G4double range,
542  const G4String& particle,
543  const G4String& material)
544 {
545  return ComputeEnergyCutFromRangeCut(range,FindParticle(particle),
546  FindMaterial(material));
547 }
548 
549 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
550 
551 inline
553  const G4String& particle,
554  const G4String& processName,
555  const G4String& material,
556  G4double cut)
557 {
558  return ComputeMeanFreePath(kinEnergy,FindParticle(particle),processName,
559  FindMaterial(material),cut);
560 }
561 
562 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
563 
564 #endif
G4double ComputeNuclearDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *)
G4double GetRangeFromRestricteDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
void PrintInverseRangeTable(const G4ParticleDefinition *)
const G4PhysicsTable * currentLambda
G4DataVector localCuts
G4double GetMeanFreePath(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, const G4Region *r=nullptr)
void FindLambdaTable(const G4ParticleDefinition *, const G4String &processName, G4double kinEnergy, G4int &proctype)
G4String currentName
const G4Material * currentMaterial
const G4Region * FindRegion(const G4String &)
void SetupMaterial(const G4Material *)
const G4ParticleDefinition * lambdaParticle
const G4ParticleDefinition * theGenericIon
G4double GetDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
void CheckMaterial(G4int Z)
G4double ComputeDEDXForCutInRange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *mat, G4double rangecut=DBL_MAX)
G4double ComputeElectronicDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *mat, G4double cut=DBL_MAX)
G4NistManager * nist
G4DynamicParticle dynParticle
const G4ParticleDefinition * FindIon(G4int Z, G4int A)
G4ionEffectiveCharge * ionEffCharge
G4AtomicShellEnumerator
G4VProcess * FindProcess(const G4ParticleDefinition *part, const G4String &processName)
G4double cutenergy[3]
G4double GetN() const
Definition: G4Element.hh:135
Float_t Z
G4VMultipleScattering * FindMscProcess(const G4ParticleDefinition *, const G4String &processName)
const XML_Char * s
Definition: expat.h:262
const G4Material * FindMaterial(const G4String &)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
TString part[npart]
Definition: Style.C:32
G4String currentProcessName
G4double GetShellIonisationCrossSectionPerAtom(const G4String &part, G4int Z, G4AtomicShellEnumerator shell, G4double kinEnergy)
const G4ParticleDefinition const G4Material *G4double range
void PrintDEDXTable(const G4ParticleDefinition *)
G4double GetKinEnergy(G4double range, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
G4double ComputeMeanFreePath(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=0.0)
const G4ParticleDefinition * FindParticle(const G4String &)
G4VEnergyLossProcess * currentProcess
G4VProcess * curProcess
double A(double temperature)
G4double ComputeTotalDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, G4double cut=DBL_MAX)
G4IonTable * ionTable
Float_t mat
G4VEmModel * currentModel
G4double ComputeGammaAttenuationLength(G4double kinEnergy, const G4Material *)
G4bool UpdateCouple(const G4Material *, G4double cut)
G4double ComputeCrossSectionPerAtom(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, G4double Z, G4double A, G4double cut=0.0)
void PrintRangeTable(const G4ParticleDefinition *)
std::vector< const G4Material * > localMaterials
G4int GetZasInt() const
Definition: G4Element.hh:132
G4VEnergyLossProcess * FindEnLossProcess(const G4ParticleDefinition *, const G4String &processName)
G4double ComputeCrossSectionPerVolume(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=0.0)
int G4int
Definition: G4Types.hh:78
G4EmCorrections * corr
G4double ComputeCrossSectionPerShell(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, G4int Z, G4int shellIdx, G4double cut=0.0)
G4double GetCrossSectionPerVolume(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, const G4Region *r=nullptr)
G4double GetCSDARange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
const G4MaterialCutsCouple * currentCouple
G4String currentMaterialName
G4String currentParticleName
G4LossTableManager * manager
G4double GetZ() const
Definition: G4Element.hh:131
G4double chargeSquare
G4bool FindEmModel(const G4ParticleDefinition *, const G4String &processName, G4double kinEnergy)
G4VEmProcess * FindDiscreteProcess(const G4ParticleDefinition *, const G4String &processName)
G4double GetRange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
static const G4double reg
std::vector< const G4MaterialCutsCouple * > localCouples
const G4MaterialCutsCouple * FindCouple(const G4Material *, const G4Region *r=nullptr)
G4EmParameters * theParameters
void SetVerbose(G4int val)
#define DBL_MAX
Definition: templates.hh:83
const G4ParticleDefinition * baseParticle
G4bool UpdateParticle(const G4ParticleDefinition *, G4double kinEnergy)
G4EmCalculator & operator=(const G4EmCalculator &right)=delete
G4VEnergyLossProcess * FindEnergyLossProcess(const G4ParticleDefinition *)
const G4ParticleDefinition * currentParticle
G4double ComputeDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=DBL_MAX)
G4double ComputeShellIonisationCrossSectionPerAtom(const G4String &part, G4int Z, G4AtomicShellEnumerator shell, G4double kinEnergy, const G4Material *mat=nullptr)
const G4Material * cutMaterial
G4double ComputeEnergyCutFromRangeCut(G4double range, const G4ParticleDefinition *, const G4Material *)
G4bool ActiveForParticle(const G4ParticleDefinition *part, G4VProcess *proc)
G4VEmModel * loweModel