Geant4 10.7.0
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 ()
 
G4double GetMaxParameter () const
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const
 
void ComputeDimensions (G4Box &box, const G4int copyNo, const G4VPhysicalVolume *physVol) const
 
- Public Member Functions inherited from G4VParameterisationBox
 G4VParameterisationBox (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
virtual ~G4VParameterisationBox ()
 
- Public Member Functions inherited from G4VDivisionParameterisation
 G4VDivisionParameterisation (EAxis axis, G4int nDiv, G4double width, G4double offset, DivisionType divType, G4VSolid *motherSolid=nullptr)
 
virtual ~G4VDivisionParameterisation ()
 
virtual G4VSolidComputeSolid (const G4int, G4VPhysicalVolume *)
 
virtual void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const =0
 
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 ()
 
virtual ~G4VPVParameterisation ()
 
virtual void ComputeTransformation (const G4int, G4VPhysicalVolume *) const =0
 
virtual G4VSolidComputeSolid (const G4int, G4VPhysicalVolume *)
 
virtual G4MaterialComputeMaterial (const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=nullptr)
 
virtual G4bool IsNested () const
 
virtual G4VVolumeMaterialScannerGetMaterialScanner ()
 
virtual void ComputeDimensions (G4Box &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Tubs &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Trd &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Trap &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Cons &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Sphere &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Orb &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Ellipsoid &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Torus &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Para &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Polycone &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Polyhedra &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Hype &, const G4int, const G4VPhysicalVolume *) const
 

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)
 
virtual G4double GetMaxParameter () const =0
 
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 68 of file G4ParameterisationBox.cc.

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

◆ ~G4ParameterisationBoxX()

G4ParameterisationBoxX::~G4ParameterisationBoxX ( )

Definition at line 98 of file G4ParameterisationBox.cc.

99{
100}

Member Function Documentation

◆ ComputeDimensions()

void G4ParameterisationBoxX::ComputeDimensions ( G4Box box,
const G4int  copyNo,
const G4VPhysicalVolume physVol 
) const
virtual

Reimplemented from G4VPVParameterisation.

Definition at line 146 of file G4ParameterisationBox.cc.

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

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 111 of file G4ParameterisationBox.cc.

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

◆ GetMaxParameter()

G4double G4ParameterisationBoxX::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 103 of file G4ParameterisationBox.cc.

104{
105 G4Box* msol = (G4Box*)(fmotherSolid);
106 return 2*msol->GetXHalfLength();
107}

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