Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Functions | Variables
emcalculator Namespace Reference

Functions

def CalculatePhotonCrossSection
 
def CalculateDEDX
 

Variables

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

Function Documentation

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.

Referenced by CalculatePhotonCrossSection().

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

References CalculateDEDX().

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

Variable Documentation

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

Definition at line 17 of file emcalculator.py.