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

#include <G4ParameterisationBox.hh>

+ Inheritance diagram for G4ParameterisationBoxY:

Public Member Functions

 G4ParameterisationBoxY (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
 ~G4ParameterisationBoxY ()
 
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 113 of file G4ParameterisationBox.hh.

Constructor & Destructor Documentation

◆ G4ParameterisationBoxY()

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

Definition at line 171 of file G4ParameterisationBox.cc.

175 : G4VParameterisationBox( axis, nDiv, width, offset, msolid, divType )
176{
178 SetType( "DivisionBoxY" );
179
180 G4Box* mbox = (G4Box*)(fmotherSolid);
181 if( divType == DivWIDTH )
182 {
183 fnDiv = CalculateNDiv( 2*mbox->GetYHalfLength(), width, offset );
184 }
185 else if( divType == DivNDIV )
186 {
187 fwidth = CalculateWidth( 2*mbox->GetYHalfLength(), nDiv, offset );
188 }
189
190#ifdef G4DIVDEBUG
191 if( verbose >= 1 )
192 {
193 G4cout << " G4ParameterisationBoxY - no divisions " << fnDiv << " = "
194 << nDiv << ". Offset " << foffset << " = " << offset
195 << ". Width " << fwidth << " = " << width << G4endl;
196 }
197#endif
198}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Definition: G4Box.hh:56
G4double GetYHalfLength() const
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const

◆ ~G4ParameterisationBoxY()

G4ParameterisationBoxY::~G4ParameterisationBoxY ( )

Definition at line 201 of file G4ParameterisationBox.cc.

202{
203}

Member Function Documentation

◆ ComputeDimensions()

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

Reimplemented from G4VPVParameterisation.

Definition at line 248 of file G4ParameterisationBox.cc.

251{
252 G4Box* msol = (G4Box*)(fmotherSolid);
253
254 G4double pDx = msol->GetXHalfLength();
255 G4double pDy = fwidth/2. - fhgap;
256 G4double pDz = msol->GetZHalfLength();
257
258 box.SetXHalfLength( pDx );
259 box.SetYHalfLength( pDy );
260 box.SetZHalfLength( pDz );
261
262#ifdef G4DIVDEBUG
263 if( verbose >= 2 )
264 {
265 G4cout << " G4ParameterisationBoxY::ComputeDimensions()" << G4endl
266 << " pDx: " << pDz << G4endl;
267 box.DumpInfo();
268 }
269#endif
270}
double G4double
Definition: G4Types.hh:83
G4double GetZHalfLength() const
void SetZHalfLength(G4double dz)
Definition: G4Box.cc:172
G4double GetXHalfLength() const
void SetYHalfLength(G4double dy)
Definition: G4Box.cc:149
void SetXHalfLength(G4double dx)
Definition: G4Box.cc:125
void DumpInfo() const

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 214 of file G4ParameterisationBox.cc.

216{
217 G4Box* msol = (G4Box*)(fmotherSolid);
218 G4double mdy = msol->GetYHalfLength();
219
220 //----- translation
221 G4ThreeVector origin(0.,0.,0.);
222 G4double posi = -mdy + foffset + (copyNo+0.5)*fwidth;
223 if( faxis == kYAxis )
224 {
225 origin.setY( posi );
226 }
227 else
228 {
229 std::ostringstream message;
230 message << "Only axes along Y are allowed ! Axis: " << faxis;
231 G4Exception("G4ParameterisationBoxY::ComputeTransformation()",
232 "GeomDiv0002", FatalException, message);
233 }
234#ifdef G4DIVDEBUG
235 if( verbose >= 2 )
236 {
237 G4cout << std::setprecision(8) << " G4ParameterisationBoxY: "
238 << copyNo << G4endl
239 << " Position " << origin << " Axis " << faxis << G4endl;
240 }
241#endif
242 //----- set translation
243 physVol->SetTranslation( origin );
244}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
void SetTranslation(const G4ThreeVector &v)
@ kYAxis
Definition: geomdefs.hh:56

◆ GetMaxParameter()

G4double G4ParameterisationBoxY::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 206 of file G4ParameterisationBox.cc.

207{
208 G4Box* msol = (G4Box*)(fmotherSolid);
209 return 2*msol->GetYHalfLength();
210}

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