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

212 : G4VParameterisationCons( axis, nDiv, width, offset, msolid, divType )
213{
215 SetType( "DivisionConsPhi" );
216
217 auto msol = (G4Cons*)(fmotherSolid);
218 if( divType == DivWIDTH )
219 {
220 fnDiv = CalculateNDiv( msol->GetDeltaPhiAngle(), width, offset );
221 }
222 else if( divType == DivNDIV )
223 {
224 fwidth = CalculateWidth( msol->GetDeltaPhiAngle(), nDiv, offset );
225 }
226
227#ifdef G4DIVDEBUG
228 if( verbose >= 1 )
229 {
230 G4cout << " G4ParameterisationConsPhi no divisions " << fnDiv << " = "
231 << nDiv << G4endl
232 << " Offset " << foffset << " = " << offset << G4endl
233 << " Width " << fwidth << " = " << width << G4endl;
234 }
235#endif
236}
#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
G4VParameterisationCons(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)

◆ ~G4ParameterisationConsPhi()

G4ParameterisationConsPhi::~G4ParameterisationConsPhi ( )
overridedefault

Member Function Documentation

◆ ComputeDimensions()

void G4ParameterisationConsPhi::ComputeDimensions ( G4Cons & tubs,
const G4int copyNo,
const G4VPhysicalVolume * physVol ) const
overridevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 277 of file G4ParameterisationCons.cc.

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

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 250 of file G4ParameterisationCons.cc.

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

◆ GetMaxParameter()

G4double G4ParameterisationConsPhi::GetMaxParameter ( ) const
overridevirtual

Implements G4VDivisionParameterisation.

Definition at line 242 of file G4ParameterisationCons.cc.

243{
244 auto msol = (G4Cons*)(fmotherSolid);
245 return msol->GetDeltaPhiAngle();
246}

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