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

#include <G4ParameterisationTubs.hh>

+ Inheritance diagram for G4ParameterisationTubsZ:

Public Member Functions

 G4ParameterisationTubsZ (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
 
 ~G4ParameterisationTubsZ ()
 
G4double GetMaxParameter () const
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const
 
void ComputeDimensions (G4Tubs &tubs, const G4int copyNo, const G4VPhysicalVolume *physVol) const
 
- Public Member Functions inherited from G4VParameterisationTubs
 G4VParameterisationTubs (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
virtual ~G4VParameterisationTubs ()
 
- 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 158 of file G4ParameterisationTubs.hh.

Constructor & Destructor Documentation

◆ G4ParameterisationTubsZ()

G4ParameterisationTubsZ::G4ParameterisationTubsZ ( EAxis  axis,
G4int  nCopies,
G4double  offset,
G4double  step,
G4VSolid motherSolid,
DivisionType  divType 
)

Definition at line 288 of file G4ParameterisationTubs.cc.

292 : G4VParameterisationTubs( axis, nDiv, width, offset, msolid, divType )
293{
295 SetType( "DivisionTubsZ" );
296
297 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
298 if( divType == DivWIDTH )
299 {
300 fnDiv = CalculateNDiv( 2*msol->GetZHalfLength(), width, offset );
301 }
302 else if( divType == DivNDIV )
303 {
304 fwidth = CalculateWidth( 2*msol->GetZHalfLength(), nDiv, offset );
305 }
306
307#ifdef G4DIVDEBUG
308 if( verbose >= 1 )
309 {
310 G4cout << " G4ParameterisationTubsZ: # divisions " << fnDiv << " = "
311 << nDiv << G4endl
312 << " Offset " << foffset << " = " << offset << G4endl
313 << " Width " << fwidth << " = " << width << G4endl;
314 }
315#endif
316}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Definition: G4Tubs.hh:75
G4double GetZHalfLength() const
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const

◆ ~G4ParameterisationTubsZ()

G4ParameterisationTubsZ::~G4ParameterisationTubsZ ( )

Definition at line 319 of file G4ParameterisationTubs.cc.

320{
321}

Member Function Documentation

◆ ComputeDimensions()

void G4ParameterisationTubsZ::ComputeDimensions ( G4Tubs tubs,
const G4int  copyNo,
const G4VPhysicalVolume physVol 
) const
virtual

Reimplemented from G4VPVParameterisation.

Definition at line 369 of file G4ParameterisationTubs.cc.

372{
373 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
374
375 G4double pRMin = msol->GetInnerRadius();
376 G4double pRMax = msol->GetOuterRadius();
377 // G4double pDz = msol->GetZHalfLength() / GetNoDiv();
378 G4double pDz = fwidth/2. - fhgap;
379 G4double pSPhi = msol->GetStartPhiAngle();
380 G4double pDPhi = msol->GetDeltaPhiAngle();
381
382 tubs.SetInnerRadius( pRMin );
383 tubs.SetOuterRadius( pRMax );
384 tubs.SetZHalfLength( pDz );
385 tubs.SetStartPhiAngle( pSPhi, false );
386 tubs.SetDeltaPhiAngle( pDPhi );
387
388#ifdef G4DIVDEBUG
389 if( verbose >= 2 )
390 {
391 G4cout << " G4ParameterisationTubsZ::ComputeDimensions()" << G4endl
392 << " pDz: " << pDz << G4endl;
393 tubs.DumpInfo();
394 }
395#endif
396}
double G4double
Definition: G4Types.hh:83
void SetDeltaPhiAngle(G4double newDPhi)
void SetStartPhiAngle(G4double newSPhi, G4bool trig=true)
G4double GetInnerRadius() const
G4double GetOuterRadius() const
G4double GetStartPhiAngle() const
void SetInnerRadius(G4double newRMin)
void SetOuterRadius(G4double newRMax)
G4double GetDeltaPhiAngle() const
void SetZHalfLength(G4double newDz)
void DumpInfo() const

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 332 of file G4ParameterisationTubs.cc.

334{
335 //----- set translation: along Z axis
336 G4Tubs* motherTubs = (G4Tubs*)(fmotherSolid);
337 G4double posi = - motherTubs->GetZHalfLength() + OffsetZ()
338 + fwidth/2 + copyNo*fwidth;
339 G4ThreeVector origin(0.,0.,posi);
340 physVol->SetTranslation( origin );
341
342 //----- calculate rotation matrix: unit
343
344#ifdef G4DIVDEBUG
345 if( verbose >= 2 )
346 {
347 G4cout << " G4ParameterisationTubsZ::ComputeTransformation()" << G4endl
348 << " Position: " << posi << " - copyNo: " << copyNo << G4endl
349 << " foffset " << foffset/deg << " - fwidth " << fwidth/deg
350 << G4endl;
351 }
352#endif
353
354 ChangeRotMatrix( physVol );
355
356#ifdef G4DIVDEBUG
357 if( verbose >= 2 )
358 {
359 G4cout << std::setprecision(8) << " G4ParameterisationTubsZ " << copyNo
360 << G4endl
361 << " Position: " << origin << " - Width: " << fwidth
362 << " - Axis: " << faxis << G4endl;
363 }
364#endif
365}
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.0) const
void SetTranslation(const G4ThreeVector &v)

◆ GetMaxParameter()

G4double G4ParameterisationTubsZ::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 324 of file G4ParameterisationTubs.cc.

325{
326 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
327 return 2*msol->GetZHalfLength();
328}

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