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

#include <G4ParameterisationTubs.hh>

+ Inheritance diagram for G4ParameterisationTubsPhi:

Public Member Functions

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

Constructor & Destructor Documentation

◆ G4ParameterisationTubsPhi()

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

Definition at line 175 of file G4ParameterisationTubs.cc.

179 : G4VParameterisationTubs( axis, nDiv, width, offset, msolid, divType )
180{
182 SetType( "DivisionTubsPhi" );
183
184 auto msol = (G4Tubs*)(fmotherSolid);
185 if( divType == DivWIDTH )
186 {
187 fnDiv = CalculateNDiv( msol->GetDeltaPhiAngle(), width, offset );
188 }
189 else if( divType == DivNDIV )
190 {
191 fwidth = CalculateWidth( msol->GetDeltaPhiAngle(), nDiv, offset );
192 }
193
194#ifdef G4DIVDEBUG
195 if( verbose >= 1 )
196 {
197 G4cout << " G4ParameterisationTubsPhi no divisions " << fnDiv << " = "
198 << nDiv << G4endl
199 << " Offset " << foffset << " = " << offset << G4endl
200 << " Width " << fwidth << " = " << width << G4endl;
201 }
202#endif
203}
#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)

◆ ~G4ParameterisationTubsPhi()

G4ParameterisationTubsPhi::~G4ParameterisationTubsPhi ( )
overridedefault

Member Function Documentation

◆ ComputeDimensions()

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

Reimplemented from G4VPVParameterisation.

Definition at line 252 of file G4ParameterisationTubs.cc.

255{
256 auto msol = (G4Tubs*)(fmotherSolid);
257
258 G4double pRMin = msol->GetInnerRadius();
259 G4double pRMax = msol->GetOuterRadius();
260 G4double pDz = msol->GetZHalfLength();
261 //----- already rotated in 'ComputeTransformation'
262 G4double pSPhi = msol->GetStartPhiAngle() + fhgap;
263 G4double pDPhi = fwidth - 2.*fhgap;
264
265 tubs.SetInnerRadius( pRMin );
266 tubs.SetOuterRadius( pRMax );
267 tubs.SetZHalfLength( pDz );
268 tubs.SetStartPhiAngle( pSPhi, false );
269 tubs.SetDeltaPhiAngle( pDPhi );
270
271#ifdef G4DIVDEBUG
272 if( verbose >= 2 )
273 {
274 G4cout << " G4ParameterisationTubsPhi::ComputeDimensions" << G4endl
275 << " pSPhi: " << pSPhi << " - pDPhi: " << pDPhi << G4endl;
276 tubs.DumpInfo();
277 }
278#endif
279}
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 G4ParameterisationTubsPhi::ComputeTransformation ( const G4int copyNo,
G4VPhysicalVolume * physVol ) const
overridevirtual

Implements G4VDivisionParameterisation.

Definition at line 217 of file G4ParameterisationTubs.cc.

219{
220 //----- translation
221 G4ThreeVector origin(0.,0.,0.);
222 //----- set translation
223 physVol->SetTranslation( origin );
224
225 //----- calculate rotation matrix (so that all volumes point to the centre)
226 G4double posi = foffset + copyNo*fwidth;
227
228#ifdef G4DIVDEBUG
229 if( verbose >= 2 )
230 {
231 G4cout << " G4ParameterisationTubsPhi - position: " << posi/CLHEP::deg << G4endl
232 << " copyNo: " << copyNo << " - foffset: " << foffset/CLHEP::deg
233 << " - fwidth: " << fwidth/CLHEP::deg << G4endl;
234 }
235#endif
236
237 ChangeRotMatrix( physVol, -posi );
238
239#ifdef G4DIVDEBUG
240 if( verbose >= 2 )
241 {
242 G4cout << std::setprecision(8) << " G4ParameterisationTubsPhi " << copyNo
243 << G4endl
244 << " Position: " << origin << " - Width: " << fwidth
245 << " - Axis: " << faxis << G4endl;
246 }
247#endif
248}
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.0) const
void SetTranslation(const G4ThreeVector &v)

◆ GetMaxParameter()

G4double G4ParameterisationTubsPhi::GetMaxParameter ( ) const
overridevirtual

Implements G4VDivisionParameterisation.

Definition at line 209 of file G4ParameterisationTubs.cc.

210{
211 auto msol = (G4Tubs*)(fmotherSolid);
212 return msol->GetDeltaPhiAngle();
213}

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