Geant4-11
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Static Protected Attributes | Private Attributes
G4VPhysicalVolume Class Referenceabstract

#include <G4VPhysicalVolume.hh>

Inheritance diagram for G4VPhysicalVolume:
G4PVPlacement G4PVReplica G4VExternalPhysicalVolume G4PVDivision G4PVParameterised G4ReplicatedSlice

Public Member Functions

virtual G4bool CheckOverlaps (G4int res=1000, G4double tol=0., G4bool verbose=true, G4int errMax=1)
 
EVolume DeduceVolumeType () const
 
 G4VPhysicalVolume (__void__ &)
 
 G4VPhysicalVolume (const G4VPhysicalVolume &)=delete
 
 G4VPhysicalVolume (G4RotationMatrix *pRot, const G4ThreeVector &tlate, const G4String &pName, G4LogicalVolume *pLogical, G4VPhysicalVolume *pMother)
 
virtual G4int GetCopyNo () const =0
 
const G4RotationMatrixGetFrameRotation () const
 
G4ThreeVector GetFrameTranslation () const
 
G4int GetInstanceID () const
 
G4LogicalVolumeGetLogicalVolume () const
 
G4LogicalVolumeGetMotherLogical () const
 
virtual G4int GetMultiplicity () const
 
const G4StringGetName () const
 
G4RotationMatrixGetObjectRotation () const
 
G4RotationMatrix GetObjectRotationValue () const
 
G4ThreeVector GetObjectTranslation () const
 
virtual G4VPVParameterisationGetParameterisation () const =0
 
virtual G4int GetRegularStructureId () const =0
 
virtual void GetReplicationData (EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const =0
 
G4RotationMatrixGetRotation ()
 
const G4RotationMatrixGetRotation () const
 
const G4ThreeVector GetTranslation () const
 
virtual G4bool IsMany () const =0
 
virtual G4bool IsParameterised () const =0
 
virtual G4bool IsRegularStructure () const =0
 
virtual G4bool IsReplicated () const =0
 
G4VPhysicalVolumeoperator= (const G4VPhysicalVolume &)=delete
 
G4bool operator== (const G4VPhysicalVolume &p) const
 
virtual void SetCopyNo (G4int CopyNo)=0
 
void SetLogicalVolume (G4LogicalVolume *pLogical)
 
void SetMotherLogical (G4LogicalVolume *pMother)
 
void SetName (const G4String &pName)
 
void SetRotation (G4RotationMatrix *)
 
void SetTranslation (const G4ThreeVector &v)
 
virtual EVolume VolumeType () const =0
 
virtual ~G4VPhysicalVolume ()
 

Static Public Member Functions

static void Clean ()
 
static const G4PVManagerGetSubInstanceManager ()
 

Protected Member Functions

void InitialiseWorker (G4VPhysicalVolume *pMasterObject, G4RotationMatrix *pRot, const G4ThreeVector &tlate)
 
void TerminateWorker (G4VPhysicalVolume *pMasterObject)
 

Protected Attributes

G4int instanceID
 

Static Protected Attributes

static G4GEOM_DLL G4PVManager subInstanceManager
 

Private Attributes

G4LogicalVolumeflmother = nullptr
 
G4LogicalVolumeflogical = nullptr
 
G4String fname
 
G4PVDatapvdata = nullptr
 

Detailed Description

Definition at line 78 of file G4VPhysicalVolume.hh.

Constructor & Destructor Documentation

◆ G4VPhysicalVolume() [1/3]

G4VPhysicalVolume::G4VPhysicalVolume ( G4RotationMatrix pRot,
const G4ThreeVector tlate,
const G4String pName,
G4LogicalVolume pLogical,
G4VPhysicalVolume pMother 
)

Definition at line 54 of file G4VPhysicalVolume.cc.

59 : flogical(pLogical), fname(pName)
60{
62
63 this->SetRotation( pRot ); // G4MT_rot = pRot;
64 this->SetTranslation( tlate ); // G4MT_trans = tlate;
65
66 // Initialize 'Shadow' data structure - for use by object persistency
67 pvdata = new G4PVData();
68 pvdata->frot = pRot;
69 pvdata->tx = tlate.x();
70 pvdata->ty = tlate.y();
71 pvdata->tz = tlate.z();
72
74}
double z() const
double x() const
double y() const
G4int CreateSubInstance()
G4RotationMatrix * frot
static void Register(G4VPhysicalVolume *pSolid)
static G4GEOM_DLL G4PVManager subInstanceManager
G4LogicalVolume * flogical
void SetTranslation(const G4ThreeVector &v)
void SetRotation(G4RotationMatrix *)

References G4GeomSplitter< T >::CreateSubInstance(), G4PVData::frot, instanceID, pvdata, G4PhysicalVolumeStore::Register(), SetRotation(), SetTranslation(), subInstanceManager, G4PVData::tx, G4PVData::ty, G4PVData::tz, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ ~G4VPhysicalVolume()

G4VPhysicalVolume::~G4VPhysicalVolume ( )
virtual

Definition at line 91 of file G4VPhysicalVolume.cc.

92{
93 delete pvdata;
95}
static void DeRegister(G4VPhysicalVolume *pSolid)

References G4PhysicalVolumeStore::DeRegister(), and pvdata.

◆ G4VPhysicalVolume() [2/3]

G4VPhysicalVolume::G4VPhysicalVolume ( const G4VPhysicalVolume )
delete

◆ G4VPhysicalVolume() [3/3]

G4VPhysicalVolume::G4VPhysicalVolume ( __void__ &  )

Member Function Documentation

◆ CheckOverlaps()

G4bool G4VPhysicalVolume::CheckOverlaps ( G4int  res = 1000,
G4double  tol = 0.,
G4bool  verbose = true,
G4int  errMax = 1 
)
virtual

◆ Clean()

void G4VPhysicalVolume::Clean ( )
static

◆ DeduceVolumeType()

EVolume G4VPhysicalVolume::DeduceVolumeType ( ) const
inline

◆ GetCopyNo()

virtual G4int G4VPhysicalVolume::GetCopyNo ( ) const
pure virtual

◆ GetFrameRotation()

const G4RotationMatrix * G4VPhysicalVolume::GetFrameRotation ( ) const

◆ GetFrameTranslation()

G4ThreeVector G4VPhysicalVolume::GetFrameTranslation ( ) const

Definition at line 213 of file G4VPhysicalVolume.cc.

214{
216}
CLHEP::Hep3Vector G4ThreeVector
#define G4MT_ty
#define G4MT_tz
#define G4MT_tx

References G4MT_tx, G4MT_ty, and G4MT_tz.

◆ GetInstanceID()

G4int G4VPhysicalVolume::GetInstanceID ( ) const
inline

◆ GetLogicalVolume()

G4LogicalVolume * G4VPhysicalVolume::GetLogicalVolume ( ) const
inline

Referenced by G4AdjointCrossSurfChecker::AddanExtSurfaceOfAvolume(), G4VSceneHandler::AddCompound(), G4VtkSceneHandler::AddCompound(), G4LogicalVolume::AddDaughter(), G4OpenGLStoredSceneHandler::AddPrimitivePreambleInternal(), G4GMocrenFileSceneHandler::AddSolid(), G4VtkSceneHandler::AddSolid(), G4ParallelWorldProcess::AtRestDoIt(), G4ParallelWorldScoringProcess::AtRestDoIt(), G4FastSimulationManagerProcess::AtRestGetPhysicalInteractionLength(), G4ReplicaNavigation::BackLocate(), G4Region::BelongsTo(), G4PhantomParameterisation::BuildContainerSolid(), G4SmartVoxelHeader::BuildNodes(), G4GeometryManager::BuildOptimisations(), G4Track::CalculateVelocityForOpticalPhoton(), G4PVDivision::CheckAndSetParameters(), G4PVPlacement::CheckOverlaps(), G4GeometryWorkspace::CloneReplicaSolid(), G4AdjointPrimaryGenerator::ComputeAccumulatedDepthVectorAlongBackRay(), G4VPVParameterisation::ComputeMaterial(), G4Navigator::ComputeSafety(), G4ITNavigator1::ComputeSafety(), G4ITNavigator2::ComputeSafety(), G4NormalNavigation::ComputeSafety(), G4VoxelNavigation::ComputeSafety(), G4ReplicaNavigation::ComputeSafety(), G4ParameterisedNavigation::ComputeSafety(), G4VoxelSafety::ComputeSafety(), G4PolarizedAnnihilation::ComputeSaturationFactor(), G4PolarizedCompton::ComputeSaturationFactor(), G4PolarizedIonisation::ComputeSaturationFactor(), G4VNestedParameterisation::ComputeSolid(), G4VPVParameterisation::ComputeSolid(), G4PhantomParameterisation::ComputeSolid(), G4ParameterisedNavigation::ComputeStep(), G4RegularNavigation::ComputeStep(), G4VoxelNavigation::ComputeStep(), G4ReplicaNavigation::ComputeStep(), G4NormalNavigation::ComputeStep(), G4Navigator::ComputeStep(), G4ITNavigator1::ComputeStep(), G4ITNavigator2::ComputeStep(), G4PropagatorInField::ComputeStep(), G4RegularNavigation::ComputeStepSkippingEqualMaterials(), G4tgbVolume::ConstructG4Volumes(), G4TheRayTracer::CreateBitMap(), G4AdjointCrossSurfChecker::CrossingAnInterfaceBetweenTwoVolumes(), G4Radioactivation::DecayIt(), G4RadioactiveDecay::DecayIt(), G4AdjointPosOnPhysVolGenerator::DefinePhysicalVolume(), G4RunManagerKernel::DefineWorldVolume(), G4GeometryManager::DeleteOptimisations(), G4LogicalVolumeModel::DescribeYourselfTo(), G4PhysicalVolumeModel::DescribeYourselfTo(), G4VFieldModel::DescribeYourselfTo(), G4GeometryWorkspace::DestroyWorkspace(), G4GlobalFastSimulationManager::DisplayRegion(), G4GDMLWriteStructure::DivisionvolWrite(), G4VisManager::Draw(), G4TrajectoryDrawByOriginVolume::Draw(), G4tgbGeometryDumper::DumpPhysVol(), G4tgbGeometryDumper::DumpPVParameterised(), G4tgbGeometryDumper::DumpPVPlacement(), G4TrajectoryOriginVolumeFilter::Evaluate(), export_G4VPhysicalVolume(), G4PropagatorInField::FindAndSetFieldManager(), G4VReadOutGeometry::FindROTouchable(), G4FastTrack::FRecordsAffineTransformation(), G4Mesh::G4Mesh(), G4PVDivision::G4PVDivision(), G4PVParameterised::G4PVParameterised(), G4PVPlacement::G4PVPlacement(), G4PVReplica::G4PVReplica(), G4ReplicatedSlice::G4ReplicatedSlice(), G4VExternalPhysicalVolume::G4VExternalPhysicalVolume(), G4GDMLRead::GeneratePhysvolName(), G4RTPrimaryGeneratorAction::GeneratePrimaries(), G4Navigator::GetGlobalExitNormal(), G4Navigator::GetLocalExitNormal(), G4ITNavigator1::GetLocalExitNormal(), G4VIntersectionLocator::GetLocalSurfaceNormal(), G4LogicalVolume::GetMass(), G4Channeling::GetMatData(), G4Channeling::GetMeanFreePath(), G4VXTRenergyLoss::GetMeanFreePath(), G4NeutrinoElectronProcess::GetMeanFreePath(), G4ElNeutrinoNucleusProcess::GetMeanFreePath(), G4MuNeutrinoNucleusProcess::GetMeanFreePath(), G4VTransitionRadiation::GetMeanFreePath(), G4Navigator::GetMotherToDaughterTransform(), G4ITNavigator1::GetMotherToDaughterTransform(), G4ITNavigator2::GetMotherToDaughterTransform(), G4TransportationManager::GetParallelWorld(), G4ITTransportationManager::GetParallelWorld(), G4tgbGeometryDumper::GetTopPhysVol(), G4AdjointCrossSurfChecker::GoingInOrOutOfaVolumeByExtSurface(), G4GeometryWorkspace::InitialisePhysicalVolumes(), G4BOptnForceCommonTruncatedExp::Initialize(), G4IStore::IsInWorld(), G4WeightWindowStore::IsInWorld(), G4VEnergyLossProcess::IsRegionForCubcutProcessor(), G4ParameterisedNavigation::LevelLocate(), G4RegularNavigation::LevelLocate(), G4LatticeManager::LoadLattice(), G4ITNavigator1::LocateGlobalPointAndSetup(), G4ITNavigator2::LocateGlobalPointAndSetup(), G4Navigator::LocateGlobalPointAndSetup(), G4Navigator::LocateGlobalPointWithinVolume(), G4ITNavigator1::LocateGlobalPointWithinVolume(), G4ITNavigator2::LocateGlobalPointWithinVolume(), G4FastSimHitMaker::make(), GFlashHitMaker::make(), G4GDMLWriteParamvol::ParametersWrite(), G4GDMLWriteParamvol::ParamvolAlgorithmWrite(), G4GDMLWriteParamvol::ParamvolWrite(), G4GDMLWriteStructure::PhysvolWrite(), G4ParallelWorldProcess::PostStepDoIt(), G4ParallelWorldScoringProcess::PostStepDoIt(), G4ScoreSplittingProcess::PostStepDoIt(), G4Channeling::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4MicroElecSurface::PostStepDoIt(), G4ForwardXrayTR::PostStepDoIt(), G4VXTRenergyLoss::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4ITTransportation::PostStepDoIt(), G4VTransitionRadiation::PostStepDoIt(), G4CoupledTransportation::PostStepDoIt(), G4Transportation::PostStepDoIt(), G4MaxTimeCuts::PostStepGetPhysicalInteractionLength(), G4MinEkineCuts::PostStepGetPhysicalInteractionLength(), G4BiasingProcessInterface::PostStepGetPhysicalInteractionLength(), G4FastSimulationManagerProcess::PostStepGetPhysicalInteractionLength(), G4StepLimiter::PostStepGetPhysicalInteractionLength(), G4UserSpecialCuts::PostStepGetPhysicalInteractionLength(), G4AdjointForcedInteractionForGamma::PostStepGetPhysicalInteractionLength(), G4LowECapture::PostStepGetPhysicalInteractionLength(), G4NavigationLogger::PreComputeStepLog(), G4PSFlatSurfaceCurrent::ProcessHits(), G4PSFlatSurfaceFlux::ProcessHits(), G4PSSphereSurfaceFlux::ProcessHits(), G4PSVolumeFlux::ProcessHits(), G4ErrorFreeTrajState::PropagateErrorIoni(), G4ErrorFreeTrajState::PropagateErrorMSC(), G4ITNavigator2::RecheckDistanceToCurrentBoundary(), G4ReflectionFactory::ReflectPVDivision(), G4ReflectionFactory::ReflectPVPlacement(), G4ReflectionFactory::ReflectPVReplica(), G4LatticeManager::RegisterLattice(), G4GDMLWriteStructure::ReplicavolWrite(), G4PropagatorInField::ReportLoopingParticle(), G4NavigationLogger::ReportOutsideMother(), G4NavigationLogger::ReportVolumeAndIntersection(), G4ASCIITreeSceneHandler::RequestPrimitives(), G4VoxelSafety::SafetyForVoxelHeader(), G4VoxelSafety::SafetyForVoxelNode(), G4PolarizedAnnihilationModel::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), G4PolarizedIonisationModel::SampleSecondaries(), G4DNAIRT_geometries::Sampling(), G4ProductionCutsTable::ScanAndSetCouple(), G4Region::ScanVolumeTree(), G4LogicalVolume::SetFieldManager(), G4ITStepProcessor::SetInitialStep(), G4SteppingManager::SetInitialStep(), G4VVisCommandGeometrySet::SetLVVisAtts(), G4VisCommandsTouchable::SetNewValue(), G4RTPrimaryGeneratorAction::SetUp(), G4ScoringProbe::SetupGeometry(), G4ScoringBox::SetupGeometry(), G4ScoringCylinder::SetupGeometry(), G4Navigator::SetupHierarchy(), G4ITNavigator1::SetupHierarchy(), G4ITNavigator2::SetupHierarchy(), G4GlobalFastSimulationManager::ShowSetup(), G4SteppingManager::Stepping(), G4ParallelWorldProcess::SwitchMaterial(), G4GeomTestVolume::TestOverlapInTree(), G4GeomTestVolume::TestRecursiveOverlap(), G4LogicalVolume::TotalVolumeEntities(), G4GDMLWriteStructure::TraverseVolumeTree(), G4Channeling::UpdateParameters(), and G4MSSteppingAction::UserSteppingAction().

◆ GetMotherLogical()

G4LogicalVolume * G4VPhysicalVolume::GetMotherLogical ( ) const
inline

◆ GetMultiplicity()

G4int G4VPhysicalVolume::GetMultiplicity ( ) const
virtual

◆ GetName()

const G4String & G4VPhysicalVolume::GetName ( ) const
inline

Referenced by G4TransportationManager::ActivateNavigator(), G4VSceneHandler::AddCompound(), G4LogicalVolume::AddDaughter(), G4HepRepFileSceneHandler::AddHepRepInstance(), G4GDMLWrite::AddModule(), G4GMocrenFileSceneHandler::AddPrimitive(), G4GMocrenFileSceneHandler::AddSolid(), G4VtkSceneHandler::AddSolid(), G4DNABrownianTransportation::AlongStepGetPhysicalInteractionLength(), G4CoupledTransportation::AlongStepGetPhysicalInteractionLength(), G4GDMLWriteStructure::BorderSurfaceCache(), G4SmartVoxelHeader::BuildNodes(), G4ReplicatedSlice::CheckAndSetParameters(), G4PVDivision::CheckAndSetParameters(), G4PVReplica::CheckOnlyDaughter(), G4PVParameterised::CheckOverlaps(), G4PVPlacement::CheckOverlaps(), G4Navigator::CheckOverlapsIterative(), G4Navigator::ComputeSafety(), G4ITNavigator1::ComputeSafety(), G4ITNavigator2::ComputeSafety(), G4VoxelNavigation::ComputeSafety(), G4VoxelSafety::ComputeSafety(), G4PolarizedAnnihilation::ComputeSaturationFactor(), G4PolarizedCompton::ComputeSaturationFactor(), G4PolarizedIonisation::ComputeSaturationFactor(), G4ParameterisedNavigation::ComputeStep(), G4ReplicaNavigation::ComputeStep(), G4Navigator::ComputeStep(), G4ITNavigator1::ComputeStep(), G4ITNavigator2::ComputeStep(), G4PropagatorInField::ComputeStep(), G4tgbPlaceParamCircle::ComputeTransformation(), G4tgbPlaceParamLinear::ComputeTransformation(), G4tgbPlaceParamSquare::ComputeTransformation(), G4ImportanceConfigurator::Configure(), G4WeightCutOffConfigurator::Configure(), G4WeightWindowConfigurator::Configure(), G4tgbDetectorConstruction::Construct(), G4tgbDetectorBuilder::ConstructDetector(), G4tgbVolume::ConstructG4PhysVol(), G4FastSimulationPhysics::ConstructProcess(), G4Qt3DSceneHandler::CreateNewNode(), G4PathFinder::CreateTouchableHandle(), G4ITPathFinder::CreateTouchableHandle(), G4AdjointCrossSurfChecker::CrossingAnInterfaceBetweenTwoVolumes(), G4TransportationManager::DeActivateNavigator(), G4DNAMolecularDissociation::DecayIt(), G4PhysicalVolumeStore::DeRegister(), G4TransportationManager::DeRegisterNavigator(), G4TransportationManager::DeRegisterWorld(), G4ITTransportationManager::DeRegisterWorld(), G4PhysicalVolumeModel::DescribeAndDescend(), G4LogicalVolumeModel::DescribeYourselfTo(), G4GDMLWriteStructure::DivisionvolWrite(), G4TrajectoryDrawByOriginVolume::Draw(), G4tgbVolumeMgr::DumpG4PhysVolLeaf(), G4LogicalBorderSurface::DumpInfo(), G4tgbGeometryDumper::DumpPVParameterised(), G4tgbGeometryDumper::DumpPVPlacement(), G4tgbGeometryDumper::DumpPVReplica(), G4RunManagerKernel::DumpRegion(), G4HadronicProcess::DumpState(), G4MuonicAtomDecay::DumpState(), G4tgbVolumeMgr::DumpSummary(), G4ExceptionHandler::DumpTrackInfo(), G4ASCIITreeSceneHandler::EndModeling(), G4TrajectoryOriginVolumeFilter::Evaluate(), export_G4VPhysicalVolume(), G4FastTrack::FRecordsAffineTransformation(), G4FastSimulationManagerProcess::G4FastSimulationManagerProcess(), G4IStore::G4IStore(), G4PhysicalVolumeModel::G4PhysicalVolumeModel(), G4PVParameterised::G4PVParameterised(), G4PVReplica::G4PVReplica(), G4Navigator::GetGlobalExitNormal(), G4PSDoseDeposit3D::GetIndex(), G4PSEnergyDeposit3D::GetIndex(), G4LatticeManager::GetLattice(), G4Navigator::GetLocalExitNormal(), G4ITNavigator1::GetLocalExitNormal(), G4ModelingParameters::PVPointerCopyNo::GetName(), G4TransportationManager::GetNavigator(), G4ITTransportationManager::GetNavigator(), G4tgbVolumeMgr::GetTopPhysVol(), G4AdjointCrossSurfChecker::GoingInOrOutOfaVolume(), G4AdjointCrossSurfChecker::GoingInOrOutOfaVolumeByExtSurface(), G4SPSPosDistribution::IsSourceConfined(), G4LatticeManager::LoadLattice(), G4ITMultiNavigator::LocateGlobalPointAndSetup(), G4ITNavigator1::LocateGlobalPointAndSetup(), G4ITNavigator2::LocateGlobalPointAndSetup(), G4MultiNavigator::LocateGlobalPointAndSetup(), G4Navigator::LocateGlobalPointAndSetup(), operator<<(), G4GDMLWriteParamvol::ParametersWrite(), Path(), G4GDMLWriteStructure::PhysvolWrite(), G4UCNAbsorption::PostStepDoIt(), G4UCNBoundaryProcess::PostStepDoIt(), G4UCNMultiScattering::PostStepDoIt(), G4MicroElecSurface::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4CoupledTransportation::PostStepDoIt(), G4ITSteppingVerbose::PostStepVerbose(), G4NavigationLogger::PreComputeStepLog(), G4MultiNavigator::PrepareNavigators(), G4ITMultiNavigator::PrepareNavigators(), G4ITSteppingVerbose::PreStepVerbose(), G4MultiNavigator::PrintLimited(), G4PathFinder::PrintLimited(), G4ITMultiNavigator::PrintLimited(), G4ITPathFinder::PrintLimited(), G4Navigator::PrintState(), G4ITNavigator1::PrintState(), G4PropagatorInField::printStatus(), G4PhysicalVolumeMassScene::ProcessVolume(), G4PhysicalVolumesSearchScene::ProcessVolume(), G4ReflectionFactory::ReflectPVDivision(), G4ReflectionFactory::ReflectPVParameterised(), G4ReflectionFactory::ReflectPVPlacement(), G4ReflectionFactory::ReflectPVReplica(), G4PhysicalVolumeStore::Register(), G4tgbVolumeMgr::RegisterMe(), G4PropagatorInField::ReportLoopingParticle(), G4PropagatorInField::ReportStuckParticle(), G4NavigationLogger::ReportVolumeAndIntersection(), G4ASCIITreeSceneHandler::RequestPrimitives(), G4VoxelSafety::SafetyForVoxelHeader(), G4Region::ScanVolumeTree(), G4VisCommandsTouchable::SetNewValue(), G4WeightCutOffProcess::SetParallelWorld(), G4WeightWindowProcess::SetParallelWorld(), G4ParallelWorldProcess::SetParallelWorld(), G4ParallelWorldScoringProcess::SetParallelWorld(), G4IStore::SetParallelWorldVolume(), G4IStore::SetWorldVolume(), G4WeightWindowStore::SetWorldVolume(), G4FastSimulationManagerProcess::SetWorldVolume(), G4GlobalFastSimulationManager::ShowSetup(), G4ITSteppingVerbose::ShowStep(), G4SteppingVerbose::ShowStep(), G4SteppingVerboseWithUnits::ShowStep(), G4ITSteppingVerbose::StepInfo(), G4SteppingVerbose::StepInfo(), G4SteppingVerboseWithUnits::StepInfo(), G4ITSteppingVerbose::StepInfoForLeadingTrack(), G4GDMLRead::StripNames(), G4ErrorGeomVolumeTarget::TargetReached(), G4GeomTestVolume::TestOverlapInTree(), G4ITSteppingVerbose::TrackingEnded(), G4SteppingVerbose::TrackingStarted(), G4SteppingVerboseWithUnits::TrackingStarted(), G4ITSteppingVerbose::TrackingStarted(), G4ParallelWorldScoringProcess::Verbose(), G4ScoreSplittingProcess::Verbose(), G4ITSteppingVerbose::VerboseTrack(), G4SteppingVerbose::VerboseTrack(), G4SteppingVerboseWithUnits::VerboseTrack(), G4VScoringMesh::WorkerConstruct(), and G4RunManagerKernel::WorkerUpdateWorldVolume().

◆ GetObjectRotation()

G4RotationMatrix * G4VPhysicalVolume::GetObjectRotation ( ) const

Definition at line 175 of file G4VPhysicalVolume.cc.

176{
177 static G4RotationMatrix aRotM;
178 static G4RotationMatrix IdentityRM;
179
180 G4RotationMatrix* retval = &IdentityRM;
181
182 // Insure against frot being a null pointer
183 if(this->GetRotation())
184 {
185 aRotM = GetRotation()->inverse();
186 retval= &aRotM;
187 }
188 return retval;
189}
HepRotation inverse() const
const G4RotationMatrix * GetRotation() const

References GetRotation(), and CLHEP::HepRotation::inverse().

◆ GetObjectRotationValue()

G4RotationMatrix G4VPhysicalVolume::GetObjectRotationValue ( ) const

Definition at line 191 of file G4VPhysicalVolume.cc.

192{
193 G4RotationMatrix aRotM; // Initialised to identity
194
195 // Insure against G4MT_rot being a null pointer
196 if(G4MT_rot)
197 {
198 aRotM = G4MT_rot->inverse();
199 }
200 return aRotM;
201}

References G4MT_rot.

Referenced by G4PhysicalVolumeModel::CreateCurrentAttValues(), export_G4VPhysicalVolume(), G4GDMLWriteParamvol::ParametersWrite(), and G4ReflectionFactory::ReflectPVPlacement().

◆ GetObjectTranslation()

G4ThreeVector G4VPhysicalVolume::GetObjectTranslation ( ) const

◆ GetParameterisation()

virtual G4VPVParameterisation * G4VPhysicalVolume::GetParameterisation ( ) const
pure virtual

◆ GetRegularStructureId()

virtual G4int G4VPhysicalVolume::GetRegularStructureId ( ) const
pure virtual

◆ GetReplicationData()

virtual void G4VPhysicalVolume::GetReplicationData ( EAxis axis,
G4int nReplicas,
G4double width,
G4double offset,
G4bool consuming 
) const
pure virtual

◆ GetRotation() [1/2]

G4RotationMatrix * G4VPhysicalVolume::GetRotation ( )

Definition at line 165 of file G4VPhysicalVolume.cc.

166{
167 return G4MT_rot;
168}

References G4MT_rot.

◆ GetRotation() [2/2]

const G4RotationMatrix * G4VPhysicalVolume::GetRotation ( ) const

◆ GetSubInstanceManager()

const G4PVManager & G4VPhysicalVolume::GetSubInstanceManager ( )
static

Definition at line 140 of file G4VPhysicalVolume.cc.

141{
142 return subInstanceManager;
143}

References subInstanceManager.

Referenced by G4GeometryWorkspace::G4GeometryWorkspace().

◆ GetTranslation()

const G4ThreeVector G4VPhysicalVolume::GetTranslation ( ) const

◆ InitialiseWorker()

void G4VPhysicalVolume::InitialiseWorker ( G4VPhysicalVolume pMasterObject,
G4RotationMatrix pRot,
const G4ThreeVector tlate 
)
protected

Definition at line 111 of file G4VPhysicalVolume.cc.

115{
117
118 this->SetRotation( pRot ); // G4MT_rot = pRot;
119 this->SetTranslation( tlate ); // G4MT_trans = tlate;
120 // G4PhysicalVolumeStore::Register(this);
121}
void SlaveCopySubInstanceArray()

References SetRotation(), SetTranslation(), G4GeomSplitter< T >::SlaveCopySubInstanceArray(), and subInstanceManager.

Referenced by G4PVReplica::InitialiseWorker().

◆ IsMany()

virtual G4bool G4VPhysicalVolume::IsMany ( ) const
pure virtual

◆ IsParameterised()

virtual G4bool G4VPhysicalVolume::IsParameterised ( ) const
pure virtual

◆ IsRegularStructure()

virtual G4bool G4VPhysicalVolume::IsRegularStructure ( ) const
pure virtual

◆ IsReplicated()

virtual G4bool G4VPhysicalVolume::IsReplicated ( ) const
pure virtual

◆ operator=()

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

◆ operator==()

G4bool G4VPhysicalVolume::operator== ( const G4VPhysicalVolume p) const
inline

◆ SetCopyNo()

virtual void G4VPhysicalVolume::SetCopyNo ( G4int  CopyNo)
pure virtual

◆ SetLogicalVolume()

void G4VPhysicalVolume::SetLogicalVolume ( G4LogicalVolume pLogical)
inline

◆ SetMotherLogical()

void G4VPhysicalVolume::SetMotherLogical ( G4LogicalVolume pMother)
inline

◆ SetName()

void G4VPhysicalVolume::SetName ( const G4String pName)

◆ SetRotation()

void G4VPhysicalVolume::SetRotation ( G4RotationMatrix pRot)

◆ SetTranslation()

void G4VPhysicalVolume::SetTranslation ( const G4ThreeVector v)

Definition at line 155 of file G4VPhysicalVolume.cc.

156{
157 G4MT_tx=vec.x(); G4MT_ty=vec.y(); G4MT_tz=vec.z();
158}

References G4MT_tx, G4MT_ty, G4MT_tz, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by G4ParameterisationBoxX::ComputeTransformation(), G4ParameterisationBoxY::ComputeTransformation(), G4ParameterisationBoxZ::ComputeTransformation(), G4ParameterisationConsRho::ComputeTransformation(), G4ParameterisationConsPhi::ComputeTransformation(), G4ParameterisationConsZ::ComputeTransformation(), G4ParameterisationParaX::ComputeTransformation(), G4ParameterisationParaY::ComputeTransformation(), G4ParameterisationParaZ::ComputeTransformation(), G4ParameterisationPolyconeRho::ComputeTransformation(), G4ParameterisationPolyconePhi::ComputeTransformation(), G4ParameterisationPolyconeZ::ComputeTransformation(), G4ParameterisationPolyhedraRho::ComputeTransformation(), G4ParameterisationPolyhedraPhi::ComputeTransformation(), G4ParameterisationPolyhedraZ::ComputeTransformation(), G4ParameterisationTrdX::ComputeTransformation(), G4ParameterisationTrdY::ComputeTransformation(), G4ParameterisationTrdZ::ComputeTransformation(), G4ParameterisationTubsRho::ComputeTransformation(), G4ParameterisationTubsPhi::ComputeTransformation(), G4ParameterisationTubsZ::ComputeTransformation(), G4tgbPlaceParamCircle::ComputeTransformation(), G4tgbPlaceParamLinear::ComputeTransformation(), G4tgbPlaceParamSquare::ComputeTransformation(), G4ReplicaNavigation::ComputeTransformation(), G4PartialPhantomParameterisation::ComputeTransformation(), G4PhantomParameterisation::ComputeTransformation(), G4GDMLParameterisation::ComputeTransformation(), export_G4VPhysicalVolume(), G4VPhysicalVolume(), and InitialiseWorker().

◆ TerminateWorker()

void G4VPhysicalVolume::TerminateWorker ( G4VPhysicalVolume pMasterObject)
protected

Definition at line 134 of file G4VPhysicalVolume.cc.

135{
136}

◆ VolumeType()

virtual EVolume G4VPhysicalVolume::VolumeType ( ) const
pure virtual

Field Documentation

◆ flmother

G4LogicalVolume* G4VPhysicalVolume::flmother = nullptr
private

Definition at line 246 of file G4VPhysicalVolume.hh.

◆ flogical

G4LogicalVolume* G4VPhysicalVolume::flogical = nullptr
private

Definition at line 242 of file G4VPhysicalVolume.hh.

◆ fname

G4String G4VPhysicalVolume::fname
private

Definition at line 245 of file G4VPhysicalVolume.hh.

Referenced by SetName().

◆ instanceID

G4int G4VPhysicalVolume::instanceID
protected

Definition at line 234 of file G4VPhysicalVolume.hh.

Referenced by G4VPhysicalVolume().

◆ pvdata

G4PVData* G4VPhysicalVolume::pvdata = nullptr
private

Definition at line 248 of file G4VPhysicalVolume.hh.

Referenced by G4VPhysicalVolume(), and ~G4VPhysicalVolume().

◆ subInstanceManager

G4PVManager G4VPhysicalVolume::subInstanceManager
staticprotected

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