Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ParameterisationPolyhedraPhi Class Reference

#include <G4ParameterisationPolyhedra.hh>

Inheritance diagram for G4ParameterisationPolyhedraPhi:
Collaboration 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 *)
 
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 G4MaterialComputeMaterial (const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=0)
 
virtual G4bool IsNested () const
 
virtual G4VVolumeMaterialScannerGetMaterialScanner ()
 

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
 
void CheckOffset (G4double maxPar)
 
void CheckNDivAndWidth (G4double maxPar)
 
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 G4ThreadLocal
G4RotationMatrix
fRot = 0
 
static const G4int verbose = 5
 

Detailed Description

Definition at line 137 of file G4ParameterisationPolyhedra.hh.

Constructor & Destructor Documentation

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

Definition at line 285 of file G4ParameterisationPolyhedra.cc.

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

Here is the call graph for this function:

G4ParameterisationPolyhedraPhi::~G4ParameterisationPolyhedraPhi ( )

Definition at line 315 of file G4ParameterisationPolyhedra.cc.

316 {
317 }

Member Function Documentation

void G4ParameterisationPolyhedraPhi::CheckParametersValidity ( )
virtual

Reimplemented from G4VDivisionParameterisation.

Definition at line 327 of file G4ParameterisationPolyhedra.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

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

Reimplemented from G4VPVParameterisation.

Definition at line 408 of file G4ParameterisationPolyhedra.cc.

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

Here is the call graph for this function:

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

Implements G4VDivisionParameterisation.

Definition at line 372 of file G4ParameterisationPolyhedra.cc.

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

Here is the call graph for this function:

G4double G4ParameterisationPolyhedraPhi::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 320 of file G4ParameterisationPolyhedra.cc.

321 {
323  return msol->GetEndPhi() - msol->GetStartPhi();
324 }
G4double GetEndPhi() const
G4double GetStartPhi() const

Here is the call graph for this function:


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