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

#include <G4LogicalVolume.hh>

+ Inheritance diagram for G4LogicalVolume:

Public Member Functions

 G4LogicalVolume (G4VSolid *pSolid, G4Material *pMaterial, const G4String &name, G4FieldManager *pFieldMgr=nullptr, G4VSensitiveDetector *pSDetector=nullptr, G4UserLimits *pULimits=nullptr, G4bool optimise=true)
 
virtual ~G4LogicalVolume ()
 
 G4LogicalVolume (const G4LogicalVolume &)=delete
 
G4LogicalVolumeoperator= (const G4LogicalVolume &)=delete
 
const G4StringGetName () const
 
void SetName (const G4String &pName)
 
std::size_t GetNoDaughters () const
 
G4VPhysicalVolumeGetDaughter (const std::size_t 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
 
EVolume CharacteriseDaughters () const
 
EVolume DeduceDaughtersType () 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=nullptr)
 
void ResetMass ()
 
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__ &)
 
virtual G4bool IsExtended () const
 
G4FieldManagerGetMasterFieldManager () const
 
G4VSensitiveDetectorGetMasterSensitiveDetector () const
 
G4VSolidGetMasterSolid () const
 
G4int GetInstanceID () const
 
void Lock ()
 
void InitialiseWorker (G4LogicalVolume *ptrMasterObject, G4VSolid *pSolid, G4VSensitiveDetector *pSDetector)
 
void TerminateWorker (G4LogicalVolume *ptrMasterObject)
 
void AssignFieldManager (G4FieldManager *fldMgr)
 
G4bool ChangeDaughtersType (EVolume atype)
 

Static Public Member Functions

static const G4LVManagerGetSubInstanceManager ()
 
static void Clean ()
 
static G4VSolidGetSolid (G4LVData &instLVdata)
 
static void SetSolid (G4LVData &instLVdata, G4VSolid *pSolid)
 

Detailed Description

Definition at line 182 of file G4LogicalVolume.hh.

Constructor & Destructor Documentation

◆ G4LogicalVolume() [1/3]

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

Definition at line 68 of file G4LogicalVolume.cc.

75 : fDaughters(0,(G4VPhysicalVolume*)nullptr), fDaughtersVolumeType(kNormal),
76 fOptimise(optimise)
77{
78 // Initialize 'Shadow'/master pointers - for use in copying to workers
79 //
80 fSolid = pSolid;
81 fSensitiveDetector = pSDetector;
82 fFieldManager = pFieldMgr;
83
84 instanceID = subInstanceManager.CreateSubInstance();
85 AssignFieldManager(pFieldMgr);
86
87 G4MT_mass = 0.;
88 G4MT_ccouple = nullptr;
89
90 SetSolid(pSolid);
91 SetMaterial(pMaterial);
92 SetName(name);
93 SetSensitiveDetector(pSDetector);
94 SetUserLimits(pULimits);
95
96 // Initialize 'Shadow' data structure - for use by object persistency
97 //
98 lvdata = new G4LVData();
99 lvdata->fSolid = pSolid;
100 lvdata->fMaterial = pMaterial;
101
102 //
103 // Add to store
104 //
106}
#define G4MT_mass
#define G4MT_ccouple
G4int CreateSubInstance()
G4Material * fMaterial
G4VSolid * fSolid
static void Register(G4LogicalVolume *pVolume)
void SetName(const G4String &pName)
void SetUserLimits(G4UserLimits *pULimits)
void AssignFieldManager(G4FieldManager *fldMgr)
void SetMaterial(G4Material *pMaterial)
void SetSolid(G4VSolid *pSolid)
void SetSensitiveDetector(G4VSensitiveDetector *pSDetector)
@ kNormal
Definition: geomdefs.hh:84

◆ ~G4LogicalVolume()

G4LogicalVolume::~G4LogicalVolume ( )
virtual

Definition at line 134 of file G4LogicalVolume.cc.

135{
136 if (!fLock && fRootRegion) // De-register root region first if not locked
137 { // and flagged as root logical-volume
138 fRegion->RemoveRootLogicalVolume(this, true);
139 }
140 delete lvdata;
142}
static void DeRegister(G4LogicalVolume *pVolume)
void RemoveRootLogicalVolume(G4LogicalVolume *lv, G4bool scan=true)
Definition: G4Region.cc:328

◆ G4LogicalVolume() [2/3]

G4LogicalVolume::G4LogicalVolume ( const G4LogicalVolume )
delete

◆ G4LogicalVolume() [3/3]

G4LogicalVolume::G4LogicalVolume ( __void__ &  )

Definition at line 113 of file G4LogicalVolume.cc.

114 : fDaughters(0,(G4VPhysicalVolume*)nullptr), fName("")
115{
116 instanceID = subInstanceManager.CreateSubInstance();
117
118 SetSensitiveDetector(nullptr); // G4MT_sdetector = nullptr;
119 SetFieldManager(nullptr, false); // G4MT_fmanager = nullptr;
120
121 G4MT_mass = 0.;
122 G4MT_ccouple = nullptr;
123
124 // Add to store
125 //
127}
void SetFieldManager(G4FieldManager *pFieldMgr, G4bool forceToAllDaughters)

Member Function Documentation

◆ AddDaughter()

void G4LogicalVolume::AddDaughter ( G4VPhysicalVolume p)

Definition at line 281 of file G4LogicalVolume.cc.

282{
283 EVolume daughterType = pNewDaughter->VolumeType();
284
285 // The type of the navigation needed is determined by the first daughter
286 //
287 if( fDaughters.empty() )
288 {
289 fDaughtersVolumeType = daughterType;
290 }
291 else
292 {
293 // Check consistency of detector description
294
295 // 1. A replica or parameterised volume can have only one daughter
296 //
297 if( fDaughters[0]->IsReplicated() )
298 {
299 std::ostringstream message;
300 message << "ERROR - Attempt to place a volume in a mother volume"
301 << G4endl
302 << " already containing a replicated volume." << G4endl
303 << " A volume can either contain several placements" << G4endl
304 << " or a unique replica or parameterised volume !" << G4endl
305 << " Mother logical volume: " << GetName() << G4endl
306 << " Placing volume: " << pNewDaughter->GetName()
307 << G4endl;
308 G4Exception("G4LogicalVolume::AddDaughter()", "GeomMgt0002",
309 FatalException, message,
310 "Replica or parameterised volume must be the only daughter!");
311 }
312 else
313 {
314 // 2. Ensure that Placement and External physical volumes do not mix
315 //
316 if( daughterType != fDaughtersVolumeType )
317 {
318 std::ostringstream message;
319 message << "ERROR - Attempt to place a volume in a mother volume"
320 << G4endl
321 << " already containing a different type of volume." << G4endl
322 << " A volume can either contain" << G4endl
323 << " - one or more placements, OR" << G4endl
324 << " - one or more 'external' type physical volumes." << G4endl
325 << " Mother logical volume: " << GetName() << G4endl
326 << " Volume being placed: " << pNewDaughter->GetName()
327 << G4endl;
328 G4Exception("G4LogicalVolume::AddDaughter()", "GeomMgt0002",
329 FatalException, message,
330 "Cannot mix placements and external physical volumes !");
331 }
332 }
333 }
334
335 // Invalidate previous calculation of mass - if any - for all threads
336 //
337 G4MT_mass = 0.;
338 fDaughters.push_back(pNewDaughter);
339
340 G4LogicalVolume* pDaughterLogical = pNewDaughter->GetLogicalVolume();
341
342 // Propagate the Field Manager, if the daughter has no field Manager
343 //
344 G4FieldManager* pDaughterFieldManager = pDaughterLogical->GetFieldManager();
345
346 // Avoid propagating the fieldManager pointer if null
347 // and daughter's one is null as well...
348 //
349 if( (G4MT_fmanager != nullptr ) && (pDaughterFieldManager == nullptr) )
350 {
351 pDaughterLogical->SetFieldManager(G4MT_fmanager, false);
352 }
353 if (fRegion)
354 {
356 fRegion->RegionModified(true);
357 }
358}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
#define G4MT_fmanager
#define G4endl
Definition: G4ios.hh:57
void PropagateRegion()
G4FieldManager * GetFieldManager() const
const G4String & GetName() const
void RegionModified(G4bool flag)
EVolume
Definition: geomdefs.hh:83

Referenced by G4PVDivision::G4PVDivision(), G4PVParameterised::G4PVParameterised(), G4PVPlacement::G4PVPlacement(), G4PVReplica::G4PVReplica(), and G4VExternalPhysicalVolume::G4VExternalPhysicalVolume().

◆ AssignFieldManager()

void G4LogicalVolume::AssignFieldManager ( G4FieldManager fldMgr)

Definition at line 240 of file G4LogicalVolume.cc.

241{
242 G4MT_fmanager= fldMgr;
243 if(G4Threading::IsMasterThread()) { fFieldManager = fldMgr; }
244}
G4bool IsMasterThread()
Definition: G4Threading.cc:124

Referenced by G4LogicalVolume(), InitialiseWorker(), and SetFieldManager().

◆ ChangeDaughtersType()

G4bool G4LogicalVolume::ChangeDaughtersType ( EVolume  atype)

Definition at line 654 of file G4LogicalVolume.cc.

655{
656 G4bool works = false;
657 if( aType == kExternal )
658 {
659 // It is the responsibility of External Navigator to handle types selected
660 //
661 fDaughtersVolumeType = aType;
662 works = true;
663 }
664 else
665 {
666 EVolume expectedVType = DeduceDaughtersType();
667 works = (expectedVType == aType);
668 if ( works )
669 {
670 fDaughtersVolumeType = aType;
671 }
672 }
673 return works;
674}
bool G4bool
Definition: G4Types.hh:86
EVolume DeduceDaughtersType() const
@ kExternal
Definition: geomdefs.hh:87

◆ CharacteriseDaughters()

EVolume G4LogicalVolume::CharacteriseDaughters ( ) const
inline

◆ Clean()

void G4LogicalVolume::Clean ( )
static

Definition at line 198 of file G4LogicalVolume.cc.

199{
200 subInstanceManager.FreeSlave();
201}

Referenced by G4LogicalVolumeStore::~G4LogicalVolumeStore().

◆ ClearDaughters()

void G4LogicalVolume::ClearDaughters ( )

Definition at line 385 of file G4LogicalVolume.cc.

386{
387 fDaughters.erase(fDaughters.cbegin(), fDaughters.cend());
388 if (fRegion != nullptr)
389 {
390 fRegion->RegionModified(true);
391 }
392 G4MT_mass = 0.;
393}

◆ DeduceDaughtersType()

EVolume G4LogicalVolume::DeduceDaughtersType ( ) const
inline

Referenced by ChangeDaughtersType().

◆ GetBiasWeight()

G4double G4LogicalVolume::GetBiasWeight ( ) const
inline

◆ GetDaughter()

◆ GetFastSimulationManager()

◆ GetFieldManager()

◆ GetInstanceID()

G4int G4LogicalVolume::GetInstanceID ( ) const
inline

◆ GetMass()

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

Definition at line 562 of file G4LogicalVolume.cc.

565{
566 // Return the cached non-zero value, if not forced
567 //
568 if ( (G4MT_mass) && (!forced) ) { return G4MT_mass; }
569
570 // Global density and computed mass associated to the logical
571 // volume without considering its daughters
572 //
573 G4Material* logMaterial = parMaterial ? parMaterial : GetMaterial();
574 if (logMaterial == nullptr)
575 {
576 std::ostringstream message;
577 message << "No material associated to the logical volume: "
578 << fName << " !" << G4endl
579 << "Sorry, cannot compute the mass ...";
580 G4Exception("G4LogicalVolume::GetMass()", "GeomMgt0002",
581 FatalException, message);
582 return 0.0;
583 }
584 if ( GetSolid() == nullptr )
585 {
586 std::ostringstream message;
587 message << "No solid is associated to the logical volume: "
588 << fName << " !" << G4endl
589 << "Sorry, cannot compute the mass ...";
590 G4Exception("G4LogicalVolume::GetMass()", "GeomMgt0002",
591 FatalException, message);
592 return 0.0;
593 }
594 G4double globalDensity = logMaterial->GetDensity();
595 G4double motherMass = GetSolid()->GetCubicVolume() * globalDensity;
596 G4double massSum = motherMass;
597
598 // For each daughter in the tree, subtract the mass occupied
599 // and if required by the propagate flag, add the real daughter's
600 // one computed recursively
601
602 for (auto itDau = fDaughters.cbegin(); itDau != fDaughters.cend(); ++itDau)
603 {
604 G4VPhysicalVolume* physDaughter = (*itDau);
605 G4LogicalVolume* logDaughter = physDaughter->GetLogicalVolume();
606 G4double subMass = 0.0;
607 G4VSolid* daughterSolid = nullptr;
608 G4Material* daughterMaterial = nullptr;
609
610 // Compute the mass to subtract and to add for each daughter
611 // considering its multiplicity (i.e. replicated or not) and
612 // eventually its parameterisation (by solid and/or by material)
613 //
614 for (auto i=0; i<physDaughter->GetMultiplicity(); ++i)
615 {
616 G4VPVParameterisation* physParam = physDaughter->GetParameterisation();
617 if (physParam)
618 {
619 daughterSolid = physParam->ComputeSolid(i, physDaughter);
620 daughterSolid->ComputeDimensions(physParam, i, physDaughter);
621 daughterMaterial = physParam->ComputeMaterial(i, physDaughter);
622 }
623 else
624 {
625 daughterSolid = logDaughter->GetSolid();
626 daughterMaterial = logDaughter->GetMaterial();
627 }
628 subMass = daughterSolid->GetCubicVolume() * globalDensity;
629
630 // Subtract the daughter's portion for the mass and, if required,
631 // add the real daughter's mass computed recursively
632 //
633 massSum -= subMass;
634 if (propagate)
635 {
636 massSum += logDaughter->GetMass(true, true, daughterMaterial);
637 }
638 }
639 }
640 G4MT_mass = massSum;
641 return massSum;
642}
double G4double
Definition: G4Types.hh:83
G4VSolid * GetSolid() const
G4double GetMass(G4bool forced=false, G4bool propagate=true, G4Material *parMaterial=nullptr)
G4Material * GetMaterial() const
G4double GetDensity() const
Definition: G4Material.hh:175
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
virtual G4Material * ComputeMaterial(const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=nullptr)
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

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

◆ GetMasterFieldManager()

G4FieldManager * G4LogicalVolume::GetMasterFieldManager ( ) const
inline

◆ GetMasterSensitiveDetector()

G4VSensitiveDetector * G4LogicalVolume::GetMasterSensitiveDetector ( ) const
inline

◆ GetMasterSolid()

G4VSolid * G4LogicalVolume::GetMasterSolid ( ) const
inline

◆ GetMaterial()

◆ GetMaterialCutsCouple()

◆ GetName()

const G4String & G4LogicalVolume::GetName ( ) const
inline

Referenced by G4AdjointCrossSurfChecker::AddanExtSurfaceOfAvolume(), AddDaughter(), G4OpenInventorSceneHandler::AddPrimitive(), G4VtkSceneHandler::AddSolid(), G4VBiasingOperator::AttachTo(), G4SmartVoxelHeader::BuildNodes(), G4SmartVoxelHeader::BuildReplicaVoxels(), G4SmartVoxelHeader::BuildVoxelsWithinLimits(), G4PVParameterised::CheckOverlaps(), G4PVPlacement::CheckOverlaps(), checkVol(), G4tgbVolume::ConstructG4LogVol(), G4tgbVolume::ConstructG4PhysVol(), G4tgbVolume::ConstructG4Volumes(), G4PhysicalVolumeModel::CreateCurrentAttValues(), G3Division::CreatePVReplica(), G4AdjointCrossSurfChecker::CrossingAnInterfaceBetweenTwoVolumes(), G4Radioactivation::DecayIt(), G4RadioactiveDecay::DecayIt(), G4RunManagerKernel::DefineWorldVolume(), G4LogicalVolumeStore::DeRegister(), G4ReflectionFactory::Divide(), G4GDMLReadStructure::DivisionvolRead(), G4GDMLWriteStructure::DivisionvolWrite(), G4TrajectoryDrawByOriginVolume::Draw(), G4tgbVolumeMgr::DumpG4LogVolLeaf(), G4tgbGeometryDumper::DumpLogVol(), G4tgbGeometryDumper::DumpPhysVol(), G4tgbGeometryDumper::DumpPVParameterised(), G4tgbGeometryDumper::DumpPVPlacement(), G4tgbGeometryDumper::DumpPVReplica(), G4TrajectoryOriginVolumeFilter::Evaluate(), G4BuildGeom(), G4PVReplica::G4PVReplica(), G4GDMLRead::GeneratePhysvolName(), G4Navigator::GetLocalExitNormal(), G4ITNavigator1::GetLocalExitNormal(), G4ITNavigator2::GetLocalExitNormal(), G4tgbVolumeMgr::GetTopLogVol(), G4tgbVolumeMgr::GetTopPhysVol(), G4GDMLReadStructure::GetWorldVolume(), G4AdjointCrossSurfChecker::GoingInOrOutOfaVolumeByExtSurface(), G4GDMLWriteParamvol::ParamvolAlgorithmWrite(), G4GDMLReadParamvol::ParamvolRead(), G4GDMLWriteParamvol::ParamvolWrite(), G4GDMLReadStructure::PhysvolRead(), G4GDMLWriteStructure::PhysvolWrite(), G4ReflectionFactory::Place(), G4LogicalVolumeStore::Register(), G4tgbVolumeMgr::RegisterMe(), G4RunManager::ReOptimize(), G4GDMLReadStructure::ReplicaRead(), G4ReflectionFactory::Replicate(), G4GDMLWriteStructure::ReplicavolWrite(), G4ASCIITreeSceneHandler::RequestPrimitives(), G4VoxelSafety::SafetyForVoxelHeader(), G4PolarizedAnnihilationModel::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), G4Region::ScanVolumeTree(), G4RadioactiveDecay::SelectAllVolumes(), G4VVisCommandGeometrySet::Set(), G4VVisCommandGeometrySet::SetLVVisAtts(), G4VisCommandGeometryList::SetNewValue(), G4VisCommandGeometryRestore::SetNewValue(), G4VUserDetectorConstruction::SetSensitiveDetector(), G4GDMLWriteSetup::SetupWrite(), G4PolarizationManager::SetVolumePolarization(), G4GDMLWriteStructure::SkinSurfaceCache(), G4GDMLRead::StripNames(), and G4GDMLWriteStructure::TraverseVolumeTree().

◆ GetNoDaughters()

◆ GetRegion()

◆ GetSensitiveDetector()

◆ GetSmartless()

G4double G4LogicalVolume::GetSmartless ( ) const
inline

◆ GetSolid() [1/2]

G4VSolid * G4LogicalVolume::GetSolid ( ) const

Definition at line 413 of file G4LogicalVolume.cc.

414{
415 return this->GetSolid( subInstanceManager.offset[instanceID] );
416}
static G4GEOM_DLL G4ThreadLocal T * offset

Referenced by G4AdjointCrossSurfChecker::AddanExtSurfaceOfAvolume(), G4VSceneHandler::AddCompound(), G4PseudoScene::AddCompound(), G4GMocrenFileSceneHandler::AddPrimitive(), G4GMocrenFileSceneHandler::AddSolid(), G4ReplicaNavigation::BackLocate(), G4PhantomParameterisation::BuildContainerSolid(), G4SmartVoxelHeader::BuildNodes(), G4SmartVoxelHeader::BuildReplicaVoxels(), G4SmartVoxelHeader::BuildVoxelsWithinLimits(), G4PVParameterised::CheckOverlaps(), G4PVPlacement::CheckOverlaps(), G4GeometryWorkspace::CloneReplicaSolid(), 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(), G4ITNavigator1::ComputeStep(), G4ITNavigator2::ComputeStep(), G4RegularNavigation::ComputeStepSkippingEqualMaterials(), G4tgbVolume::ConstructG4PhysVol(), G4TheRayTracer::CreateBitMap(), G4PhysicalVolumeModel::CreateCurrentAttValues(), G3Division::CreatePVReplica(), G4AdjointPosOnPhysVolGenerator::DefinePhysicalVolume(), G4LogicalVolumeModel::DescribeYourselfTo(), G4VFieldModel::DescribeYourselfTo(), G4VisManager::Draw(), G4tgbGeometryDumper::DumpLogVol(), G4Mesh::G4Mesh(), G4RTPrimaryGeneratorAction::GeneratePrimaries(), G4Navigator::GetGlobalExitNormal(), G4Navigator::GetLocalExitNormal(), G4ITNavigator1::GetLocalExitNormal(), G4ITNavigator2::GetLocalExitNormal(), GetMass(), G4TransportationManager::GetParallelWorld(), G4ITTransportationManager::GetParallelWorld(), GetSolid(), G4BOptnForceCommonTruncatedExp::Initialize(), G4ITNavigator1::LocateGlobalPointAndSetup(), G4Navigator::LocateGlobalPointAndSetup(), G4GDMLWriteParamvol::ParametersWrite(), G4NeutrinoElectronProcess::PostStepDoIt(), G4VXTRenergyLoss::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4TauNeutrinoNucleusProcess::PostStepDoIt(), G4NavigationLogger::PreComputeStepLog(), G4PSFlatSurfaceCurrent::ProcessHits(), G4PSFlatSurfaceFlux::ProcessHits(), G4PSSphereSurfaceFlux::ProcessHits(), G4PSVolumeFlux::ProcessHits(), G4ITNavigator2::RecheckDistanceToCurrentBoundary(), G4NavigationLogger::ReportOutsideMother(), G4NavigationLogger::ReportVolumeAndIntersection(), G4VoxelSafety::SafetyForVoxelHeader(), G4VoxelSafety::SafetyForVoxelNode(), G4VisCommandsTouchable::SetNewValue(), G4RTPrimaryGeneratorAction::SetUp(), G4GeomTestVolume::TestOverlapInTree(), and G4GDMLWriteStructure::TraverseVolumeTree().

◆ GetSolid() [2/2]

G4VSolid * G4LogicalVolume::GetSolid ( G4LVData instLVdata)
static

Definition at line 408 of file G4LogicalVolume.cc.

409{
410 return instLVdata.fSolid;
411}

◆ GetSubInstanceManager()

const G4LVManager & G4LogicalVolume::GetSubInstanceManager ( )
static

Definition at line 222 of file G4LogicalVolume.cc.

223{
224 return subInstanceManager;
225}

Referenced by G4GeometryWorkspace::G4GeometryWorkspace().

◆ GetUserLimits()

◆ GetVisAttributes()

◆ GetVoxelHeader()

◆ InitialiseWorker()

void G4LogicalVolume::InitialiseWorker ( G4LogicalVolume ptrMasterObject,
G4VSolid pSolid,
G4VSensitiveDetector pSDetector 
)

Definition at line 164 of file G4LogicalVolume.cc.

168{
169 subInstanceManager.SlaveCopySubInstanceArray();
170
171 SetSolid(pSolid);
172 SetSensitiveDetector(pSDetector); // How this object is available now ?
173 AssignFieldManager(fFieldManager);
174 // Should be set - but a per-thread copy is not available yet
175 // Must not call SetFieldManager(), which propagates FieldMgr
176
177#ifdef CLONE_FIELD_MGR
178 // Create a field FieldManager by cloning
179 //
180 G4FieldManager workerFldMgr = fFieldManager->GetWorkerClone(G4bool* created);
181 if( created || (GetFieldManager() != workerFldMgr) )
182 {
183 SetFieldManager(fFieldManager, false); // which propagates FieldMgr
184 }
185 else
186 {
187 // Field manager existed and is equal to current one
188 //
189 AssignFieldManager(workerFldMgr);
190 }
191#endif
192}
void SlaveCopySubInstanceArray()

Referenced by G4GeometryWorkspace::CloneReplicaSolid(), and G4GeometryWorkspace::InitialisePhysicalVolumes().

◆ IsAncestor()

G4bool G4LogicalVolume::IsAncestor ( const G4VPhysicalVolume p) const

Definition at line 510 of file G4LogicalVolume.cc.

511{
512 G4bool isDaughter = IsDaughter(aVolume);
513 if (!isDaughter)
514 {
515 for (auto itDau = fDaughters.cbegin(); itDau != fDaughters.cend(); ++itDau)
516 {
517 isDaughter = (*itDau)->GetLogicalVolume()->IsAncestor(aVolume);
518 if (isDaughter) break;
519 }
520 }
521 return isDaughter;
522}
G4bool IsDaughter(const G4VPhysicalVolume *p) const

◆ IsDaughter()

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

Referenced by IsAncestor().

◆ IsExtended()

G4bool G4LogicalVolume::IsExtended ( ) const
virtual

Reimplemented in G4LogicalCrystalVolume.

Definition at line 250 of file G4LogicalVolume.cc.

251{
252 return false;
253}

◆ IsRegion()

G4bool G4LogicalVolume::IsRegion ( ) const
inline

◆ IsRootRegion()

G4bool G4LogicalVolume::IsRootRegion ( ) const
inline

◆ IsToOptimise()

G4bool G4LogicalVolume::IsToOptimise ( ) const
inline

◆ Lock()

void G4LogicalVolume::Lock ( )
inline

◆ operator=()

G4LogicalVolume & G4LogicalVolume::operator= ( const G4LogicalVolume )
delete

◆ operator==()

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

◆ PropagateRegion()

void G4LogicalVolume::PropagateRegion ( )
inline

Referenced by AddDaughter().

◆ RemoveDaughter()

void G4LogicalVolume::RemoveDaughter ( const G4VPhysicalVolume p)

Definition at line 364 of file G4LogicalVolume.cc.

365{
366 for (auto i=fDaughters.cbegin(); i!=fDaughters.cend(); ++i )
367 {
368 if (**i==*p)
369 {
370 fDaughters.erase(i);
371 break;
372 }
373 }
374 if (fRegion)
375 {
376 fRegion->RegionModified(true);
377 }
378 G4MT_mass = 0.;
379}

Referenced by G4PhysicalVolumeStore::DeRegister().

◆ ResetMass()

void G4LogicalVolume::ResetMass ( )

Definition at line 399 of file G4LogicalVolume.cc.

400{
401 G4MT_mass= 0.0;
402}

Referenced by SetSolid().

◆ SetBiasWeight()

void G4LogicalVolume::SetBiasWeight ( G4double  w)
inline

◆ SetFieldManager()

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

Definition at line 260 of file G4LogicalVolume.cc.

262{
263 AssignFieldManager(pNewFieldMgr);
264
265 auto NoDaughters = GetNoDaughters();
266 while ( (NoDaughters--)>0 )
267 {
268 G4LogicalVolume* DaughterLogVol;
269 DaughterLogVol = GetDaughter(NoDaughters)->GetLogicalVolume();
270 if ( forceAllDaughters || (DaughterLogVol->GetFieldManager() == nullptr) )
271 {
272 DaughterLogVol->SetFieldManager(pNewFieldMgr, forceAllDaughters);
273 }
274 }
275}
std::size_t GetNoDaughters() const
G4VPhysicalVolume * GetDaughter(const std::size_t i) const

Referenced by AddDaughter(), G4VUserDetectorConstruction::CloneF(), G4LogicalVolume(), InitialiseWorker(), SetFieldManager(), and G4WorkerThread::UpdateGeometryAndPhysicsVectorFromMaster().

◆ SetMaterial()

void G4LogicalVolume::SetMaterial ( G4Material pMaterial)

Definition at line 448 of file G4LogicalVolume.cc.

449{
450 G4MT_material = pMaterial;
451 G4MT_mass = 0.0;
452}

Referenced by G4LogicalVolume(), and G4ScoreSplittingProcess::PostStepDoIt().

◆ SetMaterialCutsCouple()

void G4LogicalVolume::SetMaterialCutsCouple ( G4MaterialCutsCouple cuts)

Definition at line 497 of file G4LogicalVolume.cc.

498{
499 G4MT_ccouple = cuts;
500}

◆ SetName()

void G4LogicalVolume::SetName ( const G4String pName)

Definition at line 148 of file G4LogicalVolume.cc.

149{
150 fName = pName;
152}
void SetMapValid(G4bool val)
static G4LogicalVolumeStore * GetInstance()

Referenced by G4LogicalVolume(), and G4GDMLRead::StripNames().

◆ SetOptimisation()

void G4LogicalVolume::SetOptimisation ( G4bool  optim)
inline

◆ SetRegion()

void G4LogicalVolume::SetRegion ( G4Region reg)
inline

◆ SetRegionRootFlag()

void G4LogicalVolume::SetRegionRootFlag ( G4bool  rreg)
inline

◆ SetSensitiveDetector()

◆ SetSmartless()

void G4LogicalVolume::SetSmartless ( G4double  s)
inline

◆ SetSolid() [1/2]

void G4LogicalVolume::SetSolid ( G4LVData instLVdata,
G4VSolid pSolid 
)
static

Definition at line 429 of file G4LogicalVolume.cc.

430{
431 instLVdata.fSolid = pSolid;
432 instLVdata.fMass = 0.0;
433}
G4double fMass

◆ SetSolid() [2/2]

◆ SetUserLimits()

void G4LogicalVolume::SetUserLimits ( G4UserLimits pULimits)
inline

Referenced by G4LogicalVolume().

◆ SetVisAttributes() [1/2]

void G4LogicalVolume::SetVisAttributes ( const G4VisAttributes VA)

Definition at line 680 of file G4LogicalVolume.cc.

681{
682 if (G4Threading::IsWorkerThread()) return;
683 fVisAttributes = std::make_shared<const G4VisAttributes>(VA);
684}
G4bool IsWorkerThread()
Definition: G4Threading.cc:123

◆ SetVisAttributes() [2/2]

◆ SetVoxelHeader()

void G4LogicalVolume::SetVoxelHeader ( G4SmartVoxelHeader pVoxel)
inline

◆ TerminateWorker()

void G4LogicalVolume::TerminateWorker ( G4LogicalVolume ptrMasterObject)

Definition at line 211 of file G4LogicalVolume.cc.

213{
214}

Referenced by G4GeometryWorkspace::DestroyWorkspace().

◆ TotalVolumeEntities()

G4int G4LogicalVolume::TotalVolumeEntities ( ) const

Definition at line 531 of file G4LogicalVolume.cc.

532{
533 G4int vols = 1;
534 for (auto itDau = fDaughters.cbegin(); itDau != fDaughters.cend(); ++itDau)
535 {
536 G4VPhysicalVolume* physDaughter = (*itDau);
537 vols += physDaughter->GetMultiplicity()
538 *physDaughter->GetLogicalVolume()->TotalVolumeEntities();
539 }
540 return vols;
541}
int G4int
Definition: G4Types.hh:85
G4int TotalVolumeEntities() const

Referenced by TotalVolumeEntities().

◆ UpdateMaterial()

void G4LogicalVolume::UpdateMaterial ( G4Material pMaterial)

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