Geant4  10.02.p03
G4LogicalVolume Class Reference

#include <G4LogicalVolume.hh>

Collaboration diagram for G4LogicalVolume:

Public Member Functions

 G4LogicalVolume (G4VSolid *pSolid, G4Material *pMaterial, const G4String &name, G4FieldManager *pFieldMgr=0, G4VSensitiveDetector *pSDetector=0, G4UserLimits *pULimits=0, G4bool optimise=true)
 
 ~G4LogicalVolume ()
 
const G4StringGetName () 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)
 

Private Types

typedef std::vector< G4VPhysicalVolume * > G4PhysicalVolumeList
 

Private Member Functions

 G4LogicalVolume (const G4LogicalVolume &)
 
G4LogicalVolumeoperator= (const G4LogicalVolume &)
 

Private Attributes

G4PhysicalVolumeList fDaughters
 
G4String fName
 
G4UserLimitsfUserLimits
 
G4SmartVoxelHeaderfVoxel
 
G4bool fOptimise
 
G4bool fRootRegion
 
G4bool fLock
 
G4double fSmartless
 
const G4VisAttributesfVisAttributes
 
G4RegionfRegion
 
G4double fBiasWeight
 
G4int instanceID
 
G4VSolidfSolid
 
G4VSensitiveDetectorfSensitiveDetector
 
G4FieldManagerfFieldManager
 
G4LVDatalvdata
 

Static Private Attributes

static G4GEOM_DLL G4LVManager subInstanceManager
 

Detailed Description

Definition at line 189 of file G4LogicalVolume.hh.

Member Typedef Documentation

◆ G4PhysicalVolumeList

Definition at line 191 of file G4LogicalVolume.hh.

Constructor & Destructor Documentation

◆ G4LogicalVolume() [1/3]

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

Definition at line 77 of file G4LogicalVolume.cc.

84  : fDaughters(0,(G4VPhysicalVolume*)0),
85  fVoxel(0), fOptimise(optimise), fRootRegion(false), fLock(false),
87 {
88  // Initialize 'Shadow'/master pointers - for use in copying to workers
89  fSolid = pSolid;
90  fSensitiveDetector = pSDetector;
91  fFieldManager = pFieldMgr;
92 
94  AssignFieldManager(pFieldMgr); // G4MT_fmanager = pFieldMgr;
95 
96  // fMasterFieldMgr= pFieldMgr;
97  G4MT_mass = 0.;
98  G4MT_ccouple = 0;
99 
100  SetSolid(pSolid);
101  SetMaterial(pMaterial);
102  SetName(name);
103  SetSensitiveDetector(pSDetector);
104  SetUserLimits(pULimits);
105 
106  // Initialize 'Shadow' data structure - for use by object persistency
107  lvdata = new G4LVData();
108  lvdata->fSolid = pSolid;
109  lvdata->fMaterial = pMaterial;
110 
111  //
112  // Add to store
113  //
115 }
G4SmartVoxelHeader * fVoxel
G4FieldManager * fFieldManager
static G4GEOM_DLL G4LVManager subInstanceManager
void SetUserLimits(G4UserLimits *pULimits)
#define G4MT_ccouple
void SetSolid(G4VSolid *pSolid)
G4int CreateSubInstance()
void AssignFieldManager(G4FieldManager *fldMgr)
#define G4MT_mass
G4Material * fMaterial
void SetName(const G4String &pName)
G4VSensitiveDetector * fSensitiveDetector
void SetMaterial(G4Material *pMaterial)
static void Register(G4LogicalVolume *pVolume)
void SetSensitiveDetector(G4VSensitiveDetector *pSDetector)
G4VSolid * fSolid
G4PhysicalVolumeList fDaughters
const G4VisAttributes * fVisAttributes
Here is the call graph for this function:

◆ ~G4LogicalVolume()

G4LogicalVolume::~G4LogicalVolume ( )

Definition at line 147 of file G4LogicalVolume.cc.

148 {
149  if (!fLock && fRootRegion) // De-register root region first if not locked
150  { // and flagged as root logical-volume
151  fRegion->RemoveRootLogicalVolume(this, true);
152  }
153  delete lvdata;
155 }
static void DeRegister(G4LogicalVolume *pVolume)
void RemoveRootLogicalVolume(G4LogicalVolume *lv, G4bool scan=true)
Definition: G4Region.cc:319
Here is the call graph for this function:

◆ G4LogicalVolume() [2/3]

G4LogicalVolume::G4LogicalVolume ( __void__ &  )

Definition at line 122 of file G4LogicalVolume.cc.

123  : fDaughters(0,(G4VPhysicalVolume*)0),
124  fName(""), fUserLimits(0),
125  fVoxel(0), fOptimise(true), fRootRegion(false), fLock(false),
126  fSmartless(2.), fVisAttributes(0), fRegion(0), fBiasWeight(1.),
128 {
130 
131  SetSensitiveDetector(0); // G4MT_sdetector = 0;
132  SetFieldManager(0, false); // G4MT_fmanager = 0;
133 
134  G4MT_mass = 0.;
135  G4MT_ccouple = 0;
136 
137  // Add to store
138  //
140 }
G4SmartVoxelHeader * fVoxel
G4FieldManager * fFieldManager
static G4GEOM_DLL G4LVManager subInstanceManager
#define G4MT_ccouple
G4int CreateSubInstance()
void SetFieldManager(G4FieldManager *pFieldMgr, G4bool forceToAllDaughters)
#define G4MT_mass
G4VSensitiveDetector * fSensitiveDetector
static void Register(G4LogicalVolume *pVolume)
G4UserLimits * fUserLimits
void SetSensitiveDetector(G4VSensitiveDetector *pSDetector)
G4PhysicalVolumeList fDaughters
const G4VisAttributes * fVisAttributes
Here is the call graph for this function:

◆ G4LogicalVolume() [3/3]

G4LogicalVolume::G4LogicalVolume ( const G4LogicalVolume )
private

Member Function Documentation

◆ AddDaughter()

void G4LogicalVolume::AddDaughter ( G4VPhysicalVolume p)

Definition at line 263 of file G4LogicalVolume.cc.

264 {
265  if( !fDaughters.empty() && fDaughters[0]->IsReplicated() )
266  {
267  std::ostringstream message;
268  message << "ERROR - Attempt to place a volume in a mother volume" << G4endl
269  << " already containing a replicated volume." << G4endl
270  << " A volume can either contain several placements" << G4endl
271  << " or a unique replica or parameterised volume !" << G4endl
272  << " Mother logical volume: " << GetName() << G4endl
273  << " Placing volume: " << pNewDaughter->GetName() << G4endl;
274  G4Exception("G4LogicalVolume::AddDaughter()", "GeomMgt0002",
275  FatalException, message,
276  "Replica or parameterised volume must be the only daughter !");
277  }
278 
279  // Invalidate previous calculation of mass - if any - for all threads
280  G4MT_mass = 0.;
281  // SignalVolumeChange(); // fVolumeChanged= true;
282  fDaughters.push_back(pNewDaughter);
283 
284  G4LogicalVolume* pDaughterLogical = pNewDaughter->GetLogicalVolume();
285 
286  // Propagate the Field Manager, if the daughter has no field Manager.
287  //
288  G4FieldManager* pDaughterFieldManager = pDaughterLogical->GetFieldManager();
289 
290  if( pDaughterFieldManager == 0 )
291  {
292  pDaughterLogical->SetFieldManager(G4MT_fmanager, false);
293  }
294  if (fRegion)
295  {
296  PropagateRegion();
297  fRegion->RegionModified(true);
298  }
299 }
void RegionModified(G4bool flag)
#define G4MT_fmanager
void SetFieldManager(G4FieldManager *pFieldMgr, G4bool forceToAllDaughters)
const G4String & GetName() const
G4FieldManager * GetFieldManager() const
#define G4MT_mass
void PropagateRegion()
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
G4PhysicalVolumeList fDaughters
Here is the call graph for this function:
Here is the caller graph for this function:

◆ AssignFieldManager()

void G4LogicalVolume::AssignFieldManager ( G4FieldManager fldMgr)

Definition at line 230 of file G4LogicalVolume.cc.

231 {
232  G4MT_fmanager= fldMgr;
234 }
G4FieldManager * fFieldManager
#define G4MT_fmanager
G4bool IsMasterThread()
Definition: G4Threading.cc:136
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CharacteriseDaughters()

EVolume G4LogicalVolume::CharacteriseDaughters ( ) const
inline
Here is the caller graph for this function:

◆ ClearDaughters()

void G4LogicalVolume::ClearDaughters ( )

Definition at line 327 of file G4LogicalVolume.cc.

328 {
329  fDaughters.erase(fDaughters.begin(), fDaughters.end());
330  if (fRegion)
331  {
332  fRegion->RegionModified(true);
333  }
334  G4MT_mass = 0.;
335 }
void RegionModified(G4bool flag)
#define G4MT_mass
G4PhysicalVolumeList fDaughters
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetBiasWeight()

G4double G4LogicalVolume::GetBiasWeight ( ) const
inline
Here is the caller graph for this function:

◆ GetDaughter()

G4VPhysicalVolume* G4LogicalVolume::GetDaughter ( const G4int  i) const
inline
Here is the caller graph for this function:

◆ GetFastSimulationManager()

G4FastSimulationManager* G4LogicalVolume::GetFastSimulationManager ( ) const
inline
Here is the caller graph for this function:

◆ GetFieldManager()

G4FieldManager * G4LogicalVolume::GetFieldManager ( ) const

Definition at line 221 of file G4LogicalVolume.cc.

222 {
223  return G4MT_fmanager;
224 }
#define G4MT_fmanager
Here is the caller graph for this function:

◆ GetInstanceID()

G4int G4LogicalVolume::GetInstanceID ( ) const
inline

◆ GetMass()

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

Definition at line 512 of file G4LogicalVolume.cc.

515 {
516  // Return the cached non-zero value, if not forced
517  //
518  if ( (G4MT_mass) && (!forced) ) return G4MT_mass;
519 
520  // Global density and computed mass associated to the logical
521  // volume without considering its daughters
522  //
523  G4Material* logMaterial = parMaterial ? parMaterial : GetMaterial(); // G4MT_material;
524  if (!logMaterial)
525  {
526  std::ostringstream message;
527  message << "No material associated to the logical volume: " << fName << " !"
528  << G4endl
529  << "Sorry, cannot compute the mass ...";
530  G4Exception("G4LogicalVolume::GetMass()", "GeomMgt0002",
531  FatalException, message);
532  return 0;
533  }
534  if (! GetSolid() ) // !G4MT_solid)
535  {
536  std::ostringstream message;
537  message << "No solid is associated to the logical volume: " << fName << " !"
538  << G4endl
539  << "Sorry, cannot compute the mass ...";
540  G4Exception("G4LogicalVolume::GetMass()", "GeomMgt0002",
541  FatalException, message);
542  return 0;
543  }
544  G4double globalDensity = logMaterial->GetDensity();
545  G4double motherMass= GetSolid()->GetCubicVolume() * globalDensity;
546 
547  // G4MT_mass =
548  // SetMass( motherMmass );
549  G4double massSum= motherMass;
550 
551  // For each daughter in the tree, subtract the mass occupied
552  // and if required by the propagate flag, add the real daughter's
553  // one computed recursively
554 
555  for (G4PhysicalVolumeList::const_iterator itDau = fDaughters.begin();
556  itDau != fDaughters.end(); itDau++)
557  {
558  G4VPhysicalVolume* physDaughter = (*itDau);
559  G4LogicalVolume* logDaughter = physDaughter->GetLogicalVolume();
560  G4double subMass=0.;
561  G4VSolid* daughterSolid = 0;
562  G4Material* daughterMaterial = 0;
563 
564  // Compute the mass to subtract and to add for each daughter
565  // considering its multiplicity (i.e. replicated or not) and
566  // eventually its parameterisation (by solid and/or by material)
567  //
568  for (G4int i=0; i<physDaughter->GetMultiplicity(); i++)
569  {
571  physParam = physDaughter->GetParameterisation();
572  if (physParam)
573  {
574  daughterSolid = physParam->ComputeSolid(i, physDaughter);
575  daughterSolid->ComputeDimensions(physParam, i, physDaughter);
576  daughterMaterial = physParam->ComputeMaterial(i, physDaughter);
577  }
578  else
579  {
580  daughterSolid = logDaughter->GetSolid();
581  daughterMaterial = logDaughter->GetMaterial();
582  }
583  subMass = daughterSolid->GetCubicVolume() * globalDensity;
584 
585  // Subtract the daughter's portion for the mass and, if required,
586  // add the real daughter's mass computed recursively
587  //
588  massSum -= subMass;
589  if (propagate)
590  {
591  massSum += logDaughter->GetMass(true, true, daughterMaterial);
592  }
593  }
594  }
595  G4MT_mass= massSum;
596  return massSum;
597 }
virtual G4Material * ComputeMaterial(const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=0)
virtual G4int GetMultiplicity() const
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4VSolid.cc:137
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
G4Material * GetMaterial() const
virtual G4double GetCubicVolume()
Definition: G4VSolid.cc:188
G4double GetDensity() const
Definition: G4Material.hh:180
int G4int
Definition: G4Types.hh:78
virtual G4VPVParameterisation * GetParameterisation() const =0
#define G4MT_mass
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
G4double GetMass(G4bool forced=false, G4bool propagate=true, G4Material *parMaterial=0)
double G4double
Definition: G4Types.hh:76
G4LogicalVolume * GetLogicalVolume() const
G4VSolid * GetSolid() const
G4PhysicalVolumeList fDaughters
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetMasterFieldManager()

G4FieldManager* G4LogicalVolume::GetMasterFieldManager ( ) const
inline
Here is the caller graph for this function:

◆ GetMasterSensitiveDetector()

G4VSensitiveDetector* G4LogicalVolume::GetMasterSensitiveDetector ( ) const
inline
Here is the caller graph for this function:

◆ GetMasterSolid()

G4VSolid* G4LogicalVolume::GetMasterSolid ( ) const
inline
Here is the caller graph for this function:

◆ GetMaterial()

G4Material * G4LogicalVolume::GetMaterial ( ) const

Definition at line 387 of file G4LogicalVolume.cc.

388 {
389  return G4MT_material;
390 }
#define G4MT_material

◆ GetMaterialCutsCouple()

const G4MaterialCutsCouple * G4LogicalVolume::GetMaterialCutsCouple ( ) const

Definition at line 436 of file G4LogicalVolume.cc.

437 {
438  return G4MT_ccouple;
439 }
#define G4MT_ccouple
Here is the caller graph for this function:

◆ GetName()

const G4String& G4LogicalVolume::GetName ( ) const
inline

◆ GetNoDaughters()

G4int G4LogicalVolume::GetNoDaughters ( ) const
inline
Here is the caller graph for this function:

◆ GetRegion()

G4Region* G4LogicalVolume::GetRegion ( ) const
inline
Here is the caller graph for this function:

◆ GetSensitiveDetector()

G4VSensitiveDetector * G4LogicalVolume::GetSensitiveDetector ( ) const

Definition at line 417 of file G4LogicalVolume.cc.

418 {
419  return G4MT_sdetector;
420 }
#define G4MT_sdetector
Here is the caller graph for this function:

◆ GetSmartless()

G4double G4LogicalVolume::GetSmartless ( ) const
inline
Here is the caller graph for this function:

◆ GetSolid() [1/2]

G4VSolid * G4LogicalVolume::GetSolid ( ) const

Definition at line 355 of file G4LogicalVolume.cc.

356 {
357  // return G4MT_solid;
358  // return ((subInstanceManager.offset[instanceID]).fSolid);
359  return this->GetSolid( subInstanceManager.offset[instanceID] );
360 }
static G4GEOM_DLL G4ThreadLocal T * offset
static G4GEOM_DLL G4LVManager subInstanceManager
G4VSolid * GetSolid() const

◆ GetSolid() [2/2]

G4VSolid * G4LogicalVolume::GetSolid ( G4LVData instLVdata)
static

Definition at line 350 of file G4LogicalVolume.cc.

351 {
352  return instLVdata.fSolid;
353 }
G4VSolid * fSolid

◆ GetSubInstanceManager()

const G4LVManager & G4LogicalVolume::GetSubInstanceManager ( )
static

Definition at line 212 of file G4LogicalVolume.cc.

213 {
214  return subInstanceManager;
215 }
static G4GEOM_DLL G4LVManager subInstanceManager
Here is the caller graph for this function:

◆ GetUserLimits()

G4UserLimits* G4LogicalVolume::GetUserLimits ( ) const
inline
Here is the caller graph for this function:

◆ GetVisAttributes()

const G4VisAttributes* G4LogicalVolume::GetVisAttributes ( ) const
inline
Here is the caller graph for this function:

◆ GetVoxelHeader()

G4SmartVoxelHeader* G4LogicalVolume::GetVoxelHeader ( ) const
inline
Here is the caller graph for this function:

◆ InitialiseWorker()

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

Definition at line 168 of file G4LogicalVolume.cc.

171 {
173 
174  SetSolid(pSolid);
175  SetSensitiveDetector(pSDetector); // How this object is available now ?
176  AssignFieldManager(fFieldManager); // Should be set - but a per-thread copy is not available yet
177  // G4MT_fmanager= fFieldManager;
178  // Must not call SetFieldManager(fFieldManager, false); which propagates FieldMgr
179 
180 #ifdef CLONE_FIELD_MGR
181  // Create a field FieldManager by cloning
182  G4FieldManager workerFldMgr= fFieldManager->GetWorkerClone(G4bool* created);
183  if( created || (GetFieldManager()!=workerFldMgr) )
184  {
185  SetFieldManager(fFieldManager, false); // which propagates FieldMgr
186  }else{
187  // Field manager existed and is equal to current one
188  AssignFieldManager(workerFldMgr);
189  }
190 #endif
191 }
G4FieldManager * fFieldManager
static G4GEOM_DLL G4LVManager subInstanceManager
void SetSolid(G4VSolid *pSolid)
void SetFieldManager(G4FieldManager *pFieldMgr, G4bool forceToAllDaughters)
void AssignFieldManager(G4FieldManager *fldMgr)
bool G4bool
Definition: G4Types.hh:79
G4FieldManager * GetFieldManager() const
void SlaveCopySubInstanceArray()
void SetSensitiveDetector(G4VSensitiveDetector *pSDetector)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ IsAncestor()

G4bool G4LogicalVolume::IsAncestor ( const G4VPhysicalVolume p) const

Definition at line 458 of file G4LogicalVolume.cc.

459 {
460  G4bool isDaughter = IsDaughter(aVolume);
461  if (!isDaughter)
462  {
463  for (G4PhysicalVolumeList::const_iterator itDau = fDaughters.begin();
464  itDau != fDaughters.end(); itDau++)
465  {
466  isDaughter = (*itDau)->GetLogicalVolume()->IsAncestor(aVolume);
467  if (isDaughter) break;
468  }
469  }
470  return isDaughter;
471 }
G4bool IsDaughter(const G4VPhysicalVolume *p) const
bool G4bool
Definition: G4Types.hh:79
G4PhysicalVolumeList fDaughters
Here is the call graph for this function:
Here is the caller graph for this function:

◆ IsDaughter()

G4bool G4LogicalVolume::IsDaughter ( const G4VPhysicalVolume p) const
inline
Here is the caller graph for this function:

◆ IsRegion()

G4bool G4LogicalVolume::IsRegion ( ) const
inline
Here is the caller graph for this function:

◆ IsRootRegion()

G4bool G4LogicalVolume::IsRootRegion ( ) const
inline
Here is the caller graph for this function:

◆ IsToOptimise()

G4bool G4LogicalVolume::IsToOptimise ( ) const
inline
Here is the caller graph for this function:

◆ Lock()

void G4LogicalVolume::Lock ( )
inline

◆ operator=()

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

◆ operator==()

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

◆ PropagateRegion()

void G4LogicalVolume::PropagateRegion ( )
inline
Here is the caller graph for this function:

◆ RemoveDaughter()

void G4LogicalVolume::RemoveDaughter ( const G4VPhysicalVolume p)

Definition at line 305 of file G4LogicalVolume.cc.

306 {
307  G4PhysicalVolumeList::iterator i;
308  for ( i=fDaughters.begin(); i!=fDaughters.end(); ++i )
309  {
310  if (**i==*p)
311  {
312  fDaughters.erase(i);
313  break;
314  }
315  }
316  if (fRegion)
317  {
318  fRegion->RegionModified(true);
319  }
320  G4MT_mass = 0.;
321 }
void RegionModified(G4bool flag)
#define G4MT_mass
G4PhysicalVolumeList fDaughters
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ResetMass()

void G4LogicalVolume::ResetMass ( )

Definition at line 341 of file G4LogicalVolume.cc.

342 {
343  G4MT_mass= 0.0;
344 }
#define G4MT_mass
Here is the caller graph for this function:

◆ SetBiasWeight()

void G4LogicalVolume::SetBiasWeight ( G4double  w)
inline
Here is the caller graph for this function:

◆ SetFieldManager()

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

Definition at line 241 of file G4LogicalVolume.cc.

243 {
244  // G4MT_fmanager = pNewFieldMgr;
245  AssignFieldManager(pNewFieldMgr);
246 
247  G4int NoDaughters = GetNoDaughters();
248  while ( (NoDaughters--)>0 )
249  {
250  G4LogicalVolume* DaughterLogVol;
251  DaughterLogVol = GetDaughter(NoDaughters)->GetLogicalVolume();
252  if ( forceAllDaughters || (DaughterLogVol->GetFieldManager() == 0) )
253  {
254  DaughterLogVol->SetFieldManager(pNewFieldMgr, forceAllDaughters);
255  }
256  }
257 }
G4int GetNoDaughters() const
int G4int
Definition: G4Types.hh:78
void SetFieldManager(G4FieldManager *pFieldMgr, G4bool forceToAllDaughters)
void AssignFieldManager(G4FieldManager *fldMgr)
G4FieldManager * GetFieldManager() const
G4VPhysicalVolume * GetDaughter(const G4int i) const
G4LogicalVolume * GetLogicalVolume() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetMaterial()

void G4LogicalVolume::SetMaterial ( G4Material pMaterial)

Definition at line 396 of file G4LogicalVolume.cc.

397 {
398  G4MT_material=pMaterial;
399  G4MT_mass = 0.;
400 }
#define G4MT_mass
#define G4MT_material
Here is the caller graph for this function:

◆ SetMaterialCutsCouple()

void G4LogicalVolume::SetMaterialCutsCouple ( G4MaterialCutsCouple cuts)

Definition at line 445 of file G4LogicalVolume.cc.

446 {
447  G4MT_ccouple = cuts;
448 }
#define G4MT_ccouple
Here is the caller graph for this function:

◆ SetName()

void G4LogicalVolume::SetName ( const G4String pName)
inline
Here is the caller graph for this function:

◆ SetOptimisation()

void G4LogicalVolume::SetOptimisation ( G4bool  optim)
inline
Here is the caller graph for this function:

◆ SetRegion()

void G4LogicalVolume::SetRegion ( G4Region reg)
inline
Here is the caller graph for this function:

◆ SetRegionRootFlag()

void G4LogicalVolume::SetRegionRootFlag ( G4bool  rreg)
inline
Here is the caller graph for this function:

◆ SetSensitiveDetector()

void G4LogicalVolume::SetSensitiveDetector ( G4VSensitiveDetector pSDetector)

Definition at line 426 of file G4LogicalVolume.cc.

427 {
428  G4MT_sdetector = pSDetector;
430 }
#define G4MT_sdetector
G4bool IsMasterThread()
Definition: G4Threading.cc:136
G4VSensitiveDetector * fSensitiveDetector
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetSmartless()

void G4LogicalVolume::SetSmartless ( G4double  s)
inline
Here is the caller graph for this function:

◆ SetSolid() [1/2]

void G4LogicalVolume::SetSolid ( G4VSolid pSolid)

Definition at line 366 of file G4LogicalVolume.cc.

367 {
368 
369  // ((subInstanceManager.offset[instanceID]).fSolid) = pSolid;
370  G4MT_solid=pSolid;
371  // G4MT_mass = 0.;
372  this->ResetMass();
373 }
#define G4MT_solid
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetSolid() [2/2]

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

Definition at line 375 of file G4LogicalVolume.cc.

376 {
377  instLVdata.fSolid = pSolid;
378  // G4MT_solid=pSolid;
379  instLVdata.fMass= 0;
380  // A fast way to reset the mass ... ie G4MT_mass = 0.;
381 }
G4double fMass
G4VSolid * fSolid

◆ SetUserLimits()

void G4LogicalVolume::SetUserLimits ( G4UserLimits pULimits)
inline
Here is the caller graph for this function:

◆ SetVisAttributes() [1/2]

void G4LogicalVolume::SetVisAttributes ( const G4VisAttributes pVA)
inline

◆ SetVisAttributes() [2/2]

void G4LogicalVolume::SetVisAttributes ( const G4VisAttributes VA)

Definition at line 599 of file G4LogicalVolume.cc.

600 {
602 }
const G4VisAttributes * fVisAttributes

◆ SetVoxelHeader()

void G4LogicalVolume::SetVoxelHeader ( G4SmartVoxelHeader pVoxel)
inline
Here is the caller graph for this function:

◆ TerminateWorker()

void G4LogicalVolume::TerminateWorker ( G4LogicalVolume ptrMasterObject)

Definition at line 202 of file G4LogicalVolume.cc.

203 {
204 }
Here is the caller graph for this function:

◆ TotalVolumeEntities()

G4int G4LogicalVolume::TotalVolumeEntities ( ) const

Definition at line 480 of file G4LogicalVolume.cc.

481 {
482  G4int vols = 1;
483  for (G4PhysicalVolumeList::const_iterator itDau = fDaughters.begin();
484  itDau != fDaughters.end(); itDau++)
485  {
486  G4VPhysicalVolume* physDaughter = (*itDau);
487  vols += physDaughter->GetMultiplicity()
488  *physDaughter->GetLogicalVolume()->TotalVolumeEntities();
489  }
490  return vols;
491 }
virtual G4int GetMultiplicity() const
G4int TotalVolumeEntities() const
int G4int
Definition: G4Types.hh:78
G4LogicalVolume * GetLogicalVolume() const
G4PhysicalVolumeList fDaughters
Here is the call graph for this function:
Here is the caller graph for this function:

◆ UpdateMaterial()

void G4LogicalVolume::UpdateMaterial ( G4Material pMaterial)

Definition at line 406 of file G4LogicalVolume.cc.

407 {
408  G4MT_material=pMaterial;
409  if(fRegion) { G4MT_ccouple = fRegion->FindCouple(pMaterial); }
410  G4MT_mass = 0.;
411 }
#define G4MT_ccouple
#define G4MT_mass
G4MaterialCutsCouple * FindCouple(G4Material *mat)
#define G4MT_material
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ fBiasWeight

G4double G4LogicalVolume::fBiasWeight
private

Definition at line 419 of file G4LogicalVolume.hh.

◆ fDaughters

G4PhysicalVolumeList G4LogicalVolume::fDaughters
private

Definition at line 396 of file G4LogicalVolume.hh.

◆ fFieldManager

G4FieldManager* G4LogicalVolume::fFieldManager
private

Definition at line 433 of file G4LogicalVolume.hh.

◆ fLock

G4bool G4LogicalVolume::fLock
private

Definition at line 410 of file G4LogicalVolume.hh.

◆ fName

G4String G4LogicalVolume::fName
private

Definition at line 398 of file G4LogicalVolume.hh.

◆ fOptimise

G4bool G4LogicalVolume::fOptimise
private

Definition at line 406 of file G4LogicalVolume.hh.

◆ fRegion

G4Region* G4LogicalVolume::fRegion
private

Definition at line 417 of file G4LogicalVolume.hh.

◆ fRootRegion

G4bool G4LogicalVolume::fRootRegion
private

Definition at line 408 of file G4LogicalVolume.hh.

◆ fSensitiveDetector

G4VSensitiveDetector* G4LogicalVolume::fSensitiveDetector
private

Definition at line 432 of file G4LogicalVolume.hh.

◆ fSmartless

G4double G4LogicalVolume::fSmartless
private

Definition at line 412 of file G4LogicalVolume.hh.

◆ fSolid

G4VSolid* G4LogicalVolume::fSolid
private

Definition at line 431 of file G4LogicalVolume.hh.

◆ fUserLimits

G4UserLimits* G4LogicalVolume::fUserLimits
private

Definition at line 402 of file G4LogicalVolume.hh.

◆ fVisAttributes

const G4VisAttributes* G4LogicalVolume::fVisAttributes
private

Definition at line 415 of file G4LogicalVolume.hh.

◆ fVoxel

G4SmartVoxelHeader* G4LogicalVolume::fVoxel
private

Definition at line 404 of file G4LogicalVolume.hh.

◆ instanceID

G4int G4LogicalVolume::instanceID
private

Definition at line 422 of file G4LogicalVolume.hh.

◆ lvdata

G4LVData* G4LogicalVolume::lvdata
private

Definition at line 434 of file G4LogicalVolume.hh.

◆ subInstanceManager

G4LVManager G4LogicalVolume::subInstanceManager
staticprivate

Definition at line 424 of file G4LogicalVolume.hh.


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