Geant4 11.2.2
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 () override
 
void CheckParametersValidity () override
 
G4double GetMaxParameter () const override
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const override
 
void ComputeDimensions (G4Polyhedra &phedra, const G4int copyNo, const G4VPhysicalVolume *physVol) const override
 
- Public Member Functions inherited from G4VParameterisationPolyhedra
 G4VParameterisationPolyhedra (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
 ~G4VParameterisationPolyhedra () override
 
- Public Member Functions inherited from G4VDivisionParameterisation
 G4VDivisionParameterisation (EAxis axis, G4int nDiv, G4double width, G4double offset, DivisionType divType, G4VSolid *motherSolid=nullptr)
 
 ~G4VDivisionParameterisation () override
 
G4VSolidComputeSolid (const G4int, G4VPhysicalVolume *) override
 
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 G4MaterialComputeMaterial (const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=nullptr)
 
virtual G4bool IsNested () const
 
virtual G4VVolumeMaterialScannerGetMaterialScanner ()
 

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
 
void CheckOffset (G4double maxPar)
 
void CheckNDivAndWidth (G4double maxPar)
 
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 79 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 133 of file G4ParameterisationPolyhedra.cc.

137 : G4VParameterisationPolyhedra( axis, nDiv, width, offset, msolid, divType )
138{
139
141 SetType( "DivisionPolyhedraRho" );
142
143 auto msol = (G4Polyhedra*)(fmotherSolid);
144 G4PolyhedraHistorical* original_pars = msol->GetOriginalParameters();
145
146 if( divType == DivWIDTH )
147 {
148 fnDiv = CalculateNDiv( original_pars->Rmax[0]
149 - original_pars->Rmin[0], width, offset );
150 }
151 else if( divType == DivNDIV )
152 {
153 fwidth = CalculateWidth( original_pars->Rmax[0]
154 - original_pars->Rmin[0], nDiv, offset );
155 }
156
157#ifdef G4DIVDEBUG
158 if( verbose >= 1 )
159 {
160 G4cout << " G4ParameterisationPolyhedraRho - # divisions " << fnDiv
161 << " = " << nDiv << G4endl
162 << " Offset " << foffset << " = " << offset << G4endl
163 << " Width " << fwidth << " = " << width << G4endl;
164 }
165#endif
166}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const
G4VParameterisationPolyhedra(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)

◆ ~G4ParameterisationPolyhedraRho()

G4ParameterisationPolyhedraRho::~G4ParameterisationPolyhedraRho ( )
overridedefault

Member Function Documentation

◆ CheckParametersValidity()

void G4ParameterisationPolyhedraRho::CheckParametersValidity ( )
overridevirtual

Reimplemented from G4VDivisionParameterisation.

Definition at line 172 of file G4ParameterisationPolyhedra.cc.

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

Referenced by G4ParameterisationPolyhedraRho().

◆ ComputeDimensions()

void G4ParameterisationPolyhedraRho::ComputeDimensions ( G4Polyhedra & phedra,
const G4int copyNo,
const G4VPhysicalVolume * physVol ) const
overridevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 246 of file G4ParameterisationPolyhedra.cc.

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

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 210 of file G4ParameterisationPolyhedra.cc.

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

◆ GetMaxParameter()

G4double G4ParameterisationPolyhedraRho::GetMaxParameter ( ) const
overridevirtual

Implements G4VDivisionParameterisation.

Definition at line 201 of file G4ParameterisationPolyhedra.cc.

202{
203 auto msol = (G4Polyhedra*)(fmotherSolid);
204 G4PolyhedraHistorical* original_pars = msol->GetOriginalParameters();
205 return original_pars->Rmax[0] - original_pars->Rmin[0];
206}

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