Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ParameterisationTubsPhi Class Reference

#include <G4ParameterisationTubs.hh>

Inheritance diagram for G4ParameterisationTubsPhi:
Collaboration diagram for G4ParameterisationTubsPhi:

Public Member Functions

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

Detailed Description

Definition at line 120 of file G4ParameterisationTubs.hh.

Constructor & Destructor Documentation

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

Definition at line 183 of file G4ParameterisationTubs.cc.

186  : G4VParameterisationTubs( axis, nDiv, width, offset, msolid, divType )
187 {
189  SetType( "DivisionTubsPhi" );
190 
191  G4Tubs* msol = (G4Tubs*)(fmotherSolid);
192  if( divType == DivWIDTH )
193  {
194  fnDiv = CalculateNDiv( msol->GetDeltaPhiAngle(), width, offset );
195  }
196  else if( divType == DivNDIV )
197  {
198  fwidth = CalculateWidth( msol->GetDeltaPhiAngle(), nDiv, offset );
199  }
200 
201 #ifdef G4DIVDEBUG
202  if( verbose >= 1 )
203  {
204  G4cout << " G4ParameterisationTubsPhi no divisions " << fnDiv << " = "
205  << nDiv << G4endl
206  << " Offset " << foffset << " = " << offset << G4endl
207  << " Width " << fwidth << " = " << width << G4endl;
208  }
209 #endif
210 }
void SetType(const G4String &type)
Definition: G4Tubs.hh:85
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
G4double GetDeltaPhiAngle() const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

G4ParameterisationTubsPhi::~G4ParameterisationTubsPhi ( )

Definition at line 213 of file G4ParameterisationTubs.cc.

214 {
215 }

Member Function Documentation

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

Reimplemented from G4VPVParameterisation.

Definition at line 262 of file G4ParameterisationTubs.cc.

264 {
265  G4Tubs* msol = (G4Tubs*)(fmotherSolid);
266 
267  G4double pRMin = msol->GetInnerRadius();
268  G4double pRMax = msol->GetOuterRadius();
269  G4double pDz = msol->GetZHalfLength();
270  //----- already rotated in 'ComputeTransformation'
271  G4double pSPhi = msol->GetStartPhiAngle() + fhgap;
272  G4double pDPhi = fwidth - 2.*fhgap;
273 
274  tubs.SetInnerRadius( pRMin );
275  tubs.SetOuterRadius( pRMax );
276  tubs.SetZHalfLength( pDz );
277  tubs.SetStartPhiAngle( pSPhi, false );
278  tubs.SetDeltaPhiAngle( pDPhi );
279 
280 #ifdef G4DIVDEBUG
281  if( verbose >= 2 )
282  {
283  G4cout << " G4ParameterisationTubsPhi::ComputeDimensions" << G4endl
284  << " pSPhi: " << pSPhi << " - pDPhi: " << pDPhi << G4endl;
285  tubs.DumpInfo();
286  }
287 #endif
288 }
Definition: G4Tubs.hh:85
void SetStartPhiAngle(G4double newSPhi, G4bool trig=true)
void DumpInfo() const
void SetDeltaPhiAngle(G4double newDPhi)
G4GLOB_DLL std::ostream G4cout
G4double GetStartPhiAngle() const
G4double GetInnerRadius() const
void SetInnerRadius(G4double newRMin)
void SetOuterRadius(G4double newRMax)
G4double GetZHalfLength() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double GetOuterRadius() const
void SetZHalfLength(G4double newDz)

Here is the call graph for this function:

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

Implements G4VDivisionParameterisation.

Definition at line 227 of file G4ParameterisationTubs.cc.

228 {
229  //----- translation
230  G4ThreeVector origin(0.,0.,0.);
231  //----- set translation
232  physVol->SetTranslation( origin );
233 
234  //----- calculate rotation matrix (so that all volumes point to the centre)
235  G4double posi = foffset + copyNo*fwidth;
236 
237 #ifdef G4DIVDEBUG
238  if( verbose >= 2 )
239  {
240  G4cout << " G4ParameterisationTubsPhi - position: " << posi/deg << G4endl
241  << " copyNo: " << copyNo << " - foffset: " << foffset/deg
242  << " - fwidth: " << fwidth/deg << G4endl;
243  }
244 #endif
245 
246  ChangeRotMatrix( physVol, -posi );
247 
248 #ifdef G4DIVDEBUG
249  if( verbose >= 2 )
250  {
251  G4cout << std::setprecision(8) << " G4ParameterisationTubsPhi " << copyNo
252  << G4endl
253  << " Position: " << origin << " - Width: " << fwidth
254  << " - Axis: " << faxis << G4endl;
255  }
256 #endif
257 }
G4GLOB_DLL std::ostream G4cout
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.) const
void SetTranslation(const G4ThreeVector &v)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static constexpr double deg
Definition: G4SIunits.hh:152

Here is the call graph for this function:

G4double G4ParameterisationTubsPhi::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 218 of file G4ParameterisationTubs.cc.

219 {
220  G4Tubs* msol = (G4Tubs*)(fmotherSolid);
221  return msol->GetDeltaPhiAngle();
222 }
Definition: G4Tubs.hh:85
G4double GetDeltaPhiAngle() const

Here is the call graph for this function:


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