32 #include <boost/python.hpp> 33 #include "G4Version.hh" 34 #include "G4LogicalVolume.hh" 35 #include "G4Material.hh" 36 #include "G4VSolid.hh" 37 #include "G4FieldManager.hh" 38 #include "G4VSensitiveDetector.hh" 39 #include "G4UserLimits.hh" 40 #include "G4SmartVoxelHeader.hh" 41 #include "G4MaterialCutsCouple.hh" 42 #include "G4FastSimulationManager.hh" 43 #include "G4VisAttributes.hh" 53 = &G4LogicalVolume::SetVisAttributes;
56 = &G4LogicalVolume::SetVisAttributes;
58 G4VSolid*(G4LogicalVolume::*
f1_GetSolid)()
const = &G4LogicalVolume::GetSolid;
61 = &G4LogicalVolume::SetSolid;
74 class_<G4LogicalVolume, G4LogicalVolume*, boost::noncopyable>
75 (
"G4LogicalVolume",
"logical volume class", no_init)
77 .def(init<G4VSolid*, G4Material*, const G4String& >())
78 .def(init<G4VSolid*, G4Material*,
const G4String&,
80 .def(init<G4VSolid*, G4Material*,
const G4String&,
81 G4FieldManager*, G4VSensitiveDetector* >())
82 .def(init<G4VSolid*, G4Material*,
const G4String&,
83 G4FieldManager*, G4VSensitiveDetector*,
85 .def(init<G4VSolid*, G4Material*,
const G4String&,
86 G4FieldManager*, G4VSensitiveDetector*,
87 G4UserLimits*, G4bool >())
89 .def(
"GetName", &G4LogicalVolume::GetName,
90 return_value_policy<reference_existing_object>())
91 .def(
"SetName", &G4LogicalVolume::SetName)
93 .def(
"GetNoDaughters", &G4LogicalVolume::GetNoDaughters)
94 .def(
"GetDaughter", &G4LogicalVolume::GetDaughter,
95 return_internal_reference<>())
96 .def(
"AddDaughter", &G4LogicalVolume::AddDaughter)
97 .def(
"IsDaughter", &G4LogicalVolume::IsDaughter)
98 .def(
"IsAncestor", &G4LogicalVolume::IsAncestor)
99 .def(
"RemoveDaughter", &G4LogicalVolume::RemoveDaughter)
100 .def(
"ClearDaughters", &G4LogicalVolume::ClearDaughters)
101 .def(
"TotalVolumeEntities", &G4LogicalVolume::TotalVolumeEntities)
104 return_internal_reference<>())
106 .def(
"GetMaterial", &G4LogicalVolume::GetMaterial,
107 return_internal_reference<>())
108 .def(
"SetMaterial", &G4LogicalVolume::SetMaterial)
109 .def(
"UpdateMaterial", &G4LogicalVolume::UpdateMaterial)
111 .def(
"GetMass", &G4LogicalVolume::GetMass, f_GetMass())
112 .def(
"GetFieldManager", &G4LogicalVolume::GetFieldManager,
113 return_internal_reference<>())
114 .def(
"SetFieldManager", &G4LogicalVolume::SetFieldManager)
115 .def(
"GetSensitiveDetector", &G4LogicalVolume::GetSensitiveDetector,
116 return_internal_reference<>())
117 .def(
"SetSensitiveDetector", &G4LogicalVolume::SetSensitiveDetector)
118 .def(
"GetUserLimits", &G4LogicalVolume::GetUserLimits,
119 return_internal_reference<>())
120 .def(
"SetUserLimits", &G4LogicalVolume::SetUserLimits)
122 .def(
"GetVoxelHeader", &G4LogicalVolume::GetVoxelHeader,
123 return_internal_reference<>())
124 .def(
"SetVoxelHeader", &G4LogicalVolume::SetVoxelHeader)
125 .def(
"GetSmartless", &G4LogicalVolume::GetSmartless)
126 .def(
"SetSmartless", &G4LogicalVolume::SetSmartless)
127 .def(
"IsToOptimise", &G4LogicalVolume::IsToOptimise)
128 .def(
"SetOptimisation", &G4LogicalVolume::SetOptimisation)
130 .def(
"IsRootRegion", &G4LogicalVolume::IsRootRegion)
131 .def(
"SetRegionRootFlag", &G4LogicalVolume::SetRegionRootFlag)
132 .def(
"IsRegion", &G4LogicalVolume::IsRegion)
133 .def(
"SetRegion", &G4LogicalVolume::SetRegion)
134 .def(
"GetRegion", &G4LogicalVolume::GetRegion,
135 return_internal_reference<>())
136 .def(
"PropagateRegion", &G4LogicalVolume::PropagateRegion)
137 .def(
"GetMaterialCutsCouple", &G4LogicalVolume::GetMaterialCutsCouple,
138 return_internal_reference<>())
139 .def(
"SetMaterialCutsCouple", &G4LogicalVolume::SetMaterialCutsCouple)
141 .def(
"GetVisAttributes", &G4LogicalVolume::GetVisAttributes,
142 return_internal_reference<>())
146 .def(
"GetFastSimulationManager",
147 &G4LogicalVolume::GetFastSimulationManager,
148 return_internal_reference<>())
150 .def(
"SetBiasWeight", &G4LogicalVolume::SetBiasWeight)
151 .def(
"GetBiasWeight", &G4LogicalVolume::GetBiasWeight)
void(G4LogicalVolume::* f2_SetVisAttributes)(const G4VisAttributes &)
BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_CreateTubeVolume, CreateTubeVolume, 4, 6) BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(f_CreateConeVolume
G4VSolid *(G4LogicalVolume::* f1_GetSolid)() const
void(G4LogicalVolume::* f1_SetSolid)(G4VSolid *)
void(G4LogicalVolume::* f1_SetVisAttributes)(const G4VisAttributes *)
void export_G4LogicalVolume()