Geant4 10.7.0
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 ()
 
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 ()
 
virtual ~G4VPVParameterisation ()
 
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 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 318 of file G4ParameterisationCons.cc.

322 : G4VParameterisationCons( axis, nDiv, width, offset, msolid, divType )
323{
325 SetType( "DivisionConsZ" );
326
327 G4Cons* msol = (G4Cons*)(fmotherSolid);
328 if( divType == DivWIDTH )
329 {
330 fnDiv = CalculateNDiv( 2*msol->GetZHalfLength(), width, offset );
331 }
332 else if( divType == DivNDIV )
333 {
334 fwidth = CalculateWidth( 2*msol->GetZHalfLength(), nDiv, offset );
335 }
336
337#ifdef G4DIVDEBUG
338 if( verbose >= 1 )
339 {
340 G4cout << " G4ParameterisationConsZ: # divisions " << fnDiv << " = "
341 << nDiv << G4endl
342 << " Offset " << foffset << " = " << offset << G4endl
343 << " Width " << fwidth << " = " << width << G4endl
344 << " - Axis: " << faxis << G4endl;
345 }
346#endif
347}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Definition: G4Cons.hh:78
G4double GetZHalfLength() const
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const

◆ ~G4ParameterisationConsZ()

G4ParameterisationConsZ::~G4ParameterisationConsZ ( )

Definition at line 350 of file G4ParameterisationCons.cc.

351{
352}

Member Function Documentation

◆ ComputeDimensions()

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

Reimplemented from G4VPVParameterisation.

Definition at line 391 of file G4ParameterisationCons.cc.

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

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 363 of file G4ParameterisationCons.cc.

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

◆ GetMaxParameter()

G4double G4ParameterisationConsZ::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 355 of file G4ParameterisationCons.cc.

356{
357 G4Cons* msol = (G4Cons*)(fmotherSolid);
358 return 2*msol->GetZHalfLength();
359}

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