Geant4-11
Functions
pymodG4processes.cc File Reference
#include <boost/python.hpp>

Go to the source code of this file.

Functions

 BOOST_PYTHON_MODULE (G4processes)
 
void export_G4CrossSectionHandler ()
 
void export_G4EmCalculator ()
 
void export_G4ProcessManager ()
 
void export_G4ProcessTable ()
 
void export_G4ProcessType ()
 
void export_G4ProcVector ()
 
void export_G4ProductionCutsTable ()
 
void export_G4VCrossSectionHandler ()
 
void export_G4VProcess ()
 

Function Documentation

◆ BOOST_PYTHON_MODULE()

BOOST_PYTHON_MODULE ( G4processes  )

◆ export_G4CrossSectionHandler()

void export_G4CrossSectionHandler ( )

Definition at line 39 of file pyG4CrossSectionHandler.cc.

40{
41 class_<G4CrossSectionHandler, bases<G4VCrossSectionHandler>,
42 boost::noncopyable>
43 ("G4CrossSectionHandler", "cross section handler")
44 .def(init<>())
45 ;
46}

Referenced by BOOST_PYTHON_MODULE().

◆ export_G4EmCalculator()

void export_G4EmCalculator ( )

Definition at line 213 of file pyG4EmCalculator.cc.

214{
215 class_<G4EmCalculator, boost::noncopyable>
216 ("G4EmCalculator", "Provide access to dE/dx and cross section")
217 // ---
218 .def("GetDEDX", f1_GetDEDX, f_GetDEDX())
219 .def("GetDEDX", f2_GetDEDX, f_GetDEDX())
220 .def("GetRange", f1_GetRange, f_GetRange())
221 .def("GetRange", f2_GetDEDX, f_GetRange())
222 .def("GetKinEnergy", f1_GetKinEnergy, f_GetKinEnergy())
223 .def("GetKinEnergy", f2_GetKinEnergy, f_GetKinEnergy())
224 .def("GetCrossSectionPerVolume",
225 f1_GetCrossSectionPerVolume, f_GetCrossSectionPerVolume())
226 .def("GetCrossSectionPerVolume",
227 f2_GetCrossSectionPerVolume, f_GetCrossSectionPerVolume())
228 .def("GetMeanFreePath", f1_GetMeanFreePath, f_GetMeanFreePath())
229 .def("GetMeanFreePath", f2_GetMeanFreePath, f_GetMeanFreePath())
230 // ---
231 .def("PrintDEDXTable", &G4EmCalculator::PrintDEDXTable)
232 .def("PrintRangeTable", &G4EmCalculator::PrintRangeTable)
233 .def("PrintInverseRangeTable", &G4EmCalculator::PrintInverseRangeTable)
234 // ---
235 .def("ComputeDEDX", f1_ComputeDEDX, f_ComputeDEDX())
236 .def("ComputeDEDX", f2_ComputeDEDX, f_ComputeDEDX())
237 .def("ComputeNuclearDEDX", f1_ComputeNuclearDEDX)
238 .def("ComputeNuclearDEDX", f2_ComputeNuclearDEDX)
239 .def("ComputeElectronicDEDX", f1_ComputeElectronicDEDX,
240 f_ComputeElectronicDEDX())
241 .def("ComputeDEDX", f2_ComputeElectronicDEDX,
242 f_ComputeElectronicDEDX())
243 .def("ComputeTotalDEDX", f1_ComputeTotalDEDX, f_ComputeTotalDEDX())
244 .def("ComputeTotalDEDX", f2_ComputeTotalDEDX, f_ComputeTotalDEDX())
245 // ---
246 .def("ComputeCrossSectionPerVolume",
247 f1_ComputeCrossSectionPerVolume, f_ComputeCrossSectionPerVolume())
248 .def("ComputeCrossSectionPerVolume",
249 f2_ComputeCrossSectionPerVolume, f_ComputeCrossSectionPerVolume())
250 .def("ComputeCrossSectionPerAtom",
251 f1_ComputeCrossSectionPerAtom, f_ComputeCrossSectionPerAtom())
252 .def("ComputeCrossSectionPerAtom",
253 f2_ComputeCrossSectionPerAtom, g_ComputeCrossSectionPerAtom())
254 .def("ComputeEnergyCutFromRangeCut", f1_ComputeEnergyCutFromRangeCut)
255 .def("ComputeEnergyCutFromRangeCut", f2_ComputeEnergyCutFromRangeCut)
256 // ---
257 .def("ComputeMeanFreePath",
258 f1_ComputeMeanFreePath, f_ComputeMeanFreePath())
259 .def("ComputeMeanFreePath",
260 f2_ComputeMeanFreePath, f_ComputeMeanFreePath())
261 // ---
262 .def("FindParticle", &G4EmCalculator::FindParticle,
263 return_value_policy<reference_existing_object>())
264 .def("FindMaterial", &G4EmCalculator::FindMaterial,
265 return_value_policy<reference_existing_object>())
266 .def("FindRegion", &G4EmCalculator::FindRegion,
267 return_value_policy<reference_existing_object>())
268 .def("FindCouple", &G4EmCalculator::FindCouple,
269 f_FindCouple()[return_value_policy<reference_existing_object>()])
270 // ---
271 .def("SetVerbose", &G4EmCalculator::SetVerbose)
272 ;
273}
const G4ParticleDefinition * FindParticle(const G4String &)
void SetVerbose(G4int val)
void PrintDEDXTable(const G4ParticleDefinition *)
const G4Region * FindRegion(const G4String &)
const G4MaterialCutsCouple * FindCouple(const G4Material *, const G4Region *r=nullptr)
void PrintRangeTable(const G4ParticleDefinition *)
void PrintInverseRangeTable(const G4ParticleDefinition *)
const G4Material * FindMaterial(const G4String &)
G4double(G4EmCalculator::* f2_ComputeMeanFreePath)(G4double, const G4String &, const G4String &, const G4String &, G4double)
G4double(G4EmCalculator::* f2_ComputeEnergyCutFromRangeCut)(G4double range, const G4String &, const G4String &)
G4double(G4EmCalculator::* f1_GetDEDX)(G4double, const G4ParticleDefinition *, const G4Material *, const G4Region *)
G4double(G4EmCalculator::* f1_ComputeNuclearDEDX)(G4double, const G4ParticleDefinition *, const G4Material *)
G4double(G4EmCalculator::* f2_GetCrossSectionPerVolume)(G4double, const G4String &, const G4String &, const G4String &, const G4String &)
G4double(G4EmCalculator::* f2_GetDEDX)(G4double, const G4String &, const G4String &, const G4String &)
G4double(G4EmCalculator::* f1_GetRange)(G4double, const G4ParticleDefinition *, const G4Material *, const G4Region *)
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_GetKinEnergy)(G4double, const G4String &, const G4String &, const G4String &)
BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_ComputeTotalDEDX, ComputeTotalDEDX, 3, 4) G4double(G4EmCalculator G4double(G4EmCalculator::* f2_ComputeCrossSectionPerVolume)(G4double, const G4String &, const G4String &, const G4String &, 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)
G4double(G4EmCalculator::* f1_ComputeElectronicDEDX)(G4double, const G4ParticleDefinition *, const G4Material *, G4double)
G4double(G4EmCalculator::* f1_GetCrossSectionPerVolume)(G4double, const G4ParticleDefinition *, const G4String &, const G4Material *, const G4Region *)
G4double(G4EmCalculator::* f2_ComputeNuclearDEDX)(G4double, const G4String &, const G4String &)
BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_ComputeCrossSectionPerVolume, ComputeCrossSectionPerVolume, 4, 5) G4double(G4EmCalculator G4double(G4EmCalculator::* f2_ComputeCrossSectionPerAtom)(G4double, const G4String &, const G4String &, const G4Element *, G4double)
G4double(G4EmCalculator::* f1_ComputeMeanFreePath)(G4double, const G4ParticleDefinition *, const G4String &, const G4Material *, G4double)
G4double(G4EmCalculator::* f1_ComputeEnergyCutFromRangeCut)(G4double, const G4ParticleDefinition *, const G4Material *)
G4double(G4EmCalculator::* f2_ComputeElectronicDEDX)(G4double, const G4String &, const G4String &, G4double)
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_GetKinEnergy)(G4double, const G4ParticleDefinition *, const G4Material *, const G4Region *)

References pyG4EmCalculator::f1_ComputeElectronicDEDX, pyG4EmCalculator::f1_ComputeEnergyCutFromRangeCut, pyG4EmCalculator::f1_ComputeMeanFreePath, pyG4EmCalculator::f1_ComputeNuclearDEDX, pyG4EmCalculator::f1_GetCrossSectionPerVolume, pyG4EmCalculator::f1_GetDEDX, pyG4EmCalculator::f1_GetKinEnergy, pyG4EmCalculator::f1_GetRange, pyG4EmCalculator::f2_ComputeCrossSectionPerAtom, pyG4EmCalculator::f2_ComputeCrossSectionPerVolume, pyG4EmCalculator::f2_ComputeDEDX, pyG4EmCalculator::f2_ComputeElectronicDEDX, pyG4EmCalculator::f2_ComputeEnergyCutFromRangeCut, pyG4EmCalculator::f2_ComputeMeanFreePath, pyG4EmCalculator::f2_ComputeNuclearDEDX, pyG4EmCalculator::f2_ComputeTotalDEDX, pyG4EmCalculator::f2_GetCrossSectionPerVolume, pyG4EmCalculator::f2_GetDEDX, pyG4EmCalculator::f2_GetKinEnergy, pyG4EmCalculator::f2_GetMeanFreePath, G4EmCalculator::FindCouple(), G4EmCalculator::FindMaterial(), G4EmCalculator::FindParticle(), G4EmCalculator::FindRegion(), G4EmCalculator::PrintDEDXTable(), G4EmCalculator::PrintInverseRangeTable(), G4EmCalculator::PrintRangeTable(), and G4EmCalculator::SetVerbose().

Referenced by BOOST_PYTHON_MODULE().

◆ export_G4ProcessManager()

void export_G4ProcessManager ( )

Definition at line 170 of file pyG4ProcessManager.cc.

171{
172 class_<G4ProcessManager, G4ProcessManager*, boost::noncopyable>
173 ("G4ProcessManager", "process manager class", no_init)
174 // ---
175 .def("GetProcessList", f_GetProcessList)
176 .def("GetProcessListLength", &G4ProcessManager::GetProcessListLength)
177 .def("GetProcessIndex", &G4ProcessManager::GetProcessIndex)
178 .def("GetProcessVector", f_GetProcessVector,
179 g_GetProcessVector())
180 .def("GetAtRestProcessVector", f_GetAtRestProcessVector,
181 g_GetAtRestProcessVector())
182 .def("GetAlongStepProcessVector", f_GetAlongStepProcessVector,
183 g_GetAlongStepProcessVector())
184 .def("GetPostStepProcessVector", f_GetPostStepProcessVector,
185 g_GetPostStepProcessVector())
186 .def("GetProcessVectorIndex",
188 f_GetProcessVectorIndex())
189 .def("GetAtRestIndex", &G4ProcessManager::GetAtRestIndex,
190 f_GetAtRestIndex())
191 .def("GetAlongStepIndex", &G4ProcessManager::GetAlongStepIndex,
192 f_GetAlongStepIndex())
193 .def("GetPostStepIndex", &G4ProcessManager::GetPostStepIndex,
194 f_GetPostStepIndex())
195 // ----
196 .def("AddProcess", &G4ProcessManager::AddProcess,
197 f_AddProcess())
198 .def("AddRestProcess", &G4ProcessManager::AddRestProcess,
199 f_AddRestProcess())
200 .def("AddDiscreteProcess", &G4ProcessManager::AddDiscreteProcess,
201 f_AddDiscreteProcess())
202 .def("AddContinuousProcess", &G4ProcessManager::AddContinuousProcess,
203 f_AddContinuousProcess())
204 // ---
205 .def("GetProcessOrdering", &G4ProcessManager::GetProcessOrdering)
206 .def("SetProcessOrdering", &G4ProcessManager::SetProcessOrdering,
207 f_SetProcessOrdering())
208 .def("SetProcessOrderingToFirst",
210 .def("SetProcessOrderingToLast",
212 // ---
213 .def("RemoveProcess", f1_RemoveProcess,
214 return_value_policy<reference_existing_object>())
215 .def("RemoveProcess", f2_RemoveProcess,
216 return_value_policy<reference_existing_object>())
217 // ---
218 .def("SetProcessActivation", f1_SetProcessActivation,
219 return_value_policy<reference_existing_object>())
220 .def("SetProcessActivation", f2_SetProcessActivation,
221 return_value_policy<reference_existing_object>())
222 .def("GetProcessActivation", f1_GetProcessActivation)
223 .def("GetProcessActivation", f2_GetProcessActivation)
224 // ---
225 .def("GetParticleType", &G4ProcessManager::GetParticleType,
226 return_internal_reference<>())
227 .def("SetParticleType", &G4ProcessManager::SetParticleType)
228 .def("DumpInfo", &G4ProcessManager::DumpInfo)
229 .def("SetVerboseLevel", &G4ProcessManager::SetVerboseLevel)
230 .def("GetVerboseLevel", &G4ProcessManager::GetVerboseLevel)
231 ;
232
233 // enums...
234 enum_<G4ProcessVectorTypeIndex>("G4ProcessVectorTypeIndex")
235 .value("typeGPIL", typeGPIL)
236 .value("typeGPIL", typeDoIt)
237 ;
238
239 enum_<G4ProcessVectorDoItIndex>("G4ProcessVectorDoItIndex")
240 .value("idxAll", idxAll)
241 .value("idxAtRest", idxAtRest)
242 .value("idxAlongStep", idxAlongStep)
243 .value("idxPostStep", idxPostStep)
244 ;
245
246 enum_<G4ProcessVectorOrdering>("G4ProcessVectorOrdering")
247 .value("ordInActive", ordInActive)
248 .value("ordDefault", ordDefault)
249 .value("ordLast", ordLast)
250 ;
251}
@ typeGPIL
@ typeDoIt
@ ordInActive
@ ordDefault
@ ordLast
@ idxPostStep
@ idxAtRest
@ idxAll
@ idxAlongStep
void SetProcessOrdering(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt, G4int ordDoIt=ordDefault)
G4int AddRestProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4ParticleDefinition * GetParticleType() const
void SetParticleType(const G4ParticleDefinition *)
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4int GetAlongStepIndex(G4VProcess *aProcess, G4ProcessVectorTypeIndex typ=typeGPIL) const
G4int AddContinuousProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4int GetProcessListLength() const
G4int GetAtRestIndex(G4VProcess *aProcess, G4ProcessVectorTypeIndex typ=typeGPIL) const
void SetVerboseLevel(G4int value)
G4int GetProcessVectorIndex(G4VProcess *aProcess, G4ProcessVectorDoItIndex idx, G4ProcessVectorTypeIndex typ=typeGPIL) const
G4int GetProcessOrdering(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt)
G4int GetPostStepIndex(G4VProcess *aProcess, G4ProcessVectorTypeIndex typ=typeGPIL) const
G4int GetProcessIndex(G4VProcess *) const
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
void SetProcessOrderingToLast(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt)
G4int GetVerboseLevel() const
void SetProcessOrderingToFirst(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt)
list f_GetProcessList(const G4ProcessManager *procMgr)
G4VProcess *(G4ProcessManager::* f2_SetProcessActivation)(G4int, G4bool)
G4bool(G4ProcessManager::* f1_GetProcessActivation)(G4VProcess *) const
G4bool(G4ProcessManager::* f2_GetProcessActivation)(G4int) const
list f_GetProcessVector(const G4ProcessManager *procMgr, G4ProcessVectorDoItIndex idx, G4ProcessVectorTypeIndex typ=typeGPIL)
G4VProcess *(G4ProcessManager::* f1_SetProcessActivation)(G4VProcess *, G4bool)
G4VProcess *(G4ProcessManager::* f2_RemoveProcess)(G4int)

References G4ProcessManager::AddContinuousProcess(), G4ProcessManager::AddDiscreteProcess(), G4ProcessManager::AddProcess(), G4ProcessManager::AddRestProcess(), G4ProcessManager::DumpInfo(), pyG4ProcessManager::f1_GetProcessActivation, pyG4ProcessManager::f1_SetProcessActivation, pyG4ProcessManager::f2_GetProcessActivation, pyG4ProcessManager::f2_RemoveProcess, pyG4ProcessManager::f2_SetProcessActivation, pyG4ProcessManager::f_GetProcessList(), pyG4ProcessManager::f_GetProcessVector(), G4ProcessManager::GetAlongStepIndex(), G4ProcessManager::GetAtRestIndex(), G4ProcessManager::GetParticleType(), G4ProcessManager::GetPostStepIndex(), G4ProcessManager::GetProcessIndex(), G4ProcessManager::GetProcessListLength(), G4ProcessManager::GetProcessOrdering(), G4ProcessManager::GetProcessVectorIndex(), G4ProcessManager::GetVerboseLevel(), idxAll, idxAlongStep, idxAtRest, idxPostStep, ordDefault, ordInActive, ordLast, G4ProcessManager::SetParticleType(), G4ProcessManager::SetProcessOrdering(), G4ProcessManager::SetProcessOrderingToFirst(), G4ProcessManager::SetProcessOrderingToLast(), G4ProcessManager::SetVerboseLevel(), typeDoIt, and typeGPIL.

Referenced by BOOST_PYTHON_MODULE().

◆ export_G4ProcessTable()

void export_G4ProcessTable ( )

Definition at line 144 of file pyG4ProcessTable.cc.

145{
146 class_<G4ProcessTable, boost::noncopyable>
147 ("G4ProcessTable", "process table", no_init)
148 // ---
149 .def("GetProcessTable", &G4ProcessTable::GetProcessTable,
150 return_value_policy<reference_existing_object>())
151 .staticmethod("GetProcessTable")
152 .def("Length", &G4ProcessTable::Length)
153 //.def("Insert", &G4ProcessTable::Insert) // protected
154 //.def("Remove", &G4ProcessTable::Remove) // protected
155 // ---
156 .def("FindProcess", f1_FindProcess,
157 return_value_policy<reference_existing_object>())
158 .def("FindProcess", f2_FindProcess,
159 return_value_policy<reference_existing_object>())
160 .def("FindProcess", f3_FindProcess,
161 return_value_policy<reference_existing_object>())
162 .def("FindProcess", f3_FindProcess,
163 return_value_policy<reference_existing_object>())
164 // ---
165 .def("FindProcesses", f1_FindProcesses)
166 .def("FindProcesses", f2_FindProcesses)
167 .def("FindProcesses", f3_FindProcesses)
168 .def("FindProcesses", f4_FindProcesses)
169 // ---
170 .def("SetProcessActivation", f1_SetProcessActivation)
171 .def("SetProcessActivation", f2_SetProcessActivation)
172 .def("SetProcessActivation", f3_SetProcessActivation)
173 .def("SetProcessActivation", f4_SetProcessActivation)
174 .def("SetProcessActivation", f5_SetProcessActivation)
175 .def("SetProcessActivation", f6_SetProcessActivation)
176 .def("SetProcessActivation", f7_SetProcessActivation)
177 .def("SetProcessActivation", f8_SetProcessActivation)
178 // ---
179 .def("GetNameList", &G4ProcessTable::GetNameList,
180 return_internal_reference<>())
181 .def("DumpInfo", &G4ProcessTable::DumpInfo, f_DumpInfo())
182 .def("SetVerboseLevel", &G4ProcessTable::SetVerboseLevel)
183 .def("GetVerboseLevel", &G4ProcessTable::GetVerboseLevel)
184 ;
185}
void DumpInfo(G4VProcess *process, const G4ParticleDefinition *particle=nullptr)
static G4ProcessTable * GetProcessTable()
G4int Length() const
G4ProcNameVector * GetNameList()
G4int GetVerboseLevel() const
void SetVerboseLevel(G4int value)
list f3_FindProcesses(G4ProcessTable *procTable, const G4String &pname)
G4VProcess *(G4ProcessTable::* f1_FindProcess)(const G4String &, const G4String &) const
void(G4ProcessTable::* f7_SetProcessActivation)(G4ProcessType, const G4ParticleDefinition *, G4bool)
list f1_FindProcesses(G4ProcessTable *procTable)
void(G4ProcessTable::* f2_SetProcessActivation)(const G4String &, const G4String &, G4bool)
list f2_FindProcesses(G4ProcessTable *procTable, const G4ProcessManager *procManager)
void(G4ProcessTable::* f1_SetProcessActivation)(const G4String &, G4bool)
void(G4ProcessTable::* f8_SetProcessActivation)(G4ProcessType, G4ProcessManager *, G4bool)
void(G4ProcessTable::* f5_SetProcessActivation)(G4ProcessType, G4bool)
G4VProcess *(G4ProcessTable::* f3_FindProcess)(const G4String &, const G4ProcessManager *) const
G4VProcess *(G4ProcessTable::* f2_FindProcess)(const G4String &, const G4ParticleDefinition *) const
void(G4ProcessTable::* f3_SetProcessActivation)(const G4String &, const G4ParticleDefinition *, G4bool)
void(G4ProcessTable::* f4_SetProcessActivation)(const G4String &, G4ProcessManager *, G4bool)
list f4_FindProcesses(G4ProcessTable *procTable, G4ProcessType ptype)
void(G4ProcessTable::* f6_SetProcessActivation)(G4ProcessType, const G4String &, G4bool)

References G4ProcessTable::DumpInfo(), pyG4ProcessTable::f1_FindProcess, pyG4ProcessTable::f1_FindProcesses(), pyG4ProcessTable::f1_SetProcessActivation, pyG4ProcessTable::f2_FindProcess, pyG4ProcessTable::f2_FindProcesses(), pyG4ProcessTable::f2_SetProcessActivation, pyG4ProcessTable::f3_FindProcess, pyG4ProcessTable::f3_FindProcesses(), pyG4ProcessTable::f3_SetProcessActivation, pyG4ProcessTable::f4_FindProcesses(), pyG4ProcessTable::f4_SetProcessActivation, pyG4ProcessTable::f5_SetProcessActivation, pyG4ProcessTable::f6_SetProcessActivation, pyG4ProcessTable::f7_SetProcessActivation, pyG4ProcessTable::f8_SetProcessActivation, G4ProcessTable::GetNameList(), G4ProcessTable::GetProcessTable(), G4ProcessTable::GetVerboseLevel(), G4ProcessTable::Length(), and G4ProcessTable::SetVerboseLevel().

Referenced by BOOST_PYTHON_MODULE().

◆ export_G4ProcessType()

void export_G4ProcessType ( )

Definition at line 39 of file pyG4ProcessType.cc.

40{
41 enum_<G4ProcessType>("G4ProcessType")
42 .value("fNotDefined", fNotDefined)
43 .value("fTransportation", fTransportation)
44 .value("fElectromagnetic", fElectromagnetic)
45 .value("fOptical", fOptical)
46 .value("fHadronic", fHadronic)
47 .value("fPhotolepton_hadron", fPhotolepton_hadron)
48 .value("fDecay", fDecay)
49 .value("fGeneral", fGeneral)
50 .value("fParameterisation", fParameterisation)
51 .value("fUserDefined", fUserDefined)
52 ;
53}
@ fOptical
@ fParameterisation
@ fGeneral
@ fDecay
@ fElectromagnetic
@ fHadronic
@ fUserDefined
@ fTransportation
@ fPhotolepton_hadron
@ fNotDefined

References fDecay, fElectromagnetic, fGeneral, fHadronic, fNotDefined, fOptical, fParameterisation, fPhotolepton_hadron, fTransportation, and fUserDefined.

Referenced by BOOST_PYTHON_MODULE().

◆ export_G4ProcVector()

void export_G4ProcVector ( )

Definition at line 48 of file pyG4ProcVector.cc.

49{
50 class_<G4ProcVector> ("G4ProcVector", "process vector")
51 .def(vector_indexing_suite<G4ProcVector>())
52 ;
53}

Referenced by BOOST_PYTHON_MODULE().

◆ export_G4ProductionCutsTable()

void export_G4ProductionCutsTable ( )

Definition at line 41 of file pyG4ProductionCutsTable.cc.

42{
43 class_<G4ProductionCutsTable, boost::noncopyable>
44 ("G4ProductionCutsTable", "production cuts table", no_init)
45 .def("GetProductionCutsTable",
47 return_value_policy<reference_existing_object>())
48 .staticmethod("GetProductionCutsTable")
49
50 // internally used methods are limmitted to be exposed...
51
52 // ---
53 .def("GetLowEdgeEnergy", &G4ProductionCutsTable::GetLowEdgeEnergy)
54 .def("GetHighEdgeEnergy", &G4ProductionCutsTable::GetHighEdgeEnergy)
55 .def("SetEnergyRange", &G4ProductionCutsTable::SetEnergyRange)
56 .def("DumpCouples", &G4ProductionCutsTable::DumpCouples)
57 .def("IsModified", &G4ProductionCutsTable::IsModified)
58 // ---
59#if G4VERSION_NUMBER >= 830
60 .def("ConvertRangeToEnergy", &G4ProductionCutsTable::ConvertRangeToEnergy)
61#endif
62 // ---
63 .def("SetVerboseLevel", &G4ProductionCutsTable::SetVerboseLevel)
64 .def("GetVerboseLevel", &G4ProductionCutsTable::GetVerboseLevel)
65 ;
66}
G4double GetLowEdgeEnergy() const
void SetVerboseLevel(G4int value)
G4double GetHighEdgeEnergy() const
void SetEnergyRange(G4double lowedge, G4double highedge)
static G4ProductionCutsTable * GetProductionCutsTable()
G4double ConvertRangeToEnergy(const G4ParticleDefinition *particle, const G4Material *material, G4double range)

References G4ProductionCutsTable::ConvertRangeToEnergy(), G4ProductionCutsTable::DumpCouples(), G4ProductionCutsTable::GetHighEdgeEnergy(), G4ProductionCutsTable::GetLowEdgeEnergy(), G4ProductionCutsTable::GetProductionCutsTable(), G4ProductionCutsTable::GetVerboseLevel(), G4ProductionCutsTable::IsModified(), G4ProductionCutsTable::SetEnergyRange(), and G4ProductionCutsTable::SetVerboseLevel().

Referenced by BOOST_PYTHON_MODULE().

◆ export_G4VCrossSectionHandler()

void export_G4VCrossSectionHandler ( )

Definition at line 60 of file pyG4VCrossSectionHandler.cc.

61{
62 class_<G4VCrossSectionHandler, boost::noncopyable>
63 ("G4VCrossSectionHandler", "cross section handler", no_init)
64 // ---
65 .def("Initialise", &G4VCrossSectionHandler::Initialise,
66 f_Initialise())
67 .def("SelectRandomElement", &G4VCrossSectionHandler::SelectRandomElement,
68 return_value_policy<reference_existing_object>())
69 .def("SelectRandomShell", &G4VCrossSectionHandler::SelectRandomShell)
70 .def("FindValue", f1_FindValue)
71 .def("FindValue", f2_FindValue)
72 .def("ValueForMaterial", &G4VCrossSectionHandler::ValueForMaterial)
73 .def("LoadData", &G4VCrossSectionHandler::LoadData)
74 .def("LoadShellData", &G4VCrossSectionHandler::LoadShellData)
75 .def("PrintData", &G4VCrossSectionHandler::PrintData)
76 .def("Clear", &G4VCrossSectionHandler::Clear)
77 ;
78}
G4double ValueForMaterial(const G4Material *material, G4double e) const
void LoadShellData(const G4String &dataFile)
void LoadData(const G4String &dataFile)
G4int SelectRandomShell(G4int Z, G4double e) const
const G4Element * SelectRandomElement(const G4MaterialCutsCouple *material, G4double e) const
void Initialise(G4VDataSetAlgorithm *interpolation=nullptr, G4double minE=250 *CLHEP::eV, G4double maxE=100 *CLHEP::GeV, G4int numberOfBins=200, G4double unitE=CLHEP::MeV, G4double unitData=CLHEP::barn, G4int minZ=1, G4int maxZ=99)
G4double(G4VCrossSectionHandler::* f1_FindValue)(G4int, G4double) const
G4double(G4VCrossSectionHandler::* f2_FindValue)(G4int, G4double, G4int) const

References G4VCrossSectionHandler::Clear(), pyG4VCrossSectionHandler::f1_FindValue, pyG4VCrossSectionHandler::f2_FindValue, G4VCrossSectionHandler::Initialise(), G4VCrossSectionHandler::LoadData(), G4VCrossSectionHandler::LoadShellData(), G4VCrossSectionHandler::PrintData(), G4VCrossSectionHandler::SelectRandomElement(), G4VCrossSectionHandler::SelectRandomShell(), and G4VCrossSectionHandler::ValueForMaterial().

Referenced by BOOST_PYTHON_MODULE().

◆ export_G4VProcess()

void export_G4VProcess ( )

Definition at line 50 of file pyG4VProcess.cc.

51{
52 class_<G4VProcess, G4VProcess*, boost::noncopyable>
53 ("G4VProcess", "base class for process", no_init)
54 // ---
55 // Note that only limited methods are exposed.
56 .def("SetPILfactor", &G4VProcess::SetPILfactor)
57 .def("GetPILfactor", &G4VProcess::GetPILfactor)
58 .def("IsApplicable", &G4VProcess::IsApplicable)
59 .def("BuildPhysicsTable", &G4VProcess::BuildPhysicsTable)
60 .def("PreparePhysicsTable", &G4VProcess::PreparePhysicsTable)
61 .def("StorePhysicsTable", &G4VProcess::StorePhysicsTable)
62 .def("RetrievePhysicsTable", &G4VProcess::RetrievePhysicsTable)
63 .def("GetPhysicsTableFileName", &G4VProcess::GetPhysicsTableFileName,
64 f_GetPhysicsTableFileName()
65 [return_value_policy<return_by_value>()])
66 .def("GetProcessName", &G4VProcess::GetProcessName,
67 return_value_policy<return_by_value>())
68 .def("GetProcessType", &G4VProcess::GetProcessType)
69 .def("DumpInfo", &G4VProcess::DumpInfo)
70 .def("SetVerboseLevel", &G4VProcess::SetVerboseLevel)
71 .def("GetVerboseLevel", &G4VProcess::GetVerboseLevel)
72 ;
73}
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:412
G4int GetVerboseLevel() const
Definition: G4VProcess.hh:418
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &, G4bool)
Definition: G4VProcess.hh:211
void SetPILfactor(G4double value)
Definition: G4VProcess.hh:449
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:388
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:182
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:187
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &, G4bool)
Definition: G4VProcess.hh:206
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:194
virtual void DumpInfo() const
Definition: G4VProcess.cc:167
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
Definition: G4VProcess.cc:182
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
G4double GetPILfactor() const
Definition: G4VProcess.hh:455

References G4VProcess::BuildPhysicsTable(), G4VProcess::DumpInfo(), G4VProcess::GetPhysicsTableFileName(), G4VProcess::GetPILfactor(), G4VProcess::GetProcessName(), G4VProcess::GetProcessType(), G4VProcess::GetVerboseLevel(), G4VProcess::IsApplicable(), G4VProcess::PreparePhysicsTable(), G4VProcess::RetrievePhysicsTable(), G4VProcess::SetPILfactor(), G4VProcess::SetVerboseLevel(), and G4VProcess::StorePhysicsTable().

Referenced by BOOST_PYTHON_MODULE().