Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LogicalVolume Class Reference

#include <G4LogicalVolume.hh>

Public Member Functions

 G4LogicalVolume (G4VSolid *pSolid, G4Material *pMaterial, const G4String &name, G4FieldManager *pFieldMgr=0, G4VSensitiveDetector *pSDetector=0, G4UserLimits *pULimits=0, G4bool optimise=true)
 
 ~G4LogicalVolume ()
 
G4String GetName () const
 
void SetName (const G4String &pName)
 
G4int GetNoDaughters () const
 
G4VPhysicalVolumeGetDaughter (const G4int i) const
 
void AddDaughter (G4VPhysicalVolume *p)
 
G4bool IsDaughter (const G4VPhysicalVolume *p) const
 
G4bool IsAncestor (const G4VPhysicalVolume *p) const
 
void RemoveDaughter (const G4VPhysicalVolume *p)
 
void ClearDaughters ()
 
G4int TotalVolumeEntities () const
 
G4VSolidGetSolid () const
 
void SetSolid (G4VSolid *pSolid)
 
G4MaterialGetMaterial () const
 
void SetMaterial (G4Material *pMaterial)
 
void UpdateMaterial (G4Material *pMaterial)
 
G4double GetMass (G4bool forced=false, G4bool propagate=true, G4Material *parMaterial=0)
 
G4FieldManagerGetFieldManager () const
 
void SetFieldManager (G4FieldManager *pFieldMgr, G4bool forceToAllDaughters)
 
G4VSensitiveDetectorGetSensitiveDetector () const
 
void SetSensitiveDetector (G4VSensitiveDetector *pSDetector)
 
G4UserLimitsGetUserLimits () const
 
void SetUserLimits (G4UserLimits *pULimits)
 
G4SmartVoxelHeaderGetVoxelHeader () const
 
void SetVoxelHeader (G4SmartVoxelHeader *pVoxel)
 
G4double GetSmartless () const
 
void SetSmartless (G4double s)
 
G4bool IsToOptimise () const
 
void SetOptimisation (G4bool optim)
 
G4bool IsRootRegion () const
 
void SetRegionRootFlag (G4bool rreg)
 
G4bool IsRegion () const
 
void SetRegion (G4Region *reg)
 
G4RegionGetRegion () const
 
void PropagateRegion ()
 
const G4MaterialCutsCoupleGetMaterialCutsCouple () const
 
void SetMaterialCutsCouple (G4MaterialCutsCouple *cuts)
 
G4bool operator== (const G4LogicalVolume &lv) const
 
const G4VisAttributesGetVisAttributes () const
 
void SetVisAttributes (const G4VisAttributes *pVA)
 
void SetVisAttributes (const G4VisAttributes &VA)
 
G4FastSimulationManagerGetFastSimulationManager () const
 
void SetBiasWeight (G4double w)
 
G4double GetBiasWeight () const
 
 G4LogicalVolume (__void__ &)
 
void Lock ()
 

Detailed Description

Definition at line 133 of file G4LogicalVolume.hh.

Constructor & Destructor Documentation

◆ G4LogicalVolume() [1/2]

G4LogicalVolume::G4LogicalVolume ( G4VSolid pSolid,
G4Material pMaterial,
const G4String name,
G4FieldManager pFieldMgr = 0,
G4VSensitiveDetector pSDetector = 0,
G4UserLimits pULimits = 0,
G4bool  optimise = true 
)

Definition at line 56 of file G4LogicalVolume.cc.

63 : fDaughters(0,(G4VPhysicalVolume*)0), fFieldManager(pFieldMgr),
64 fVoxel(0), fOptimise(optimise), fRootRegion(false), fLock(false),
65 fSmartless(2.), fMass(0.), fVisAttributes(0), fRegion(0), fCutsCouple(0)
66{
67 SetSolid(pSolid);
68 SetMaterial(pMaterial);
69 SetName(name);
70 SetSensitiveDetector(pSDetector);
71 SetUserLimits(pULimits);
72 //
73 // Add to store
74 //
76}
static void Register(G4LogicalVolume *pVolume)
void SetName(const G4String &pName)
void SetUserLimits(G4UserLimits *pULimits)
void SetMaterial(G4Material *pMaterial)
void SetSolid(G4VSolid *pSolid)
void SetSensitiveDetector(G4VSensitiveDetector *pSDetector)

◆ ~G4LogicalVolume()

G4LogicalVolume::~G4LogicalVolume ( )

Definition at line 99 of file G4LogicalVolume.cc.

100{
101 if (!fLock && fRootRegion) // De-register root region first if not locked
102 { // and flagged as root logical-volume
103 fRegion->RemoveRootLogicalVolume(this, true);
104 }
106}
static void DeRegister(G4LogicalVolume *pVolume)
void RemoveRootLogicalVolume(G4LogicalVolume *lv, G4bool scan=true)
Definition: G4Region.cc:256

◆ G4LogicalVolume() [2/2]

G4LogicalVolume::G4LogicalVolume ( __void__ &  )

Definition at line 83 of file G4LogicalVolume.cc.

84 : fDaughters(0,(G4VPhysicalVolume*)0), fFieldManager(0),
85 fMaterial(0), fName(""), fSensitiveDetector(0), fSolid(0), fUserLimits(0),
86 fVoxel(0), fOptimise(true), fRootRegion(false), fLock(false), fSmartless(2.),
87 fMass(0.), fVisAttributes(0), fRegion(0), fCutsCouple(0), fBiasWeight(0.)
88{
89 // Add to store
90 //
92}

Member Function Documentation

◆ AddDaughter()

void G4LogicalVolume::AddDaughter ( G4VPhysicalVolume p)
inline

◆ ClearDaughters()

void G4LogicalVolume::ClearDaughters ( )
inline

◆ GetBiasWeight()

G4double G4LogicalVolume::GetBiasWeight ( ) const
inline

◆ GetDaughter()

◆ GetFastSimulationManager()

◆ GetFieldManager()

G4FieldManager * G4LogicalVolume::GetFieldManager ( ) const
inline

◆ GetMass()

G4double G4LogicalVolume::GetMass ( G4bool  forced = false,
G4bool  propagate = true,
G4Material parMaterial = 0 
)

Definition at line 193 of file G4LogicalVolume.cc.

196{
197 // Return the cached non-zero value, if not forced
198 //
199 if ( (fMass) && (!forced) ) return fMass;
200
201 // Global density and computed mass associated to the logical
202 // volume without considering its daughters
203 //
204 G4Material* logMaterial = parMaterial ? parMaterial : fMaterial;
205 if (!logMaterial)
206 {
207 std::ostringstream message;
208 message << "No material associated to the logical volume: " << fName << " !"
209 << G4endl
210 << "Sorry, cannot compute the mass ...";
211 G4Exception("G4LogicalVolume::GetMass()", "GeomMgt0002",
212 FatalException, message);
213 return 0;
214 }
215 if (!fSolid)
216 {
217 std::ostringstream message;
218 message << "No solid is associated to the logical volume: " << fName << " !"
219 << G4endl
220 << "Sorry, cannot compute the mass ...";
221 G4Exception("G4LogicalVolume::GetMass()", "GeomMgt0002",
222 FatalException, message);
223 return 0;
224 }
225 G4double globalDensity = logMaterial->GetDensity();
226 fMass = fSolid->GetCubicVolume() * globalDensity;
227
228 // For each daughter in the tree, subtract the mass occupied
229 // and if required by the propagate flag, add the real daughter's
230 // one computed recursively
231
232 for (G4PhysicalVolumeList::const_iterator itDau = fDaughters.begin();
233 itDau != fDaughters.end(); itDau++)
234 {
235 G4VPhysicalVolume* physDaughter = (*itDau);
236 G4LogicalVolume* logDaughter = physDaughter->GetLogicalVolume();
237 G4double subMass=0.;
238 G4VSolid* daughterSolid = 0;
239 G4Material* daughterMaterial = 0;
240
241 // Compute the mass to subtract and to add for each daughter
242 // considering its multiplicity (i.e. replicated or not) and
243 // eventually its parameterisation (by solid and/or by material)
244 //
245 for (G4int i=0; i<physDaughter->GetMultiplicity(); i++)
246 {
248 physParam = physDaughter->GetParameterisation();
249 if (physParam)
250 {
251 daughterSolid = physParam->ComputeSolid(i, physDaughter);
252 daughterSolid->ComputeDimensions(physParam, i, physDaughter);
253 daughterMaterial = physParam->ComputeMaterial(i, physDaughter);
254 }
255 else
256 {
257 daughterSolid = logDaughter->GetSolid();
258 daughterMaterial = logDaughter->GetMaterial();
259 }
260 subMass = daughterSolid->GetCubicVolume() * globalDensity;
261
262 // Subtract the daughter's portion for the mass and, if required,
263 // add the real daughter's mass computed recursively
264 //
265 fMass -= subMass;
266 if (propagate)
267 {
268 fMass += logDaughter->GetMass(true, true, daughterMaterial);
269 }
270 }
271 }
272
273 return fMass;
274}
@ FatalException
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4VSolid * GetSolid() const
G4double GetMass(G4bool forced=false, G4bool propagate=true, G4Material *parMaterial=0)
G4Material * GetMaterial() const
G4double GetDensity() const
Definition: G4Material.hh:179
virtual G4Material * ComputeMaterial(const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=0)
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
G4LogicalVolume * GetLogicalVolume() const
virtual G4int GetMultiplicity() const
virtual G4VPVParameterisation * GetParameterisation() const =0
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4VSolid.cc:137
virtual G4double GetCubicVolume()
Definition: G4VSolid.cc:188
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Referenced by GetMass(), and G4ASCIITreeSceneHandler::RequestPrimitives().

◆ GetMaterial()

◆ GetMaterialCutsCouple()

◆ GetName()

G4String G4LogicalVolume::GetName ( ) const
inline

Referenced by G4AdjointCrossSurfChecker::AddanExtSurfaceOfAvolume(), G4HepRepSceneHandler::AddSolid(), G4SmartVoxelHeader::BuildNodes(), G4SmartVoxelHeader::BuildReplicaVoxels(), G4SmartVoxelHeader::BuildVoxelsWithinLimits(), G4PVParameterised::CheckOverlaps(), G4PVPlacement::CheckOverlaps(), checkVol(), G4tgbVolume::ConstructG4LogVol(), G4tgbVolume::ConstructG4PhysVol(), G4tgbVolume::ConstructG4Volumes(), G4PhysicalVolumeModel::CreateCurrentAttValues(), G3Division::CreatePVReplica(), G4AdjointCrossSurfChecker::CrossingAnInterfaceBetweenTwoVolumes(), G4RadioactiveDecay::DecayIt(), G4RunManagerKernel::DefineWorldVolume(), G4RadioactiveDecay::DeselectAVolume(), G4ReflectionFactory::Divide(), G4GDMLReadStructure::DivisionvolRead(), G4GDMLWriteStructure::DivisionvolWrite(), G4TrajectoryDrawByOriginVolume::Draw(), G4tgbVolumeMgr::DumpG4LogVolLeaf(), G4LogicalSkinSurface::DumpInfo(), G4tgbGeometryDumper::DumpLogVol(), G4tgbGeometryDumper::DumpPhysVol(), G4tgbGeometryDumper::DumpPVParameterised(), G4tgbGeometryDumper::DumpPVPlacement(), G4tgbGeometryDumper::DumpPVReplica(), G4TrajectoryOriginVolumeFilter::Evaluate(), G4BuildGeom(), G4PVReplica::G4PVReplica(), G4GDMLRead::GeneratePhysvolName(), G4Navigator::GetLocalExitNormal(), G4ITNavigator::GetLocalExitNormal(), G4PolarizedCompton::GetMeanFreePath(), G4eplusPolarizedAnnihilation::GetMeanFreePath(), G4tgbVolumeMgr::GetTopLogVol(), G4tgbVolumeMgr::GetTopPhysVol(), G4GDMLReadStructure::GetWorldVolume(), G4AdjointCrossSurfChecker::GoingInOrOutOfaVolumeByExtSurface(), G4GDMLWriteParamvol::ParamvolAlgorithmWrite(), G4GDMLReadParamvol::ParamvolRead(), G4GDMLWriteParamvol::ParamvolWrite(), G4GDMLReadStructure::PhysvolRead(), G4GDMLWriteStructure::PhysvolWrite(), G4ReflectionFactory::Place(), G4eplusPolarizedAnnihilation::PostStepGetPhysicalInteractionLength(), G4PolarizedCompton::PostStepGetPhysicalInteractionLength(), G4tgbVolumeMgr::RegisterMe(), G4GDMLReadStructure::ReplicaRead(), G4ReflectionFactory::Replicate(), G4GDMLWriteStructure::ReplicavolWrite(), G4ASCIITreeSceneHandler::RequestPrimitives(), G4VoxelSafety::SafetyForVoxelHeader(), G4PolarizedComptonModel::SampleSecondaries(), G4Region::ScanVolumeTree(), G4RadioactiveDecay::SelectAllVolumes(), G4RadioactiveDecay::SelectAVolume(), G4VVisCommandGeometrySet::Set(), G4VVisCommandGeometrySet::SetLVVisAtts(), G4VisCommandGeometryList::SetNewValue(), G4VisCommandGeometryRestore::SetNewValue(), G4GDMLWriteSetup::SetupWrite(), G4PolarizationManager::SetVolumePolarization(), G4GDMLWriteStructure::SkinSurfaceCache(), G4GDMLRead::StripNames(), and G4GDMLWriteStructure::TraverseVolumeTree().

◆ GetNoDaughters()

◆ GetRegion()

◆ GetSensitiveDetector()

◆ GetSmartless()

G4double G4LogicalVolume::GetSmartless ( ) const
inline

◆ GetSolid()

G4VSolid * G4LogicalVolume::GetSolid ( ) const
inline

Referenced by G4AdjointCrossSurfChecker::AddanExtSurfaceOfAvolume(), G4GMocrenFileSceneHandler::AddPrimitive(), G4ReplicaNavigation::BackLocate(), G4PhantomParameterisation::BuildContainerSolid(), G4SmartVoxelHeader::BuildNodes(), G4SmartVoxelHeader::BuildReplicaVoxels(), G4SmartVoxelHeader::BuildVoxelsWithinLimits(), G4PVParameterised::CheckOverlaps(), G4PVPlacement::CheckOverlaps(), G4NormalNavigation::ComputeSafety(), G4VoxelNavigation::ComputeSafety(), G4ReplicaNavigation::ComputeSafety(), G4ParameterisedNavigation::ComputeSafety(), G4VoxelSafety::ComputeSafety(), G4VNestedParameterisation::ComputeSolid(), G4VPVParameterisation::ComputeSolid(), G4PhantomParameterisation::ComputeSolid(), G4ParameterisedNavigation::ComputeStep(), G4VoxelNavigation::ComputeStep(), G4ReplicaNavigation::ComputeStep(), G4NormalNavigation::ComputeStep(), G4Navigator::ComputeStep(), G4ITNavigator::ComputeStep(), G4RegularNavigation::ComputeStepSkippingEqualMaterials(), G4PSCellFlux::ComputeVolume(), G4PSDoseDeposit::ComputeVolume(), G4PSPassageCellFlux::ComputeVolume(), G4tgbVolume::ConstructG4PhysVol(), G4TheRayTracer::CreateBitMap(), G4PhysicalVolumeModel::CreateCurrentAttValues(), G3Division::CreatePVReplica(), G4AdjointPosOnPhysVolGenerator::DefinePhysicalVolume(), G4VisManager::Draw(), G4tgbGeometryDumper::DumpLogVol(), G4Navigator::GetLocalExitNormal(), G4ITNavigator::GetLocalExitNormal(), GetMass(), G4TransportationManager::GetParallelWorld(), G4Navigator::LocateGlobalPointAndSetup(), G4ITNavigator::LocateGlobalPointAndSetup(), G4GDMLWriteParamvol::ParametersWrite(), G4VXTRenergyLoss::PostStepDoIt(), G4NavigationLogger::PreComputeStepLog(), G4PSCylinderSurfaceCurrent::ProcessHits(), G4PSCylinderSurfaceFlux::ProcessHits(), G4PSFlatSurfaceCurrent::ProcessHits(), G4PSFlatSurfaceFlux::ProcessHits(), G4PSSphereSurfaceCurrent::ProcessHits(), G4PSSphereSurfaceFlux::ProcessHits(), G4VoxelSafety::SafetyForVoxelHeader(), G4VoxelSafety::SafetyForVoxelNode(), G4GeomTestVolume::TestOneLine(), and G4GDMLWriteStructure::TraverseVolumeTree().

◆ GetUserLimits()

◆ GetVisAttributes()

◆ GetVoxelHeader()

◆ IsAncestor()

G4bool G4LogicalVolume::IsAncestor ( const G4VPhysicalVolume p) const

Definition at line 139 of file G4LogicalVolume.cc.

140{
141 G4bool isDaughter = IsDaughter(aVolume);
142 if (!isDaughter)
143 {
144 for (G4PhysicalVolumeList::const_iterator itDau = fDaughters.begin();
145 itDau != fDaughters.end(); itDau++)
146 {
147 isDaughter = (*itDau)->GetLogicalVolume()->IsAncestor(aVolume);
148 if (isDaughter) break;
149 }
150 }
151 return isDaughter;
152}
bool G4bool
Definition: G4Types.hh:67
G4bool IsDaughter(const G4VPhysicalVolume *p) const

◆ IsDaughter()

G4bool G4LogicalVolume::IsDaughter ( const G4VPhysicalVolume p) const
inline

Referenced by IsAncestor().

◆ IsRegion()

G4bool G4LogicalVolume::IsRegion ( ) const
inline

◆ IsRootRegion()

G4bool G4LogicalVolume::IsRootRegion ( ) const
inline

◆ IsToOptimise()

G4bool G4LogicalVolume::IsToOptimise ( ) const
inline

◆ Lock()

void G4LogicalVolume::Lock ( )
inline

◆ operator==()

G4bool G4LogicalVolume::operator== ( const G4LogicalVolume lv) const

◆ PropagateRegion()

void G4LogicalVolume::PropagateRegion ( )
inline

◆ RemoveDaughter()

void G4LogicalVolume::RemoveDaughter ( const G4VPhysicalVolume p)
inline

◆ SetBiasWeight()

void G4LogicalVolume::SetBiasWeight ( G4double  w)
inline

◆ SetFieldManager()

void G4LogicalVolume::SetFieldManager ( G4FieldManager pFieldMgr,
G4bool  forceToAllDaughters 
)

Definition at line 113 of file G4LogicalVolume.cc.

115{
116 fFieldManager = pNewFieldMgr;
117
118 G4int NoDaughters = GetNoDaughters();
119 while ( (NoDaughters--)>0 )
120 {
121 G4LogicalVolume* DaughterLogVol;
122 DaughterLogVol = GetDaughter(NoDaughters)->GetLogicalVolume();
123 if ( forceAllDaughters || (DaughterLogVol->GetFieldManager() == 0) )
124 {
125 DaughterLogVol->SetFieldManager(pNewFieldMgr, forceAllDaughters);
126 }
127 }
128}
G4int GetNoDaughters() const
void SetFieldManager(G4FieldManager *pFieldMgr, G4bool forceToAllDaughters)
G4VPhysicalVolume * GetDaughter(const G4int i) const
G4FieldManager * GetFieldManager() const

Referenced by SetFieldManager().

◆ SetMaterial()

void G4LogicalVolume::SetMaterial ( G4Material pMaterial)
inline

◆ SetMaterialCutsCouple()

void G4LogicalVolume::SetMaterialCutsCouple ( G4MaterialCutsCouple cuts)
inline

◆ SetName()

void G4LogicalVolume::SetName ( const G4String pName)
inline

◆ SetOptimisation()

void G4LogicalVolume::SetOptimisation ( G4bool  optim)
inline

◆ SetRegion()

void G4LogicalVolume::SetRegion ( G4Region reg)
inline

◆ SetRegionRootFlag()

void G4LogicalVolume::SetRegionRootFlag ( G4bool  rreg)
inline

◆ SetSensitiveDetector()

void G4LogicalVolume::SetSensitiveDetector ( G4VSensitiveDetector pSDetector)
inline

Referenced by G4LogicalVolume().

◆ SetSmartless()

void G4LogicalVolume::SetSmartless ( G4double  s)
inline

◆ SetSolid()

◆ SetUserLimits()

void G4LogicalVolume::SetUserLimits ( G4UserLimits pULimits)
inline

Referenced by G4LogicalVolume().

◆ SetVisAttributes() [1/2]

void G4LogicalVolume::SetVisAttributes ( const G4VisAttributes VA)

Definition at line 276 of file G4LogicalVolume.cc.

277{
278 fVisAttributes = new G4VisAttributes(VA);
279}

◆ SetVisAttributes() [2/2]

◆ SetVoxelHeader()

void G4LogicalVolume::SetVoxelHeader ( G4SmartVoxelHeader pVoxel)
inline

◆ TotalVolumeEntities()

G4int G4LogicalVolume::TotalVolumeEntities ( ) const

Definition at line 161 of file G4LogicalVolume.cc.

162{
163 G4int vols = 1;
164 for (G4PhysicalVolumeList::const_iterator itDau = fDaughters.begin();
165 itDau != fDaughters.end(); itDau++)
166 {
167 G4VPhysicalVolume* physDaughter = (*itDau);
168 vols += physDaughter->GetMultiplicity()
169 *physDaughter->GetLogicalVolume()->TotalVolumeEntities();
170 }
171 return vols;
172}
G4int TotalVolumeEntities() const

Referenced by TotalVolumeEntities().

◆ UpdateMaterial()


The documentation for this class was generated from the following files: