Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ParameterisationPolyconePhi Class Reference

#include <G4ParameterisationPolycone.hh>

Inheritance diagram for G4ParameterisationPolyconePhi:
Collaboration diagram for G4ParameterisationPolyconePhi:

Public Member Functions

 G4ParameterisationPolyconePhi (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
 
 ~G4ParameterisationPolyconePhi ()
 
G4double GetMaxParameter () const
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const
 
void ComputeDimensions (G4Polycone &pcone, const G4int copyNo, const G4VPhysicalVolume *physVol) const
 
- Public Member Functions inherited from G4VParameterisationPolycone
 G4VParameterisationPolycone (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
virtual ~G4VParameterisationPolycone ()
 
- 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
 
virtual void CheckParametersValidity ()
 
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 131 of file G4ParameterisationPolycone.hh.

Constructor & Destructor Documentation

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

Definition at line 246 of file G4ParameterisationPolycone.cc.

249  : G4VParameterisationPolycone( axis, nDiv, width, offset, msolid, divType )
250 {
252  SetType( "DivisionPolyconePhi" );
253 
254  G4Polycone* msol = (G4Polycone*)(fmotherSolid);
255  G4double deltaPhi = msol->GetEndPhi() - msol->GetStartPhi();
256 
257  if( divType == DivWIDTH )
258  {
259  fnDiv = CalculateNDiv( deltaPhi, width, offset );
260  }
261  else if( divType == DivNDIV )
262  {
263  fwidth = CalculateWidth( deltaPhi, nDiv, offset );
264  }
265 
266 #ifdef G4DIVDEBUG
267  if( verbose >= 1 )
268  {
269  G4cout << " G4ParameterisationPolyconePhi - # divisions " << fnDiv
270  << " = " << nDiv << G4endl
271  << " Offset " << foffset/deg << " = " << offset/deg << G4endl
272  << " Width " << fwidth/deg << " = " << width/deg << G4endl;
273  }
274 #endif
275 }
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4double GetEndPhi() const
G4GLOB_DLL std::ostream G4cout
G4VParameterisationPolycone(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
G4double GetStartPhi() const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const
#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:

G4ParameterisationPolyconePhi::~G4ParameterisationPolyconePhi ( )

Definition at line 278 of file G4ParameterisationPolycone.cc.

279 {
280 }

Member Function Documentation

void G4ParameterisationPolyconePhi::ComputeDimensions ( G4Polycone pcone,
const G4int  copyNo,
const G4VPhysicalVolume physVol 
) const
virtual

Reimplemented from G4VPVParameterisation.

Definition at line 328 of file G4ParameterisationPolycone.cc.

330 {
331  G4Polycone* msol = (G4Polycone*)(fmotherSolid);
332 
333  G4PolyconeHistorical* origparamMother = msol->GetOriginalParameters();
334  G4PolyconeHistorical origparam( *origparamMother );
335  origparam.Start_angle = origparamMother->Start_angle;
336  origparam.Opening_angle = fwidth;
337 
338  pcone.SetOriginalParameters(&origparam); // copy values & transfer pointers
339  pcone.Reset(); // reset to new solid parameters
340 
341 #ifdef G4DIVDEBUG
342  if( verbose >= 2 )
343  {
344  G4cout << "G4ParameterisationPolyconePhi::ComputeDimensions():" << G4endl;
345  pcone.DumpInfo();
346  }
347 #endif
348 }
void DumpInfo() const
G4GLOB_DLL std::ostream G4cout
void SetOriginalParameters(G4PolyconeHistorical *pars)
G4PolyconeHistorical * GetOriginalParameters() const
G4bool Reset()
Definition: G4Polycone.cc:447
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

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

Implements G4VDivisionParameterisation.

Definition at line 292 of file G4ParameterisationPolycone.cc.

293 {
294  //----- translation
295  G4ThreeVector origin(0.,0.,0.);
296  //----- set translation
297  physVol->SetTranslation( origin );
298 
299  //----- calculate rotation matrix (so that all volumes point to the centre)
300  G4double posi = foffset + copyNo*fwidth;
301 
302 #ifdef G4DIVDEBUG
303  if( verbose >= 2 )
304  {
305  G4cout << " G4ParameterisationPolyconePhi - position: " << posi/deg
306  << G4endl
307  << " copyNo: " << copyNo << " - foffset: " << foffset/deg
308  << " - fwidth: " << fwidth/deg << G4endl;
309  }
310 #endif
311 
312  ChangeRotMatrix( physVol, -posi );
313 
314 #ifdef G4DIVDEBUG
315  if( verbose >= 2 )
316  {
317  G4cout << std::setprecision(8) << " G4ParameterisationPolyconePhi "
318  << copyNo << G4endl
319  << " Position: " << origin << " - Width: " << fwidth
320  << " - Axis: " << faxis << G4endl;
321  }
322 #endif
323 }
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 G4ParameterisationPolyconePhi::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 283 of file G4ParameterisationPolycone.cc.

284 {
285  G4Polycone* msol = (G4Polycone*)(fmotherSolid);
286  return msol->GetEndPhi() - msol->GetStartPhi();
287 }
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: