Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParameterisationTubsRho Class Reference

#include <G4ParameterisationTubs.hh>

+ Inheritance 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 *)
 
virtual void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const =0
 
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 void ComputeTransformation (const G4int, G4VPhysicalVolume *) const =0
 
virtual G4VSolidComputeSolid (const G4int, G4VPhysicalVolume *)
 
virtual G4MaterialComputeMaterial (const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=0)
 
virtual G4bool IsNested () const
 
virtual G4VVolumeMaterialScannerGetMaterialScanner ()
 
virtual void ComputeDimensions (G4Box &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Tubs &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Trd &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Trap &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Cons &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Sphere &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Orb &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Torus &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Para &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Polycone &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Polyhedra &, const G4int, const G4VPhysicalVolume *) const
 
virtual void ComputeDimensions (G4Hype &, const G4int, const G4VPhysicalVolume *) const
 

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)
 
virtual G4double GetMaxParameter () const =0
 
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 G4int verbose = 5
 

Detailed Description

Definition at line 74 of file G4ParameterisationTubs.hh.

Constructor & Destructor Documentation

◆ G4ParameterisationTubsRho()

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

Definition at line 71 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 {
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}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
Definition: G4Tubs.hh:77
G4double GetInnerRadius() const
G4double GetOuterRadius() const
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const

◆ ~G4ParameterisationTubsRho()

G4ParameterisationTubsRho::~G4ParameterisationTubsRho ( )

Definition at line 105 of file G4ParameterisationTubs.cc.

106{
107}

Member Function Documentation

◆ ComputeDimensions()

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

Reimplemented from G4VPVParameterisation.

Definition at line 152 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}
double G4double
Definition: G4Types.hh:64
G4double GetZHalfLength() const
void SetDeltaPhiAngle(G4double newDPhi)
void SetStartPhiAngle(G4double newSPhi, G4bool trig=true)
G4double GetStartPhiAngle() const
void SetInnerRadius(G4double newRMin)
void SetOuterRadius(G4double newRMax)
G4double GetDeltaPhiAngle() const
void SetZHalfLength(G4double newDz)
void DumpInfo() const

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 119 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}
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.) const
void SetTranslation(const G4ThreeVector &v)

◆ GetMaxParameter()

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}

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