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

#include <G4ParameterisationPolyhedra.hh>

+ Inheritance diagram for G4ParameterisationPolyhedraPhi:

Public Member Functions

 G4ParameterisationPolyhedraPhi (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
 
 ~G4ParameterisationPolyhedraPhi () 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 130 of file G4ParameterisationPolyhedra.hh.

Constructor & Destructor Documentation

◆ G4ParameterisationPolyhedraPhi()

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

Definition at line 279 of file G4ParameterisationPolyhedra.cc.

283 : G4VParameterisationPolyhedra( axis, nDiv, width, offset, msolid, divType )
284{
286 SetType( "DivisionPolyhedraPhi" );
287
288 auto msol = (G4Polyhedra*)(fmotherSolid);
289 G4double deltaPhi = msol->GetEndPhi() - msol->GetStartPhi();
290
291 if( divType == DivWIDTH )
292 {
293 fnDiv = msol->GetNumSide();
294 }
295
296 fwidth = CalculateWidth( deltaPhi, fnDiv, 0.0 );
297
298#ifdef G4DIVDEBUG
299 if( verbose >= 1 )
300 {
301 G4cout << " G4ParameterisationPolyhedraPhi - # divisions " << fnDiv
302 << " = " << nDiv << G4endl
303 << " Offset " << foffset << " = " << offset << G4endl
304 << " Width " << fwidth << " = " << width << G4endl;
305 }
306#endif
307}
double G4double
Definition G4Types.hh:83
#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
G4VParameterisationPolyhedra(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)

◆ ~G4ParameterisationPolyhedraPhi()

G4ParameterisationPolyhedraPhi::~G4ParameterisationPolyhedraPhi ( )
overridedefault

Member Function Documentation

◆ CheckParametersValidity()

void G4ParameterisationPolyhedraPhi::CheckParametersValidity ( )
overridevirtual

Reimplemented from G4VDivisionParameterisation.

Definition at line 320 of file G4ParameterisationPolyhedra.cc.

321{
323
324 auto msol = (G4Polyhedra*)(fmotherSolid);
325
327 {
328 std::ostringstream message;
329 message << "In solid " << msol->GetName() << G4endl
330 << " Division along PHI will be done splitting "
331 << "in the defined numSide." << G4endl
332 << "WIDTH will not be used !";
333 G4Exception("G4ParameterisationPolyhedraPhi::CheckParametersValidity()",
334 "GeomDiv1001", JustWarning, message);
335 }
336 if( foffset != 0. )
337 {
338 std::ostringstream message;
339 message << "In solid " << msol->GetName() << G4endl
340 << "Division along PHI will be done splitting "
341 << "in the defined numSide." << G4endl
342 << "OFFSET will not be used !";
343 G4Exception("G4ParameterisationPolyhedraPhi::CheckParametersValidity()",
344 "GeomDiv1001", JustWarning, message);
345 }
346
347 G4PolyhedraHistorical* origparamMother = msol->GetOriginalParameters();
348
349 if( origparamMother->numSide != fnDiv && fDivisionType != DivWIDTH)
350 {
351 std::ostringstream message;
352 message << "Configuration not supported." << G4endl
353 << "Division along PHI will be done splitting in the defined"
354 << G4endl
355 << "numSide, i.e, the number of division would be :"
356 << origparamMother->numSide << " instead of " << fnDiv << " !";
357 G4Exception("G4ParameterisationPolyhedraPhi::CheckParametersValidity()",
358 "GeomDiv0001", FatalException, message);
359 }
360}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)

Referenced by G4ParameterisationPolyhedraPhi().

◆ ComputeDimensions()

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

Reimplemented from G4VPVParameterisation.

Definition at line 400 of file G4ParameterisationPolyhedra.cc.

403{
404 auto msol = (G4Polyhedra*)(fmotherSolid);
405
406 G4PolyhedraHistorical* origparamMother = msol->GetOriginalParameters();
407 G4PolyhedraHistorical origparam( *origparamMother );
408
409 origparam.numSide = 1;
410 origparam.Start_angle = origparamMother->Start_angle;
411 origparam.Opening_angle = fwidth;
412
413 phedra.SetOriginalParameters(&origparam); // copy values & transfer pointers
414 phedra.Reset(); // reset to new solid parameters
415
416#ifdef G4DIVDEBUG
417 if( verbose >= 2 )
418 {
419 G4cout << "G4ParameterisationPolyhedraPhi::ComputeDimensions():" << G4endl;
420 phedra.DumpInfo();
421 }
422#endif
423}
void SetOriginalParameters(G4PolyhedraHistorical *pars)
G4bool Reset()
void DumpInfo() const

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 364 of file G4ParameterisationPolyhedra.cc.

366{
367 //----- translation
368 G4ThreeVector origin(0.,0.,0.);
369 //----- set translation
370 physVol->SetTranslation( origin );
371
372 //----- calculate rotation matrix (so that all volumes point to the centre)
373 G4double posi = copyNo*fwidth;
374
375#ifdef G4DIVDEBUG
376 if( verbose >= 2 )
377 {
378 G4cout << " G4ParameterisationPolyhedraPhi - position: " << posi/CLHEP::deg
379 << G4endl
380 << " copyNo: " << copyNo
381 << " - fwidth: " << fwidth/CLHEP::deg << G4endl;
382 }
383#endif
384
385 ChangeRotMatrix( physVol, -posi );
386
387#ifdef G4DIVDEBUG
388 if( verbose >= 2 )
389 {
390 G4cout << std::setprecision(8) << " G4ParameterisationPolyhedraPhi "
391 << copyNo << G4endl
392 << " Position: " << origin << " - Width: " << fwidth
393 << " - Axis: " << faxis << G4endl;
394 }
395#endif
396}
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.0) const
void SetTranslation(const G4ThreeVector &v)

◆ GetMaxParameter()

G4double G4ParameterisationPolyhedraPhi::GetMaxParameter ( ) const
overridevirtual

Implements G4VDivisionParameterisation.

Definition at line 313 of file G4ParameterisationPolyhedra.cc.

314{
315 auto msol = (G4Polyhedra*)(fmotherSolid);
316 return msol->GetEndPhi() - msol->GetStartPhi();
317}

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