Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Static Public Member Functions
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
 
EVolume CharacteriseDaughters () 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)
 
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__ &)
 
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)
 

Static Public Member Functions

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

Detailed Description

Definition at line 187 of file G4LogicalVolume.hh.

Constructor & Destructor Documentation

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

Definition at line 138 of file G4LogicalVolume.cc.

References AssignFieldManager(), G4GeomSplitter< T >::CreateSubInstance(), G4LogicalVolumeStore::Register(), SetMaterial(), SetName(), SetSensitiveDetector(), SetSolid(), and SetUserLimits().

145  : fDaughters(0,(G4VPhysicalVolume*)0),
146  fVoxel(0), fOptimise(optimise), fRootRegion(false), fLock(false),
147  fSmartless(2.), fVisAttributes(0), fRegion(0)
148 {
149  // Initialize 'Shadow'/master pointers - for use in copying to workers
150  fSolid = pSolid;
151  fSensitiveDetector = pSDetector;
152  fFieldManager = pFieldMgr;
153 
154  instanceID = subInstanceManager.CreateSubInstance();
155  AssignFieldManager(pFieldMgr); // G4MT_fmanager = pFieldMgr;
156 
157  // fMasterFieldMgr= pFieldMgr;
158  G4MT_mass = 0.;
159  G4MT_ccouple = 0;
160 
161  SetSolid(pSolid);
162  SetMaterial(pMaterial);
163  SetName(name);
164  SetSensitiveDetector(pSDetector);
165  SetUserLimits(pULimits);
166  //
167  // Add to store
168  //
170 }
void SetUserLimits(G4UserLimits *pULimits)
void SetSolid(G4VSolid *pSolid)
G4int CreateSubInstance()
void AssignFieldManager(G4FieldManager *fldMgr)
void SetName(const G4String &pName)
void SetMaterial(G4Material *pMaterial)
static void Register(G4LogicalVolume *pVolume)
void SetSensitiveDetector(G4VSensitiveDetector *pSDetector)
G4LogicalVolume::~G4LogicalVolume ( )

Definition at line 206 of file G4LogicalVolume.cc.

References G4LogicalVolumeStore::DeRegister(), and G4Region::RemoveRootLogicalVolume().

207 {
208  if (!fLock && fRootRegion) // De-register root region first if not locked
209  { // and flagged as root logical-volume
210  fRegion->RemoveRootLogicalVolume(this, true);
211  }
213 }
static void DeRegister(G4LogicalVolume *pVolume)
void RemoveRootLogicalVolume(G4LogicalVolume *lv, G4bool scan=true)
Definition: G4Region.cc:283
G4LogicalVolume::G4LogicalVolume ( __void__ &  )

Definition at line 177 of file G4LogicalVolume.cc.

References G4GeomSplitter< T >::CreateSubInstance(), G4LogicalVolumeStore::Register(), SetFieldManager(), SetMaterial(), SetSensitiveDetector(), and SetSolid().

178  : fDaughters(0,(G4VPhysicalVolume*)0),
179  fName(""), fUserLimits(0),
180  fVoxel(0), fOptimise(true), fRootRegion(false), fLock(false),
181  fSmartless(2.), fVisAttributes(0), fRegion(0), fBiasWeight(0.),
182  fSolid(0), fSensitiveDetector(0), fFieldManager(0)
183 {
184  instanceID = subInstanceManager.CreateSubInstance();
185 
186  // G4MT_solid = 0,
187  SetSolid(0);
188  SetSensitiveDetector(0); // G4MT_sdetector = 0;
189  SetFieldManager(0, false); // G4MT_fmanager = 0;
190  SetMaterial(0); // G4MT_material = 0;
191  // this->SetMass(0); //
192  G4MT_mass = 0.;
193  // this->SetCutsCouple(0);
194  G4MT_ccouple = 0;
195 
196  // Add to store
197  //
199 }
void SetSolid(G4VSolid *pSolid)
G4int CreateSubInstance()
void SetFieldManager(G4FieldManager *pFieldMgr, G4bool forceToAllDaughters)
void SetMaterial(G4Material *pMaterial)
static void Register(G4LogicalVolume *pVolume)
void SetSensitiveDetector(G4VSensitiveDetector *pSDetector)

Member Function Documentation

void G4LogicalVolume::AddDaughter ( G4VPhysicalVolume p)
inline
void G4LogicalVolume::AssignFieldManager ( G4FieldManager fldMgr)
inline
EVolume G4LogicalVolume::CharacteriseDaughters ( ) const
inline
void G4LogicalVolume::ClearDaughters ( )
inline

Referenced by export_G4LogicalVolume().

G4double G4LogicalVolume::GetBiasWeight ( ) const
inline

Referenced by export_G4LogicalVolume().

G4VPhysicalVolume* G4LogicalVolume::GetDaughter ( const G4int  i) const
inline
G4FastSimulationManager* G4LogicalVolume::GetFastSimulationManager ( ) const
inline
G4FieldManager* G4LogicalVolume::GetFieldManager ( ) const
inline
G4int G4LogicalVolume::GetInstanceID ( ) const
inline
G4double G4LogicalVolume::GetMass ( G4bool  forced = false,
G4bool  propagate = true,
G4Material parMaterial = 0 
)

Definition at line 301 of file G4LogicalVolume.cc.

References G4VSolid::ComputeDimensions(), G4VPVParameterisation::ComputeMaterial(), G4VPVParameterisation::ComputeSolid(), FatalException, G4endl, G4Exception(), G4VSolid::GetCubicVolume(), G4Material::GetDensity(), G4VPhysicalVolume::GetLogicalVolume(), GetMass(), GetMaterial(), G4VPhysicalVolume::GetMultiplicity(), G4VPhysicalVolume::GetParameterisation(), and GetSolid().

Referenced by B1RunAction::EndOfRunAction(), B1ConRunAction::EndOfRunAction(), export_G4LogicalVolume(), GetMass(), and G4ASCIITreeSceneHandler::RequestPrimitives().

304 {
305  // Return the cached non-zero value, if not forced
306  //
307  if ( (G4MT_mass) && (!forced) ) return G4MT_mass;
308 
309  // Global density and computed mass associated to the logical
310  // volume without considering its daughters
311  //
312  G4Material* logMaterial = parMaterial ? parMaterial : GetMaterial(); // G4MT_material;
313  if (!logMaterial)
314  {
315  std::ostringstream message;
316  message << "No material associated to the logical volume: " << fName << " !"
317  << G4endl
318  << "Sorry, cannot compute the mass ...";
319  G4Exception("G4LogicalVolume::GetMass()", "GeomMgt0002",
320  FatalException, message);
321  return 0;
322  }
323  if (! GetSolid() ) // !G4MT_solid)
324  {
325  std::ostringstream message;
326  message << "No solid is associated to the logical volume: " << fName << " !"
327  << G4endl
328  << "Sorry, cannot compute the mass ...";
329  G4Exception("G4LogicalVolume::GetMass()", "GeomMgt0002",
330  FatalException, message);
331  return 0;
332  }
333  G4double globalDensity = logMaterial->GetDensity();
334  G4double motherMass= GetSolid()->GetCubicVolume() * globalDensity;
335 
336  // G4MT_mass =
337  // SetMass( motherMmass );
338  G4double massSum= motherMass;
339 
340  // For each daughter in the tree, subtract the mass occupied
341  // and if required by the propagate flag, add the real daughter's
342  // one computed recursively
343 
344  for (G4PhysicalVolumeList::const_iterator itDau = fDaughters.begin();
345  itDau != fDaughters.end(); itDau++)
346  {
347  G4VPhysicalVolume* physDaughter = (*itDau);
348  G4LogicalVolume* logDaughter = physDaughter->GetLogicalVolume();
349  G4double subMass=0.;
350  G4VSolid* daughterSolid = 0;
351  G4Material* daughterMaterial = 0;
352 
353  // Compute the mass to subtract and to add for each daughter
354  // considering its multiplicity (i.e. replicated or not) and
355  // eventually its parameterisation (by solid and/or by material)
356  //
357  for (G4int i=0; i<physDaughter->GetMultiplicity(); i++)
358  {
360  physParam = physDaughter->GetParameterisation();
361  if (physParam)
362  {
363  daughterSolid = physParam->ComputeSolid(i, physDaughter);
364  daughterSolid->ComputeDimensions(physParam, i, physDaughter);
365  daughterMaterial = physParam->ComputeMaterial(i, physDaughter);
366  }
367  else
368  {
369  daughterSolid = logDaughter->GetSolid();
370  daughterMaterial = logDaughter->GetMaterial();
371  }
372  subMass = daughterSolid->GetCubicVolume() * globalDensity;
373 
374  // Subtract the daughter's portion for the mass and, if required,
375  // add the real daughter's mass computed recursively
376  //
377  massSum -= subMass;
378  if (propagate)
379  {
380  massSum += logDaughter->GetMass(true, true, daughterMaterial);
381  }
382  }
383  }
384  G4MT_mass= massSum;
385  return massSum;
386 }
virtual G4Material * ComputeMaterial(const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=0)
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4VSolid.cc:137
G4Material * GetMaterial() const
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
virtual G4double GetCubicVolume()
Definition: G4VSolid.cc:188
G4double GetDensity() const
Definition: G4Material.hh:178
int G4int
Definition: G4Types.hh:78
virtual G4VPVParameterisation * GetParameterisation() const =0
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4LogicalVolume * GetLogicalVolume() const
virtual G4int GetMultiplicity() const
#define G4endl
Definition: G4ios.hh:61
G4double GetMass(G4bool forced=false, G4bool propagate=true, G4Material *parMaterial=0)
double G4double
Definition: G4Types.hh:76
G4VSolid * GetSolid() const
G4FieldManager* G4LogicalVolume::GetMasterFieldManager ( ) const
inline
G4VSensitiveDetector* G4LogicalVolume::GetMasterSensitiveDetector ( ) const
inline
G4VSolid* G4LogicalVolume::GetMasterSolid ( ) const
inline
G4Material* G4LogicalVolume::GetMaterial ( ) const
inline

Referenced by CML2ExpVoxels::add(), G4Track::CalculateVelocityForOpticalPhoton(), G4VPVParameterisation::ComputeMaterial(), CellParameterisation::ComputeMaterial(), G4MIRDLiver::Construct(), G4MIRDLeftLeg::Construct(), G4MIRDLeftLegBone::Construct(), G4MIRDRightClavicle::Construct(), G4MIRDRightKidney::Construct(), G4MIRDRightLeg::Construct(), G4MIRDRightLegBone::Construct(), G4MIRDLeftOvary::Construct(), G4MIRDRightScapula::Construct(), G4MIRDLeftAdrenal::Construct(), G4MIRDLeftTeste::Construct(), G4MIRDSkull::Construct(), G4MIRDSmallIntestine::Construct(), G4MIRDLeftScapula::Construct(), G4MIRDThyroid::Construct(), G4MIRDLowerLargeIntestine::Construct(), G4MIRDUpperLargeIntestine::Construct(), G4MIRDLeftClavicle::Construct(), G4MIRDMaleGenitalia::Construct(), G4MIRDPelvis::Construct(), G4MIRDRightAdrenal::Construct(), G4MIRDLeftKidney::Construct(), G4MIRDRightLung::Construct(), G4MIRDRightOvary::Construct(), G4MIRDLeftLung::Construct(), G4MIRDRightTeste::Construct(), G4MIRDBrain::Construct(), G4MIRDSpleen::Construct(), G4MIRDStomach::Construct(), G4MIRDThymus::Construct(), G4MIRDHeart::Construct(), G4MIRDRightBreast::Construct(), G4MIRDTrunk::Construct(), G4MIRDMiddleLowerSpine::Construct(), G4MIRDUpperSpine::Construct(), G4MIRDUterus::Construct(), G4MIRDPancreas::Construct(), G4MIRDRibCage::Construct(), G4MIRDLeftBreast::Construct(), G4MIRDRightArmBone::Construct(), G4MIRDHead::Construct(), G4MIRDUrinaryBladder::Construct(), G4MIRDLeftArmBone::Construct(), G4tgbGeometryDumper::DumpLogVol(), export_G4LogicalVolume(), GetMass(), PhysicsList::GetRange(), G4LatticeManager::LoadLattice(), G4Transportation::PostStepDoIt(), G4CoupledTransportation::PostStepDoIt(), G4MonopoleTransportation::PostStepDoIt(), G4VXTRenergyLoss::PostStepDoIt(), G4ITTransportation::PostStepDoIt(), XrayFluoMercuryDetectorConstruction::PrintApparateParameters(), XrayFluoPlaneDetectorConstruction::PrintApparateParameters(), XrayFluoDetectorConstruction::PrintApparateParameters(), CML2SDWithVoxels::ProcessHits(), G4LatticeManager::RegisterLattice(), G4Region::ScanVolumeTree(), G4GDMLWriteStructure::TraverseVolumeTree(), and SteppingAction::UserSteppingAction().

const G4MaterialCutsCouple* G4LogicalVolume::GetMaterialCutsCouple ( ) const
inline
G4String G4LogicalVolume::GetName ( ) const
inline

Referenced by G4AdjointCrossSurfChecker::AddanExtSurfaceOfAvolume(), G4HepRepSceneHandler::AddSolid(), G4VBiasingOperator::AttachTo(), G4SmartVoxelHeader::BuildNodes(), G4SmartVoxelHeader::BuildReplicaVoxels(), G4SmartVoxelHeader::BuildVoxelsWithinLimits(), G4PVParameterised::CheckOverlaps(), G4PVPlacement::CheckOverlaps(), checkVol(), GB01DetectorConstruction::Construct(), G3toG4DetectorConstruction::Construct(), G4tgbVolume::ConstructG4LogVol(), G4tgbVolume::ConstructG4PhysVol(), G4tgbVolume::ConstructG4Volumes(), CCalG4Hall::constructIn(), CCalG4Hcal::constructIn(), CCalG4Ecal::constructIn(), CCalG4Hcal::constructScintillatorLayer(), GB02DetectorConstruction::ConstructSDandField(), RE01CalorimeterHit::CreateAttValues(), B5HodoscopeHit::CreateAttValues(), B5EmCalorimeterHit::CreateAttValues(), 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(), export_G4LogicalVolume(), G4BuildGeom(), G4PVReplica::G4PVReplica(), G4GDMLRead::GeneratePhysvolName(), G4Navigator::GetLocalExitNormal(), G4ITNavigator::GetLocalExitNormal(), G4eplusPolarizedAnnihilation::GetMeanFreePath(), G4PolarizedCompton::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(), G4HumanPhantomSD::ProcessHits(), G4tgbVolumeMgr::RegisterMe(), CCalSensitiveDetectors::registerVolume(), G4RunManager::ReOptimize(), 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(), DicomDetectorConstruction::SetScorer(), CCalSensitiveDetectors::setSensitive(), G4GDMLWriteSetup::SetupWrite(), CCalG4Able::setVisType(), G4PolarizationManager::SetVolumePolarization(), G4GDMLWriteStructure::SkinSurfaceCache(), G4GDMLRead::StripNames(), and G4GDMLWriteStructure::TraverseVolumeTree().

G4int G4LogicalVolume::GetNoDaughters ( ) const
inline
G4Region* G4LogicalVolume::GetRegion ( ) const
inline
G4VSensitiveDetector* G4LogicalVolume::GetSensitiveDetector ( ) const
inline
G4double G4LogicalVolume::GetSmartless ( ) const
inline
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(), G4GeometryWorkspace::CloneParameterisedSolids(), G4GeometryWorkspace::CloneReplicaSolid(), G4VoxelSafety::ComputeSafety(), G4ParameterisedNavigation::ComputeSafety(), G4NormalNavigation::ComputeSafety(), G4VoxelNavigation::ComputeSafety(), G4ReplicaNavigation::ComputeSafety(), G4PhantomParameterisation::ComputeSolid(), G4VPVParameterisation::ComputeSolid(), G4VNestedParameterisation::ComputeSolid(), G4ParameterisedNavigation::ComputeStep(), G4NormalNavigation::ComputeStep(), G4Navigator::ComputeStep(), G4VoxelNavigation::ComputeStep(), G4ReplicaNavigation::ComputeStep(), G4ITNavigator::ComputeStep(), G4RegularNavigation::ComputeStepSkippingEqualMaterials(), G4PSDoseDeposit::ComputeVolume(), G4PSPassageCellFlux::ComputeVolume(), G4PSCellFlux::ComputeVolume(), G4MIRDLiver::Construct(), G4MIRDLeftTeste::Construct(), G4MIRDLowerLargeIntestine::Construct(), G4MIRDMaleGenitalia::Construct(), G4MIRDPelvis::Construct(), G4MIRDLeftAdrenal::Construct(), G4MIRDRightAdrenal::Construct(), G4MIRDRightScapula::Construct(), G4MIRDRightClavicle::Construct(), G4MIRDSkull::Construct(), G4MIRDRightKidney::Construct(), G4MIRDRightLeg::Construct(), G4MIRDRightLegBone::Construct(), G4MIRDLeftClavicle::Construct(), G4MIRDLeftKidney::Construct(), G4MIRDSmallIntestine::Construct(), G4MIRDLeftLeg::Construct(), G4MIRDLeftLegBone::Construct(), G4MIRDThyroid::Construct(), G4MIRDUpperLargeIntestine::Construct(), G4MIRDLeftOvary::Construct(), G4MIRDLeftScapula::Construct(), G4MIRDBrain::Construct(), G4MIRDRightBreast::Construct(), G4MIRDRightTeste::Construct(), G4MIRDRightOvary::Construct(), G4MIRDPancreas::Construct(), G4MIRDMiddleLowerSpine::Construct(), G4MIRDRibCage::Construct(), G4MIRDHeart::Construct(), G4MIRDRightLung::Construct(), G4MIRDLeftBreast::Construct(), G4MIRDSpleen::Construct(), G4MIRDStomach::Construct(), G4MIRDThymus::Construct(), G4MIRDTrunk::Construct(), G4MIRDLeftLung::Construct(), G4MIRDUterus::Construct(), G4MIRDUpperSpine::Construct(), G4MIRDHead::Construct(), G4MIRDRightArmBone::Construct(), G4MIRDLeftArmBone::Construct(), G4MIRDUrinaryBladder::Construct(), G4tgbVolume::ConstructG4PhysVol(), G4TheRayTracer::CreateBitMap(), G4PhysicalVolumeModel::CreateCurrentAttValues(), G3Division::CreatePVReplica(), G4AdjointPosOnPhysVolGenerator::DefinePhysicalVolume(), G4VisManager::Draw(), G4tgbGeometryDumper::DumpLogVol(), F04SimpleSolenoid::F04SimpleSolenoid(), G4RTPrimaryGeneratorAction::GeneratePrimaries(), B2PrimaryGeneratorAction::GeneratePrimaries(), B4PrimaryGeneratorAction::GeneratePrimaries(), B1PrimaryGeneratorAction::GeneratePrimaries(), G4Navigator::GetLocalExitNormal(), G4ITNavigator::GetLocalExitNormal(), GetMass(), G4TransportationManager::GetParallelWorld(), G4BOptnForceCommonTruncatedExp::Initialize(), G4Navigator::LocateGlobalPointAndSetup(), G4ITNavigator::LocateGlobalPointAndSetup(), G4GDMLWriteParamvol::ParametersWrite(), G4VXTRenergyLoss::PostStepDoIt(), G4NavigationLogger::PreComputeStepLog(), G4PSSphereSurfaceCurrent::ProcessHits(), G4PSCylinderSurfaceFlux::ProcessHits(), G4PSCylinderSurfaceCurrent::ProcessHits(), G4PSFlatSurfaceCurrent::ProcessHits(), G4PSSphereSurfaceFlux::ProcessHits(), G4PSFlatSurfaceFlux::ProcessHits(), G4VoxelSafety::SafetyForVoxelHeader(), G4VoxelSafety::SafetyForVoxelNode(), G4RTPrimaryGeneratorAction::SetUp(), CexmcPhysicsList< BasePhysics, StudiedPhysics, ProductionModel >::SetupConstructionHook(), and G4GDMLWriteStructure::TraverseVolumeTree().

static G4VSolid* G4LogicalVolume::GetSolid ( G4LVData instLVdata)
inlinestatic
const G4LVManager & G4LogicalVolume::GetSubInstanceManager ( )
static

Definition at line 127 of file G4LogicalVolume.cc.

Referenced by G4GeometryWorkspace::G4GeometryWorkspace().

128 {
129  return subInstanceManager;
130 }
G4UserLimits* G4LogicalVolume::GetUserLimits ( ) const
inline
const G4VisAttributes* G4LogicalVolume::GetVisAttributes ( ) const
inline
G4SmartVoxelHeader* G4LogicalVolume::GetVoxelHeader ( ) const
inline
void G4LogicalVolume::InitialiseWorker ( G4LogicalVolume ptrMasterObject,
G4VSolid pSolid,
G4VSensitiveDetector pSDetector 
)

Definition at line 83 of file G4LogicalVolume.cc.

References AssignFieldManager(), GetFieldManager(), SetFieldManager(), SetSensitiveDetector(), SetSolid(), and G4GeomSplitter< T >::SlaveCopySubInstanceArray().

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

86 {
87  subInstanceManager.SlaveCopySubInstanceArray();
88 
89  SetSolid(pSolid);
90  SetSensitiveDetector(pSDetector); // How this object is available now ?
91  AssignFieldManager(fFieldManager); // Should be set - but a per-thread copy is not available yet
92  // G4MT_fmanager= fFieldManager;
93  // Must not call SetFieldManager(fFieldManager, false); which propagates FieldMgr
94 
95 #ifdef CLONE_FIELD_MGR
96  // Create a field FieldManager by cloning
97  G4FieldManager workerFldMgr= fFieldManager->GetWorkerClone(G4bool* created);
98  if( created || (GetFieldManager()!=workerFldMgr) )
99  {
100  SetFieldManager(fFieldManager, false); // which propagates FieldMgr
101  }else{
102  // Field manager existed and is equal to current one
103  AssignFieldManager(workerFldMgr);
104  }
105 #endif
106 }
void SetSolid(G4VSolid *pSolid)
void SetFieldManager(G4FieldManager *pFieldMgr, G4bool forceToAllDaughters)
void AssignFieldManager(G4FieldManager *fldMgr)
G4FieldManager * GetFieldManager() const
bool G4bool
Definition: G4Types.hh:79
void SlaveCopySubInstanceArray()
void SetSensitiveDetector(G4VSensitiveDetector *pSDetector)
G4bool G4LogicalVolume::IsAncestor ( const G4VPhysicalVolume p) const

Definition at line 247 of file G4LogicalVolume.cc.

References IsDaughter().

Referenced by export_G4LogicalVolume().

248 {
249  G4bool isDaughter = IsDaughter(aVolume);
250  if (!isDaughter)
251  {
252  for (G4PhysicalVolumeList::const_iterator itDau = fDaughters.begin();
253  itDau != fDaughters.end(); itDau++)
254  {
255  isDaughter = (*itDau)->GetLogicalVolume()->IsAncestor(aVolume);
256  if (isDaughter) break;
257  }
258  }
259  return isDaughter;
260 }
bool G4bool
Definition: G4Types.hh:79
G4bool IsDaughter(const G4VPhysicalVolume *p) const
G4bool G4LogicalVolume::IsDaughter ( const G4VPhysicalVolume p) const
inline
G4bool G4LogicalVolume::IsRegion ( ) const
inline

Referenced by export_G4LogicalVolume().

G4bool G4LogicalVolume::IsRootRegion ( ) const
inline
G4bool G4LogicalVolume::IsToOptimise ( ) const
inline

Referenced by export_G4LogicalVolume().

void G4LogicalVolume::Lock ( )
inline
G4bool G4LogicalVolume::operator== ( const G4LogicalVolume lv) const
void G4LogicalVolume::PropagateRegion ( )
inline

Referenced by export_G4LogicalVolume().

void G4LogicalVolume::RemoveDaughter ( const G4VPhysicalVolume p)
inline
void G4LogicalVolume::ResetMass ( )
void G4LogicalVolume::SetBiasWeight ( G4double  w)
inline

Referenced by export_G4LogicalVolume().

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

Definition at line 220 of file G4LogicalVolume.cc.

References AssignFieldManager(), GetDaughter(), GetFieldManager(), G4VPhysicalVolume::GetLogicalVolume(), GetNoDaughters(), and SetFieldManager().

Referenced by G4VUserDetectorConstruction::CloneF(), B5DetectorConstruction::ConstructSDandField(), F03DetectorConstruction::ConstructSDandField(), export_G4LogicalVolume(), G4LogicalVolume(), InitialiseWorker(), SetFieldManager(), and G4WorkerThread::UpdateGeometryAndPhysicsVectorFromMaster().

222 {
223  // G4MT_fmanager = pNewFieldMgr;
224  AssignFieldManager(pNewFieldMgr);
225 
226  G4int NoDaughters = GetNoDaughters();
227  while ( (NoDaughters--)>0 )
228  {
229  G4LogicalVolume* DaughterLogVol;
230  DaughterLogVol = GetDaughter(NoDaughters)->GetLogicalVolume();
231  if ( forceAllDaughters || (DaughterLogVol->GetFieldManager() == 0) )
232  {
233  DaughterLogVol->SetFieldManager(pNewFieldMgr, forceAllDaughters);
234  }
235  }
236 }
G4VPhysicalVolume * GetDaughter(const G4int i) const
int G4int
Definition: G4Types.hh:78
void SetFieldManager(G4FieldManager *pFieldMgr, G4bool forceToAllDaughters)
void AssignFieldManager(G4FieldManager *fldMgr)
G4FieldManager * GetFieldManager() const
G4int GetNoDaughters() const
G4LogicalVolume * GetLogicalVolume() const
void G4LogicalVolume::SetMaterial ( G4Material pMaterial)
inline

Referenced by export_G4LogicalVolume(), G4LogicalVolume(), G4ScoreSplittingProcess::PostStepDoIt(), DetectorConstruction::SetAbsMaterial(), RE06DetectorConstruction::SetAbsorberMaterial(), DetectorConstruction::SetAbsorberMaterial(), Em10DetectorConstruction::SetAbsorberMaterial(), F03DetectorConstruction::SetAbsorberMaterial(), F01DetectorConstruction::SetAbsorberMaterial(), F02DetectorConstruction::SetAbsorberMaterial(), ExG4DetectorConstruction02::SetBoxMaterial(), B2bDetectorConstruction::SetChamberMaterial(), ExN02DetectorConstruction::setChamberMaterial(), ExP01DetectorConstruction::setChamberMaterial(), DetectorConstruction::SetContainerMaterial(), GammaRayTelDetectorConstruction::SetConverterMaterial(), exrdmDetectorConstruction::SetDetectorMaterial(), DetectorConstruction::SetEcalMaterial(), RE06DetectorConstruction::SetGapMaterial(), DetectorConstruction::SetGapMaterial(), DetectorConstruction::SetGasMaterial(), ExG4DetectorConstruction01::SetMaterial(), XrayFluoMercuryDetectorConstruction::SetMercuryMaterial(), XrayFluoPlaneDetectorConstruction::SetPlaneMaterial(), XrayFluoDetectorConstruction::SetSampleMaterial(), DetectorConstruction::SetTallyMaterial(), DetectorConstruction::SetTarget1Material(), DetectorConstruction::SetTarget2Material(), B2bDetectorConstruction::SetTargetMaterial(), B2aDetectorConstruction::SetTargetMaterial(), DetectorConstruction::SetTargetMaterial(), ExN02DetectorConstruction::setTargetMaterial(), ExP01DetectorConstruction::setTargetMaterial(), exrdmDetectorConstruction::SetTargetMaterial(), ExG4DetectorConstruction02::SetWorldMaterial(), DetectorConstruction::SetWorldMaterial(), F03DetectorConstruction::SetWorldMaterial(), F01DetectorConstruction::SetWorldMaterial(), F02DetectorConstruction::SetWorldMaterial(), and Em10DetectorConstruction::SetWorldMaterial().

void G4LogicalVolume::SetMaterialCutsCouple ( G4MaterialCutsCouple cuts)
inline

Referenced by export_G4LogicalVolume().

void G4LogicalVolume::SetName ( const G4String pName)
inline
void G4LogicalVolume::SetOptimisation ( G4bool  optim)
inline

Referenced by export_G4LogicalVolume().

void G4LogicalVolume::SetRegion ( G4Region reg)
inline
void G4LogicalVolume::SetRegionRootFlag ( G4bool  rreg)
inline
void G4LogicalVolume::SetSensitiveDetector ( G4VSensitiveDetector pSDetector)
inline
void G4LogicalVolume::SetSmartless ( G4double  s)
inline
void G4LogicalVolume::SetSolid ( G4VSolid pSolid)
inline
static void G4LogicalVolume::SetSolid ( G4LVData instLVdata,
G4VSolid pSolid 
)
inlinestatic
void G4LogicalVolume::SetUserLimits ( G4UserLimits pULimits)
inline
void G4LogicalVolume::SetVisAttributes ( const G4VisAttributes pVA)
inline

Referenced by HadrontherapyModulator::BuildModulator(), DicomPhantomParameterisationColour::ComputeMaterial(), CellParameterisation::ComputeMaterial(), Par01DetectorConstruction::Construct(), ExN04DetectorConstruction::Construct(), exGPSGeometryConstruction::Construct(), G4MIRDLiver::Construct(), G4MIRDRightScapula::Construct(), G4MIRDSkull::Construct(), G4MIRDSmallIntestine::Construct(), G4MIRDThyroid::Construct(), G4MIRDLeftAdrenal::Construct(), G4MIRDUpperLargeIntestine::Construct(), G4MIRDLeftKidney::Construct(), G4MIRDLeftClavicle::Construct(), G4MIRDLeftLeg::Construct(), G4MIRDLeftLegBone::Construct(), G4MIRDLeftOvary::Construct(), G4MIRDLeftScapula::Construct(), G4MIRDLeftTeste::Construct(), B03DetectorConstruction::Construct(), G4MIRDLowerLargeIntestine::Construct(), G4MIRDMaleGenitalia::Construct(), G4MIRDPelvis::Construct(), G4MIRDRightAdrenal::Construct(), G4MIRDRightClavicle::Construct(), G4MIRDRightKidney::Construct(), G4MIRDRightLeg::Construct(), G4MIRDRightLegBone::Construct(), G4MIRDRightTeste::Construct(), G4MIRDBrain::Construct(), G4MIRDStomach::Construct(), G4MIRDSpleen::Construct(), G4MIRDHeart::Construct(), G4MIRDTrunk::Construct(), G4MIRDUpperSpine::Construct(), G4MIRDLeftBreast::Construct(), G4MIRDUterus::Construct(), RE05DetectorConstruction::Construct(), RE01DetectorConstruction::Construct(), G4MIRDThymus::Construct(), G4MIRDLeftLung::Construct(), G4MIRDMiddleLowerSpine::Construct(), G4MIRDPancreas::Construct(), G4MIRDRibCage::Construct(), G4MIRDRightBreast::Construct(), G4MIRDRightOvary::Construct(), G4MIRDRightLung::Construct(), G4MIRDHead::Construct(), G4MIRDLeftArmBone::Construct(), B02DetectorConstruction::Construct(), G4MIRDUrinaryBladder::Construct(), G4MIRDRightArmBone::Construct(), B01DetectorConstruction::Construct(), B3DetectorConstruction::Construct(), ExGflashDetectorConstruction::Construct(), FCALCryostatVolumes::Construct(), FCALEMModule::Construct(), FCALHadModule::Construct(), eRositaDetectorConstruction::Construct(), B5DetectorConstruction::Construct(), exrdmDetectorConstruction::Construct(), ExN02DetectorConstruction::Construct(), ExP01DetectorConstruction::Construct(), FCALTestbeamSetup::Construct(), G02DetectorConstruction::Construct(), ExErrorDetectorConstruction::Construct(), UltraDetectorConstruction::Construct(), DMXDetectorConstruction::Construct(), CML2Ph_BoxInBox::Construct(), RE02DetectorConstruction::Construct(), G4tgbVolume::ConstructG4LogVol(), CML2WorldConstruction::create(), ElectronBenchmarkDetector::CreateExitWindow(), ElectronBenchmarkDetector::CreateHeliumBag(), ElectronBenchmarkDetector::CreateMonitor(), ElectronBenchmarkDetector::CreatePrimaryFoil(), ElectronBenchmarkDetector::CreateScorer(), ElectronBenchmarkDetector::CreateWorld(), F04ElementField::F04ElementField(), G4BuildGeom(), G4GDMLReadStructure::GetWorldVolume(), F04ElementField::SetColor(), G4VVisCommandGeometrySet::SetLVVisAtts(), G4VisCommandGeometryRestore::SetNewValue(), CCalG4Able::setVisType(), and G03ColorReader::VolumeRead().

void G4LogicalVolume::SetVisAttributes ( const G4VisAttributes VA)

Definition at line 388 of file G4LogicalVolume.cc.

389 {
390  fVisAttributes = new G4VisAttributes(VA);
391 }
void G4LogicalVolume::SetVoxelHeader ( G4SmartVoxelHeader pVoxel)
inline
void G4LogicalVolume::TerminateWorker ( G4LogicalVolume ptrMasterObject)

Definition at line 117 of file G4LogicalVolume.cc.

Referenced by G4GeometryWorkspace::DestroyWorkspace().

118 {
119 }
G4int G4LogicalVolume::TotalVolumeEntities ( ) const

Definition at line 269 of file G4LogicalVolume.cc.

References G4VPhysicalVolume::GetLogicalVolume(), G4VPhysicalVolume::GetMultiplicity(), and TotalVolumeEntities().

Referenced by export_G4LogicalVolume(), and TotalVolumeEntities().

270 {
271  G4int vols = 1;
272  for (G4PhysicalVolumeList::const_iterator itDau = fDaughters.begin();
273  itDau != fDaughters.end(); itDau++)
274  {
275  G4VPhysicalVolume* physDaughter = (*itDau);
276  vols += physDaughter->GetMultiplicity()
277  *physDaughter->GetLogicalVolume()->TotalVolumeEntities();
278  }
279  return vols;
280 }
int G4int
Definition: G4Types.hh:78
G4int TotalVolumeEntities() const
G4LogicalVolume * GetLogicalVolume() const
virtual G4int GetMultiplicity() const
void G4LogicalVolume::UpdateMaterial ( G4Material pMaterial)
inline

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