Geant4  10.01.p01
pyG4EmCalculator.cc
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 // $Id: pyG4EmCalculator.cc 76884 2013-11-18 12:54:03Z gcosmo $
27 // ====================================================================
28 // pyG4EmCalculator.cc
29 //
30 // 2006 Q
31 // ====================================================================
32 #include <boost/python.hpp>
33 #include "G4Version.hh"
34 #include "G4EmCalculator.hh"
35 #include "G4Region.hh"
36 #include "G4Material.hh"
37 #include "G4Element.hh"
38 #include "G4ParticleDefinition.hh"
39 #include "G4MaterialCutsCouple.hh"
40 
41 using namespace boost::python;
42 
43 // ====================================================================
44 // thin wrappers
45 // ====================================================================
46 namespace pyG4EmCalculator {
47 
48 // GetDEDX
52 
54  (G4double, const G4String&, const G4String&, const G4String&)
56 
57 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_GetDEDX, GetDEDX, 3, 4)
58 
59 // GetRange
62  = &G4EmCalculator::GetRange;
63 
65  (G4double, const G4String&, const G4String&, const G4String&)
66  = &G4EmCalculator::GetRange;
67 
68 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_GetRange, GetRange, 3, 4)
69 
70 // GetKinEnergy
73  = &G4EmCalculator::GetKinEnergy;
74 
76  (G4double, const G4String&, const G4String&, const G4String&)
77  = &G4EmCalculator::GetKinEnergy;
78 
79 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_GetKinEnergy, GetKinEnergy, 3, 4)
80 
81 // GetCrossSectionPerVolume
84  const G4String&, const G4Material*, const G4Region*)
85  = &G4EmCalculator::GetCrossSectionPerVolume;
86 
88  (G4double, const G4String&, const G4String&,
89  const G4String&, const G4String&)
90  = &G4EmCalculator::GetCrossSectionPerVolume;
91 
92 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_GetCrossSectionPerVolume,
93  GetCrossSectionPerVolume, 4, 5)
94 
95 // GetMeanFreePath
96 G4double (G4EmCalculator::*f1_GetMeanFreePath)
98  const G4String&, const G4Material*, const G4Region*)
99  = &G4EmCalculator::GetMeanFreePath;
100 
102  (G4double, const G4String&, const G4String&,
103  const G4String&, const G4String&)
104  = &G4EmCalculator::GetMeanFreePath;
105 
107  GetMeanFreePath, 4, 5)
108 
109 // ComputeDEDX
110 G4double (G4EmCalculator::*f1_ComputeDEDX)
111  (G4double, const G4ParticleDefinition*,
112  const G4String&, const G4Material*, G4double)
113  = &G4EmCalculator::ComputeDEDX;
114 
115 G4double (G4EmCalculator::*f2_ComputeDEDX)
116  (G4double, const G4String&, const G4String&, const G4String&, G4double)
117  = &G4EmCalculator::ComputeDEDX;
118 
119 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_ComputeDEDX, ComputeDEDX, 4, 5)
120 
121 // ComputeNuclearDEDX
123  (G4double, const G4ParticleDefinition*, const G4Material*)
124  = &G4EmCalculator::ComputeNuclearDEDX;
125 
127  (G4double, const G4String&, const G4String&)
128  = &G4EmCalculator::ComputeNuclearDEDX;
129 
130 // ComputeElectronicDEDX
132  (G4double, const G4ParticleDefinition*, const G4Material*, G4double)
133  = &G4EmCalculator::ComputeElectronicDEDX;
134 
136  (G4double, const G4String&, const G4String&, G4double)
137  = &G4EmCalculator::ComputeElectronicDEDX;
138 
139 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_ComputeElectronicDEDX,
140  ComputeElectronicDEDX, 3, 4)
141 
142 // ComputeTotalDEDX
143 G4double (G4EmCalculator::*f1_ComputeTotalDEDX)
144  (G4double, const G4ParticleDefinition*, const G4Material*, G4double)
145  = &G4EmCalculator::ComputeTotalDEDX;
146 
148  (G4double, const G4String&, const G4String&, G4double)
149  = &G4EmCalculator::ComputeTotalDEDX;
150 
151 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_ComputeTotalDEDX, ComputeTotalDEDX,
152  3, 4)
153 // ComputeCrossSectionPerVolume
154 G4double (G4EmCalculator::*f1_ComputeCrossSectionPerVolume)
155  (G4double, const G4ParticleDefinition*,
156  const G4String&, const G4Material*, G4double)
157  = &G4EmCalculator::ComputeCrossSectionPerVolume;
158 
160  (G4double, const G4String&, const G4String&, const G4String&, G4double)
161  = &G4EmCalculator::ComputeCrossSectionPerVolume;
162 
163 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_ComputeCrossSectionPerVolume,
164  ComputeCrossSectionPerVolume, 4, 5)
165 
166 // ComputeCrossSectionPerAtom
167 G4double (G4EmCalculator::*f1_ComputeCrossSectionPerAtom)
168  (G4double, const G4ParticleDefinition*, const G4String&,
169  G4double, G4double, G4double)
171 
173  (G4double, const G4String&, const G4String&, const G4Element*, G4double)
174  = &G4EmCalculator::ComputeCrossSectionPerAtom;
175 
176 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_ComputeCrossSectionPerAtom,
177  ComputeCrossSectionPerAtom, 5, 6)
178 
179 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(g_ComputeCrossSectionPerAtom,
180  ComputeCrossSectionPerAtom, 4, 5)
181 
182 // ComputeEnergyCutFromRangeCut
184  (G4double, const G4ParticleDefinition*, const G4Material*)
185  = &G4EmCalculator::ComputeEnergyCutFromRangeCut;
186 
188  (G4double range, const G4String&, const G4String&)
189  = &G4EmCalculator::ComputeEnergyCutFromRangeCut;
190 
191 // ComputeMeanFreePath
193  (G4double, const G4ParticleDefinition*,
194  const G4String&, const G4Material*, G4double)
195  = &G4EmCalculator::ComputeMeanFreePath;
196 
198  (G4double, const G4String&, const G4String&, const G4String&, G4double)
199  = &G4EmCalculator::ComputeMeanFreePath;
200 
201 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_ComputeMeanFreePath,
202  ComputeMeanFreePath, 4, 5)
203 
204 // FindCouple
205 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_FindCouple, FindCouple, 1, 2)
206 
207 }
208 
209 using namespace pyG4EmCalculator;
210 
211 // ====================================================================
212 // module definition
213 // ====================================================================
215 {
216  class_<G4EmCalculator, boost::noncopyable>
217  ("G4EmCalculator", "Provide access to dE/dx and cross section")
218  // ---
219  .def("GetDEDX", f1_GetDEDX, f_GetDEDX())
220  .def("GetDEDX", f2_GetDEDX, f_GetDEDX())
221  .def("GetRange", f1_GetRange, f_GetRange())
222  .def("GetRange", f2_GetDEDX, f_GetRange())
223  .def("GetKinEnergy", f1_GetKinEnergy, f_GetKinEnergy())
224  .def("GetKinEnergy", f2_GetKinEnergy, f_GetKinEnergy())
225  .def("GetCrossSectionPerVolume",
226  f1_GetCrossSectionPerVolume, f_GetCrossSectionPerVolume())
227  .def("GetCrossSectionPerVolume",
228  f2_GetCrossSectionPerVolume, f_GetCrossSectionPerVolume())
229  .def("GetMeanFreePath", f1_GetMeanFreePath, f_GetMeanFreePath())
230  .def("GetMeanFreePath", f2_GetMeanFreePath, f_GetMeanFreePath())
231  // ---
232  .def("PrintDEDXTable", &G4EmCalculator::PrintDEDXTable)
233  .def("PrintRangeTable", &G4EmCalculator::PrintRangeTable)
234  .def("PrintInverseRangeTable", &G4EmCalculator::PrintInverseRangeTable)
235  // ---
236  .def("ComputeDEDX", f1_ComputeDEDX, f_ComputeDEDX())
237  .def("ComputeDEDX", f2_ComputeDEDX, f_ComputeDEDX())
238  .def("ComputeNuclearDEDX", f1_ComputeNuclearDEDX)
239  .def("ComputeNuclearDEDX", f2_ComputeNuclearDEDX)
240  .def("ComputeElectronicDEDX", f1_ComputeElectronicDEDX,
241  f_ComputeElectronicDEDX())
242  .def("ComputeDEDX", f2_ComputeElectronicDEDX,
243  f_ComputeElectronicDEDX())
244  .def("ComputeTotalDEDX", f1_ComputeTotalDEDX, f_ComputeTotalDEDX())
245  .def("ComputeTotalDEDX", f2_ComputeTotalDEDX, f_ComputeTotalDEDX())
246  // ---
247  .def("ComputeCrossSectionPerVolume",
248  f1_ComputeCrossSectionPerVolume, f_ComputeCrossSectionPerVolume())
249  .def("ComputeCrossSectionPerVolume",
250  f2_ComputeCrossSectionPerVolume, f_ComputeCrossSectionPerVolume())
251  .def("ComputeCrossSectionPerAtom",
252  f1_ComputeCrossSectionPerAtom, f_ComputeCrossSectionPerAtom())
253  .def("ComputeCrossSectionPerAtom",
254  f2_ComputeCrossSectionPerAtom, g_ComputeCrossSectionPerAtom())
255  .def("ComputeEnergyCutFromRangeCut", f1_ComputeEnergyCutFromRangeCut)
256  .def("ComputeEnergyCutFromRangeCut", f2_ComputeEnergyCutFromRangeCut)
257  // ---
258  .def("ComputeMeanFreePath",
259  f1_ComputeMeanFreePath, f_ComputeMeanFreePath())
260  .def("ComputeMeanFreePath",
261  f2_ComputeMeanFreePath, f_ComputeMeanFreePath())
262  // ---
263  .def("FindParticle", &G4EmCalculator::FindParticle,
264  return_value_policy<reference_existing_object>())
265  .def("FindMaterial", &G4EmCalculator::FindMaterial,
266  return_value_policy<reference_existing_object>())
267  .def("FindRegion", &G4EmCalculator::FindRegion,
268  return_value_policy<reference_existing_object>())
269  .def("FindCouple", &G4EmCalculator::FindCouple,
270  f_FindCouple()[return_value_policy<reference_existing_object>()])
271  // ---
272  .def("SetVerbose", &G4EmCalculator::SetVerbose)
273  ;
274 }
G4double(G4EmCalculator::* f2_GetCrossSectionPerVolume)(G4double, const G4String &, const G4String &, const G4String &, const G4String &)
const G4Material * FindMaterial(const G4String &)
G4double(G4EmCalculator::* f1_ComputeMeanFreePath)(G4double, const G4ParticleDefinition *, const G4String &, const G4Material *, G4double)
BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_GetMeanFreePath, GetMeanFreePath, 4, 5) G4double(G4EmCalculator G4double(G4EmCalculator::* f2_ComputeDEDX)(G4double, const G4String &, const G4String &, const G4String &, G4double)
void PrintRangeTable(const G4ParticleDefinition *)
void export_G4EmCalculator()
void PrintInverseRangeTable(const G4ParticleDefinition *)
G4double(G4EmCalculator::* f2_ComputeEnergyCutFromRangeCut)(G4double range, const G4String &, const G4String &)
G4double GetDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
G4double(G4EmCalculator::* f1_GetKinEnergy)(G4double, const G4ParticleDefinition *, const G4Material *, const G4Region *)
G4double(G4EmCalculator::* f1_GetDEDX)(G4double, const G4ParticleDefinition *, const G4Material *, const G4Region *)
G4double(G4EmCalculator::* f2_GetDEDX)(G4double, const G4String &, const G4String &, const G4String &)
BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_ComputeCrossSectionPerAtom, ComputeCrossSectionPerAtom, 5, 6) BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(g_ComputeCrossSectionPerAtom
const G4MaterialCutsCouple * FindCouple(const G4Material *, const G4Region *r=0)
BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_ComputeElectronicDEDX, ComputeElectronicDEDX, 3, 4) G4double(G4EmCalculator G4double(G4EmCalculator::* f2_ComputeTotalDEDX)(G4double, const G4String &, const G4String &, G4double)
G4double(G4EmCalculator::* f1_GetRange)(G4double, const G4ParticleDefinition *, const G4Material *, const G4Region *)
G4double(G4EmCalculator::* f2_ComputeElectronicDEDX)(G4double, const G4String &, const G4String &, G4double)
G4double(G4EmCalculator::* f2_ComputeMeanFreePath)(G4double, const G4String &, const G4String &, const G4String &, G4double)
G4double(G4EmCalculator::* f2_ComputeNuclearDEDX)(G4double, const G4String &, const G4String &)
G4double(G4EmCalculator::* f1_ComputeElectronicDEDX)(G4double, const G4ParticleDefinition *, const G4Material *, G4double)
BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_ComputeTotalDEDX, ComputeTotalDEDX, 3, 4) G4double(G4EmCalculator G4double(G4EmCalculator::* f2_ComputeCrossSectionPerVolume)(G4double, const G4String &, const G4String &, const G4String &, G4double)
G4double(G4EmCalculator::* f2_GetKinEnergy)(G4double, const G4String &, const G4String &, const G4String &)
const G4ParticleDefinition * FindParticle(const G4String &)
const G4Region * FindRegion(const G4String &)
G4double(G4EmCalculator::* f1_GetCrossSectionPerVolume)(G4double, const G4ParticleDefinition *, const G4String &, const G4Material *, const G4Region *)
G4double(G4EmCalculator::* f1_ComputeEnergyCutFromRangeCut)(G4double, const G4ParticleDefinition *, const G4Material *)
void SetVerbose(G4int val)
void PrintDEDXTable(const G4ParticleDefinition *)
double G4double
Definition: G4Types.hh:76
BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_ComputeCrossSectionPerVolume, ComputeCrossSectionPerVolume, 4, 5) G4double(G4EmCalculator G4double(G4EmCalculator::* f2_ComputeCrossSectionPerAtom)(G4double, const G4String &, const G4String &, const G4Element *, G4double)
BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_GetCrossSectionPerVolume, GetCrossSectionPerVolume, 4, 5) G4double(G4EmCalculator G4double(G4EmCalculator::* f2_GetMeanFreePath)(G4double, const G4String &, const G4String &, const G4String &, const G4String &)
G4double(G4EmCalculator::* f2_GetRange)(G4double, const G4String &, const G4String &, const G4String &)
G4double(G4EmCalculator::* f1_ComputeNuclearDEDX)(G4double, const G4ParticleDefinition *, const G4Material *)