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

#include <G4ReplicatedSlice.hh>

+ Inheritance diagram for G4ReplicatedSlice:

Public Member Functions

 G4ReplicatedSlice (const G4String &pName, G4LogicalVolume *pLogical, G4LogicalVolume *pMotherLogical, const EAxis pAxis, const G4int nReplicas, const G4double width, const G4double half_gap, const G4double offset)
 
 G4ReplicatedSlice (const G4String &pName, G4LogicalVolume *pLogical, G4LogicalVolume *pMotherLogical, const EAxis pAxis, const G4int nReplicas, const G4double half_gap, const G4double offset)
 
 G4ReplicatedSlice (const G4String &pName, G4LogicalVolume *pLogical, G4LogicalVolume *pMotherLogical, const EAxis pAxis, const G4double width, const G4double half_gap, const G4double offset)
 
 G4ReplicatedSlice (const G4String &pName, G4LogicalVolume *pLogical, G4VPhysicalVolume *pMotherPhysical, const EAxis pAxis, const G4int nReplicas, const G4double width, const G4double half_gap, const G4double offset)
 
 G4ReplicatedSlice (const G4String &pName, G4LogicalVolume *pLogical, G4VPhysicalVolume *pMotherPhysical, const EAxis pAxis, const G4int nReplicas, const G4double half_gap, const G4double offset)
 
 G4ReplicatedSlice (const G4String &pName, G4LogicalVolume *pLogical, G4VPhysicalVolume *pMotherPhysical, const EAxis pAxis, const G4double width, const G4double half_gap, const G4double offset)
 
virtual ~G4ReplicatedSlice ()
 
 G4ReplicatedSlice (const G4ReplicatedSlice &)=delete
 
G4ReplicatedSliceoperator= (const G4ReplicatedSlice &)=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
 
EVolume VolumeType () const final
 
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 66 of file G4ReplicatedSlice.hh.

Constructor & Destructor Documentation

◆ G4ReplicatedSlice() [1/7]

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

Definition at line 44 of file G4ReplicatedSlice.cc.

52 : G4PVReplica(pName, nDivs, pAxis, pLogical, pMotherLogical)
53{
54 CheckAndSetParameters(pAxis, nDivs, width, half_gap, offset,
55 DivNDIVandWIDTH, pMotherLogical, pLogical);
56}

◆ G4ReplicatedSlice() [2/7]

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

Definition at line 59 of file G4ReplicatedSlice.cc.

66 : G4PVReplica(pName, nDivs, pAxis, pLogical, pMotherLogical)
67{
68 CheckAndSetParameters(pAxis, nDivs, 0., half_gap, offset,
69 DivNDIV, pMotherLogical, pLogical);
70}

◆ G4ReplicatedSlice() [3/7]

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

Definition at line 73 of file G4ReplicatedSlice.cc.

80 : G4PVReplica(pName, 0, pAxis, pLogical, pMotherLogical)
81{
82 CheckAndSetParameters(pAxis, 0, width, half_gap, offset,
83 DivWIDTH, pMotherLogical, pLogical);
84}

◆ G4ReplicatedSlice() [4/7]

G4ReplicatedSlice::G4ReplicatedSlice ( const G4String pName,
G4LogicalVolume pLogical,
G4VPhysicalVolume pMotherPhysical,
const EAxis  pAxis,
const G4int  nReplicas,
const G4double  width,
const G4double  half_gap,
const G4double  offset 
)

Definition at line 87 of file G4ReplicatedSlice.cc.

95 : G4PVReplica(pName, nDivs, pAxis, pLogical,
96 pMotherPhysical ? pMotherPhysical->GetLogicalVolume() : nullptr)
97{
98 if (pMotherPhysical == nullptr)
99 {
100 std::ostringstream message;
101 message << "Invalid setup." << G4endl
102 << "NULL pointer specified as mother for volume: " << pName;
103 G4Exception("G4ReplicatedSlice::G4ReplicatedSlice()", "GeomDiv0002",
104 FatalException, message);
105 return;
106 }
107 CheckAndSetParameters(pAxis, nDivs, width, half_gap, offset,
108 DivNDIVandWIDTH, pMotherPhysical->GetLogicalVolume(), pLogical);
109}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
#define G4endl
Definition: G4ios.hh:57
G4LogicalVolume * GetLogicalVolume() const

◆ G4ReplicatedSlice() [5/7]

G4ReplicatedSlice::G4ReplicatedSlice ( const G4String pName,
G4LogicalVolume pLogical,
G4VPhysicalVolume pMotherPhysical,
const EAxis  pAxis,
const G4int  nReplicas,
const G4double  half_gap,
const G4double  offset 
)

Definition at line 112 of file G4ReplicatedSlice.cc.

119 : G4PVReplica(pName, nDivs, pAxis, pLogical,
120 pMotherPhysical ? pMotherPhysical->GetLogicalVolume() : nullptr)
121{
122 if (pMotherPhysical == nullptr)
123 {
124 std::ostringstream message;
125 message << "Invalid setup." << G4endl
126 << "NULL pointer specified as mother for volume: " << pName;
127 G4Exception("G4ReplicatedSlice::G4ReplicatedSlice()", "GeomDiv0002",
128 FatalException, message);
129 return;
130 }
131 CheckAndSetParameters(pAxis, nDivs, 0., half_gap, offset,
132 DivNDIV, pMotherPhysical->GetLogicalVolume(), pLogical);
133}

◆ G4ReplicatedSlice() [6/7]

G4ReplicatedSlice::G4ReplicatedSlice ( const G4String pName,
G4LogicalVolume pLogical,
G4VPhysicalVolume pMotherPhysical,
const EAxis  pAxis,
const G4double  width,
const G4double  half_gap,
const G4double  offset 
)

Definition at line 136 of file G4ReplicatedSlice.cc.

143 : G4PVReplica(pName, 0, pAxis, pLogical,
144 pMotherPhysical ? pMotherPhysical->GetLogicalVolume() : nullptr)
145{
146 if (pMotherPhysical == nullptr)
147 {
148 std::ostringstream message;
149 message << "Invalid setup." << G4endl
150 << "NULL pointer specified as mother for volume: " << pName;
151 G4Exception("G4ReplicatedSlice::G4ReplicatedSlice()", "GeomDiv0002",
152 FatalException, message);
153 return;
154 }
155 CheckAndSetParameters(pAxis, 0, width, half_gap, offset,
156 DivWIDTH, pMotherPhysical->GetLogicalVolume(), pLogical);
157}

◆ ~G4ReplicatedSlice()

G4ReplicatedSlice::~G4ReplicatedSlice ( )
virtual

Definition at line 279 of file G4ReplicatedSlice.cc.

280{
281 delete GetRotation();
282}
const G4RotationMatrix * GetRotation() const

◆ G4ReplicatedSlice() [7/7]

G4ReplicatedSlice::G4ReplicatedSlice ( const G4ReplicatedSlice )
delete

Member Function Documentation

◆ GetDivisionAxis()

EAxis G4ReplicatedSlice::GetDivisionAxis ( ) const

Definition at line 285 of file G4ReplicatedSlice.cc.

286{
287 return fdivAxis;
288}

◆ GetMultiplicity()

G4int G4ReplicatedSlice::GetMultiplicity ( ) const
virtual

Reimplemented from G4PVReplica.

Definition at line 309 of file G4ReplicatedSlice.cc.

310{
311 return fnReplicas;
312}

◆ GetParameterisation()

G4VPVParameterisation * G4ReplicatedSlice::GetParameterisation ( ) const
virtual

Reimplemented from G4PVReplica.

Definition at line 315 of file G4ReplicatedSlice.cc.

316{
317 return fparam;
318}
G4VDivisionParameterisation * fparam

◆ GetRegularStructureId()

G4int G4ReplicatedSlice::GetRegularStructureId ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 574 of file G4ReplicatedSlice.cc.

575{
576 return 0;
577}

◆ GetReplicationData()

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

Reimplemented from G4PVReplica.

Definition at line 327 of file G4ReplicatedSlice.cc.

332{
333 axis = faxis;
334 nDivs = fnReplicas;
335 width = fwidth;
336 offset = foffset;
337 consuming = false;
338}

◆ IsMany()

G4bool G4ReplicatedSlice::IsMany ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 297 of file G4ReplicatedSlice.cc.

298{
299 return false;
300}

◆ IsParameterised()

G4bool G4ReplicatedSlice::IsParameterised ( ) const
virtual

Reimplemented from G4PVReplica.

Definition at line 291 of file G4ReplicatedSlice.cc.

292{
293 return true;
294}

◆ IsRegularStructure()

G4bool G4ReplicatedSlice::IsRegularStructure ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 566 of file G4ReplicatedSlice.cc.

567{
568 return false;
569}

◆ IsReplicated()

G4bool G4ReplicatedSlice::IsReplicated ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 303 of file G4ReplicatedSlice.cc.

304{
305 return true;
306}

◆ operator=()

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

◆ VolumeType()

EVolume G4ReplicatedSlice::VolumeType ( ) const
finalvirtual

Reimplemented from G4PVReplica.

Definition at line 321 of file G4ReplicatedSlice.cc.

322{
323 return kParameterised;
324}
@ kParameterised
Definition: geomdefs.hh:86

Member Data Documentation

◆ faxis

EAxis G4ReplicatedSlice::faxis
protected

Definition at line 174 of file G4ReplicatedSlice.hh.

Referenced by GetReplicationData().

◆ fdivAxis

EAxis G4ReplicatedSlice::fdivAxis
protected

Definition at line 175 of file G4ReplicatedSlice.hh.

Referenced by GetDivisionAxis().

◆ fnReplicas

G4int G4ReplicatedSlice::fnReplicas = 0
protected

Definition at line 176 of file G4ReplicatedSlice.hh.

Referenced by GetMultiplicity(), and GetReplicationData().

◆ foffset

G4double G4ReplicatedSlice::foffset = 0.0
protected

Definition at line 177 of file G4ReplicatedSlice.hh.

Referenced by GetReplicationData().

◆ fparam

G4VDivisionParameterisation* G4ReplicatedSlice::fparam = nullptr
protected

Definition at line 178 of file G4ReplicatedSlice.hh.

Referenced by GetParameterisation().

◆ fwidth

G4double G4ReplicatedSlice::fwidth = 0.0
protected

Definition at line 177 of file G4ReplicatedSlice.hh.

Referenced by GetReplicationData().


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