Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4ParameterisationTrdZ Class Reference

#include <G4ParameterisationTrd.hh>

Inheritance diagram for G4ParameterisationTrdZ:
G4VParameterisationTrd G4VDivisionParameterisation G4VPVParameterisation

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 *)
 
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 G4MaterialComputeMaterial (const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=0)
 
virtual G4bool IsNested () const
 
virtual G4VVolumeMaterialScannerGetMaterialScanner ()
 

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)
 
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 G4ThreadLocal G4int verbose = 5
 

Detailed Description

Definition at line 184 of file G4ParameterisationTrd.hh.

Constructor & Destructor Documentation

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

Definition at line 421 of file G4ParameterisationTrd.cc.

References G4VDivisionParameterisation::CalculateNDiv(), G4VDivisionParameterisation::CalculateWidth(), G4VDivisionParameterisation::CheckParametersValidity(), DivNDIV, DivWIDTH, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::fnDiv, G4VDivisionParameterisation::foffset, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4Trd::GetZHalfLength(), G4VDivisionParameterisation::SetType(), G4VDivisionParameterisation::verbose, and width.

424  : G4VParameterisationTrd( axis, nDiv, width, offset, msolid, divType )
425 {
427  SetType( "DivTrdZ" );
428 
429  G4Trd* msol = (G4Trd*)(fmotherSolid);
430  if( divType == DivWIDTH )
431  {
432  fnDiv = CalculateNDiv( 2*msol->GetZHalfLength(),
433  width, offset );
434  }
435  else if( divType == DivNDIV )
436  {
437  fwidth = CalculateWidth( 2*msol->GetZHalfLength(),
438  nDiv, offset );
439  }
440 
441 #ifdef G4DIVDEBUG
442  if( verbose >= 1 )
443  {
444  G4cout << " G4ParameterisationTrdZ no divisions " << fnDiv
445  << " = " << nDiv << G4endl
446  << " Offset " << foffset << " = " << offset << G4endl
447  << " Width " << fwidth << " = " << width << G4endl;
448  }
449 #endif
450 }
void SetType(const G4String &type)
#define width
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
Definition: G4Trd.hh:71
G4double GetZHalfLength() const
G4GLOB_DLL std::ostream G4cout
G4VParameterisationTrd(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const
#define G4endl
Definition: G4ios.hh:61
G4ParameterisationTrdZ::~G4ParameterisationTrdZ ( )

Definition at line 453 of file G4ParameterisationTrd.cc.

454 {
455 }

Member Function Documentation

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

Reimplemented from G4VPVParameterisation.

Definition at line 504 of file G4ParameterisationTrd.cc.

References G4VSolid::DumpInfo(), G4VDivisionParameterisation::fhgap, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4Trd::GetXHalfLength1(), G4Trd::GetXHalfLength2(), G4Trd::GetYHalfLength1(), G4Trd::GetYHalfLength2(), G4Trd::GetZHalfLength(), G4VDivisionParameterisation::OffsetZ(), G4Trd::SetAllParameters(), and G4VDivisionParameterisation::verbose.

506 {
507  //---- The division along Z of a Trd will result a Trd
508  G4Trd* msol = (G4Trd*)(fmotherSolid);
509 
510  G4double pDx1 = msol->GetXHalfLength1();
511  G4double DDx = (msol->GetXHalfLength2() - msol->GetXHalfLength1() );
512  G4double pDy1 = msol->GetYHalfLength1();
513  G4double DDy = (msol->GetYHalfLength2() - msol->GetYHalfLength1() );
514  G4double pDz = fwidth/2. - fhgap;
515  G4double zLength = 2*msol->GetZHalfLength();
516 
517  trd.SetAllParameters( pDx1+DDx*(OffsetZ()+copyNo*fwidth+fhgap)/zLength,
518  pDx1+DDx*(OffsetZ()+(copyNo+1)*fwidth-fhgap)/zLength,
519  pDy1+DDy*(OffsetZ()+copyNo*fwidth+fhgap)/zLength,
520  pDy1+DDy*(OffsetZ()+(copyNo+1)*fwidth-fhgap)/zLength,
521  pDz );
522 
523 #ifdef G4DIVDEBUG
524  if( verbose >= 1 )
525  {
526  G4cout << " G4ParameterisationTrdZ::ComputeDimensions()"
527  << " - Mother TRD " << G4endl;
528  msol->DumpInfo();
529  G4cout << " - Parameterised TRD: "
530  << copyNo << G4endl;
531  trd.DumpInfo();
532  }
533 #endif
534 }
void SetAllParameters(G4double pdx1, G4double pdx2, G4double pdy1, G4double pdy2, G4double pdz)
Definition: G4Trd.cc:169
G4double GetYHalfLength1() const
Definition: G4Trd.hh:71
G4double GetZHalfLength() const
void DumpInfo() const
G4double GetXHalfLength2() const
G4GLOB_DLL std::ostream G4cout
G4double GetYHalfLength2() const
#define G4endl
Definition: G4ios.hh:61
G4double GetXHalfLength1() const
double G4double
Definition: G4Types.hh:76
void G4ParameterisationTrdZ::ComputeTransformation ( const G4int  copyNo,
G4VPhysicalVolume physVol 
) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 467 of file G4ParameterisationTrd.cc.

References FatalException, G4VDivisionParameterisation::faxis, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::foffset, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4Exception(), G4Trd::GetZHalfLength(), kZAxis, G4VDivisionParameterisation::OffsetZ(), G4VPhysicalVolume::SetTranslation(), CLHEP::Hep3Vector::setZ(), and G4VDivisionParameterisation::verbose.

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

Implements G4VDivisionParameterisation.

Definition at line 458 of file G4ParameterisationTrd.cc.

References G4VDivisionParameterisation::fmotherSolid, and G4Trd::GetZHalfLength().

459 {
460  G4Trd* msol = (G4Trd*)(fmotherSolid);
461  return 2*msol->GetZHalfLength();
462 }
Definition: G4Trd.hh:71
G4double GetZHalfLength() const

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