Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VPhysicalVolume Class Referenceabstract

#include <G4VPhysicalVolume.hh>

+ Inheritance diagram for G4VPhysicalVolume:

Public Member Functions

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

Static Public Member Functions

static const G4PVManagerGetSubInstanceManager ()
 
static void Clean ()
 

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
 

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
void SetTranslation(const G4ThreeVector &v)
void SetRotation(G4RotationMatrix *)

◆ ~G4VPhysicalVolume()

G4VPhysicalVolume::~G4VPhysicalVolume ( )
virtual

Definition at line 91 of file G4VPhysicalVolume.cc.

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

◆ G4VPhysicalVolume() [2/3]

G4VPhysicalVolume::G4VPhysicalVolume ( const G4VPhysicalVolume & )
delete

◆ G4VPhysicalVolume() [3/3]

G4VPhysicalVolume::G4VPhysicalVolume ( __void__ & )

Definition at line 79 of file G4VPhysicalVolume.cc.

80 : fname("")
81{
82 // Register to store
83 //
85
87}

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()

◆ GetFrameRotation()

◆ 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

◆ GetInstanceID()

G4int G4VPhysicalVolume::GetInstanceID ( ) const
inline

◆ GetLogicalVolume()

G4LogicalVolume * G4VPhysicalVolume::GetLogicalVolume ( ) const
inline

Referenced by G4AdjointCrossSurfChecker::AddanExtSurfaceOfAvolume(), G4PseudoScene::AddCompound(), G4VSceneHandler::AddCompound(), G4LogicalVolume::AddDaughter(), G4GMocrenFileSceneHandler::AddSolid(), G4VtkSceneHandler::AddSolid(), G4ParallelWorldProcess::AtRestDoIt(), G4ParallelWorldScoringProcess::AtRestDoIt(), G4FastSimulationManagerProcess::AtRestGetPhysicalInteractionLength(), G4ReplicaNavigation::BackLocate(), G4Region::BelongsTo(), G4PhantomParameterisation::BuildContainerSolid(), G4SmartVoxelHeader::BuildNodes(), G4Track::CalculateVelocityForOpticalPhoton(), G4PVPlacement::CheckOverlaps(), G4GeometryWorkspace::CloneReplicaSolid(), G4AdjointPrimaryGenerator::ComputeAccumulatedDepthVectorAlongBackRay(), G4VPVParameterisation::ComputeMaterial(), G4ITNavigator1::ComputeSafety(), G4ITNavigator2::ComputeSafety(), G4NormalNavigation::ComputeSafety(), G4ParameterisedNavigation::ComputeSafety(), G4ReplicaNavigation::ComputeSafety(), G4VoxelNavigation::ComputeSafety(), G4VoxelSafety::ComputeSafety(), G4PhantomParameterisation::ComputeSolid(), G4VNestedParameterisation::ComputeSolid(), G4VPVParameterisation::ComputeSolid(), G4ITNavigator1::ComputeStep(), G4ITNavigator2::ComputeStep(), G4Navigator::ComputeStep(), G4NormalNavigation::ComputeStep(), G4ParameterisedNavigation::ComputeStep(), G4RegularNavigation::ComputeStep(), G4ReplicaNavigation::ComputeStep(), G4VoxelNavigation::ComputeStep(), G4RegularNavigation::ComputeStepSkippingEqualMaterials(), G4tgbVolume::ConstructG4Volumes(), G4TheRayTracer::CreateBitMap(), G4AdjointCrossSurfChecker::CrossingAnInterfaceBetweenTwoVolumes(), G4Radioactivation::DecayIt(), G4RadioactiveDecay::DecayIt(), G4AdjointPosOnPhysVolGenerator::DefinePhysicalVolume(), G4RunManagerKernel::DefineWorldVolume(), G4LogicalVolumeModel::DescribeYourselfTo(), G4PhysicalVolumeModel::DescribeYourselfTo(), G4VFieldModel::DescribeYourselfTo(), G4GDMLWriteStructure::DivisionvolWrite(), G4TrajectoryDrawByOriginVolume::Draw(), G4VisManager::Draw(), G4tgbGeometryDumper::DumpPhysVol(), G4tgbGeometryDumper::DumpPVParameterised(), G4tgbGeometryDumper::DumpPVPlacement(), G4TrajectoryOriginVolumeFilter::Evaluate(), G4PropagatorInField::FindAndSetFieldManager(), G4Mesh::G4Mesh(), G4PVDivision::G4PVDivision(), G4PVParameterised::G4PVParameterised(), G4PVPlacement::G4PVPlacement(), G4PVPlacement::G4PVPlacement(), G4PVReplica::G4PVReplica(), G4ReplicatedSlice::G4ReplicatedSlice(), G4ReplicatedSlice::G4ReplicatedSlice(), G4ReplicatedSlice::G4ReplicatedSlice(), G4VExternalPhysicalVolume::G4VExternalPhysicalVolume(), G4GDMLRead::GeneratePhysvolName(), G4RTPrimaryGeneratorAction::GeneratePrimaries(), G4Navigator::GetGlobalExitNormal(), G4ITNavigator1::GetLocalExitNormal(), G4Navigator::GetLocalExitNormal(), G4LogicalVolume::GetMass(), G4Channeling::GetMeanFreePath(), G4ElNeutrinoNucleusProcess::GetMeanFreePath(), G4MuNeutrinoNucleusProcess::GetMeanFreePath(), G4NeutrinoElectronProcess::GetMeanFreePath(), G4NuVacOscProcess::GetMeanFreePath(), G4TauNeutrinoNucleusProcess::GetMeanFreePath(), G4VTransitionRadiation::GetMeanFreePath(), G4VXTRenergyLoss::GetMeanFreePath(), G4XrayReflection::GetMeanFreePath(), G4ITNavigator1::GetMotherToDaughterTransform(), G4ITNavigator2::GetMotherToDaughterTransform(), G4Navigator::GetMotherToDaughterTransform(), G4ITTransportationManager::GetParallelWorld(), G4TransportationManager::GetParallelWorld(), G4tgbGeometryDumper::GetTopPhysVol(), G4AdjointCrossSurfChecker::GoingInOrOutOfaVolumeByExtSurface(), G4BOptnForceCommonTruncatedExp::Initialize(), G4ParameterisedNavigation::LevelLocate(), G4RegularNavigation::LevelLocate(), G4LatticeManager::LoadLattice(), G4ITNavigator1::LocateGlobalPointAndSetup(), G4ITNavigator2::LocateGlobalPointAndSetup(), G4Navigator::LocateGlobalPointAndSetup(), G4ITNavigator1::LocateGlobalPointWithinVolume(), G4ITNavigator2::LocateGlobalPointWithinVolume(), G4Navigator::LocateGlobalPointWithinVolume(), G4FastSimHitMaker::make(), GFlashHitMaker::make(), G4GDMLWriteParamvol::ParametersWrite(), G4GDMLWriteParamvol::ParamvolAlgorithmWrite(), G4GDMLWriteParamvol::ParamvolWrite(), G4GDMLWriteStructure::PhysvolWrite(), G4Channeling::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4ForwardXrayTR::PostStepDoIt(), G4ITTransportation::PostStepDoIt(), G4MicroElecSurface::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4NuVacOscProcess::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4ParallelWorldProcess::PostStepDoIt(), G4ParallelWorldScoringProcess::PostStepDoIt(), G4ScoreSplittingProcess::PostStepDoIt(), G4TauNeutrinoNucleusProcess::PostStepDoIt(), G4VTransitionRadiation::PostStepDoIt(), G4VXTRenergyLoss::PostStepDoIt(), G4AdjointForcedInteractionForGamma::PostStepGetPhysicalInteractionLength(), G4BiasingProcessInterface::PostStepGetPhysicalInteractionLength(), G4FastSimulationManagerProcess::PostStepGetPhysicalInteractionLength(), G4LowECapture::PostStepGetPhysicalInteractionLength(), G4MaxTimeCuts::PostStepGetPhysicalInteractionLength(), G4MinEkineCuts::PostStepGetPhysicalInteractionLength(), G4StepLimiter::PostStepGetPhysicalInteractionLength(), G4UserSpecialCuts::PostStepGetPhysicalInteractionLength(), G4NavigationLogger::PreComputeStepLog(), G4PSFlatSurfaceCurrent::ProcessHits(), G4PSFlatSurfaceFlux::ProcessHits(), G4PSSphereSurfaceFlux::ProcessHits(), G4PSVolumeFlux::ProcessHits(), G4SafetyCalculator::QuickLocateWithinVolume(), G4ITNavigator2::RecheckDistanceToCurrentBoundary(), G4LatticeManager::RegisterLattice(), G4ParameterisedNavigation::RelocateWithinVolume(), G4VoxelNavigation::RelocateWithinVolume(), G4GDMLWriteStructure::ReplicavolWrite(), G4PropagatorInField::ReportLoopingParticle(), G4TransportationLogger::ReportLoopingTrack(), G4NavigationLogger::ReportOutsideMother(), G4NavigationLogger::ReportVolumeAndIntersection(), G4ASCIITreeSceneHandler::RequestPrimitives(), G4VoxelSafety::SafetyForVoxelHeader(), G4VoxelSafety::SafetyForVoxelNode(), G4SafetyCalculator::SafetyInCurrentVolume(), G4PolarizedAnnihilationModel::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), G4PolarizedIonisationModel::SampleSecondaries(), G4Region::ScanVolumeTree(), G4LogicalVolume::SetFieldManager(), G4ITStepProcessor::SetInitialStep(), G4SteppingManager::SetInitialStep(), G4VVisCommandGeometrySet::SetLVVisAtts(), G4VisCommandsTouchable::SetNewValue(), G4Transportation::SetTouchableInformation(), G4RTPrimaryGeneratorAction::SetUp(), G4ScoringBox::SetupGeometry(), G4ScoringCylinder::SetupGeometry(), G4ScoringProbe::SetupGeometry(), G4ITNavigator1::SetupHierarchy(), G4ITNavigator2::SetupHierarchy(), G4Navigator::SetupHierarchy(), G4GlobalFastSimulationManager::ShowSetup(), G4VSceneHandler::StandardSpecialMeshRendering(), G4SteppingManager::Stepping(), G4ParallelWorldProcess::SwitchMaterial(), G4GeomTestVolume::TestOverlapInTree(), G4GeomTestVolume::TestRecursiveOverlap(), G4LogicalVolume::TotalVolumeEntities(), G4GDMLWriteStructure::TraverseVolumeTree(), and G4MSSteppingAction::UserSteppingAction().

◆ GetMotherLogical()

◆ GetMultiplicity()

◆ GetName()

const G4String & G4VPhysicalVolume::GetName ( ) const
inline

Referenced by G4TransportationManager::ActivateNavigator(), G4VSceneHandler::AddCompound(), G4LogicalVolume::AddDaughter(), G4GDMLWrite::AddModule(), G4GMocrenFileSceneHandler::AddPrimitive(), G4GMocrenFileSceneHandler::AddSolid(), G4VtkSceneHandler::AddSolid(), G4CoupledTransportation::AlongStepGetPhysicalInteractionLength(), G4DNABrownianTransportation::AlongStepGetPhysicalInteractionLength(), G4GDMLWriteStructure::BorderSurfaceCache(), G4SmartVoxelHeader::BuildNodes(), G4PVParameterised::CheckOverlaps(), G4PVPlacement::CheckOverlaps(), G4Navigator::CheckOverlapsIterative(), G4SafetyCalculator::CompareSafetyValues(), G4ITNavigator1::ComputeSafety(), G4ITNavigator2::ComputeSafety(), G4VoxelNavigation::ComputeSafety(), G4VoxelSafety::ComputeSafety(), G4ITNavigator1::ComputeStep(), G4ITNavigator2::ComputeStep(), G4Navigator::ComputeStep(), G4ParameterisedNavigation::ComputeStep(), G4PropagatorInField::ComputeStep(), G4ReplicaNavigation::ComputeStep(), G4tgbPlaceParamCircle::ComputeTransformation(), G4tgbPlaceParamLinear::ComputeTransformation(), G4tgbPlaceParamSquare::ComputeTransformation(), G4ImportanceConfigurator::Configure(), G4WeightCutOffConfigurator::Configure(), G4WeightWindowConfigurator::Configure(), G4tgbDetectorConstruction::Construct(), G4tgbDetectorBuilder::ConstructDetector(), G4tgbVolume::ConstructG4PhysVol(), G4FastSimulationPhysics::ConstructProcess(), G4Qt3DSceneHandler::CreateNewNode(), G4ITPathFinder::CreateTouchableHandle(), G4PathFinder::CreateTouchableHandle(), G4AdjointCrossSurfChecker::CrossingAnInterfaceBetweenTwoVolumes(), G4TransportationManager::DeActivateNavigator(), G4DNAMolecularDissociation::DecayIt(), G4PhysicalVolumeStore::DeRegister(), G4TransportationManager::DeRegisterNavigator(), G4PhysicalVolumeModel::DescribeAndDescend(), G4GDMLWriteStructure::DivisionvolWrite(), G4TrajectoryDrawByOriginVolume::Draw(), G4VSceneHandler::Draw3DRectMeshAsDots(), G4VSceneHandler::Draw3DRectMeshAsSurfaces(), G4VSceneHandler::DrawTetMeshAsDots(), G4VSceneHandler::DrawTetMeshAsSurfaces(), G4tgbVolumeMgr::DumpG4PhysVolLeaf(), G4LogicalBorderSurface::DumpInfo(), G4tgbGeometryDumper::DumpPVParameterised(), G4tgbGeometryDumper::DumpPVPlacement(), G4tgbGeometryDumper::DumpPVReplica(), G4RunManagerKernel::DumpRegion(), G4HadronicProcess::DumpState(), G4tgbVolumeMgr::DumpSummary(), G4ASCIITreeSceneHandler::EndModeling(), G4TrajectoryOriginVolumeFilter::Evaluate(), G4FastSimulationManagerProcess::G4FastSimulationManagerProcess(), G4FastSimulationManagerProcess::G4FastSimulationManagerProcess(), G4FastSimulationManagerProcess::G4FastSimulationManagerProcess(), G4IStore::G4IStore(), G4PhysicalVolumeModel::G4PhysicalVolumeModel(), G4PVParameterised::G4PVParameterised(), G4PVReplica::G4PVReplica(), G4Navigator::GetGlobalExitNormal(), G4PSDoseDeposit3D::GetIndex(), G4PSEnergyDeposit3D::GetIndex(), G4LatticeManager::GetLattice(), G4ITNavigator1::GetLocalExitNormal(), G4Navigator::GetLocalExitNormal(), G4ModelingParameters::PVPointerCopyNo::GetName(), G4ITTransportationManager::GetNavigator(), G4TransportationManager::GetNavigator(), G4tgbVolumeMgr::GetTopPhysVol(), G4AdjointCrossSurfChecker::GoingInOrOutOfaVolume(), G4AdjointCrossSurfChecker::GoingInOrOutOfaVolumeByExtSurface(), G4LatticeManager::LoadLattice(), G4ITMultiNavigator::LocateGlobalPointAndSetup(), G4ITNavigator1::LocateGlobalPointAndSetup(), G4ITNavigator2::LocateGlobalPointAndSetup(), G4MultiNavigator::LocateGlobalPointAndSetup(), G4Navigator::LocateGlobalPointAndSetup(), operator<<(), operator<<(), operator<<(), G4GDMLWriteParamvol::ParametersWrite(), G4GDMLWriteStructure::PhysvolWrite(), G4CoupledTransportation::PostStepDoIt(), G4MicroElecSurface::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4UCNAbsorption::PostStepDoIt(), G4UCNBoundaryProcess::PostStepDoIt(), G4UCNMultiScattering::PostStepDoIt(), G4ITSteppingVerbose::PostStepVerbose(), G4NavigationLogger::PreComputeStepLog(), G4ITMultiNavigator::PrepareNavigators(), G4MultiNavigator::PrepareNavigators(), G4ITSteppingVerbose::PreStepVerbose(), G4ITMultiNavigator::PrintLimited(), G4ITPathFinder::PrintLimited(), G4MultiNavigator::PrintLimited(), G4PathFinder::PrintLimited(), G4ITNavigator1::PrintState(), G4Navigator::PrintState(), G4PropagatorInField::printStatus(), G4PhysicalVolumeStore::Register(), G4tgbVolumeMgr::RegisterMe(), G4PropagatorInField::ReportLoopingParticle(), G4PropagatorInField::ReportStuckParticle(), G4NavigationLogger::ReportVolumeAndIntersection(), G4ASCIITreeSceneHandler::RequestPrimitives(), G4VoxelSafety::SafetyForVoxelHeader(), G4SafetyCalculator::SafetyInCurrentVolume(), G4Region::ScanVolumeTree(), G4VisCommandsTouchable::SetNewValue(), G4ParallelWorldProcess::SetParallelWorld(), G4ParallelWorldScoringProcess::SetParallelWorld(), G4WeightCutOffProcess::SetParallelWorld(), G4WeightWindowProcess::SetParallelWorld(), G4IStore::SetParallelWorldVolume(), G4FastSimulationManagerProcess::SetWorldVolume(), G4FastSimulationManagerProcess::SetWorldVolume(), G4IStore::SetWorldVolume(), G4WeightWindowStore::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(), G4ITSteppingVerbose::TrackingStarted(), G4SteppingVerbose::TrackingStarted(), G4SteppingVerboseWithUnits::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() != nullptr)
184 {
185 aRotM = GetRotation()->inverse();
186 retval= &aRotM;
187 }
188 return retval;
189}
HepRotation inverse() const
const G4RotationMatrix * GetRotation() const

◆ 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}

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

◆ GetObjectTranslation()

◆ GetParameterisation()

◆ GetRegularStructureId()

◆ GetReplicationData()

◆ GetRotation() [1/2]

G4RotationMatrix * G4VPhysicalVolume::GetRotation ( )

Definition at line 165 of file G4VPhysicalVolume.cc.

166{
167 return G4MT_rot;
168}

◆ GetRotation() [2/2]

◆ GetSubInstanceManager()

const G4PVManager & G4VPhysicalVolume::GetSubInstanceManager ( )
static

Definition at line 140 of file G4VPhysicalVolume.cc.

141{
142 return subInstanceManager;
143}

Referenced by G4GeometryWorkspace::G4GeometryWorkspace().

◆ GetTranslation()

◆ 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()

Referenced by G4PVReplica::InitialiseWorker().

◆ IsMany()

virtual G4bool G4VPhysicalVolume::IsMany ( ) const
pure virtual

◆ IsParameterised()

◆ IsRegularStructure()

virtual G4bool G4VPhysicalVolume::IsRegularStructure ( ) const
pure virtual

◆ IsReplicated()

◆ operator=()

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

◆ operator==()

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

◆ SetCopyNo()

◆ SetLogicalVolume()

void G4VPhysicalVolume::SetLogicalVolume ( G4LogicalVolume * pLogical)
inline

◆ SetMotherLogical()

◆ SetName()

void G4VPhysicalVolume::SetName ( const G4String & pName)

◆ SetRotation()

◆ 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}

Referenced by G4ParameterisationBoxX::ComputeTransformation(), G4ParameterisationBoxY::ComputeTransformation(), G4ParameterisationBoxZ::ComputeTransformation(), G4ParameterisationConsPhi::ComputeTransformation(), G4ParameterisationConsRho::ComputeTransformation(), G4ParameterisationConsZ::ComputeTransformation(), G4ParameterisationParaX::ComputeTransformation(), G4ParameterisationParaY::ComputeTransformation(), G4ParameterisationParaZ::ComputeTransformation(), G4ParameterisationPolyconePhi::ComputeTransformation(), G4ParameterisationPolyconeRho::ComputeTransformation(), G4ParameterisationPolyconeZ::ComputeTransformation(), G4ParameterisationPolyhedraPhi::ComputeTransformation(), G4ParameterisationPolyhedraRho::ComputeTransformation(), G4ParameterisationPolyhedraZ::ComputeTransformation(), G4ParameterisationTrdX::ComputeTransformation(), G4ParameterisationTrdY::ComputeTransformation(), G4ParameterisationTrdZ::ComputeTransformation(), G4ParameterisationTubsPhi::ComputeTransformation(), G4ParameterisationTubsRho::ComputeTransformation(), G4ParameterisationTubsZ::ComputeTransformation(), G4PartialPhantomParameterisation::ComputeTransformation(), G4PhantomParameterisation::ComputeTransformation(), G4ReplicaNavigation::ComputeTransformation(), G4ReplicaNavigation::ComputeTransformation(), G4tgbPlaceParamCircle::ComputeTransformation(), G4tgbPlaceParamLinear::ComputeTransformation(), G4tgbPlaceParamSquare::ComputeTransformation(), 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

Member Data Documentation

◆ instanceID

G4int G4VPhysicalVolume::instanceID
protected

Definition at line 234 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: