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

#include <G4ParameterisationCons.hh>

+ Inheritance diagram for G4ParameterisationConsPhi:

Public Member Functions

 G4ParameterisationConsPhi (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
 
 ~G4ParameterisationConsPhi ()
 
G4double GetMaxParameter () const
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const
 
void ComputeDimensions (G4Cons &tubs, const G4int copyNo, const G4VPhysicalVolume *physVol) const
 
- Public Member Functions inherited from G4VParameterisationCons
 G4VParameterisationCons (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
virtual ~G4VParameterisationCons ()
 
- 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 114 of file G4ParameterisationCons.hh.

Constructor & Destructor Documentation

◆ G4ParameterisationConsPhi()

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

Definition at line 212 of file G4ParameterisationCons.cc.

216 : G4VParameterisationCons( axis, nDiv, width, offset, msolid, divType )
217{
219 SetType( "DivisionConsPhi" );
220
221 G4Cons* msol = (G4Cons*)(fmotherSolid);
222 if( divType == DivWIDTH )
223 {
224 fnDiv = CalculateNDiv( msol->GetDeltaPhiAngle(), width, offset );
225 }
226 else if( divType == DivNDIV )
227 {
228 fwidth = CalculateWidth( msol->GetDeltaPhiAngle(), nDiv, offset );
229 }
230
231#ifdef G4DIVDEBUG
232 if( verbose >= 1 )
233 {
234 G4cout << " G4ParameterisationConsPhi no divisions " << fnDiv << " = "
235 << nDiv << G4endl
236 << " Offset " << foffset << " = " << offset << G4endl
237 << " Width " << fwidth << " = " << width << G4endl;
238 }
239#endif
240}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Definition: G4Cons.hh:78
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

◆ ~G4ParameterisationConsPhi()

G4ParameterisationConsPhi::~G4ParameterisationConsPhi ( )

Definition at line 243 of file G4ParameterisationCons.cc.

244{
245}

Member Function Documentation

◆ ComputeDimensions()

void G4ParameterisationConsPhi::ComputeDimensions ( G4Cons tubs,
const G4int  copyNo,
const G4VPhysicalVolume physVol 
) const
virtual

Reimplemented from G4VPVParameterisation.

Definition at line 283 of file G4ParameterisationCons.cc.

286{
287 G4Cons* msol = (G4Cons*)(fmotherSolid);
288
289 G4double pRMin1 = msol->GetInnerRadiusMinusZ();
290 G4double pRMax1 = msol->GetOuterRadiusMinusZ();
291 G4double pRMin2 = msol->GetInnerRadiusPlusZ();
292 G4double pRMax2 = msol->GetOuterRadiusPlusZ();
293 G4double pDz = msol->GetZHalfLength();
294
295 //- already rotated double pSPhi = foffset + copyNo*fwidth;
296 G4double pSPhi = foffset + msol->GetStartPhiAngle() + fhgap;
297 G4double pDPhi = fwidth - 2.*fhgap;
298
299 cons.SetInnerRadiusMinusZ( pRMin1 );
300 cons.SetOuterRadiusMinusZ( pRMax1 );
301 cons.SetInnerRadiusPlusZ( pRMin2 );
302 cons.SetOuterRadiusPlusZ( pRMax2 );
303 cons.SetZHalfLength( pDz );
304 cons.SetStartPhiAngle( pSPhi, false );
305 cons.SetDeltaPhiAngle( pDPhi );
306
307#ifdef G4DIVDEBUG
308 if( verbose >= 2 )
309 {
310 G4cout << " G4ParameterisationConsPhi::ComputeDimensions" << G4endl
311 << " pSPhi: " << pSPhi << " - pDPhi: " << pDPhi << G4endl;
312 if( verbose >= 4 ) cons.DumpInfo();
313 }
314#endif
315}
double G4double
Definition: G4Types.hh:83
G4double GetOuterRadiusPlusZ() const
G4double GetStartPhiAngle() const
G4double GetInnerRadiusMinusZ() const
G4double GetInnerRadiusPlusZ() const
G4double GetOuterRadiusMinusZ() const
G4double GetZHalfLength() const

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 256 of file G4ParameterisationCons.cc.

258{
259 //----- translation
260 G4ThreeVector origin(0.,0.,0.);
261 //----- set translation
262 physVol->SetTranslation( origin );
263
264 //----- calculate rotation matrix (so that all volumes point to the centre)
265 G4double posi = foffset + copyNo*fwidth;
266
267#ifdef G4DIVDEBUG
268 if( verbose >= 2 )
269 {
270 G4cout << " G4ParameterisationConsPhi - position: " << posi/CLHEP::deg << G4endl
271 << " Origin: " << origin << " copyNo: " << copyNo
272 << " - foffset: " << foffset/CLHEP::deg
273 << " - fwidth: " << fwidth/CLHEP::deg << G4endl
274 << " - Axis: " << faxis << G4endl;
275 }
276#endif
277
278 ChangeRotMatrix( physVol, -posi );
279}
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.0) const
void SetTranslation(const G4ThreeVector &v)

◆ GetMaxParameter()

G4double G4ParameterisationConsPhi::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 248 of file G4ParameterisationCons.cc.

249{
250 G4Cons* msol = (G4Cons*)(fmotherSolid);
251 return msol->GetDeltaPhiAngle();
252}

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