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

286 : G4VParameterisationTubs( axis, nDiv, width, offset, msolid, divType )
287{
289 SetType( "DivisionTubsZ" );
290
291 auto msol = (G4Tubs*)(fmotherSolid);
292 if( divType == DivWIDTH )
293 {
294 fnDiv = CalculateNDiv( 2*msol->GetZHalfLength(), width, offset );
295 }
296 else if( divType == DivNDIV )
297 {
298 fwidth = CalculateWidth( 2*msol->GetZHalfLength(), nDiv, offset );
299 }
300
301#ifdef G4DIVDEBUG
302 if( verbose >= 1 )
303 {
304 G4cout << " G4ParameterisationTubsZ: # divisions " << fnDiv << " = "
305 << nDiv << G4endl
306 << " Offset " << foffset << " = " << offset << G4endl
307 << " Width " << fwidth << " = " << width << G4endl;
308 }
309#endif
310}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const
G4VParameterisationTubs(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)

◆ ~G4ParameterisationTubsZ()

G4ParameterisationTubsZ::~G4ParameterisationTubsZ ( )
overridedefault

Member Function Documentation

◆ ComputeDimensions()

void G4ParameterisationTubsZ::ComputeDimensions ( G4Tubs & tubs,
const G4int copyNo,
const G4VPhysicalVolume * physVol ) const
overridevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 361 of file G4ParameterisationTubs.cc.

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

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 324 of file G4ParameterisationTubs.cc.

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

◆ GetMaxParameter()

G4double G4ParameterisationTubsZ::GetMaxParameter ( ) const
overridevirtual

Implements G4VDivisionParameterisation.

Definition at line 316 of file G4ParameterisationTubs.cc.

317{
318 auto msol = (G4Tubs*)(fmotherSolid);
319 return 2*msol->GetZHalfLength();
320}

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