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

#include <G4ParameterisationCons.hh>

+ Inheritance diagram for G4ParameterisationConsZ:

Public Member Functions

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

Constructor & Destructor Documentation

◆ G4ParameterisationConsZ()

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

Definition at line 312 of file G4ParameterisationCons.cc.

316 : G4VParameterisationCons( axis, nDiv, width, offset, msolid, divType )
317{
319 SetType( "DivisionConsZ" );
320
321 auto msol = (G4Cons*)(fmotherSolid);
322 if( divType == DivWIDTH )
323 {
324 fnDiv = CalculateNDiv( 2*msol->GetZHalfLength(), width, offset );
325 }
326 else if( divType == DivNDIV )
327 {
328 fwidth = CalculateWidth( 2*msol->GetZHalfLength(), nDiv, offset );
329 }
330
331#ifdef G4DIVDEBUG
332 if( verbose >= 1 )
333 {
334 G4cout << " G4ParameterisationConsZ: # divisions " << fnDiv << " = "
335 << nDiv << G4endl
336 << " Offset " << foffset << " = " << offset << G4endl
337 << " Width " << fwidth << " = " << width << G4endl
338 << " - Axis: " << faxis << G4endl;
339 }
340#endif
341}
#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)

◆ ~G4ParameterisationConsZ()

G4ParameterisationConsZ::~G4ParameterisationConsZ ( )
overridedefault

Member Function Documentation

◆ ComputeDimensions()

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

Reimplemented from G4VPVParameterisation.

Definition at line 383 of file G4ParameterisationCons.cc.

386{
387 auto msol = (G4Cons*)(fmotherSolid);
388
389 G4double mHalfLength = msol->GetZHalfLength() - fhgap;
390 G4double aRInner = (msol->GetInnerRadiusPlusZ()
391 - msol->GetInnerRadiusMinusZ()) / (2*mHalfLength);
392 G4double bRInner = (msol->GetInnerRadiusPlusZ()
393 + msol->GetInnerRadiusMinusZ()) / 2;
394 G4double aROuter = (msol->GetOuterRadiusPlusZ()
395 - msol->GetOuterRadiusMinusZ()) / (2*mHalfLength);
396 G4double bROuter = (msol->GetOuterRadiusPlusZ()
397 + msol->GetOuterRadiusMinusZ()) / 2;
398 G4double xMinusZ = -mHalfLength + OffsetZ() + fwidth*copyNo + fhgap;
399 G4double xPlusZ = -mHalfLength + OffsetZ() + fwidth*(copyNo+1) - fhgap;
400 cons.SetInnerRadiusMinusZ( aRInner * xMinusZ + bRInner );
401 cons.SetOuterRadiusMinusZ( aROuter * xMinusZ + bROuter );
402 cons.SetInnerRadiusPlusZ( aRInner * xPlusZ + bRInner );
403 cons.SetOuterRadiusPlusZ( aROuter * xPlusZ + bROuter );
404
405 G4double pDz = fwidth / 2. - fhgap;
406 G4double pSPhi = msol->GetStartPhiAngle();
407 G4double pDPhi = msol->GetDeltaPhiAngle();
408
409 cons.SetZHalfLength( pDz );
410 cons.SetStartPhiAngle( pSPhi, false );
411 cons.SetDeltaPhiAngle( pDPhi );
412
413#ifdef G4DIVDEBUG
414 if( verbose >= 2 )
415 {
416 G4cout << " G4ParameterisationConsZ::ComputeDimensions()" << G4endl
417 << " pDz: " << pDz << G4endl;
418 if( verbose >= 4 ) cons.DumpInfo();
419 }
420#endif
421
422}
double G4double
Definition G4Types.hh:83

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 355 of file G4ParameterisationCons.cc.

357{
358 //----- set translation: along Z axis
359 auto motherCons = (G4Cons*)(GetMotherSolid());
360 G4double posi = - motherCons->GetZHalfLength() + OffsetZ()
361 + fwidth/2 + copyNo*fwidth;
362 G4ThreeVector origin(0.,0.,posi);
363 physVol->SetTranslation( origin );
364
365 //----- calculate rotation matrix: unit
366
367#ifdef G4DIVDEBUG
368 if( verbose >= 2 )
369 {
370 G4cout << " G4ParameterisationConsZ::ComputeTransformation()" << G4endl
371 << " Origin: " << origin << " - copyNo: " << copyNo << G4endl
372 << " foffset: " << foffset << " - fwidth: " << fwidth
373 << G4endl;
374 }
375#endif
376
377 ChangeRotMatrix( physVol );
378}
G4VSolid * GetMotherSolid() const
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.0) const
void SetTranslation(const G4ThreeVector &v)

◆ GetMaxParameter()

G4double G4ParameterisationConsZ::GetMaxParameter ( ) const
overridevirtual

Implements G4VDivisionParameterisation.

Definition at line 347 of file G4ParameterisationCons.cc.

348{
349 auto msol = (G4Cons*)(fmotherSolid);
350 return 2*msol->GetZHalfLength();
351}

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