Geant4 9.6.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=0)
 
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=0)
 
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 (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.) 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
 
G4double fwidth
 
G4double foffset
 
DivisionType fDivisionType
 
G4VSolidfmotherSolid
 
G4bool fReflectedSolid
 
G4bool fDeleteSolid
 
G4int theVoluFirstCopyNo
 
G4double kCarTolerance
 
G4double fhgap
 
- Static Protected Attributes inherited from G4VDivisionParameterisation
static G4int verbose = 5
 

Detailed Description

Definition at line 160 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 291 of file G4ParameterisationTubs.cc.

295 : G4VParameterisationTubs( axis, nDiv, width, offset, msolid, divType )
296{
298 SetType( "DivisionTubsZ" );
299
300 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
301 if( divType == DivWIDTH )
302 {
303 fnDiv = CalculateNDiv( 2*msol->GetZHalfLength(), width, offset );
304 }
305 else if( divType == DivNDIV )
306 {
307 fwidth = CalculateWidth( 2*msol->GetZHalfLength(), nDiv, offset );
308 }
309
310#ifdef G4DIVDEBUG
311 if( verbose >= 1 )
312 {
313 G4cout << " G4ParameterisationTubsZ: # divisions " << fnDiv << " = "
314 << nDiv << G4endl
315 << " Offset " << foffset << " = " << offset << G4endl
316 << " Width " << fwidth << " = " << width << G4endl;
317 }
318#endif
319}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
Definition: G4Tubs.hh:77
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 322 of file G4ParameterisationTubs.cc.

323{
324}

Member Function Documentation

◆ ComputeDimensions()

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

Reimplemented from G4VPVParameterisation.

Definition at line 372 of file G4ParameterisationTubs.cc.

375{
376 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
377
378 G4double pRMin = msol->GetInnerRadius();
379 G4double pRMax = msol->GetOuterRadius();
380 // G4double pDz = msol->GetZHalfLength() / GetNoDiv();
381 G4double pDz = fwidth/2. - fhgap;
382 G4double pSPhi = msol->GetStartPhiAngle();
383 G4double pDPhi = msol->GetDeltaPhiAngle();
384
385 tubs.SetInnerRadius( pRMin );
386 tubs.SetOuterRadius( pRMax );
387 tubs.SetZHalfLength( pDz );
388 tubs.SetStartPhiAngle( pSPhi, false );
389 tubs.SetDeltaPhiAngle( pDPhi );
390
391#ifdef G4DIVDEBUG
392 if( verbose >= 2 )
393 {
394 G4cout << " G4ParameterisationTubsZ::ComputeDimensions()" << G4endl
395 << " pDz: " << pDz << G4endl;
396 tubs.DumpInfo();
397 }
398#endif
399}
double G4double
Definition: G4Types.hh:64
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 335 of file G4ParameterisationTubs.cc.

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

◆ GetMaxParameter()

G4double G4ParameterisationTubsZ::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 327 of file G4ParameterisationTubs.cc.

328{
329 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
330 return 2*msol->GetZHalfLength();
331}

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