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

#include <G4ParameterisationTubs.hh>

+ Inheritance diagram for G4ParameterisationTubsRho:

Public Member Functions

 G4ParameterisationTubsRho (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
 
 ~G4ParameterisationTubsRho () 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 68 of file G4ParameterisationTubs.hh.

Constructor & Destructor Documentation

◆ G4ParameterisationTubsRho()

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

Definition at line 66 of file G4ParameterisationTubs.cc.

70 : G4VParameterisationTubs( axis, nDiv, width, offset, msolid, divType )
71{
73 SetType( "DivisionTubsRho" );
74
75 auto msol = (G4Tubs*)(fmotherSolid);
76 if( divType == DivWIDTH )
77 {
78 fnDiv = CalculateNDiv( msol->GetOuterRadius() - msol->GetInnerRadius(),
79 width, offset );
80 }
81 else if( divType == DivNDIV )
82 {
83 fwidth = CalculateWidth( msol->GetOuterRadius() - msol->GetInnerRadius(),
84 nDiv, offset );
85 }
86
87#ifdef G4DIVDEBUG
88 if( verbose >= 1 )
89 {
90 G4cout << " G4ParameterisationTubsRho - no divisions " << fnDiv << " = "
91 << nDiv << G4endl
92 << " Offset " << foffset << " = " << offset << G4endl
93 << " Width " << fwidth << " = " << width << G4endl
94 << " DivType " << divType << G4endl;
95 }
96#endif
97}
#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)

◆ ~G4ParameterisationTubsRho()

G4ParameterisationTubsRho::~G4ParameterisationTubsRho ( )
overridedefault

Member Function Documentation

◆ ComputeDimensions()

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

Reimplemented from G4VPVParameterisation.

Definition at line 145 of file G4ParameterisationTubs.cc.

148{
149 auto msol = (G4Tubs*)(fmotherSolid);
150
151 G4double pRMin = msol->GetInnerRadius() + foffset + fwidth*copyNo + fhgap;
152 G4double pRMax = msol->GetInnerRadius() + foffset + fwidth*(copyNo+1) - fhgap;
153 G4double pDz = msol->GetZHalfLength();
154 //- already rotated G4double pSR = foffset + copyNo*fwidth;
155 G4double pSPhi = msol->GetStartPhiAngle();
156 G4double pDPhi = msol->GetDeltaPhiAngle();;
157
158 tubs.SetInnerRadius( pRMin );
159 tubs.SetOuterRadius( pRMax );
160 tubs.SetZHalfLength( pDz );
161 tubs.SetStartPhiAngle( pSPhi, false );
162 tubs.SetDeltaPhiAngle( pDPhi );
163
164#ifdef G4DIVDEBUG
165 if( verbose >= 2 )
166 {
167 G4cout << " G4ParameterisationTubsRho::ComputeDimensions()" << G4endl
168 << " pRMin: " << pRMin << " - pRMax: " << pRMax << G4endl;
169 tubs.DumpInfo();
170 }
171#endif
172}
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 G4ParameterisationTubsRho::ComputeTransformation ( const G4int copyNo,
G4VPhysicalVolume * physVol ) const
overridevirtual

Implements G4VDivisionParameterisation.

Definition at line 112 of file G4ParameterisationTubs.cc.

114{
115 //----- translation
116 G4ThreeVector origin(0.,0.,0.);
117 //----- set translation
118 physVol->SetTranslation( origin );
119
120 //----- calculate rotation matrix: unit
121
122#ifdef G4DIVDEBUG
123 if( verbose >= 2 )
124 {
125 G4cout << " G4ParameterisationTubsRho " << G4endl
126 << " Offset: " << foffset/CLHEP::deg
127 << " - Width: " << fwidth/CLHEP::deg << G4endl;
128 }
129#endif
130
131 ChangeRotMatrix( physVol );
132
133#ifdef G4DIVDEBUG
134 if( verbose >= 2 )
135 {
136 G4cout << std::setprecision(8) << " G4ParameterisationTubsRho " << G4endl
137 << " Position: " << origin << " - Width: " << fwidth
138 << " - Axis " << faxis << G4endl;
139 }
140#endif
141}
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.0) const
void SetTranslation(const G4ThreeVector &v)

◆ GetMaxParameter()

G4double G4ParameterisationTubsRho::GetMaxParameter ( ) const
overridevirtual

Implements G4VDivisionParameterisation.

Definition at line 103 of file G4ParameterisationTubs.cc.

104{
105 auto msol = (G4Tubs*)(fmotherSolid);
106 return msol->GetOuterRadius() - msol->GetInnerRadius();
107}

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