Geant4 11.1.1
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 ()
 
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=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 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 283 of file G4ParameterisationPolyhedra.cc.

287 : G4VParameterisationPolyhedra( axis, nDiv, width, offset, msolid, divType )
288{
290 SetType( "DivisionPolyhedraPhi" );
291
293 G4double deltaPhi = msol->GetEndPhi() - msol->GetStartPhi();
294
295 if( divType == DivWIDTH )
296 {
297 fnDiv = msol->GetNumSide();
298 }
299
300 fwidth = CalculateWidth( deltaPhi, fnDiv, 0.0 );
301
302#ifdef G4DIVDEBUG
303 if( verbose >= 1 )
304 {
305 G4cout << " G4ParameterisationPolyhedraPhi - # divisions " << fnDiv
306 << " = " << nDiv << G4endl
307 << " Offset " << foffset << " = " << offset << G4endl
308 << " Width " << fwidth << " = " << width << G4endl;
309 }
310#endif
311}
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double GetEndPhi() const
G4int GetNumSide() const
G4double GetStartPhi() const
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const

◆ ~G4ParameterisationPolyhedraPhi()

G4ParameterisationPolyhedraPhi::~G4ParameterisationPolyhedraPhi ( )

Definition at line 314 of file G4ParameterisationPolyhedra.cc.

315{
316}

Member Function Documentation

◆ CheckParametersValidity()

void G4ParameterisationPolyhedraPhi::CheckParametersValidity ( )
virtual

Reimplemented from G4VDivisionParameterisation.

Definition at line 326 of file G4ParameterisationPolyhedra.cc.

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

Referenced by G4ParameterisationPolyhedraPhi().

◆ ComputeDimensions()

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

Reimplemented from G4VPVParameterisation.

Definition at line 406 of file G4ParameterisationPolyhedra.cc.

409{
411
412 G4PolyhedraHistorical* origparamMother = msol->GetOriginalParameters();
413 G4PolyhedraHistorical origparam( *origparamMother );
414
415 origparam.numSide = 1;
416 origparam.Start_angle = origparamMother->Start_angle;
417 origparam.Opening_angle = fwidth;
418
419 phedra.SetOriginalParameters(&origparam); // copy values & transfer pointers
420 phedra.Reset(); // reset to new solid parameters
421
422#ifdef G4DIVDEBUG
423 if( verbose >= 2 )
424 {
425 G4cout << "G4ParameterisationPolyhedraPhi::ComputeDimensions():" << G4endl;
426 phedra.DumpInfo();
427 }
428#endif
429}
void SetOriginalParameters(G4PolyhedraHistorical *pars)
G4bool Reset()
Definition: G4Polyhedra.cc:463
void DumpInfo() const

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 370 of file G4ParameterisationPolyhedra.cc.

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

◆ GetMaxParameter()

G4double G4ParameterisationPolyhedraPhi::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 319 of file G4ParameterisationPolyhedra.cc.

320{
322 return msol->GetEndPhi() - msol->GetStartPhi();
323}

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