Geant4 9.6.0
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=0)
 
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=0)
 
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 (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.) 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
 
- Protected Attributes inherited from G4VDivisionParameterisation
G4String ftype
 
EAxis faxis
 
G4int fnDiv
 
G4double fwidth
 
G4double foffset
 
DivisionType fDivisionType
 
G4VSolidfmotherSolid
 
G4bool fReflectedSolid
 
G4bool fDeleteSolid
 
G4int theVoluFirstCopyNo
 
G4double kCarTolerance
 
G4double fhgap
 
- Static Protected Attributes inherited from G4VDivisionParameterisation
static G4int verbose = 5
 

Detailed Description

Definition at line 179 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 418 of file G4ParameterisationTrd.cc.

422 : G4VParameterisationTrd( axis, nDiv, width, offset, msolid, divType )
423{
425 SetType( "DivTrdZ" );
426
427 G4Trd* msol = (G4Trd*)(fmotherSolid);
428 if( divType == DivWIDTH )
429 {
431 width, offset );
432 }
433 else if( divType == DivNDIV )
434 {
436 nDiv, offset );
437 }
438
439#ifdef G4DIVDEBUG
440 if( verbose >= 1 )
441 {
442 G4cout << " G4ParameterisationTrdZ no divisions " << fnDiv
443 << " = " << nDiv << G4endl
444 << " Offset " << foffset << " = " << offset << G4endl
445 << " Width " << fwidth << " = " << width << G4endl;
446 }
447#endif
448}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT 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 451 of file G4ParameterisationTrd.cc.

452{
453}

Member Function Documentation

◆ ComputeDimensions()

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

Reimplemented from G4VPVParameterisation.

Definition at line 500 of file G4ParameterisationTrd.cc.

503{
504 //---- The division along Z of a Trd will result a Trd
505 G4Trd* msol = (G4Trd*)(fmotherSolid);
506
507 G4double pDx1 = msol->GetXHalfLength1();
508 G4double DDx = (msol->GetXHalfLength2() - msol->GetXHalfLength1() );
509 G4double pDy1 = msol->GetYHalfLength1();
510 G4double DDy = (msol->GetYHalfLength2() - msol->GetYHalfLength1() );
511 G4double pDz = fwidth/2. - fhgap;
512 G4double zLength = 2*msol->GetZHalfLength();
513
514 trd.SetAllParameters( pDx1+DDx*(OffsetZ()+copyNo*fwidth+fhgap)/zLength,
515 pDx1+DDx*(OffsetZ()+(copyNo+1)*fwidth-fhgap)/zLength,
516 pDy1+DDy*(OffsetZ()+copyNo*fwidth+fhgap)/zLength,
517 pDy1+DDy*(OffsetZ()+(copyNo+1)*fwidth-fhgap)/zLength,
518 pDz );
519
520#ifdef G4DIVDEBUG
521 if( verbose >= 1 )
522 {
523 G4cout << " G4ParameterisationTrdZ::ComputeDimensions()"
524 << " - Mother TRD " << G4endl;
525 msol->DumpInfo();
526 G4cout << " - Parameterised TRD: "
527 << copyNo << G4endl;
528 trd.DumpInfo();
529 }
530#endif
531}
double G4double
Definition: G4Types.hh:64
void SetAllParameters(G4double pdx1, G4double pdx2, G4double pdy1, G4double pdy2, G4double pdz)
Definition: G4Trd.cc:169
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 464 of file G4ParameterisationTrd.cc.

466{
467 G4Trd* msol = (G4Trd*)(fmotherSolid );
468 G4double mdz = msol->GetZHalfLength();
469
470 //----- translation
471 G4ThreeVector origin(0.,0.,0.);
472 G4double posi = -mdz + OffsetZ() + (copyNo+0.5)*fwidth;
473 if( faxis == kZAxis )
474 {
475 origin.setZ( posi );
476 }
477 else
478 {
479 G4Exception("G4ParameterisationTrdZ::ComputeTransformation()",
480 "GeomDiv0002", FatalException,
481 "Only axes along Z are allowed !");
482 }
483
484#ifdef G4DIVDEBUG
485 if( verbose >= 1 )
486 {
487 G4cout << std::setprecision(8) << " G4ParameterisationTrdZ: "
488 << copyNo << G4endl
489 << " Position: " << origin << " - Offset: " << foffset
490 << " - Width: " << fwidth << " Axis " << faxis << G4endl;
491 }
492#endif
493
494 //----- set translation
495 physVol->SetTranslation( origin );
496}
@ FatalException
void SetTranslation(const G4ThreeVector &v)
@ kZAxis
Definition: geomdefs.hh:54
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

◆ GetMaxParameter()

G4double G4ParameterisationTrdZ::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 456 of file G4ParameterisationTrd.cc.

457{
458 G4Trd* msol = (G4Trd*)(fmotherSolid);
459 return 2*msol->GetZHalfLength();
460}

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