Geant4 11.1.1
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 ()
 
G4double GetMaxParameter () const
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const
 
void ComputeDimensions (G4Trd &trd, const G4int copyNo, const G4VPhysicalVolume *pv) const
 
- Public Member Functions inherited from G4VParameterisationTrd
 G4VParameterisationTrd (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
virtual ~G4VParameterisationTrd ()
 
- 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 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 399 of file G4ParameterisationTrd.cc.

403 : G4VParameterisationTrd( axis, nDiv, width, offset, msolid, divType )
404{
406 SetType( "DivTrdZ" );
407
408 G4Trd* msol = (G4Trd*)(fmotherSolid);
409 if( divType == DivWIDTH )
410 {
412 width, offset );
413 }
414 else if( divType == DivNDIV )
415 {
417 nDiv, offset );
418 }
419
420#ifdef G4DIVDEBUG
421 if( verbose >= 1 )
422 {
423 G4cout << " G4ParameterisationTrdZ no divisions " << fnDiv
424 << " = " << nDiv << G4endl
425 << " Offset " << foffset << " = " << offset << G4endl
426 << " Width " << fwidth << " = " << width << G4endl;
427 }
428#endif
429}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Definition: G4Trd.hh:63
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

◆ ~G4ParameterisationTrdZ()

G4ParameterisationTrdZ::~G4ParameterisationTrdZ ( )

Definition at line 432 of file G4ParameterisationTrd.cc.

433{
434}

Member Function Documentation

◆ ComputeDimensions()

void G4ParameterisationTrdZ::ComputeDimensions ( G4Trd trd,
const G4int  copyNo,
const G4VPhysicalVolume pv 
) const
virtual

Reimplemented from G4VPVParameterisation.

Definition at line 482 of file G4ParameterisationTrd.cc.

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

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 445 of file G4ParameterisationTrd.cc.

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

◆ GetMaxParameter()

G4double G4ParameterisationTrdZ::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 437 of file G4ParameterisationTrd.cc.

438{
439 G4Trd* msol = (G4Trd*)(fmotherSolid);
440 return 2*msol->GetZHalfLength();
441}

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