Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ParameterisationTrdY Class Reference

#include <G4ParameterisationTrd.hh>

Inheritance diagram for G4ParameterisationTrdY:
Collaboration diagram for G4ParameterisationTrdY:

Public Member Functions

 G4ParameterisationTrdY (EAxis axis, G4int nCopies, G4double width, G4double offset, G4VSolid *motherSolid, DivisionType divType)
 
 ~G4ParameterisationTrdY ()
 
void CheckParametersValidity ()
 
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
 
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
G4RotationMatrix
fRot = 0
 
static const G4int verbose = 5
 

Detailed Description

Definition at line 136 of file G4ParameterisationTrd.hh.

Constructor & Destructor Documentation

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

Definition at line 291 of file G4ParameterisationTrd.cc.

294  : G4VParameterisationTrd( axis, nDiv, width, offset, msolid, divType )
295 {
297  SetType( "DivisionTrdY" );
298 
299  G4Trd* msol = (G4Trd*)(fmotherSolid);
300  if( divType == DivWIDTH )
301  {
302  fnDiv = CalculateNDiv( 2*msol->GetYHalfLength1(),
303  width, offset );
304  }
305  else if( divType == DivNDIV )
306  {
308  nDiv, offset );
309  }
310 
311 #ifdef G4DIVDEBUG
312  if( verbose >= 1 )
313  {
314  G4cout << " G4ParameterisationTrdY no divisions " << fnDiv
315  << " = " << nDiv << G4endl
316  << " Offset " << foffset << " = " << offset << G4endl
317  << " width " << fwidth << " = " << width << G4endl;
318  }
319 #endif
320 }
G4double GetYHalfLength1() const
void SetType(const G4String &type)
#define width
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
Definition: G4Trd.hh:72
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

Here is the call graph for this function:

G4ParameterisationTrdY::~G4ParameterisationTrdY ( )

Definition at line 323 of file G4ParameterisationTrd.cc.

324 {
325 }

Member Function Documentation

void G4ParameterisationTrdY::CheckParametersValidity ( )
virtual

Reimplemented from G4VDivisionParameterisation.

Definition at line 397 of file G4ParameterisationTrd.cc.

398 {
400 
401  G4Trd* msol = (G4Trd*)(fmotherSolid);
402 
403  G4double mpDy1 = msol->GetYHalfLength1();
404  G4double mpDy2 = msol->GetYHalfLength2();
405 
406  if( std::fabs(mpDy1 - mpDy2) > kCarTolerance )
407  {
408  std::ostringstream message;
409  message << "Invalid solid specification. NOT supported." << G4endl
410  << "Making a division of a TRD along axis Y while" << G4endl
411  << "the Y half lengths are not equal is not (yet)" << G4endl
412  << "supported. It will result in non-equal" << G4endl
413  << "division solids.";
414  G4Exception("G4ParameterisationTrdY::CheckParametersValidity()",
415  "GeomDiv0001", FatalException, message);
416  }
417 }
G4double GetYHalfLength1() const
Definition: G4Trd.hh:72
G4double GetYHalfLength2() const
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

Here is the call graph for this function:

Here is the caller graph for this function:

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

Reimplemented from G4VPVParameterisation.

Definition at line 374 of file G4ParameterisationTrd.cc.

375 {
376  //---- The division along Y of a Trd will result a Trd, only
377  //--- if Y at -Z and +Z are equal, else use the G4Trap version
378  G4Trd* msol = (G4Trd*)(fmotherSolid);
379 
380  G4double pDx1 = msol->GetXHalfLength1();
381  G4double pDx2 = msol->GetXHalfLength2();
382  G4double pDz = msol->GetZHalfLength();
383  G4double pDy = fwidth/2. - fhgap;
384 
385  trd.SetAllParameters ( pDx1, pDx2, pDy, pDy, pDz );
386 
387 #ifdef G4DIVDEBUG
388  if( verbose >= 2 )
389  {
390  G4cout << " G4ParameterisationTrdY::ComputeDimensions():" << G4endl;
391  trd.DumpInfo();
392  }
393 #endif
394 }
void SetAllParameters(G4double pdx1, G4double pdx2, G4double pdy1, G4double pdy2, G4double pdz)
Definition: G4Trd.cc:173
Definition: G4Trd.hh:72
G4double GetZHalfLength() const
void DumpInfo() const
G4double GetXHalfLength2() const
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
G4double GetXHalfLength1() const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Implements G4VDivisionParameterisation.

Definition at line 337 of file G4ParameterisationTrd.cc.

338 {
339  G4Trd* msol = (G4Trd*)(fmotherSolid );
340  G4double mdy = msol->GetYHalfLength1();
341 
342  //----- translation
343  G4ThreeVector origin(0.,0.,0.);
344  G4double posi = -mdy + foffset + (copyNo+0.5)*fwidth;
345 
346  if( faxis == kYAxis )
347  {
348  origin.setY( posi );
349  }
350  else
351  {
352  std::ostringstream message;
353  message << "Only axes along Y are allowed ! Axis: " << faxis;
354  G4Exception("G4ParameterisationTrdY::ComputeTransformation()",
355  "GeomDiv0002", FatalException, message);
356  }
357 
358 #ifdef G4DIVDEBUG
359  if( verbose >= 2 )
360  {
361  G4cout << std::setprecision(8)
362  << " G4ParameterisationTrdY::ComputeTransformation " << copyNo
363  << " pos " << origin << " rot mat " << " axis " << faxis << G4endl;
364  }
365 #endif
366 
367  //----- set translation
368  physVol->SetTranslation( origin );
369 }
G4double GetYHalfLength1() const
Definition: G4Trd.hh:72
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

Here is the call graph for this function:

G4double G4ParameterisationTrdY::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 328 of file G4ParameterisationTrd.cc.

329 {
330  G4Trd* msol = (G4Trd*)(fmotherSolid);
331  return 2*msol->GetYHalfLength1();
332 }
G4double GetYHalfLength1() const
Definition: G4Trd.hh:72

Here is the call graph for this function:


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