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

#include <G4ParameterisationBox.hh>

+ Inheritance diagram for G4ParameterisationBoxX:

Public Member Functions

 G4ParameterisationBoxX (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
 ~G4ParameterisationBoxX () override
 
G4double GetMaxParameter () const override
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const override
 
void ComputeDimensions (G4Box &box, const G4int copyNo, const G4VPhysicalVolume *physVol) const override
 
- Public Member Functions inherited from G4VParameterisationBox
 G4VParameterisationBox (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
 ~G4VParameterisationBox () override
 
- Public Member Functions inherited from G4VDivisionParameterisation
 G4VDivisionParameterisation (EAxis axis, G4int nDiv, G4double width, G4double offset, DivisionType divType, G4VSolid *motherSolid=nullptr)
 
 ~G4VDivisionParameterisation () override
 
G4VSolidComputeSolid (const G4int, G4VPhysicalVolume *) override
 
const G4StringGetType () const
 
EAxis GetAxis () const
 
G4int GetNoDiv () const
 
G4double GetWidth () const
 
G4double GetOffset () const
 
G4VSolidGetMotherSolid () const
 
void SetType (const G4String &type)
 
G4int VolumeFirstCopyNo () const
 
void SetHalfGap (G4double hg)
 
G4double GetHalfGap () const
 
- Public Member Functions inherited from G4VPVParameterisation
 G4VPVParameterisation ()=default
 
virtual ~G4VPVParameterisation ()=default
 
virtual G4MaterialComputeMaterial (const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=nullptr)
 
virtual G4bool IsNested () const
 
virtual G4VVolumeMaterialScannerGetMaterialScanner ()
 

Additional Inherited Members

- Protected Member Functions inherited from G4VDivisionParameterisation
void ChangeRotMatrix (G4VPhysicalVolume *physVol, G4double rotZ=0.0) const
 
G4int CalculateNDiv (G4double motherDim, G4double width, G4double offset) const
 
G4double CalculateWidth (G4double motherDim, G4int nDiv, G4double offset) const
 
virtual void CheckParametersValidity ()
 
void CheckOffset (G4double maxPar)
 
void CheckNDivAndWidth (G4double maxPar)
 
G4double OffsetZ () const
 
- Protected Attributes inherited from G4VDivisionParameterisation
G4String ftype
 
EAxis faxis
 
G4int fnDiv = 0
 
G4double fwidth = 0.0
 
G4double foffset = 0.0
 
DivisionType fDivisionType
 
G4VSolidfmotherSolid = nullptr
 
G4bool fReflectedSolid = false
 
G4bool fDeleteSolid = false
 
G4int theVoluFirstCopyNo = 1
 
G4double kCarTolerance
 
G4double fhgap = 0.0
 
- Static Protected Attributes inherited from G4VDivisionParameterisation
static G4ThreadLocal G4RotationMatrixfRot = nullptr
 
static const G4int verbose = 5
 

Detailed Description

Definition at line 69 of file G4ParameterisationBox.hh.

Constructor & Destructor Documentation

◆ G4ParameterisationBoxX()

G4ParameterisationBoxX::G4ParameterisationBoxX ( EAxis axis,
G4int nCopies,
G4double offset,
G4double step,
G4VSolid * msolid,
DivisionType divType )

Definition at line 66 of file G4ParameterisationBox.cc.

70 : G4VParameterisationBox( axis, nDiv, width, offset, msolid, divType )
71{
73 SetType( "DivisionBoxX" );
74
75 auto mbox = (G4Box*)(fmotherSolid);
76 if( divType == DivWIDTH )
77 {
78 fnDiv = CalculateNDiv( 2*mbox->GetXHalfLength(), width, offset );
79 }
80 else if( divType == DivNDIV )
81 {
82 fwidth = CalculateWidth( 2*mbox->GetXHalfLength(), nDiv, offset );
83 }
84#ifdef G4DIVDEBUG
85 if( verbose >= 1 )
86 {
87 G4cout << " G4ParameterisationBoxX - no divisions "
88 << fnDiv << " = " << nDiv << G4endl
89 << " Offset " << foffset << " = " << offset << G4endl
90 << " Width " << fwidth << " = " << width << G4endl;
91 }
92#endif
93}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
Definition G4Box.hh:56
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const
G4VParameterisationBox(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)

◆ ~G4ParameterisationBoxX()

G4ParameterisationBoxX::~G4ParameterisationBoxX ( )
overridedefault

Member Function Documentation

◆ ComputeDimensions()

void G4ParameterisationBoxX::ComputeDimensions ( G4Box & box,
const G4int copyNo,
const G4VPhysicalVolume * physVol ) const
overridevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 142 of file G4ParameterisationBox.cc.

145{
146 auto msol = (G4Box*)(fmotherSolid);
147
148 G4double pDx = fwidth/2. - fhgap;
149 G4double pDy = msol->GetYHalfLength();
150 G4double pDz = msol->GetZHalfLength();
151
152 box.SetXHalfLength( pDx );
153 box.SetYHalfLength( pDy );
154 box.SetZHalfLength( pDz );
155
156#ifdef G4DIVDEBUG
157 if( verbose >= 2 )
158 {
159 G4cout << " G4ParameterisationBoxX::ComputeDimensions()" << G4endl
160 << " pDx: " << pDz << G4endl;
161 box.DumpInfo();
162 }
163#endif
164}
double G4double
Definition G4Types.hh:83
void SetZHalfLength(G4double dz)
Definition G4Box.cc:167
void SetYHalfLength(G4double dy)
Definition G4Box.cc:144
void SetXHalfLength(G4double dx)
Definition G4Box.cc:120
void DumpInfo() const

◆ ComputeTransformation()

void G4ParameterisationBoxX::ComputeTransformation ( const G4int copyNo,
G4VPhysicalVolume * physVol ) const
overridevirtual

Implements G4VDivisionParameterisation.

Definition at line 107 of file G4ParameterisationBox.cc.

109{
110 auto msol = (G4Box*)(fmotherSolid );
111 G4double mdx = msol->GetXHalfLength( );
112
113 //----- translation
114 G4ThreeVector origin(0.,0.,0.);
115 G4double posi = -mdx + foffset+(copyNo+0.5)*fwidth;
116
117 if( faxis == kXAxis )
118 {
119 origin.setX( posi );
120 }
121 else
122 {
123 std::ostringstream message;
124 message << "Only axes along X are allowed ! Axis: " << faxis;
125 G4Exception("G4ParameterisationBoxX::ComputeTransformation()",
126 "GeomDiv0002", FatalException, message);
127 }
128#ifdef G4DIVDEBUG
129 if( verbose >= 2 )
130 {
131 G4cout << std::setprecision(8) << " G4ParameterisationBoxX: "
132 << copyNo << G4endl
133 << " Position " << origin << " Axis " << faxis << G4endl;
134 }
135#endif
136 //----- set translation
137 physVol->SetTranslation( origin );
138}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
void SetTranslation(const G4ThreeVector &v)
@ kXAxis
Definition geomdefs.hh:55

◆ GetMaxParameter()

G4double G4ParameterisationBoxX::GetMaxParameter ( ) const
overridevirtual

Implements G4VDivisionParameterisation.

Definition at line 99 of file G4ParameterisationBox.cc.

100{
101 auto msol = (G4Box*)(fmotherSolid);
102 return 2*msol->GetXHalfLength();
103}

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