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

#include <G4ParameterisationTubs.hh>

Inheritance diagram for G4ParameterisationTubsZ:
G4VParameterisationTubs G4VDivisionParameterisation G4VPVParameterisation

Public Member Functions

 G4ParameterisationTubsZ (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
 
 ~G4ParameterisationTubsZ ()
 
G4double GetMaxParameter () const
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const
 
void ComputeDimensions (G4Tubs &tubs, const G4int copyNo, const G4VPhysicalVolume *physVol) const
 
- Public Member Functions inherited from G4VParameterisationTubs
 G4VParameterisationTubs (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
virtual ~G4VParameterisationTubs ()
 
- 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 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 165 of file G4ParameterisationTubs.hh.

Constructor & Destructor Documentation

G4ParameterisationTubsZ::G4ParameterisationTubsZ ( EAxis  axis,
G4int  nCopies,
G4double  offset,
G4double  step,
G4VSolid motherSolid,
DivisionType  divType 
)

Definition at line 292 of file G4ParameterisationTubs.cc.

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

295  : G4VParameterisationTubs( axis, nDiv, width, offset, msolid, divType )
296 {
298  SetType( "DivisionTubsZ" );
299 
300  G4Tubs* msol = (G4Tubs*)(fmotherSolid);
301  if( divType == DivWIDTH )
302  {
303  fnDiv = CalculateNDiv( 2*msol->GetZHalfLength(), width, offset );
304  }
305  else if( divType == DivNDIV )
306  {
307  fwidth = CalculateWidth( 2*msol->GetZHalfLength(), nDiv, offset );
308  }
309 
310 #ifdef G4DIVDEBUG
311  if( verbose >= 1 )
312  {
313  G4cout << " G4ParameterisationTubsZ: # divisions " << fnDiv << " = "
314  << nDiv << G4endl
315  << " Offset " << foffset << " = " << offset << G4endl
316  << " Width " << fwidth << " = " << width << G4endl;
317  }
318 #endif
319 }
void SetType(const G4String &type)
Definition: G4Tubs.hh:84
#define width
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4VParameterisationTubs(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
G4GLOB_DLL std::ostream G4cout
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const
G4double GetZHalfLength() const
#define G4endl
Definition: G4ios.hh:61
G4ParameterisationTubsZ::~G4ParameterisationTubsZ ( )

Definition at line 322 of file G4ParameterisationTubs.cc.

323 {
324 }

Member Function Documentation

void G4ParameterisationTubsZ::ComputeDimensions ( G4Tubs tubs,
const G4int  copyNo,
const G4VPhysicalVolume physVol 
) const
virtual

Reimplemented from G4VPVParameterisation.

Definition at line 373 of file G4ParameterisationTubs.cc.

References G4VSolid::DumpInfo(), G4VDivisionParameterisation::fhgap, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4Tubs::GetDeltaPhiAngle(), G4Tubs::GetInnerRadius(), G4Tubs::GetOuterRadius(), G4Tubs::GetStartPhiAngle(), G4Tubs::SetDeltaPhiAngle(), G4Tubs::SetInnerRadius(), G4Tubs::SetOuterRadius(), G4Tubs::SetStartPhiAngle(), G4Tubs::SetZHalfLength(), and G4VDivisionParameterisation::verbose.

375 {
376  G4Tubs* msol = (G4Tubs*)(fmotherSolid);
377 
378  G4double pRMin = msol->GetInnerRadius();
379  G4double pRMax = msol->GetOuterRadius();
380  // G4double pDz = msol->GetZHalfLength() / GetNoDiv();
381  G4double pDz = fwidth/2. - fhgap;
382  G4double pSPhi = msol->GetStartPhiAngle();
383  G4double pDPhi = msol->GetDeltaPhiAngle();
384 
385  tubs.SetInnerRadius( pRMin );
386  tubs.SetOuterRadius( pRMax );
387  tubs.SetZHalfLength( pDz );
388  tubs.SetStartPhiAngle( pSPhi, false );
389  tubs.SetDeltaPhiAngle( pDPhi );
390 
391 #ifdef G4DIVDEBUG
392  if( verbose >= 2 )
393  {
394  G4cout << " G4ParameterisationTubsZ::ComputeDimensions()" << G4endl
395  << " pDz: " << pDz << G4endl;
396  tubs.DumpInfo();
397  }
398 #endif
399 }
Definition: G4Tubs.hh:84
void SetStartPhiAngle(G4double newSPhi, G4bool trig=true)
void DumpInfo() const
void SetDeltaPhiAngle(G4double newDPhi)
G4GLOB_DLL std::ostream G4cout
G4double GetDeltaPhiAngle() const
G4double GetStartPhiAngle() const
G4double GetInnerRadius() const
void SetInnerRadius(G4double newRMin)
void SetOuterRadius(G4double newRMax)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double GetOuterRadius() const
void SetZHalfLength(G4double newDz)
void G4ParameterisationTubsZ::ComputeTransformation ( const G4int  copyNo,
G4VPhysicalVolume physVol 
) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 336 of file G4ParameterisationTubs.cc.

References G4VDivisionParameterisation::ChangeRotMatrix(), python.hepunit::deg, G4VDivisionParameterisation::faxis, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::foffset, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4Tubs::GetZHalfLength(), G4VDivisionParameterisation::OffsetZ(), G4VPhysicalVolume::SetTranslation(), and G4VDivisionParameterisation::verbose.

337 {
338  //----- set translation: along Z axis
339  G4Tubs* motherTubs = (G4Tubs*)(fmotherSolid);
340  G4double posi = - motherTubs->GetZHalfLength() + OffsetZ()
341  + fwidth/2 + copyNo*fwidth;
342  G4ThreeVector origin(0.,0.,posi);
343  physVol->SetTranslation( origin );
344 
345  //----- calculate rotation matrix: unit
346 
347 #ifdef G4DIVDEBUG
348  if( verbose >= 2 )
349  {
350  G4cout << " G4ParameterisationTubsZ::ComputeTransformation()" << G4endl
351  << " Position: " << posi << " - copyNo: " << copyNo << G4endl
352  << " foffset " << foffset/deg << " - fwidth " << fwidth/deg
353  << G4endl;
354  }
355 #endif
356 
357  ChangeRotMatrix( physVol );
358 
359 #ifdef G4DIVDEBUG
360  if( verbose >= 2 )
361  {
362  G4cout << std::setprecision(8) << " G4ParameterisationTubsZ " << copyNo
363  << G4endl
364  << " Position: " << origin << " - Width: " << fwidth
365  << " - Axis: " << faxis << G4endl;
366  }
367 #endif
368 }
Definition: G4Tubs.hh:84
G4GLOB_DLL std::ostream G4cout
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.) const
void SetTranslation(const G4ThreeVector &v)
G4double GetZHalfLength() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double G4ParameterisationTubsZ::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 327 of file G4ParameterisationTubs.cc.

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

328 {
329  G4Tubs* msol = (G4Tubs*)(fmotherSolid);
330  return 2*msol->GetZHalfLength();
331 }
Definition: G4Tubs.hh:84
G4double GetZHalfLength() const

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