32 #include <boost/python.hpp> 33 #include "G4Version.hh" 34 #include "G4FieldManager.hh" 36 #include "G4ChordFinder.hh" 37 #include "G4MagneticField.hh" 45 #if G4VERSION_NUMBER <= 801 47 inline const G4ChordFinder* G4FieldManager::GetChordFinder()
const 59 = &G4FieldManager::GetChordFinder;
61 = &G4FieldManager::GetChordFinder;
72 class_<G4FieldManager, G4FieldManager*, boost::noncopyable>
73 (
"G4FieldManager",
"field manager class")
76 .def(init<G4Field*>())
77 .def(init<G4Field*, G4ChordFinder*>())
78 .def(init<G4Field*, G4ChordFinder*, G4bool>())
79 .def(init<G4MagneticField*>())
81 .def(
"SetDetectorField", &G4FieldManager::SetDetectorField)
82 .def(
"GetDetectorField", &G4FieldManager::GetDetectorField,
83 return_internal_reference<>())
84 .def(
"DoesFieldExist", &G4FieldManager::DoesFieldExist)
85 .def(
"CreateChordFinder", &G4FieldManager::CreateChordFinder)
86 .def(
"SetChordFinder", &G4FieldManager::SetChordFinder)
88 return_internal_reference<>())
90 return_internal_reference<>())
91 .def(
"ConfigureForTrack", &G4FieldManager::ConfigureForTrack)
92 .def(
"GetDeltaIntersection", &G4FieldManager::GetDeltaIntersection)
93 .def(
"GetDeltaOneStep", &G4FieldManager::GetDeltaOneStep)
94 .def(
"SetAccuraciesWithDeltaOneStep",
95 &G4FieldManager::SetAccuraciesWithDeltaOneStep)
96 .def(
"SetDeltaOneStep", &G4FieldManager::SetDeltaOneStep)
97 .def(
"SetDeltaIntersection", &G4FieldManager::SetDeltaIntersection)
98 .def(
"GetMinimumEpsilonStep", &G4FieldManager::GetMinimumEpsilonStep)
99 .def(
"SetMinimumEpsilonStep", &G4FieldManager::SetMinimumEpsilonStep)
100 .def(
"GetMaximumEpsilonStep", &G4FieldManager::GetMaximumEpsilonStep)
101 .def(
"SetMaximumEpsilonStep", &G4FieldManager::SetMaximumEpsilonStep)
102 .def(
"DoesFieldChangeEnergy", &G4FieldManager::DoesFieldChangeEnergy)
103 .def(
"SetFieldChangesEnergy", &G4FieldManager::SetFieldChangesEnergy)
const G4ChordFinder *(G4FieldManager::* f2_GetChordFinder)() const
G4ChordFinder *(G4FieldManager::* f1_GetChordFinder)()
void export_G4FieldManager()