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

#include <G4ParameterisationPolyhedra.hh>

+ Inheritance diagram for G4ParameterisationPolyhedraRho:

Public Member Functions

 G4ParameterisationPolyhedraRho (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
 
 ~G4ParameterisationPolyhedraRho ()
 
void CheckParametersValidity ()
 
G4double GetMaxParameter () const
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const
 
void ComputeDimensions (G4Polyhedra &phedra, const G4int copyNo, const G4VPhysicalVolume *physVol) const
 
- Public Member Functions inherited from G4VParameterisationPolyhedra
 G4VParameterisationPolyhedra (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
virtual ~G4VParameterisationPolyhedra ()
 
- 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 85 of file G4ParameterisationPolyhedra.hh.

Constructor & Destructor Documentation

◆ G4ParameterisationPolyhedraRho()

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

Definition at line 128 of file G4ParameterisationPolyhedra.cc.

132 : G4VParameterisationPolyhedra( axis, nDiv, width, offset, msolid, divType )
133{
134
136 SetType( "DivisionPolyhedraRho" );
137
139 G4PolyhedraHistorical* original_pars = msol->GetOriginalParameters();
140
141 if( divType == DivWIDTH )
142 {
143 fnDiv = CalculateNDiv( original_pars->Rmax[0]
144 - original_pars->Rmin[0], width, offset );
145 }
146 else if( divType == DivNDIV )
147 {
148 fwidth = CalculateWidth( original_pars->Rmax[0]
149 - original_pars->Rmin[0], nDiv, offset );
150 }
151
152#ifdef G4DIVDEBUG
153 if( verbose >= 1 )
154 {
155 G4cout << " G4ParameterisationPolyhedraRho - # divisions " << fnDiv
156 << " = " << nDiv << G4endl
157 << " Offset " << foffset << " = " << offset << G4endl
158 << " Width " << fwidth << " = " << width << G4endl;
159 }
160#endif
161}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4PolyhedraHistorical * GetOriginalParameters() const
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const

◆ ~G4ParameterisationPolyhedraRho()

G4ParameterisationPolyhedraRho::~G4ParameterisationPolyhedraRho ( )

Definition at line 164 of file G4ParameterisationPolyhedra.cc.

165{
166}

Member Function Documentation

◆ CheckParametersValidity()

void G4ParameterisationPolyhedraRho::CheckParametersValidity ( )
virtual

Reimplemented from G4VDivisionParameterisation.

Definition at line 169 of file G4ParameterisationPolyhedra.cc.

170{
172
174
176 {
177 std::ostringstream message;
178 message << "In solid " << msol->GetName() << G4endl
179 << "Division along R will be done with a width "
180 << "different for each solid section." << G4endl
181 << "WIDTH will not be used !";
182 G4Exception("G4ParameterisationPolyhedraRho::CheckParametersValidity()",
183 "GeomDiv1001", JustWarning, message);
184 }
185 if( foffset != 0. )
186 {
187 std::ostringstream message;
188 message << "In solid " << msol->GetName() << G4endl
189 << "Division along R will be done with a width "
190 << "different for each solid section." << G4endl
191 << "OFFSET will not be used !";
192 G4Exception("G4ParameterisationPolyhedraRho::CheckParametersValidity()",
193 "GeomDiv1001", JustWarning, message);
194 }
195}
@ JustWarning
G4String GetName() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Referenced by G4ParameterisationPolyhedraRho().

◆ ComputeDimensions()

void G4ParameterisationPolyhedraRho::ComputeDimensions ( G4Polyhedra phedra,
const G4int  copyNo,
const G4VPhysicalVolume physVol 
) const
virtual

Reimplemented from G4VPVParameterisation.

Definition at line 243 of file G4ParameterisationPolyhedra.cc.

246{
248
249 G4PolyhedraHistorical* origparamMother = msol->GetOriginalParameters();
250 G4PolyhedraHistorical origparam( *origparamMother );
251 G4int nZplanes = origparamMother->Num_z_planes;
252
253 G4double width = 0.;
254 for( G4int ii = 0; ii < nZplanes; ii++ )
255 {
256 width = CalculateWidth( origparamMother->Rmax[ii]
257 - origparamMother->Rmin[ii], fnDiv, foffset );
258 origparam.Rmin[ii] = origparamMother->Rmin[ii]+foffset+width*copyNo;
259 origparam.Rmax[ii] = origparamMother->Rmin[ii]+foffset+width*(copyNo+1);
260 }
261
262 phedra.SetOriginalParameters(&origparam); // copy values & transfer pointers
263 phedra.Reset(); // reset to new solid parameters
264
265#ifdef G4DIVDEBUG
266 if( verbose >= -2 )
267 {
268 G4cout << "G4ParameterisationPolyhedraRho::ComputeDimensions()" << G4endl
269 << "-- Parametrised phedra copy-number: " << copyNo << G4endl;
270 phedra.DumpInfo();
271 }
272#endif
273}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
void SetOriginalParameters(G4PolyhedraHistorical *pars)
G4bool Reset()
Definition: G4Polyhedra.cc:465
void DumpInfo() const

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 207 of file G4ParameterisationPolyhedra.cc.

209{
210 //----- translation
211 G4ThreeVector origin(0.,0.,0.);
212
213 //----- set translation
214 physVol->SetTranslation( origin );
215
216 //----- calculate rotation matrix: unit
217
218#ifdef G4DIVDEBUG
219 if( verbose >= 2 )
220 {
221 G4cout << " G4ParameterisationPolyhedraRho " << G4endl
222 << " foffset: " << foffset/deg
223 << " - fwidth: " << fwidth/deg << G4endl;
224 }
225#endif
226
227 ChangeRotMatrix( physVol );
228
229#ifdef G4DIVDEBUG
230 if( verbose >= 2 )
231 {
232 G4cout << std::setprecision(8) << " G4ParameterisationPolyhedraRho "
233 << G4endl
234 << " Position: " << origin
235 << " - Width: " << fwidth
236 << " - Axis: " << faxis << G4endl;
237 }
238#endif
239}
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.) const
void SetTranslation(const G4ThreeVector &v)

◆ GetMaxParameter()

G4double G4ParameterisationPolyhedraRho::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 198 of file G4ParameterisationPolyhedra.cc.

199{
201 G4PolyhedraHistorical* original_pars = msol->GetOriginalParameters();
202 return original_pars->Rmax[0] - original_pars->Rmin[0];
203}

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