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

#include <G4ParameterisationTrd.hh>

+ Inheritance diagram for G4ParameterisationTrdZ:

Public Member Functions

 G4ParameterisationTrdZ (EAxis axis, G4int nCopies, G4double width, G4double offset, G4VSolid *motherSolid, DivisionType divType)
 
 ~G4ParameterisationTrdZ () override
 
G4double GetMaxParameter () const override
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const override
 
void ComputeDimensions (G4Trd &trd, const G4int copyNo, const G4VPhysicalVolume *pv) const override
 
- Public Member Functions inherited from G4VParameterisationTrd
 G4VParameterisationTrd (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
 ~G4VParameterisationTrd () 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 G4VParameterisationTrd
G4bool bDivInTrap = false
 
- 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 169 of file G4ParameterisationTrd.hh.

Constructor & Destructor Documentation

◆ G4ParameterisationTrdZ()

G4ParameterisationTrdZ::G4ParameterisationTrdZ ( EAxis axis,
G4int nCopies,
G4double width,
G4double offset,
G4VSolid * motherSolid,
DivisionType divType )

Definition at line 394 of file G4ParameterisationTrd.cc.

398 : G4VParameterisationTrd( axis, nDiv, width, offset, msolid, divType )
399{
401 SetType( "DivTrdZ" );
402
403 auto msol = (G4Trd*)(fmotherSolid);
404 if( divType == DivWIDTH )
405 {
406 fnDiv = CalculateNDiv( 2*msol->GetZHalfLength(),
407 width, offset );
408 }
409 else if( divType == DivNDIV )
410 {
411 fwidth = CalculateWidth( 2*msol->GetZHalfLength(),
412 nDiv, offset );
413 }
414
415#ifdef G4DIVDEBUG
416 if( verbose >= 1 )
417 {
418 G4cout << " G4ParameterisationTrdZ no divisions " << fnDiv
419 << " = " << nDiv << G4endl
420 << " Offset " << foffset << " = " << offset << G4endl
421 << " Width " << fwidth << " = " << width << G4endl;
422 }
423#endif
424}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
Definition G4Trd.hh:63
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const
G4VParameterisationTrd(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)

◆ ~G4ParameterisationTrdZ()

G4ParameterisationTrdZ::~G4ParameterisationTrdZ ( )
overridedefault

Member Function Documentation

◆ ComputeDimensions()

void G4ParameterisationTrdZ::ComputeDimensions ( G4Trd & trd,
const G4int copyNo,
const G4VPhysicalVolume * pv ) const
overridevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 475 of file G4ParameterisationTrd.cc.

478{
479 //---- The division along Z of a Trd will result a Trd
480 auto msol = (G4Trd*)(fmotherSolid);
481
482 G4double pDx1 = msol->GetXHalfLength1();
483 G4double DDx = (msol->GetXHalfLength2() - msol->GetXHalfLength1() );
484 G4double pDy1 = msol->GetYHalfLength1();
485 G4double DDy = (msol->GetYHalfLength2() - msol->GetYHalfLength1() );
486 G4double pDz = fwidth/2. - fhgap;
487 G4double zLength = 2*msol->GetZHalfLength();
488
489 trd.SetAllParameters( pDx1+DDx*(OffsetZ()+copyNo*fwidth+fhgap)/zLength,
490 pDx1+DDx*(OffsetZ()+(copyNo+1)*fwidth-fhgap)/zLength,
491 pDy1+DDy*(OffsetZ()+copyNo*fwidth+fhgap)/zLength,
492 pDy1+DDy*(OffsetZ()+(copyNo+1)*fwidth-fhgap)/zLength,
493 pDz );
494
495#ifdef G4DIVDEBUG
496 if( verbose >= 1 )
497 {
498 G4cout << " G4ParameterisationTrdZ::ComputeDimensions()"
499 << " - Mother TRD " << G4endl;
500 msol->DumpInfo();
501 G4cout << " - Parameterised TRD: "
502 << copyNo << G4endl;
503 trd.DumpInfo();
504 }
505#endif
506}
double G4double
Definition G4Types.hh:83
void SetAllParameters(G4double pdx1, G4double pdx2, G4double pdy1, G4double pdy2, G4double pdz)
Definition G4Trd.cc:126
void DumpInfo() const

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 438 of file G4ParameterisationTrd.cc.

440{
441 auto msol = (G4Trd*)(fmotherSolid );
442 G4double mdz = msol->GetZHalfLength();
443
444 //----- translation
445 G4ThreeVector origin(0.,0.,0.);
446 G4double posi = -mdz + OffsetZ() + (copyNo+0.5)*fwidth;
447 if( faxis == kZAxis )
448 {
449 origin.setZ( posi );
450 }
451 else
452 {
453 std::ostringstream message;
454 message << "Only axes along Z are allowed ! Axis: " << faxis;
455 G4Exception("G4ParameterisationTrdZ::ComputeTransformation()",
456 "GeomDiv0002", FatalException, message);
457 }
458
459#ifdef G4DIVDEBUG
460 if( verbose >= 1 )
461 {
462 G4cout << std::setprecision(8) << " G4ParameterisationTrdZ: "
463 << copyNo << G4endl
464 << " Position: " << origin << " - Offset: " << foffset
465 << " - Width: " << fwidth << " Axis " << faxis << G4endl;
466 }
467#endif
468
469 //----- set translation
470 physVol->SetTranslation( origin );
471}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
void SetTranslation(const G4ThreeVector &v)
@ kZAxis
Definition geomdefs.hh:57

◆ GetMaxParameter()

G4double G4ParameterisationTrdZ::GetMaxParameter ( ) const
overridevirtual

Implements G4VDivisionParameterisation.

Definition at line 430 of file G4ParameterisationTrd.cc.

431{
432 auto msol = (G4Trd*)(fmotherSolid);
433 return 2*msol->GetZHalfLength();
434}

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