Geant4 9.6.0
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=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 134 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 276 of file G4ParameterisationPolyhedra.cc.

280 : G4VParameterisationPolyhedra( axis, nDiv, width, offset, msolid, divType )
281{
283 SetType( "DivisionPolyhedraPhi" );
284
286 G4double deltaPhi = msol->GetEndPhi() - msol->GetStartPhi();
287
288 if( divType == DivWIDTH )
289 {
290 fnDiv = msol->GetNumSide();
291 }
292
293 fwidth = CalculateWidth( deltaPhi, fnDiv, 0.0 );
294
295#ifdef G4DIVDEBUG
296 if( verbose >= 1 )
297 {
298 G4cout << " G4ParameterisationPolyhedraPhi - # divisions " << fnDiv
299 << " = " << nDiv << G4endl
300 << " Offset " << foffset << " = " << offset << G4endl
301 << " Width " << fwidth << " = " << width << G4endl;
302 }
303#endif
304}
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT 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 307 of file G4ParameterisationPolyhedra.cc.

308{
309}

Member Function Documentation

◆ CheckParametersValidity()

void G4ParameterisationPolyhedraPhi::CheckParametersValidity ( )
virtual

Reimplemented from G4VDivisionParameterisation.

Definition at line 319 of file G4ParameterisationPolyhedra.cc.

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

Referenced by G4ParameterisationPolyhedraPhi().

◆ ComputeDimensions()

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

Reimplemented from G4VPVParameterisation.

Definition at line 399 of file G4ParameterisationPolyhedra.cc.

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

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 363 of file G4ParameterisationPolyhedra.cc.

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

◆ GetMaxParameter()

G4double G4ParameterisationPolyhedraPhi::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 312 of file G4ParameterisationPolyhedra.cc.

313{
315 return msol->GetEndPhi() - msol->GetStartPhi();
316}

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