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

#include <G4ParameterisationCons.hh>

+ Inheritance diagram for G4ParameterisationConsRho:

Public Member Functions

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

Constructor & Destructor Documentation

◆ G4ParameterisationConsRho()

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

Definition at line 81 of file G4ParameterisationCons.cc.

85 : G4VParameterisationCons( axis, nDiv, width, offset, msolid, divType )
86{
88 SetType( "DivisionConsRho" );
89
90 G4Cons* msol = (G4Cons*)(fmotherSolid);
91 if( msol->GetInnerRadiusPlusZ() == 0. )
92 {
93 std::ostringstream message;
94 message << "OuterRadiusMinusZ = 0" << G4endl
95 << "Width is calculated as that of OuterRadiusMinusZ !";
96 G4Exception("G4ParameterisationConsRho::G4ParameterisationConsRho()",
97 "GeomDiv1001", JustWarning, message);
98 }
99
100 if( divType == DivWIDTH )
101 {
103 - msol->GetInnerRadiusMinusZ(), width, offset );
104 }
105 else if( divType == DivNDIV )
106 {
107 G4Cons* mconsol = (G4Cons*)(msolid);
109 - mconsol->GetInnerRadiusMinusZ(), nDiv, offset );
110 }
111
112#ifdef G4DIVDEBUG
113 if( verbose >= 1 )
114 {
115 G4cout << " G4ParameterisationConsRho - no divisions " << fnDiv << " = "
116 << nDiv << G4endl
117 << " Offset " << foffset << " = " << offset
118 << " - Width " << fwidth << " = " << width << G4endl;
119 }
120#endif
121}
@ JustWarning
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
Definition: G4Cons.hh:75
G4double GetInnerRadiusMinusZ() const
G4double GetInnerRadiusPlusZ() const
G4double GetOuterRadiusMinusZ() const
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

◆ ~G4ParameterisationConsRho()

G4ParameterisationConsRho::~G4ParameterisationConsRho ( )

Definition at line 124 of file G4ParameterisationCons.cc.

125{
126}

Member Function Documentation

◆ ComputeDimensions()

void G4ParameterisationConsRho::ComputeDimensions ( G4Cons tubs,
const G4int  copyNo,
const G4VPhysicalVolume physVol 
) const
virtual

Reimplemented from G4VPVParameterisation.

Definition at line 170 of file G4ParameterisationCons.cc.

173{
174 G4Cons* msol = (G4Cons*)(fmotherSolid);
175
176 G4double pRMin1 = msol->GetInnerRadiusMinusZ() + foffset + fwidth * copyNo;
177 G4double pRMax1 = msol->GetInnerRadiusMinusZ() + foffset + fwidth * (copyNo+1);
178
179 //width at Z Plus
180 //- G4double fwidthPlus =
181 // fwidth * ( msol->GetOuterRadiusPlusZ()/ msol->GetInnerRadiusPlusZ())
182 //- / ( msol->GetOuterRadiusMinusZ() - msol->GetInnerRadiusMinusZ());
183 G4double fwidthPlus = CalculateWidth( msol->GetOuterRadiusPlusZ()
184 - msol->GetInnerRadiusPlusZ(), fnDiv, foffset );
185 G4double pRMin2 = msol->GetInnerRadiusPlusZ()
186 + foffset + fwidthPlus * copyNo;
187 G4double pRMax2 = msol->GetInnerRadiusPlusZ()
188 + foffset + fwidthPlus * (copyNo+1);
189 G4double pDz = msol->GetZHalfLength();
190
191 G4double d_half_gap = fhgap * pRMax2 / pRMax1;
192 //- already rotated double pSR = foffset + copyNo*fwidth;
193 G4double pSPhi = msol->GetStartPhiAngle();
194 G4double pDPhi = msol->GetDeltaPhiAngle();;
195
196 cons.SetInnerRadiusMinusZ( pRMin1 + fhgap );
197 cons.SetOuterRadiusMinusZ( pRMax1 - fhgap );
198 cons.SetInnerRadiusPlusZ( pRMin2 + d_half_gap );
199 cons.SetOuterRadiusPlusZ( pRMax2 - d_half_gap );
200 cons.SetZHalfLength( pDz );
201 cons.SetStartPhiAngle( pSPhi, false );
202 cons.SetDeltaPhiAngle( pDPhi );
203
204#ifdef G4DIVDEBUG
205 if( verbose >= 2 )
206 {
207 G4cout << " G4ParameterisationConsRho::ComputeDimensions()" << G4endl
208 << " pRMin: " << pRMin1 << " - pRMax: " << pRMax1 << G4endl;
209 if( verbose >= 4 ) cons.DumpInfo();
210 }
211#endif
212}
double G4double
Definition: G4Types.hh:64
G4double GetOuterRadiusPlusZ() const
G4double GetStartPhiAngle() const
G4double GetDeltaPhiAngle() const
G4double GetZHalfLength() const

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 137 of file G4ParameterisationCons.cc.

139{
140 //----- translation
141 G4ThreeVector origin(0.,0.,0.);
142 //----- set translation
143 physVol->SetTranslation( origin );
144
145 //----- calculate rotation matrix: unit
146
147#ifdef G4DIVDEBUG
148 if( verbose >= 2 )
149 {
150 G4cout << " G4ParameterisationConsRho " << G4endl
151 << " Offset: " << foffset
152 << " - Width: " << fwidth << G4endl;
153 }
154#endif
155
156 ChangeRotMatrix( physVol );
157
158#ifdef G4DIVDEBUG
159 if( verbose >= 2 )
160 {
161 G4cout << std::setprecision(8) << " G4ParameterisationConsRho" << G4endl
162 << " Position: " << origin << " - Width: " << fwidth
163 << " - Axis: " << faxis << G4endl;
164 }
165#endif
166}
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.) const
void SetTranslation(const G4ThreeVector &v)

◆ GetMaxParameter()

G4double G4ParameterisationConsRho::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 129 of file G4ParameterisationCons.cc.

130{
131 G4Cons* msol = (G4Cons*)(fmotherSolid);
132 return msol->GetOuterRadiusMinusZ() - msol->GetInnerRadiusMinusZ();
133}

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