Geant4  10.02.p03
emcalculator Namespace Reference

Functions

def CalculatePhotonCrossSection (mat, elist, verbose=0, plist=["compt", phot, conv)
 
def CalculateDEDX (part, mat, elist, verbose=0, plist=["eIoni", eBrem, muIoni, muBrems, hIoni)
 

Variables

list __all__ = [ 'CalculatePhotonCrossSection', 'CalculateDEDX' ]
 

Function Documentation

◆ CalculateDEDX()

def emcalculator.CalculateDEDX (   part,
  mat,
  elist,
  verbose = 0,
  plist = ["eIoni",
  eBrem,
  muIoni,
  muBrems,
  hIoni 
)
Calculate stopping powers for a give particle, material and
a list of energy, returing stopping power for the components of
"Ionization", "Radiation" and total one.

Arguments:
  part:     particle name
  mat:      material name
  elist:    list of energy
  verbose:  verbose level [0]
  plist:    list of process name
            (electron ionization/electron brems/
             muon ionization/muon brems/hadron ionization) [StandardEM set]
             
Keys of index:
  "ioni":   ionization
  "brems":  Bremsstrahlung
  "tot":    total

Example:
  dedx_list= CalculateDEDX(...)
  value= dedx_list[energy_index]["ioni"]    

Definition at line 89 of file emcalculator.py.

89  plist=["eIoni", "eBrem", "muIoni", "muBrems", "hIoni"]):
90  """
91  Calculate stopping powers for a give particle, material and
92  a list of energy, returing stopping power for the components of
93  "Ionization", "Radiation" and total one.
94 
95  Arguments:
96  part: particle name
97  mat: material name
98  elist: list of energy
99  verbose: verbose level [0]
100  plist: list of process name
101  (electron ionization/electron brems/
102  muon ionization/muon brems/hadron ionization) [StandardEM set]
103 
104  Keys of index:
105  "ioni": ionization
106  "brems": Bremsstrahlung
107  "tot": total
108 
109  Example:
110  dedx_list= CalculateDEDX(...)
111  value= dedx_list[energy_index]["ioni"]
112  """
113  if(verbose>0):
114  print "------------------------------------------------------"
115  print " Stopping Power (", part, ",", mat, ")"
116  print " Energy Ionization Radiation Total"
117  print " (MeV) (MeVcm2/g) (MeVcm2/g) (MeVcm2/g)"
118  print "------------------------------------------------------"
119 
120  procname_brems= ""
121  procname_ioni= ""
122  if ( part=="e+" or part=="e-" ):
123  procname_ioni= plist[0]
124  procname_brems= plist[1]
125  elif ( part=="mu+" or part=="mu-"):
126  procname_ioni= plist[2]
127  procname_brems= plist[3]
128  else:
129  procname_ioni= plist[4]
130  procname_brems= ""
131 
132  dedx_list= []
133  for ekin in elist:
134  dedx= {}
135  dedx["ioni"] \
136  = gEmCalculator.ComputeDEDX(ekin, part, procname_ioni, mat) * MeV*cm2/g
137  dedx["brems"] \
138  = gEmCalculator.ComputeDEDX(ekin, part, procname_brems, mat) * MeV*cm2/g
139  dedx["tot"]= dedx["ioni"]+ dedx["brems"]
140 
141  if(verbose>0):
142  print " %8.3e %8.3e %8.3e %8.3e" \
143  % (ekin/MeV, dedx["ioni"]/(MeV*cm2/g),
144  dedx["brems"]/(MeV*cm2/g), dedx["tot"]/(MeV*cm2/g) )
145 
146 
147  dedx_list.append((ekin, dedx))
148 
149  return dedx_list
150 
151 
if(nlines<=0)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CalculatePhotonCrossSection()

def emcalculator.CalculatePhotonCrossSection (   mat,
  elist,
  verbose = 0,
  plist = ["compt",
  phot,
  conv 
)
Calculate photon cross section for a given material and
a list of energy, returing a list of cross sections for
the components of "Copmton scattering", "rayleigh scattering",
"photoelectric effect", "pair creation" and total one.

Arguments:
  mat:      material name
  elist:    list of energy
  verbose:  verbose level [0]
  plist:    list of process name
            (compton/rayleigh/photoelectic/conversion) [StandardEM set]
  
Keys of index:
  "compt":     Compton Scattering
  "rayleigh":  Rayleigh Scattering
  "phot" :     photoelectric effect
  "conv" :     pair Creation
  "tot"  :     total

Example:
  xsec_list= CalculatePhotonCrossSection(...)
  value= xsec_list[energy_index]["compt"]

Definition at line 23 of file emcalculator.py.

23  plist=["compt", "", "phot", "conv"]):
24  """
25  Calculate photon cross section for a given material and
26  a list of energy, returing a list of cross sections for
27  the components of "Copmton scattering", "rayleigh scattering",
28  "photoelectric effect", "pair creation" and total one.
29 
30  Arguments:
31  mat: material name
32  elist: list of energy
33  verbose: verbose level [0]
34  plist: list of process name
35  (compton/rayleigh/photoelectic/conversion) [StandardEM set]
36 
37  Keys of index:
38  "compt": Compton Scattering
39  "rayleigh": Rayleigh Scattering
40  "phot" : photoelectric effect
41  "conv" : pair Creation
42  "tot" : total
43 
44  Example:
45  xsec_list= CalculatePhotonCrossSection(...)
46  value= xsec_list[energy_index]["compt"]
47  """
48  if(verbose>0):
49  print "-------------------------------------------------------------------"
50  print " Photon Cross Section (", mat, ")"
51  print "Energy Compton Raleigh Photo- Pair Total"
52  print " Scattering Scattering electric Creation"
53  print "(MeV) (cm2/g) (cm2/g) (cm2/g) (cm2/g) (cm2/g)"
54  print "-------------------------------------------------------------------"
55 
56  xsection_list= []
57  for ekin in elist:
58  xsec= {}
59  xsec["compt"] \
60  = gEmCalculator.ComputeCrossSectionPerVolume(ekin, "gamma", plist[0],
61  mat) * cm2/g
62  xsec["rayleigh"] \
63  = gEmCalculator.ComputeCrossSectionPerVolume(ekin, "gamma", plist[1],
64  mat) * cm2/g
65 
66  xsec["phot"] \
67  = gEmCalculator.ComputeCrossSectionPerVolume(ekin, "gamma", plist[2],
68  mat) * cm2/g
69  xsec["conv"] \
70  = gEmCalculator.ComputeCrossSectionPerVolume(ekin, "gamma", plist[3],
71  mat) * cm2/g
72 
73  xsec["tot"]= xsec["compt"] + xsec["rayleigh"] + xsec["phot"] + xsec["conv"]
74 
75  xsection_list.append((ekin, xsec))
76 
77  if(verbose>0):
78  print " %8.3e %8.3e %8.3e %8.3e %8.3e %8.3e" \
79  % (ekin/MeV, xsec["compt"]/(cm2/g), xsec["rayleigh"]/(cm2/g),
80  xsec["phot"]/(cm2/g), xsec["conv"]/(cm2/g), xsec["tot"]/(cm2/g))
81 
82  return xsection_list
83 
84 
85 # ==================================================================
86 # Stopping Power
87 # ==================================================================
if(nlines<=0)
Here is the call graph for this function:

Variable Documentation

◆ __all__

list emcalculator.__all__ = [ 'CalculatePhotonCrossSection', 'CalculateDEDX' ]
private

Definition at line 17 of file emcalculator.py.