Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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,v 1.7 2008-03-13 07:32:18 kmura Exp $
27 // $Name: not supported by cvs2svn $
28 // ====================================================================
29 // pyG4EmCalculator.cc
30 //
31 // 2006 Q
32 // ====================================================================
33 #include <boost/python.hpp>
34 #include "G4Version.hh"
35 #include "G4EmCalculator.hh"
36 #include "G4Region.hh"
37 #include "G4Material.hh"
38 #include "G4Element.hh"
39 #include "G4ParticleDefinition.hh"
40 #include "G4MaterialCutsCouple.hh"
41 
42 using namespace boost::python;
43 
44 // ====================================================================
45 // thin wrappers
46 // ====================================================================
47 namespace pyG4EmCalculator {
48 
49 // GetDEDX
53 
55  (G4double, const G4String&, const G4String&, const G4String&)
57 
58 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_GetDEDX, GetDEDX, 3, 4);
59 
60 // GetRange
64 
66  (G4double, const G4String&, const G4String&, const G4String&)
68 
69 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_GetRange, GetRange, 3, 4);
70 
71 // GetKinEnergy
75 
77  (G4double, const G4String&, const G4String&, const G4String&)
79 
80 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_GetKinEnergy, GetKinEnergy, 3, 4);
81 
82 // GetCrossSectionPerVolume
85  const G4String&, const G4Material*, const G4Region*)
87 
89  (G4double, const G4String&, const G4String&,
90  const G4String&, const G4String&)
92 
93 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_GetCrossSectionPerVolume,
94  GetCrossSectionPerVolume, 4, 5);
95 
96 // GetCrossSectionPerAtom
97 #if G4VERSION_NUMBER <= 801
100  const G4String&, const G4Material*, const G4Region*)
101  = &G4EmCalculator::GetCrossSectionPerAtom;
102 
104  (G4double, const G4String&, const G4String&,
105  const G4String&, const G4String&)
106  = &G4EmCalculator::GetCrossSectionPerAtom;
107 
108 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_GetCrossSectionPerAtom,
109  GetCrossSectionPerAtom, 4, 5);
110 #endif
111 
112 // GetMeanFreePath
115  const G4String&, const G4Material*, const G4Region*)
117 
119  (G4double, const G4String&, const G4String&,
120  const G4String&, const G4String&)
122 
123 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_GetMeanFreePath,
124  GetMeanFreePath, 4, 5);
125 
126 // ComputeDEDX
129  const G4String&, const G4Material*, G4double)
131 
133  (G4double, const G4String&, const G4String&, const G4String&, G4double)
135 
136 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_ComputeDEDX, ComputeDEDX, 4, 5);
137 
138 // ComputeNuclearDEDX
139 #if G4VERSION_NUMBER >=710
140 G4double (G4EmCalculator::*f1_ComputeNuclearDEDX)
141  (G4double, const G4ParticleDefinition*, const G4Material*)
143 
144 G4double (G4EmCalculator::*f2_ComputeNuclearDEDX)
145  (G4double, const G4String&, const G4String&)
147 #endif
148 
149 #if G4VERSION_NUMBER >= 810
150 // ComputeElectronicDEDX
151 G4double (G4EmCalculator::*f1_ComputeElectronicDEDX)
154 
155 G4double (G4EmCalculator::*f2_ComputeElectronicDEDX)
156  (G4double, const G4String&, const G4String&, G4double)
158 
159 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_ComputeElectronicDEDX,
160  ComputeElectronicDEDX, 3, 4);
161 
162 // ComputeTotalDEDX
163 G4double (G4EmCalculator::*f1_ComputeTotalDEDX)
166 
167 G4double (G4EmCalculator::*f2_ComputeTotalDEDX)
168  (G4double, const G4String&, const G4String&, G4double)
170 
171 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_ComputeTotalDEDX, ComputeTotalDEDX,
172  3, 4);
173 #endif
174 
175 
176 // ComputeCrossSectionPerVolume
179  const G4String&, const G4Material*, G4double)
181 
183  (G4double, const G4String&, const G4String&, const G4String&, G4double)
185 
186 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_ComputeCrossSectionPerVolume,
187  ComputeCrossSectionPerVolume, 4, 5);
188 
189 // ComputeCrossSectionPerAtom
194 
196  (G4double, const G4String&, const G4String&, const G4Element*, G4double)
198 
199 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_ComputeCrossSectionPerAtom,
200  ComputeCrossSectionPerAtom, 5, 6);
201 
202 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(g_ComputeCrossSectionPerAtom,
203  ComputeCrossSectionPerAtom, 4, 5);
204 
205 #if G4VERSION_NUMBER >= 830
206 // ComputeEnergyCutFromRangeCut
207 G4double (G4EmCalculator::*f1_ComputeEnergyCutFromRangeCut)
208  (G4double, const G4ParticleDefinition*, const G4Material*)
210 
211 G4double (G4EmCalculator::*f2_ComputeEnergyCutFromRangeCut)
212  (G4double range, const G4String&, const G4String&)
214 #endif
215 
216 
217 // ComputeMeanFreePath
220  const G4String&, const G4Material*, G4double)
222 
224  (G4double, const G4String&, const G4String&, const G4String&, G4double)
226 
227 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_ComputeMeanFreePath,
228  ComputeMeanFreePath, 4, 5);
229 
230 // FindCouple
231 BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_FindCouple, FindCouple, 1, 2);
232 
233 };
234 
235 using namespace pyG4EmCalculator;
236 
237 // ====================================================================
238 // module definition
239 // ====================================================================
241 {
242  class_<G4EmCalculator, boost::noncopyable>
243  ("G4EmCalculator", "Provide access to dE/dx and cross section")
244  // ---
245  .def("GetDEDX", f1_GetDEDX, f_GetDEDX())
246  .def("GetDEDX", f2_GetDEDX, f_GetDEDX())
247  .def("GetRange", f1_GetRange, f_GetRange())
248  .def("GetRange", f2_GetDEDX, f_GetRange())
249  .def("GetKinEnergy", f1_GetKinEnergy, f_GetKinEnergy())
250  .def("GetKinEnergy", f2_GetKinEnergy, f_GetKinEnergy())
251  .def("GetCrossSectionPerVolume",
252  f1_GetCrossSectionPerVolume, f_GetCrossSectionPerVolume())
253  .def("GetCrossSectionPerVolume",
254  f2_GetCrossSectionPerVolume, f_GetCrossSectionPerVolume())
255 #if G4VERSION_NUMBER <= 801
256  .def("GetCrossSectionPerAtom",
257  f1_GetCrossSectionPerAtom, f_GetCrossSectionPerAtom())
258  .def("GetCrossSectionPerAtom",
259  f2_GetCrossSectionPerAtom, f_GetCrossSectionPerAtom())
260 #endif
261  .def("GetMeanFreePath", f1_GetMeanFreePath, f_GetMeanFreePath())
262  .def("GetMeanFreePath", f2_GetMeanFreePath, f_GetMeanFreePath())
263  // ---
264  .def("PrintDEDXTable", &G4EmCalculator::PrintDEDXTable)
265  .def("PrintRangeTable", &G4EmCalculator::PrintRangeTable)
266  .def("PrintInverseRangeTable", &G4EmCalculator::PrintInverseRangeTable)
267  // ---
268  .def("ComputeDEDX", f1_ComputeDEDX, f_ComputeDEDX())
269  .def("ComputeDEDX", f2_ComputeDEDX, f_ComputeDEDX())
270 #if G4VERSION_NUMBER >=710
271  .def("ComputeNuclearDEDX", f1_ComputeNuclearDEDX)
272  .def("ComputeNuclearDEDX", f2_ComputeNuclearDEDX)
273 #endif
274 #if G4VERSION_NUMBER >= 810
275  .def("ComputeElectronicDEDX", f1_ComputeElectronicDEDX,
276  f_ComputeElectronicDEDX())
277  .def("ComputeDEDX", f2_ComputeElectronicDEDX,
278  f_ComputeElectronicDEDX())
279  .def("ComputeTotalDEDX", f1_ComputeTotalDEDX, f_ComputeTotalDEDX())
280  .def("ComputeTotalDEDX", f2_ComputeTotalDEDX, f_ComputeTotalDEDX())
281 #endif
282  // ---
283  .def("ComputeCrossSectionPerVolume",
284  f1_ComputeCrossSectionPerVolume, f_ComputeCrossSectionPerVolume())
285  .def("ComputeCrossSectionPerVolume",
286  f2_ComputeCrossSectionPerVolume, f_ComputeCrossSectionPerVolume())
287  .def("ComputeCrossSectionPerAtom",
288  f1_ComputeCrossSectionPerAtom, f_ComputeCrossSectionPerAtom())
289  .def("ComputeCrossSectionPerAtom",
290  f2_ComputeCrossSectionPerAtom, g_ComputeCrossSectionPerAtom())
291 #if G4VERSION_NUMBER >= 830
292  .def("ComputeEnergyCutFromRangeCut", f1_ComputeEnergyCutFromRangeCut)
293  .def("ComputeEnergyCutFromRangeCut", f2_ComputeEnergyCutFromRangeCut)
294 #endif
295  // ---
296  .def("ComputeMeanFreePath",
297  f1_ComputeMeanFreePath, f_ComputeMeanFreePath())
298  .def("ComputeMeanFreePath",
299  f2_ComputeMeanFreePath, f_ComputeMeanFreePath())
300  // ---
301  .def("FindParticle", &G4EmCalculator::FindParticle,
302  return_value_policy<reference_existing_object>())
303  .def("FindMaterial", &G4EmCalculator::FindMaterial,
304  return_value_policy<reference_existing_object>())
305  .def("FindRegion", &G4EmCalculator::FindRegion,
306  return_value_policy<reference_existing_object>())
307  .def("FindCouple", &G4EmCalculator::FindCouple,
308  f_FindCouple()[return_value_policy<reference_existing_object>()])
309  // ---
310  .def("SetVerbose", &G4EmCalculator::SetVerbose)
311  ;
312 }
313