Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PVDivision Class Reference

#include <G4PVDivision.hh>

+ Inheritance diagram for G4PVDivision:

Public Member Functions

 G4PVDivision (const G4String &pName, G4LogicalVolume *pLogical, G4LogicalVolume *pMother, const EAxis pAxis, const G4int nReplicas, const G4double width, const G4double offset)
 
 G4PVDivision (const G4String &pName, G4LogicalVolume *pLogical, G4LogicalVolume *pMotherLogical, const EAxis pAxis, const G4int nReplicas, const G4double offset)
 
 G4PVDivision (const G4String &pName, G4LogicalVolume *pLogical, G4LogicalVolume *pMotherLogical, const EAxis pAxis, const G4double width, const G4double offset)
 
 G4PVDivision (const G4String &pName, G4LogicalVolume *pLogical, G4VPhysicalVolume *pMother, const EAxis pAxis, const G4int nReplicas, const G4double width, const G4double offset)
 
virtual ~G4PVDivision ()
 
 G4PVDivision (const G4PVDivision &)=delete
 
G4PVDivisionoperator= (const G4PVDivision &)=delete
 
virtual G4bool IsMany () const
 
virtual G4bool IsReplicated () const
 
virtual G4int GetMultiplicity () const
 
virtual G4VPVParameterisationGetParameterisation () const
 
virtual void GetReplicationData (EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const
 
EAxis GetDivisionAxis () const
 
G4bool IsParameterised () const
 
virtual EVolume VolumeType () const
 
G4bool IsRegularStructure () const
 
G4int GetRegularStructureId () const
 
- Public Member Functions inherited from G4PVReplica
 G4PVReplica (const G4String &pName, G4LogicalVolume *pLogical, G4LogicalVolume *pMother, const EAxis pAxis, const G4int nReplicas, const G4double width, const G4double offset=0.)
 
 G4PVReplica (const G4String &pName, G4LogicalVolume *pLogical, G4VPhysicalVolume *pMother, const EAxis pAxis, const G4int nReplicas, const G4double width, const G4double offset=0.)
 
 G4PVReplica (__void__ &)
 
 G4PVReplica (const G4PVReplica &)=delete
 
G4PVReplicaoperator= (const G4PVReplica &)=delete
 
virtual ~G4PVReplica ()
 
virtual EVolume VolumeType () const
 
G4bool IsMany () const
 
G4bool IsReplicated () const
 
virtual G4int GetCopyNo () const
 
virtual void SetCopyNo (G4int CopyNo)
 
virtual G4bool IsParameterised () const
 
virtual G4VPVParameterisationGetParameterisation () const
 
virtual G4int GetMultiplicity () const
 
virtual void GetReplicationData (EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const
 
virtual void SetRegularStructureId (G4int code)
 
G4bool IsRegularStructure () const
 
G4int GetRegularStructureId () const
 
G4int GetInstanceID () const
 
void InitialiseWorker (G4PVReplica *pMasterObject)
 
void TerminateWorker (G4PVReplica *pMasterObject)
 
- Public Member Functions inherited from G4VPhysicalVolume
 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
 

Protected Attributes

EAxis faxis
 
EAxis fdivAxis
 
G4int fnReplicas = 0
 
G4double fwidth = 0.0
 
G4double foffset = 0.0
 
G4VDivisionParameterisationfparam = nullptr
 
- Protected Attributes inherited from G4PVReplica
EAxis faxis
 
G4int fnReplicas
 
G4double fwidth
 
G4double foffset
 
- Protected Attributes inherited from G4VPhysicalVolume
G4int instanceID
 

Additional Inherited Members

- Static Public Member Functions inherited from G4PVReplica
static const G4PVRManagerGetSubInstanceManager ()
 
- Static Public Member Functions inherited from G4VPhysicalVolume
static const G4PVManagerGetSubInstanceManager ()
 
static void Clean ()
 
- Protected Member Functions inherited from G4PVReplica
 G4PVReplica (const G4String &pName, G4int nReplicas, EAxis pAxis, G4LogicalVolume *pLogical, G4LogicalVolume *pMotherLogical)
 
- Protected Member Functions inherited from G4VPhysicalVolume
void InitialiseWorker (G4VPhysicalVolume *pMasterObject, G4RotationMatrix *pRot, const G4ThreeVector &tlate)
 
void TerminateWorker (G4VPhysicalVolume *pMasterObject)
 
- Static Protected Attributes inherited from G4VPhysicalVolume
static G4GEOM_DLL G4PVManager subInstanceManager
 

Detailed Description

Definition at line 74 of file G4PVDivision.hh.

Constructor & Destructor Documentation

◆ G4PVDivision() [1/5]

G4PVDivision::G4PVDivision ( const G4String pName,
G4LogicalVolume pLogical,
G4LogicalVolume pMother,
const EAxis  pAxis,
const G4int  nReplicas,
const G4double  width,
const G4double  offset 
)

Definition at line 44 of file G4PVDivision.cc.

51 : G4PVReplica(pName, nDivs, pAxis, pLogical, pMotherLogical)
52{
53 if (pMotherLogical == nullptr)
54 {
55 std::ostringstream message;
56 message << "Invalid setup." << G4endl
57 << "NULL pointer specified as mother for volume: " << pName;
58 G4Exception("G4PVDivision::G4PVDivision()", "GeomDiv0002",
59 FatalException, message);
60 return;
61 }
62 if (pLogical == pMotherLogical)
63 {
64 std::ostringstream message;
65 message << "Invalid setup." << G4endl
66 << "Cannot place a volume inside itself! Volume: " << pName;
67 G4Exception("G4PVDivision::G4PVDivision()", "GeomDiv0002",
68 FatalException, message);
69 }
70 pMotherLogical->AddDaughter(this);
71 SetMotherLogical(pMotherLogical);
72 SetParameterisation(pMotherLogical, pAxis, nDivs,
73 width, offset, DivNDIVandWIDTH);
74 CheckAndSetParameters (pAxis, nDivs, width, offset,
75 DivNDIVandWIDTH, pMotherLogical);
76}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
#define G4endl
Definition: G4ios.hh:57
void SetMotherLogical(G4LogicalVolume *pMother)

◆ G4PVDivision() [2/5]

G4PVDivision::G4PVDivision ( const G4String pName,
G4LogicalVolume pLogical,
G4LogicalVolume pMotherLogical,
const EAxis  pAxis,
const G4int  nReplicas,
const G4double  offset 
)

Definition at line 79 of file G4PVDivision.cc.

85 : G4PVReplica(pName, nDivs, pAxis, pLogical, pMotherLogical)
86{
87 if (pMotherLogical == nullptr)
88 {
89 std::ostringstream message;
90 message << "Invalid setup." << G4endl
91 << "NULL pointer specified as mother! Volume: " << pName;
92 G4Exception("G4PVDivision::G4PVDivision()", "GeomDiv0002",
93 FatalException, message);
94 return;
95 }
96 if (pLogical == pMotherLogical)
97 {
98 std::ostringstream message;
99 message << "Invalid setup." << G4endl
100 << "Cannot place a volume inside itself! Volume: " << pName;
101 G4Exception("G4PVDivision::G4PVDivision()", "GeomDiv0002",
102 FatalException, message);
103 }
104 pMotherLogical->AddDaughter(this);
105 SetMotherLogical(pMotherLogical);
106 SetParameterisation(pMotherLogical, pAxis, nDivs, 0., offset, DivNDIV);
107 CheckAndSetParameters (pAxis, nDivs, 0., offset, DivNDIV, pMotherLogical);
108}
void AddDaughter(G4VPhysicalVolume *p)

◆ G4PVDivision() [3/5]

G4PVDivision::G4PVDivision ( const G4String pName,
G4LogicalVolume pLogical,
G4LogicalVolume pMotherLogical,
const EAxis  pAxis,
const G4double  width,
const G4double  offset 
)

Definition at line 111 of file G4PVDivision.cc.

117 : G4PVReplica(pName, 0, pAxis, pLogical, pMotherLogical)
118{
119 if (pMotherLogical == nullptr)
120 {
121 std::ostringstream message;
122 message << "Invalid setup." << G4endl
123 << "NULL pointer specified as mother! Volume: " + pName;
124 G4Exception("G4PVDivision::G4PVDivision()", "GeomDiv0002",
125 FatalException, message);
126 return;
127 }
128 if (pLogical == pMotherLogical)
129 {
130 std::ostringstream message;
131 message << "Invalid setup." << G4endl
132 << "Cannot place a volume inside itself! Volume: "+ pName;
133 G4Exception("G4PVDivision::G4PVDivision()", "GeomDiv0002",
134 FatalException, message);
135 }
136 pMotherLogical->AddDaughter(this);
137 SetMotherLogical(pMotherLogical);
138 SetParameterisation(pMotherLogical, pAxis, 0, width, offset, DivWIDTH);
139 CheckAndSetParameters (pAxis, 0, width, offset, DivWIDTH, pMotherLogical);
140}

◆ G4PVDivision() [4/5]

G4PVDivision::G4PVDivision ( const G4String pName,
G4LogicalVolume pLogical,
G4VPhysicalVolume pMother,
const EAxis  pAxis,
const G4int  nReplicas,
const G4double  width,
const G4double  offset 
)

Definition at line 143 of file G4PVDivision.cc.

150 : G4PVReplica(pName, nDivs, pAxis, pLogical,
151 pMotherPhysical ? pMotherPhysical->GetLogicalVolume() : nullptr)
152{
153 if (pMotherPhysical == nullptr)
154 {
155 std::ostringstream message;
156 message << "Invalid setup." << G4endl
157 << "NULL pointer specified as mother for volume: " << pName;
158 G4Exception("G4PVDivision::G4PVDivision()", "GeomDiv0002",
159 FatalException, message);
160 return;
161 }
162 G4LogicalVolume* pMotherLogical = pMotherPhysical->GetLogicalVolume();
163 if (pLogical == pMotherLogical)
164 {
165 std::ostringstream message;
166 message << "Invalid setup." << G4endl
167 << "Cannot place a volume inside itself! Volume: " << pName;
168 G4Exception("G4PVDivision::G4PVDivision()", "GeomDiv0002",
169 FatalException, message);
170 }
171 pMotherLogical->AddDaughter(this);
172 SetMotherLogical(pMotherLogical);
173 SetParameterisation(pMotherLogical, pAxis, nDivs,
174 width, offset, DivNDIVandWIDTH);
175 CheckAndSetParameters (pAxis, nDivs, width, offset,
176 DivNDIVandWIDTH, pMotherLogical);
177}

◆ ~G4PVDivision()

G4PVDivision::~G4PVDivision ( )
virtual

Definition at line 272 of file G4PVDivision.cc.

273{
274}

◆ G4PVDivision() [5/5]

G4PVDivision::G4PVDivision ( const G4PVDivision )
delete

Member Function Documentation

◆ GetDivisionAxis()

EAxis G4PVDivision::GetDivisionAxis ( ) const

Definition at line 277 of file G4PVDivision.cc.

278{
279 return fdivAxis;
280}

Referenced by G4GDMLWriteStructure::DivisionvolWrite().

◆ GetMultiplicity()

G4int G4PVDivision::GetMultiplicity ( ) const
virtual

Reimplemented from G4PVReplica.

Definition at line 301 of file G4PVDivision.cc.

302{
303 return fnReplicas;
304}

◆ GetParameterisation()

G4VPVParameterisation * G4PVDivision::GetParameterisation ( ) const
virtual

Reimplemented from G4PVReplica.

Definition at line 307 of file G4PVDivision.cc.

308{
309 return fparam;
310}
G4VDivisionParameterisation * fparam

◆ GetRegularStructureId()

G4int G4PVDivision::GetRegularStructureId ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 568 of file G4PVDivision.cc.

569{
570 return 0;
571}

◆ GetReplicationData()

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

Reimplemented from G4PVReplica.

Definition at line 313 of file G4PVDivision.cc.

318{
319 axis = faxis;
320 nDivs = fnReplicas;
321 width = fwidth;
322 offset = foffset;
323 consuming = false;
324}
G4double foffset
G4double fwidth

Referenced by G4GDMLWriteStructure::DivisionvolWrite().

◆ IsMany()

G4bool G4PVDivision::IsMany ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 289 of file G4PVDivision.cc.

290{
291 return false;
292}

◆ IsParameterised()

G4bool G4PVDivision::IsParameterised ( ) const
virtual

Reimplemented from G4PVReplica.

Definition at line 283 of file G4PVDivision.cc.

284{
285 return true;
286}

◆ IsRegularStructure()

G4bool G4PVDivision::IsRegularStructure ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 560 of file G4PVDivision.cc.

561{
562 return false;
563}

◆ IsReplicated()

G4bool G4PVDivision::IsReplicated ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 295 of file G4PVDivision.cc.

296{
297 return true;
298}

◆ operator=()

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

◆ VolumeType()

EVolume G4PVDivision::VolumeType ( ) const
virtual

Reimplemented from G4PVReplica.

Definition at line 327 of file G4PVDivision.cc.

328{
329 return kParameterised;
330}
@ kParameterised
Definition: geomdefs.hh:86

Member Data Documentation

◆ faxis

EAxis G4PVDivision::faxis
protected

Definition at line 156 of file G4PVDivision.hh.

Referenced by GetReplicationData().

◆ fdivAxis

EAxis G4PVDivision::fdivAxis
protected

Definition at line 157 of file G4PVDivision.hh.

Referenced by GetDivisionAxis().

◆ fnReplicas

G4int G4PVDivision::fnReplicas = 0
protected

Definition at line 158 of file G4PVDivision.hh.

Referenced by GetMultiplicity(), and GetReplicationData().

◆ foffset

G4double G4PVDivision::foffset = 0.0
protected

Definition at line 159 of file G4PVDivision.hh.

Referenced by GetReplicationData().

◆ fparam

G4VDivisionParameterisation* G4PVDivision::fparam = nullptr
protected

Definition at line 160 of file G4PVDivision.hh.

Referenced by GetParameterisation().

◆ fwidth

G4double G4PVDivision::fwidth = 0.0
protected

Definition at line 159 of file G4PVDivision.hh.

Referenced by GetReplicationData().


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