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

#include <G4ParameterisationTubs.hh>

Inheritance diagram for G4ParameterisationTubsRho:
Collaboration diagram for G4ParameterisationTubsRho:

Public Member Functions

 G4ParameterisationTubsRho (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
 
 ~G4ParameterisationTubsRho ()
 
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 75 of file G4ParameterisationTubs.hh.

Constructor & Destructor Documentation

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

Definition at line 72 of file G4ParameterisationTubs.cc.

75  : G4VParameterisationTubs( axis, nDiv, width, offset, msolid, divType )
76 {
78  SetType( "DivisionTubsRho" );
79 
80  G4Tubs* msol = (G4Tubs*)(fmotherSolid);
81  if( divType == DivWIDTH )
82  {
83  fnDiv = CalculateNDiv( msol->GetOuterRadius() - msol->GetInnerRadius(),
84  width, offset );
85  }
86  else if( divType == DivNDIV )
87  {
89  nDiv, offset );
90  }
91 
92 #ifdef G4DIVDEBUG
93  if( verbose >= 1 )
94  {
95  G4cout << " G4ParameterisationTubsRho - no divisions " << fnDiv << " = "
96  << nDiv << G4endl
97  << " Offset " << foffset << " = " << offset << G4endl
98  << " Width " << fwidth << " = " << width << G4endl
99  << " DivType " << divType << G4endl;
100  }
101 #endif
102 }
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 GetInnerRadius() const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const
#define G4endl
Definition: G4ios.hh:61
G4double GetOuterRadius() const

Here is the call graph for this function:

G4ParameterisationTubsRho::~G4ParameterisationTubsRho ( )

Definition at line 105 of file G4ParameterisationTubs.cc.

106 {
107 }

Member Function Documentation

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

Reimplemented from G4VPVParameterisation.

Definition at line 153 of file G4ParameterisationTubs.cc.

155 {
156  G4Tubs* msol = (G4Tubs*)(fmotherSolid);
157 
158  G4double pRMin = msol->GetInnerRadius() + foffset + fwidth * copyNo + fhgap;
159  G4double pRMax = msol->GetInnerRadius() + foffset + fwidth * (copyNo+1) - fhgap;
160  G4double pDz = msol->GetZHalfLength();
161  //- already rotated G4double pSR = foffset + copyNo*fwidth;
162  G4double pSPhi = msol->GetStartPhiAngle();
163  G4double pDPhi = msol->GetDeltaPhiAngle();;
164 
165  tubs.SetInnerRadius( pRMin );
166  tubs.SetOuterRadius( pRMax );
167  tubs.SetZHalfLength( pDz );
168  tubs.SetStartPhiAngle( pSPhi, false );
169  tubs.SetDeltaPhiAngle( pDPhi );
170 
171 #ifdef G4DIVDEBUG
172  if( verbose >= 2 )
173  {
174  G4cout << " G4ParameterisationTubsRho::ComputeDimensions()" << G4endl
175  << " pRMin: " << pRMin << " - pRMax: " << pRMax << G4endl;
176  tubs.DumpInfo();
177  }
178 #endif
179 }
Definition: G4Tubs.hh:85
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)
G4double GetZHalfLength() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void SetZHalfLength(G4double newDz)

Here is the call graph for this function:

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

Implements G4VDivisionParameterisation.

Definition at line 120 of file G4ParameterisationTubs.cc.

121 {
122  //----- translation
123  G4ThreeVector origin(0.,0.,0.);
124  //----- set translation
125  physVol->SetTranslation( origin );
126 
127  //----- calculate rotation matrix: unit
128 
129 #ifdef G4DIVDEBUG
130  if( verbose >= 2 )
131  {
132  G4cout << " G4ParameterisationTubsRho " << G4endl
133  << " Offset: " << foffset/deg
134  << " - Width: " << fwidth/deg << G4endl;
135  }
136 #endif
137 
138  ChangeRotMatrix( physVol );
139 
140 #ifdef G4DIVDEBUG
141  if( verbose >= 2 )
142  {
143  G4cout << std::setprecision(8) << " G4ParameterisationTubsRho " << G4endl
144  << " Position: " << origin << " - Width: " << fwidth
145  << " - Axis " << faxis << G4endl;
146  }
147 #endif
148 }
G4GLOB_DLL std::ostream G4cout
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.) const
void SetTranslation(const G4ThreeVector &v)
#define G4endl
Definition: G4ios.hh:61
static constexpr double deg
Definition: G4SIunits.hh:152

Here is the call graph for this function:

G4double G4ParameterisationTubsRho::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 110 of file G4ParameterisationTubs.cc.

111 {
112  G4Tubs* msol = (G4Tubs*)(fmotherSolid);
113  return msol->GetOuterRadius() - msol->GetInnerRadius();
114 }
Definition: G4Tubs.hh:85
G4double GetInnerRadius() const
G4double GetOuterRadius() const

Here is the call graph for this function:


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