Geant4-11
Public Member Functions | Static Public Member Functions | Private Types | Private Attributes | Static Private Attributes
G4LogicalCrystalVolume Class Reference

#include <G4LogicalCrystalVolume.hh>

Inheritance diagram for G4LogicalCrystalVolume:
G4LogicalVolume

Public Member Functions

void AddDaughter (G4VPhysicalVolume *p)
 
void AssignFieldManager (G4FieldManager *fldMgr)
 
G4bool ChangeDaughtersType (EVolume atype)
 
EVolume CharacteriseDaughters () const
 
void ClearDaughters ()
 
EVolume DeduceDaughtersType () const
 
 G4LogicalCrystalVolume (G4VSolid *pSolid, G4ExtendedMaterial *pMaterial, const G4String &name, G4FieldManager *pFieldMgr=nullptr, G4VSensitiveDetector *pSDetector=nullptr, G4UserLimits *pULimits=nullptr, G4bool optimise=true, G4int h=0, G4int k=0, G4int l=0, G4double rot=0.0)
 
const G4ThreeVectorGetBasis (G4int i) const
 
G4double GetBiasWeight () const
 
const G4CrystalExtensionGetCrystal () const
 
G4VPhysicalVolumeGetDaughter (const G4int i) const
 
G4FastSimulationManagerGetFastSimulationManager () const
 
G4FieldManagerGetFieldManager () const
 
G4int GetInstanceID () const
 
G4double GetMass (G4bool forced=false, G4bool propagate=true, G4Material *parMaterial=nullptr)
 
G4FieldManagerGetMasterFieldManager () const
 
G4VSensitiveDetectorGetMasterSensitiveDetector () const
 
G4VSolidGetMasterSolid () const
 
G4MaterialGetMaterial () const
 
const G4MaterialCutsCoupleGetMaterialCutsCouple () const
 
const G4StringGetName () const
 
size_t GetNoDaughters () const
 
G4RegionGetRegion () const
 
G4VSensitiveDetectorGetSensitiveDetector () const
 
G4double GetSmartless () const
 
G4VSolidGetSolid () const
 
G4UserLimitsGetUserLimits () const
 
const G4VisAttributesGetVisAttributes () const
 
G4SmartVoxelHeaderGetVoxelHeader () const
 
void InitialiseWorker (G4LogicalVolume *ptrMasterObject, G4VSolid *pSolid, G4VSensitiveDetector *pSDetector)
 
G4bool IsAncestor (const G4VPhysicalVolume *p) const
 
G4bool IsDaughter (const G4VPhysicalVolume *p) const
 
G4bool IsExtended () const
 
G4bool IsRegion () const
 
G4bool IsRootRegion () const
 
G4bool IsToOptimise () const
 
void Lock ()
 
G4bool operator== (const G4LogicalVolume &lv) const
 
void PropagateRegion ()
 
void RemoveDaughter (const G4VPhysicalVolume *p)
 
void ResetMass ()
 
const G4ThreeVectorRotateToLattice (G4ThreeVector &dir) const
 
const G4ThreeVectorRotateToSolid (G4ThreeVector &dir) const
 
void SetBiasWeight (G4double w)
 
void SetFieldManager (G4FieldManager *pFieldMgr, G4bool forceToAllDaughters)
 
void SetMaterial (G4Material *pMaterial)
 
void SetMaterialCutsCouple (G4MaterialCutsCouple *cuts)
 
void SetMillerOrientation (G4int h, G4int k, G4int l, G4double rot=0.0)
 
void SetName (const G4String &pName)
 
void SetOptimisation (G4bool optim)
 
void SetRegion (G4Region *reg)
 
void SetRegionRootFlag (G4bool rreg)
 
void SetSensitiveDetector (G4VSensitiveDetector *pSDetector)
 
void SetSmartless (G4double s)
 
void SetSolid (G4VSolid *pSolid)
 
void SetUserLimits (G4UserLimits *pULimits)
 
void SetVerbose (G4int aInt)
 
void SetVisAttributes (const G4VisAttributes &VA)
 
void SetVisAttributes (const G4VisAttributes *pVA)
 
void SetVoxelHeader (G4SmartVoxelHeader *pVoxel)
 
void TerminateWorker (G4LogicalVolume *ptrMasterObject)
 
G4int TotalVolumeEntities () const
 
void UpdateMaterial (G4Material *pMaterial)
 
 ~G4LogicalCrystalVolume ()
 

Static Public Member Functions

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

Private Types

using G4PhysicalVolumeList = std::vector< G4VPhysicalVolume * >
 

Private Attributes

G4double fBiasWeight = 1.0
 
G4PhysicalVolumeList fDaughters
 
EVolume fDaughtersVolumeType
 
G4FieldManagerfFieldManager = nullptr
 
G4RotationMatrix fInverse
 
G4bool fLock = false
 
G4String fName
 
G4bool fOptimise = true
 
G4RotationMatrix fOrient
 
G4RegionfRegion = nullptr
 
G4bool fRootRegion = false
 
G4double fRot = 0.0
 
G4VSensitiveDetectorfSensitiveDetector = nullptr
 
G4double fSmartless = 2.0
 
G4VSolidfSolid = nullptr
 
G4UserLimitsfUserLimits = nullptr
 
const G4VisAttributesfVisAttributes = nullptr
 
G4SmartVoxelHeaderfVoxel = nullptr
 
G4int hMiller = 1
 
G4int instanceID
 
G4int kMiller = 1
 
G4int lMiller = 0
 
G4LVDatalvdata = nullptr
 
G4int verboseLevel = 0
 

Static Private Attributes

static std::vector< G4LogicalVolume * > fLCVvec
 
static G4GEOM_DLL G4LVManager subInstanceManager
 

Detailed Description

Definition at line 46 of file G4LogicalCrystalVolume.hh.

Member Typedef Documentation

◆ G4PhysicalVolumeList

using G4LogicalVolume::G4PhysicalVolumeList = std::vector<G4VPhysicalVolume *>
privateinherited

Definition at line 392 of file G4LogicalVolume.hh.

Constructor & Destructor Documentation

◆ G4LogicalCrystalVolume()

G4LogicalCrystalVolume::G4LogicalCrystalVolume ( G4VSolid pSolid,
G4ExtendedMaterial pMaterial,
const G4String name,
G4FieldManager pFieldMgr = nullptr,
G4VSensitiveDetector pSDetector = nullptr,
G4UserLimits pULimits = nullptr,
G4bool  optimise = true,
G4int  h = 0,
G4int  k = 0,
G4int  l = 0,
G4double  rot = 0.0 
)

Definition at line 40 of file G4LogicalCrystalVolume.cc.

46: G4LogicalVolume(pSolid,pMaterial,name,pFieldMgr,pSDetector,pULimits,optimise)
47{
48 SetMillerOrientation(h, k, l, rot);
49 fLCVvec.push_back(this);
50}
static std::vector< G4LogicalVolume * > fLCVvec
void SetMillerOrientation(G4int h, G4int k, G4int l, G4double rot=0.0)
G4LogicalVolume(G4VSolid *pSolid, G4Material *pMaterial, const G4String &name, G4FieldManager *pFieldMgr=nullptr, G4VSensitiveDetector *pSDetector=nullptr, G4UserLimits *pULimits=nullptr, G4bool optimise=true)
const char * name(G4int ptype)

References fLCVvec, and SetMillerOrientation().

◆ ~G4LogicalCrystalVolume()

G4LogicalCrystalVolume::~G4LogicalCrystalVolume ( )

Definition at line 54 of file G4LogicalCrystalVolume.cc.

55{
56 fLCVvec.erase( std::remove(fLCVvec.begin(),fLCVvec.end(), this ),
57 fLCVvec.end() );
58}

References fLCVvec.

Member Function Documentation

◆ AddDaughter()

void G4LogicalVolume::AddDaughter ( G4VPhysicalVolume p)
inherited

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:35
#define G4MT_fmanager
#define G4MT_mass
#define G4endl
Definition: G4ios.hh:57
void SetFieldManager(G4FieldManager *pFieldMgr, G4bool forceToAllDaughters)
void PropagateRegion()
EVolume fDaughtersVolumeType
G4FieldManager * GetFieldManager() const
const G4String & GetName() const
G4PhysicalVolumeList fDaughters
void RegionModified(G4bool flag)
EVolume
Definition: geomdefs.hh:83

References FatalException, G4LogicalVolume::fDaughters, G4LogicalVolume::fDaughtersVolumeType, G4LogicalVolume::fRegion, G4endl, G4Exception(), G4MT_fmanager, G4MT_mass, G4LogicalVolume::GetFieldManager(), G4VPhysicalVolume::GetLogicalVolume(), G4LogicalVolume::GetName(), G4VPhysicalVolume::GetName(), G4LogicalVolume::PropagateRegion(), G4Region::RegionModified(), G4LogicalVolume::SetFieldManager(), and G4VPhysicalVolume::VolumeType().

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

◆ AssignFieldManager()

void G4LogicalVolume::AssignFieldManager ( G4FieldManager fldMgr)
inherited

◆ ChangeDaughtersType()

G4bool G4LogicalVolume::ChangeDaughtersType ( EVolume  atype)
inherited

Definition at line 659 of file G4LogicalVolume.cc.

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

References G4LogicalVolume::DeduceDaughtersType(), G4LogicalVolume::fDaughtersVolumeType, and kExternal.

◆ CharacteriseDaughters()

EVolume G4LogicalVolume::CharacteriseDaughters ( ) const
inlineinherited

◆ Clean()

void G4LogicalVolume::Clean ( )
staticinherited

◆ ClearDaughters()

void G4LogicalVolume::ClearDaughters ( )
inherited

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}

References G4LogicalVolume::fDaughters, G4LogicalVolume::fRegion, G4MT_mass, and G4Region::RegionModified().

Referenced by export_G4LogicalVolume().

◆ DeduceDaughtersType()

EVolume G4LogicalVolume::DeduceDaughtersType ( ) const
inlineinherited

◆ GetBasis()

const G4ThreeVector & G4LogicalCrystalVolume::GetBasis ( G4int  i) const

Definition at line 78 of file G4LogicalCrystalVolume.cc.

79{
80 return GetCrystal()->GetUnitCell()->GetBasis(i);
81}
G4CrystalUnitCell * GetUnitCell() const
const G4ThreeVector & GetBasis(G4int idx) const
const G4CrystalExtension * GetCrystal() const

References G4CrystalUnitCell::GetBasis(), GetCrystal(), and G4CrystalExtension::GetUnitCell().

Referenced by SetMillerOrientation().

◆ GetBiasWeight()

G4double G4LogicalVolume::GetBiasWeight ( ) const
inlineinherited

◆ GetCrystal()

const G4CrystalExtension * G4LogicalCrystalVolume::GetCrystal ( ) const

Definition at line 69 of file G4LogicalCrystalVolume.cc.

70{
71 return dynamic_cast<G4CrystalExtension*>(
72 dynamic_cast<G4ExtendedMaterial*>(GetMaterial() )
73 ->RetrieveExtension("crystal"));
74}
G4Material * GetMaterial() const

References G4LogicalVolume::GetMaterial().

Referenced by GetBasis().

◆ GetDaughter()

G4VPhysicalVolume * G4LogicalVolume::GetDaughter ( const G4int  i) const
inlineinherited

◆ GetFastSimulationManager()

G4FastSimulationManager * G4LogicalVolume::GetFastSimulationManager ( ) const
inlineinherited

◆ GetFieldManager()

G4FieldManager * G4LogicalVolume::GetFieldManager ( ) const
inherited

◆ GetInstanceID()

G4int G4LogicalVolume::GetInstanceID ( ) const
inlineinherited

◆ GetMass()

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

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)
G4double GetDensity() const
Definition: G4Material.hh:176
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

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

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

◆ GetMasterFieldManager()

G4FieldManager * G4LogicalVolume::GetMasterFieldManager ( ) const
inlineinherited

◆ GetMasterSensitiveDetector()

G4VSensitiveDetector * G4LogicalVolume::GetMasterSensitiveDetector ( ) const
inlineinherited

◆ GetMasterSolid()

G4VSolid * G4LogicalVolume::GetMasterSolid ( ) const
inlineinherited

◆ GetMaterial()

G4Material * G4LogicalVolume::GetMaterial ( ) const
inherited

◆ GetMaterialCutsCouple()

const G4MaterialCutsCouple * G4LogicalVolume::GetMaterialCutsCouple ( ) const
inherited

◆ GetName()

const G4String & G4LogicalVolume::GetName ( ) const
inlineinherited

Referenced by G4AdjointCrossSurfChecker::AddanExtSurfaceOfAvolume(), G4LogicalVolume::AddDaughter(), G4HepRepFileSceneHandler::AddHepRepInstance(), G4OpenInventorSceneHandler::AddPrimitive(), G4VtkSceneHandler::AddSolid(), G4VBiasingOperator::AttachTo(), G4SmartVoxelHeader::BuildNodes(), G4GeometryManager::BuildOptimisations(), G4SmartVoxelHeader::BuildReplicaVoxels(), G4SmartVoxelHeader::BuildVoxelsWithinLimits(), G4PVReplica::CheckOnlyDaughter(), G4PVParameterised::CheckOverlaps(), G4PVPlacement::CheckOverlaps(), checkVol(), G4PolarizedAnnihilation::ComputeSaturationFactor(), G4PolarizedCompton::ComputeSaturationFactor(), G4PolarizedIonisation::ComputeSaturationFactor(), G4tgbVolume::ConstructG4LogVol(), G4tgbVolume::ConstructG4PhysVol(), G4tgbVolume::ConstructG4Volumes(), G4PhysicalVolumeModel::CreateCurrentAttValues(), G3Division::CreatePVReplica(), G4ReflectionFactory::CreateReflectedLV(), 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(), export_G4LogicalVolume(), G4GDMLParser::ExportRegions(), G4BuildGeom(), G4PVReplica::G4PVReplica(), G4GDMLRead::GeneratePhysvolName(), G4Navigator::GetLocalExitNormal(), G4ITNavigator1::GetLocalExitNormal(), G4ITNavigator2::GetLocalExitNormal(), G4tgbGeometryDumper::GetPVChildren(), G4tgbVolumeMgr::GetTopLogVol(), G4tgbVolumeMgr::GetTopPhysVol(), G4GDMLReadStructure::GetWorldVolume(), G4AdjointCrossSurfChecker::GoingInOrOutOfaVolumeByExtSurface(), G4GDMLWriteParamvol::ParamvolAlgorithmWrite(), G4GDMLReadParamvol::ParamvolRead(), G4GDMLWriteParamvol::ParamvolWrite(), G4GDMLReadStructure::PhysvolRead(), G4GDMLWriteStructure::PhysvolWrite(), G4ReflectionFactory::Place(), G4ReflectionFactory::ReflectDaughters(), G4ReflectionFactory::ReflectPVDivision(), G4ReflectionFactory::ReflectPVPlacement(), G4ReflectionFactory::ReflectPVReplica(), 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()

size_t G4LogicalVolume::GetNoDaughters ( ) const
inlineinherited

◆ GetRegion()

G4Region * G4LogicalVolume::GetRegion ( ) const
inlineinherited

◆ GetSensitiveDetector()

G4VSensitiveDetector * G4LogicalVolume::GetSensitiveDetector ( ) const
inherited

◆ GetSmartless()

G4double G4LogicalVolume::GetSmartless ( ) const
inlineinherited

◆ GetSolid() [1/2]

G4VSolid * G4LogicalVolume::GetSolid ( ) const
inherited

Definition at line 413 of file G4LogicalVolume.cc.

414{
416}
static G4GEOM_DLL G4ThreadLocal T * offset

References G4LogicalVolume::GetSolid(), G4LogicalVolume::instanceID, G4GeomSplitter< T >::offset, and G4LogicalVolume::subInstanceManager.

Referenced by G4AdjointCrossSurfChecker::AddanExtSurfaceOfAvolume(), G4VSceneHandler::AddCompound(), G4HepRepFileSceneHandler::AddHepRepInstance(), G4GMocrenFileSceneHandler::AddPrimitive(), G4OpenGLStoredSceneHandler::AddPrimitivePreambleInternal(), G4GMocrenFileSceneHandler::AddSolid(), G4ReplicaNavigation::BackLocate(), G4PhantomParameterisation::BuildContainerSolid(), G4SmartVoxelHeader::BuildNodes(), G4SmartVoxelHeader::BuildReplicaVoxels(), G4SmartVoxelHeader::BuildVoxelsWithinLimits(), G4ReplicatedSlice::CheckAndSetParameters(), G4PVDivision::CheckAndSetParameters(), 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(), G4DrawVoxels::ComputeVoxelPolyhedra(), G4tgbVolume::ConstructG4PhysVol(), G4TheRayTracer::CreateBitMap(), G4PhysicalVolumeModel::CreateCurrentAttValues(), G3Division::CreatePVReplica(), G4ReflectionFactory::CreateReflectedLV(), G4AdjointPosOnPhysVolGenerator::DefinePhysicalVolume(), G4LogicalVolumeModel::DescribeYourselfTo(), G4VFieldModel::DescribeYourselfTo(), G4VisManager::Draw(), G4tgbGeometryDumper::DumpLogVol(), G4FastTrack::FRecordsAffineTransformation(), G4RTPrimaryGeneratorAction::GeneratePrimaries(), G4Navigator::GetGlobalExitNormal(), G4Navigator::GetLocalExitNormal(), G4ITNavigator1::GetLocalExitNormal(), G4ITNavigator2::GetLocalExitNormal(), G4VIntersectionLocator::GetLocalSurfaceNormal(), G4LogicalVolume::GetMass(), G4TransportationManager::GetParallelWorld(), G4ITTransportationManager::GetParallelWorld(), G4LogicalVolume::GetSolid(), G4BOptnForceCommonTruncatedExp::Initialize(), G4ITNavigator1::LocateGlobalPointAndSetup(), G4Navigator::LocateGlobalPointAndSetup(), G4GDMLWriteParamvol::ParametersWrite(), G4NeutrinoElectronProcess::PostStepDoIt(), G4VXTRenergyLoss::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4NavigationLogger::PreComputeStepLog(), G4PSFlatSurfaceCurrent::ProcessHits(), G4PSFlatSurfaceFlux::ProcessHits(), G4PSSphereSurfaceFlux::ProcessHits(), G4PSVolumeFlux::ProcessHits(), G4ITNavigator2::RecheckDistanceToCurrentBoundary(), G4NavigationLogger::ReportOutsideMother(), G4NavigationLogger::ReportVolumeAndIntersection(), G4VoxelSafety::SafetyForVoxelHeader(), G4VoxelSafety::SafetyForVoxelNode(), G4VisCommandsTouchable::SetNewValue(), G4ReplicatedSlice::SetParameterisation(), G4PVDivision::SetParameterisation(), G4RTPrimaryGeneratorAction::SetUp(), G4GeomTestVolume::TestOverlapInTree(), and G4GDMLWriteStructure::TraverseVolumeTree().

◆ GetSolid() [2/2]

G4VSolid * G4LogicalVolume::GetSolid ( G4LVData instLVdata)
staticinherited

Definition at line 408 of file G4LogicalVolume.cc.

409{
410 return instLVdata.fSolid;
411}
G4VSolid * fSolid

References G4LVData::fSolid.

◆ GetSubInstanceManager()

const G4LVManager & G4LogicalVolume::GetSubInstanceManager ( )
staticinherited

Definition at line 222 of file G4LogicalVolume.cc.

223{
224 return subInstanceManager;
225}

References G4LogicalVolume::subInstanceManager.

Referenced by G4GeometryWorkspace::G4GeometryWorkspace().

◆ GetUserLimits()

G4UserLimits * G4LogicalVolume::GetUserLimits ( ) const
inlineinherited

◆ GetVisAttributes()

const G4VisAttributes * G4LogicalVolume::GetVisAttributes ( ) const
inlineinherited

◆ GetVoxelHeader()

G4SmartVoxelHeader * G4LogicalVolume::GetVoxelHeader ( ) const
inlineinherited

◆ InitialiseWorker()

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

Definition at line 164 of file G4LogicalVolume.cc.

168{
170
171 SetSolid(pSolid);
172 SetSensitiveDetector(pSDetector); // How this object is available now ?
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()
void AssignFieldManager(G4FieldManager *fldMgr)
void SetSolid(G4VSolid *pSolid)
void SetSensitiveDetector(G4VSensitiveDetector *pSDetector)

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

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

◆ IsAncestor()

G4bool G4LogicalVolume::IsAncestor ( const G4VPhysicalVolume p) const
inherited

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

References G4LogicalVolume::fDaughters, and G4LogicalVolume::IsDaughter().

Referenced by export_G4LogicalVolume(), G4IStore::IsInWorld(), and G4WeightWindowStore::IsInWorld().

◆ IsDaughter()

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

◆ IsExtended()

G4bool G4LogicalCrystalVolume::IsExtended ( ) const
inlinevirtual

Reimplemented from G4LogicalVolume.

Definition at line 64 of file G4LogicalCrystalVolume.hh.

64{ return true; }

◆ IsLattice()

G4bool G4LogicalCrystalVolume::IsLattice ( G4LogicalVolume aLV)
static

Definition at line 62 of file G4LogicalCrystalVolume.cc.

63{
64 return std::find(fLCVvec.begin(), fLCVvec.end(), aLV) != fLCVvec.end();
65}

References fLCVvec.

Referenced by G4Channeling::GetMeanFreePath(), and G4Channeling::PostStepDoIt().

◆ IsRegion()

G4bool G4LogicalVolume::IsRegion ( ) const
inlineinherited

◆ IsRootRegion()

G4bool G4LogicalVolume::IsRootRegion ( ) const
inlineinherited

◆ IsToOptimise()

G4bool G4LogicalVolume::IsToOptimise ( ) const
inlineinherited

◆ Lock()

void G4LogicalVolume::Lock ( )
inlineinherited

◆ operator==()

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

◆ PropagateRegion()

void G4LogicalVolume::PropagateRegion ( )
inlineinherited

◆ RemoveDaughter()

void G4LogicalVolume::RemoveDaughter ( const G4VPhysicalVolume p)
inherited

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}

References G4LogicalVolume::fDaughters, G4LogicalVolume::fRegion, G4MT_mass, and G4Region::RegionModified().

Referenced by G4PhysicalVolumeStore::DeRegister(), and export_G4LogicalVolume().

◆ ResetMass()

void G4LogicalVolume::ResetMass ( )
inherited

Definition at line 399 of file G4LogicalVolume.cc.

400{
401 G4MT_mass= 0.0;
402}

References G4MT_mass.

Referenced by G4LogicalVolume::SetSolid().

◆ RotateToLattice()

const G4ThreeVector & G4LogicalCrystalVolume::RotateToLattice ( G4ThreeVector dir) const

Definition at line 122 of file G4LogicalCrystalVolume.cc.

123{
124 return dir.transform(fOrient);
125}
Hep3Vector & transform(const HepRotation &)
Definition: ThreeVectorR.cc:20

References fOrient, and CLHEP::Hep3Vector::transform().

Referenced by G4Channeling::UpdateParameters().

◆ RotateToSolid()

const G4ThreeVector & G4LogicalCrystalVolume::RotateToSolid ( G4ThreeVector dir) const

Definition at line 128 of file G4LogicalCrystalVolume.cc.

129{
130 return dir.transform(fInverse);
131}

References fInverse, and CLHEP::Hep3Vector::transform().

Referenced by G4Channeling::PostStepDoIt().

◆ SetBiasWeight()

void G4LogicalVolume::SetBiasWeight ( G4double  w)
inlineinherited

◆ SetFieldManager()

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

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}
size_t GetNoDaughters() const
G4VPhysicalVolume * GetDaughter(const G4int i) const

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

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

◆ SetMaterial()

void G4LogicalVolume::SetMaterial ( G4Material pMaterial)
inherited

Definition at line 448 of file G4LogicalVolume.cc.

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

References G4MT_mass, and G4MT_material.

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

◆ SetMaterialCutsCouple()

void G4LogicalVolume::SetMaterialCutsCouple ( G4MaterialCutsCouple cuts)
inherited

Definition at line 497 of file G4LogicalVolume.cc.

498{
499 G4MT_ccouple = cuts;
500}

References G4MT_ccouple.

Referenced by export_G4LogicalVolume(), and G4ProductionCutsTable::ScanAndSetCouple().

◆ SetMillerOrientation()

void G4LogicalCrystalVolume::SetMillerOrientation ( G4int  h,
G4int  k,
G4int  l,
G4double  rot = 0.0 
)

Definition at line 85 of file G4LogicalCrystalVolume.cc.

89{
90 // Align Miller normal vector (hkl) with +Z axis, and rotation about axis
91
92 if (verboseLevel)
93 {
94 G4cout << "G4LatticePhysical::SetMillerOrientation(" << h << " "
95 << k << " " << l << ", " << rot/CLHEP::deg << " deg)" << G4endl;
96 }
97
98 hMiller = h;
99 kMiller = k;
100 lMiller = l;
101 fRot = rot;
102
103 G4ThreeVector norm = (h*GetBasis(0)+k*GetBasis(1)+l*GetBasis(2)).unit();
104
105 if (verboseLevel>1) G4cout << " norm = " << norm << G4endl;
106
107 // Aligns geometry +Z axis with lattice (hkl) normal
109 fOrient.rotateZ(rot).rotateY(norm.theta()).rotateZ(norm.phi());
111
112 if (verboseLevel>1) G4cout << " fOrient = " << fOrient << G4endl;
113
114 // FIXME: Is this equivalent to (phi,theta,rot) Euler angles???
115}
G4GLOB_DLL std::ostream G4cout
double phi() const
double theta() const
HepRotation inverse() const
static DLL_API const HepRotation IDENTITY
Definition: Rotation.h:366
HepRotation & rotateZ(double delta)
Definition: Rotation.cc:87
HepRotation & rotateY(double delta)
Definition: Rotation.cc:74
const G4ThreeVector & GetBasis(G4int i) const
static constexpr double deg

References CLHEP::deg, fInverse, fOrient, fRot, G4cout, G4endl, GetBasis(), hMiller, CLHEP::HepRotation::IDENTITY, CLHEP::HepRotation::inverse(), kMiller, lMiller, CLHEP::Hep3Vector::phi(), CLHEP::HepRotation::rotateY(), CLHEP::HepRotation::rotateZ(), CLHEP::Hep3Vector::theta(), and verboseLevel.

Referenced by G4LogicalCrystalVolume().

◆ SetName()

void G4LogicalVolume::SetName ( const G4String pName)
inherited

◆ SetOptimisation()

void G4LogicalVolume::SetOptimisation ( G4bool  optim)
inlineinherited

Referenced by export_G4LogicalVolume().

◆ SetRegion()

void G4LogicalVolume::SetRegion ( G4Region reg)
inlineinherited

◆ SetRegionRootFlag()

void G4LogicalVolume::SetRegionRootFlag ( G4bool  rreg)
inlineinherited

◆ SetSensitiveDetector()

void G4LogicalVolume::SetSensitiveDetector ( G4VSensitiveDetector pSDetector)
inherited

◆ SetSmartless()

void G4LogicalVolume::SetSmartless ( G4double  s)
inlineinherited

Referenced by export_G4LogicalVolume().

◆ SetSolid() [1/2]

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

Definition at line 429 of file G4LogicalVolume.cc.

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

References G4LVData::fMass, and G4LVData::fSolid.

◆ SetSolid() [2/2]

void G4LogicalVolume::SetSolid ( G4VSolid pSolid)
inherited

◆ SetUserLimits()

void G4LogicalVolume::SetUserLimits ( G4UserLimits pULimits)
inlineinherited

◆ SetVerbose()

void G4LogicalCrystalVolume::SetVerbose ( G4int  aInt)
inline

Definition at line 81 of file G4LogicalCrystalVolume.hh.

81{ verboseLevel = aInt; }

References verboseLevel.

◆ SetVisAttributes() [1/2]

void G4LogicalVolume::SetVisAttributes ( const G4VisAttributes VA)
inherited

Definition at line 644 of file G4LogicalVolume.cc.

645{
647}
const G4VisAttributes * fVisAttributes

References G4LogicalVolume::fVisAttributes.

◆ SetVisAttributes() [2/2]

void G4LogicalVolume::SetVisAttributes ( const G4VisAttributes pVA)
inlineinherited

◆ SetVoxelHeader()

void G4LogicalVolume::SetVoxelHeader ( G4SmartVoxelHeader pVoxel)
inlineinherited

◆ TerminateWorker()

void G4LogicalVolume::TerminateWorker ( G4LogicalVolume ptrMasterObject)
inherited

Definition at line 211 of file G4LogicalVolume.cc.

213{
214}

Referenced by G4GeometryWorkspace::DestroyWorkspace().

◆ TotalVolumeEntities()

G4int G4LogicalVolume::TotalVolumeEntities ( ) const
inherited

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

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

Referenced by export_G4LogicalVolume(), and G4LogicalVolume::TotalVolumeEntities().

◆ UpdateMaterial()

void G4LogicalVolume::UpdateMaterial ( G4Material pMaterial)
inherited

Field Documentation

◆ fBiasWeight

G4double G4LogicalVolume::fBiasWeight = 1.0
privateinherited

Definition at line 412 of file G4LogicalVolume.hh.

◆ fDaughters

G4PhysicalVolumeList G4LogicalVolume::fDaughters
privateinherited

◆ fDaughtersVolumeType

EVolume G4LogicalVolume::fDaughtersVolumeType
privateinherited

◆ fFieldManager

G4FieldManager* G4LogicalVolume::fFieldManager = nullptr
privateinherited

◆ fInverse

G4RotationMatrix G4LogicalCrystalVolume::fInverse
private

Definition at line 88 of file G4LogicalCrystalVolume.hh.

Referenced by RotateToSolid(), and SetMillerOrientation().

◆ fLCVvec

std::vector< G4LogicalVolume * > G4LogicalCrystalVolume::fLCVvec
staticprivate

◆ fLock

G4bool G4LogicalVolume::fLock = false
privateinherited

Definition at line 432 of file G4LogicalVolume.hh.

Referenced by G4LogicalVolume::~G4LogicalVolume().

◆ fName

G4String G4LogicalVolume::fName
privateinherited

Definition at line 399 of file G4LogicalVolume.hh.

Referenced by G4LogicalVolume::GetMass(), and G4LogicalVolume::SetName().

◆ fOptimise

G4bool G4LogicalVolume::fOptimise = true
privateinherited

Definition at line 428 of file G4LogicalVolume.hh.

◆ fOrient

G4RotationMatrix G4LogicalCrystalVolume::fOrient
private

Definition at line 87 of file G4LogicalCrystalVolume.hh.

Referenced by RotateToLattice(), and SetMillerOrientation().

◆ fRegion

G4Region* G4LogicalVolume::fRegion = nullptr
privateinherited

◆ fRootRegion

G4bool G4LogicalVolume::fRootRegion = false
privateinherited

Definition at line 430 of file G4LogicalVolume.hh.

Referenced by G4LogicalVolume::~G4LogicalVolume().

◆ fRot

G4double G4LogicalCrystalVolume::fRot = 0.0
private

Definition at line 90 of file G4LogicalCrystalVolume.hh.

Referenced by SetMillerOrientation().

◆ fSensitiveDetector

G4VSensitiveDetector* G4LogicalVolume::fSensitiveDetector = nullptr
privateinherited

◆ fSmartless

G4double G4LogicalVolume::fSmartless = 2.0
privateinherited

Definition at line 405 of file G4LogicalVolume.hh.

◆ fSolid

G4VSolid* G4LogicalVolume::fSolid = nullptr
privateinherited

Definition at line 419 of file G4LogicalVolume.hh.

Referenced by G4LogicalVolume::G4LogicalVolume().

◆ fUserLimits

G4UserLimits* G4LogicalVolume::fUserLimits = nullptr
privateinherited

Definition at line 401 of file G4LogicalVolume.hh.

◆ fVisAttributes

const G4VisAttributes* G4LogicalVolume::fVisAttributes = nullptr
privateinherited

Definition at line 408 of file G4LogicalVolume.hh.

Referenced by G4LogicalVolume::SetVisAttributes().

◆ fVoxel

G4SmartVoxelHeader* G4LogicalVolume::fVoxel = nullptr
privateinherited

Definition at line 403 of file G4LogicalVolume.hh.

◆ hMiller

G4int G4LogicalCrystalVolume::hMiller = 1
private

Definition at line 89 of file G4LogicalCrystalVolume.hh.

Referenced by SetMillerOrientation().

◆ instanceID

G4int G4LogicalVolume::instanceID
privateinherited

◆ kMiller

G4int G4LogicalCrystalVolume::kMiller = 1
private

Definition at line 89 of file G4LogicalCrystalVolume.hh.

Referenced by SetMillerOrientation().

◆ lMiller

G4int G4LogicalCrystalVolume::lMiller = 0
private

Definition at line 89 of file G4LogicalCrystalVolume.hh.

Referenced by SetMillerOrientation().

◆ lvdata

G4LVData* G4LogicalVolume::lvdata = nullptr
privateinherited

◆ subInstanceManager

G4LVManager G4LogicalVolume::subInstanceManager
staticprivateinherited

◆ verboseLevel

G4int G4LogicalCrystalVolume::verboseLevel = 0
private

Definition at line 92 of file G4LogicalCrystalVolume.hh.

Referenced by SetMillerOrientation(), and SetVerbose().


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