Geant4 11.1.1
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=nullptr)
 
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 ()=default
 
virtual ~G4VPVParameterisation ()=default
 
virtual void ComputeTransformation (const G4int, G4VPhysicalVolume *) const =0
 
virtual G4VSolidComputeSolid (const G4int, G4VPhysicalVolume *)
 
virtual G4MaterialComputeMaterial (const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=nullptr)
 
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 (G4Ellipsoid &, 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.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 = 0
 
G4double fwidth = 0.0
 
G4double foffset = 0.0
 
DivisionType fDivisionType
 
G4VSolidfmotherSolid = nullptr
 
G4bool fReflectedSolid = false
 
G4bool fDeleteSolid = false
 
G4int theVoluFirstCopyNo = 1
 
G4double kCarTolerance
 
G4double fhgap = 0.0
 
- Static Protected Attributes inherited from G4VDivisionParameterisation
static G4ThreadLocal G4RotationMatrixfRot = nullptr
 
static const G4int verbose = 5
 

Detailed Description

Definition at line 68 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 68 of file G4ParameterisationTubs.cc.

72 : G4VParameterisationTubs( axis, nDiv, width, offset, msolid, divType )
73{
75 SetType( "DivisionTubsRho" );
76
77 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
78 if( divType == DivWIDTH )
79 {
81 width, offset );
82 }
83 else if( divType == DivNDIV )
84 {
86 nDiv, offset );
87 }
88
89#ifdef G4DIVDEBUG
90 if( verbose >= 1 )
91 {
92 G4cout << " G4ParameterisationTubsRho - no divisions " << fnDiv << " = "
93 << nDiv << G4endl
94 << " Offset " << foffset << " = " << offset << G4endl
95 << " Width " << fwidth << " = " << width << G4endl
96 << " DivType " << divType << G4endl;
97 }
98#endif
99}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Definition: G4Tubs.hh:75
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 102 of file G4ParameterisationTubs.cc.

103{
104}

Member Function Documentation

◆ ComputeDimensions()

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

Reimplemented from G4VPVParameterisation.

Definition at line 149 of file G4ParameterisationTubs.cc.

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

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

◆ GetMaxParameter()

G4double G4ParameterisationTubsRho::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 107 of file G4ParameterisationTubs.cc.

108{
109 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
110 return msol->GetOuterRadius() - msol->GetInnerRadius();
111}

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