Geant4 11.1.1
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 ()
 
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 ()=default
 
virtual ~G4VPVParameterisation ()=default
 
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 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 179 of file G4ParameterisationTubs.cc.

183 : G4VParameterisationTubs( axis, nDiv, width, offset, msolid, divType )
184{
186 SetType( "DivisionTubsPhi" );
187
188 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
189 if( divType == DivWIDTH )
190 {
191 fnDiv = CalculateNDiv( msol->GetDeltaPhiAngle(), width, offset );
192 }
193 else if( divType == DivNDIV )
194 {
195 fwidth = CalculateWidth( msol->GetDeltaPhiAngle(), nDiv, offset );
196 }
197
198#ifdef G4DIVDEBUG
199 if( verbose >= 1 )
200 {
201 G4cout << " G4ParameterisationTubsPhi no divisions " << fnDiv << " = "
202 << nDiv << G4endl
203 << " Offset " << foffset << " = " << offset << G4endl
204 << " Width " << fwidth << " = " << width << G4endl;
205 }
206#endif
207}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Definition: G4Tubs.hh:75
G4double GetDeltaPhiAngle() const
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const

◆ ~G4ParameterisationTubsPhi()

G4ParameterisationTubsPhi::~G4ParameterisationTubsPhi ( )

Definition at line 210 of file G4ParameterisationTubs.cc.

211{
212}

Member Function Documentation

◆ ComputeDimensions()

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

Reimplemented from G4VPVParameterisation.

Definition at line 258 of file G4ParameterisationTubs.cc.

261{
262 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
263
264 G4double pRMin = msol->GetInnerRadius();
265 G4double pRMax = msol->GetOuterRadius();
266 G4double pDz = msol->GetZHalfLength();
267 //----- already rotated in 'ComputeTransformation'
268 G4double pSPhi = msol->GetStartPhiAngle() + fhgap;
269 G4double pDPhi = fwidth - 2.*fhgap;
270
271 tubs.SetInnerRadius( pRMin );
272 tubs.SetOuterRadius( pRMax );
273 tubs.SetZHalfLength( pDz );
274 tubs.SetStartPhiAngle( pSPhi, false );
275 tubs.SetDeltaPhiAngle( pDPhi );
276
277#ifdef G4DIVDEBUG
278 if( verbose >= 2 )
279 {
280 G4cout << " G4ParameterisationTubsPhi::ComputeDimensions" << G4endl
281 << " pSPhi: " << pSPhi << " - pDPhi: " << pDPhi << G4endl;
282 tubs.DumpInfo();
283 }
284#endif
285}
double G4double
Definition: G4Types.hh:83
G4double GetZHalfLength() const
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)
void SetZHalfLength(G4double newDz)
void DumpInfo() const

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 223 of file G4ParameterisationTubs.cc.

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

◆ GetMaxParameter()

G4double G4ParameterisationTubsPhi::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 215 of file G4ParameterisationTubs.cc.

216{
217 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
218 return msol->GetDeltaPhiAngle();
219}

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