ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
emcalculator.py
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file emcalculator.py
1 """
2 # ==================================================================
3 # Python module (Python3)
4 #
5 # Calculation of photon cross section and stopping power for
6 # chared particles
7 #
8 # Q, 2005
9 # ==================================================================
10 """
11 from Geant4 import *
12 
13 # ==================================================================
14 # public symbols
15 # ==================================================================
16 __all__ = [ 'CalculatePhotonCrossSection', 'CalculateDEDX' ]
17 
18 # ==================================================================
19 # Photon Cross Section
20 # ==================================================================
21 def CalculatePhotonCrossSection(mat, elist, verbose=0,
22  plist=["compt", "", "phot", "conv"]):
23  """
24  Calculate photon cross section for a given material and
25  a list of energy, returing a list of cross sections for
26  the components of "Copmton scattering", "rayleigh scattering",
27  "photoelectric effect", "pair creation" and total one.
28 
29  Arguments:
30  mat: material name
31  elist: list of energy
32  verbose: verbose level [0]
33  plist: list of process name
34  (compton/rayleigh/photoelectic/conversion) [StandardEM set]
35 
36  Keys of index:
37  "compt": Compton Scattering
38  "rayleigh": Rayleigh Scattering
39  "phot" : photoelectric effect
40  "conv" : pair Creation
41  "tot" : total
42 
43  Example:
44  xsec_list= CalculatePhotonCrossSection(...)
45  value= xsec_list[energy_index]["compt"]
46  """
47  if(verbose>0):
48  print("-------------------------------------------------------------------")
49  print(" Photon Cross Section (", mat, ")")
50  print("Energy Compton Raleigh Photo- Pair Total")
51  print(" Scattering Scattering electric Creation")
52  print("(MeV) (cm2/g) (cm2/g) (cm2/g) (cm2/g) (cm2/g)")
53  print("-------------------------------------------------------------------")
54 
55  xsection_list= []
56  for ekin in elist:
57  xsec= {}
58  xsec["compt"] \
59  = gEmCalculator.ComputeCrossSectionPerVolume(ekin, "gamma", plist[0],
60  mat) * cm2/g
61  xsec["rayleigh"] \
62  = gEmCalculator.ComputeCrossSectionPerVolume(ekin, "gamma", plist[1],
63  mat) * cm2/g
64 
65  xsec["phot"] \
66  = gEmCalculator.ComputeCrossSectionPerVolume(ekin, "gamma", plist[2],
67  mat) * cm2/g
68  xsec["conv"] \
69  = gEmCalculator.ComputeCrossSectionPerVolume(ekin, "gamma", plist[3],
70  mat) * cm2/g
71 
72  xsec["tot"]= xsec["compt"] + xsec["rayleigh"] + xsec["phot"] + xsec["conv"]
73 
74  xsection_list.append((ekin, xsec))
75 
76  if(verbose>0):
77  print(" %8.3e %8.3e %8.3e %8.3e %8.3e %8.3e" \
78  % (ekin/MeV, xsec["compt"]/(cm2/g), xsec["rayleigh"]/(cm2/g),
79  xsec["phot"]/(cm2/g), xsec["conv"]/(cm2/g), xsec["tot"]/(cm2/g)))
80 
81  return xsection_list
82 
83 
84 # ==================================================================
85 # Stopping Power
86 # ==================================================================
87 def CalculateDEDX(part, mat, elist, verbose=0,
88  plist=["eIoni", "eBrem", "muIoni", "muBrems", "hIoni"]):
89  """
90  Calculate stopping powers for a give particle, material and
91  a list of energy, returing stopping power for the components of
92  "Ionization", "Radiation" and total one.
93 
94  Arguments:
95  part: particle name
96  mat: material name
97  elist: list of energy
98  verbose: verbose level [0]
99  plist: list of process name
100  (electron ionization/electron brems/
101  muon ionization/muon brems/hadron ionization) [StandardEM set]
102 
103  Keys of index:
104  "ioni": ionization
105  "brems": Bremsstrahlung
106  "tot": total
107 
108  Example:
109  dedx_list= CalculateDEDX(...)
110  value= dedx_list[energy_index]["ioni"]
111  """
112  if(verbose>0):
113  print("------------------------------------------------------")
114  print(" Stopping Power (", part, ",", mat, ")")
115  print(" Energy Ionization Radiation Total")
116  print(" (MeV) (MeVcm2/g) (MeVcm2/g) (MeVcm2/g)")
117  print("------------------------------------------------------")
118 
119  procname_brems= ""
120  procname_ioni= ""
121  if ( part=="e+" or part=="e-" ):
122  procname_ioni= plist[0]
123  procname_brems= plist[1]
124  elif ( part=="mu+" or part=="mu-"):
125  procname_ioni= plist[2]
126  procname_brems= plist[3]
127  else:
128  procname_ioni= plist[4]
129  procname_brems= ""
130 
131  dedx_list= []
132  for ekin in elist:
133  dedx= {}
134  dedx["ioni"] \
135  = gEmCalculator.ComputeDEDX(ekin, part, procname_ioni, mat) * MeV*cm2/g
136  dedx["brems"] \
137  = gEmCalculator.ComputeDEDX(ekin, part, procname_brems, mat) * MeV*cm2/g
138  dedx["tot"]= dedx["ioni"]+ dedx["brems"]
139 
140  if(verbose>0):
141  print(" %8.3e %8.3e %8.3e %8.3e" \
142  % (ekin/MeV, dedx["ioni"]/(MeV*cm2/g),
143  dedx["brems"]/(MeV*cm2/g), dedx["tot"]/(MeV*cm2/g) ))
144 
145 
146  dedx_list.append((ekin, dedx))
147 
148  return dedx_list
149