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

#include <G4ParameterisationCons.hh>

+ Inheritance diagram for G4ParameterisationConsPhi:

Public Member Functions

 G4ParameterisationConsPhi (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
 
 ~G4ParameterisationConsPhi ()
 
G4double GetMaxParameter () const
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const
 
void ComputeDimensions (G4Cons &tubs, const G4int copyNo, const G4VPhysicalVolume *physVol) const
 
- Public Member Functions inherited from G4VParameterisationCons
 G4VParameterisationCons (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
virtual ~G4VParameterisationCons ()
 
- 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 118 of file G4ParameterisationCons.hh.

Constructor & Destructor Documentation

◆ G4ParameterisationConsPhi()

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

Definition at line 215 of file G4ParameterisationCons.cc.

219 : G4VParameterisationCons( axis, nDiv, width, offset, msolid, divType )
220{
222 SetType( "DivisionConsPhi" );
223
224 G4Cons* msol = (G4Cons*)(fmotherSolid);
225 if( divType == DivWIDTH )
226 {
227 fnDiv = CalculateNDiv( msol->GetDeltaPhiAngle(), width, offset );
228 }
229 else if( divType == DivNDIV )
230 {
231 fwidth = CalculateWidth( msol->GetDeltaPhiAngle(), nDiv, offset );
232 }
233
234#ifdef G4DIVDEBUG
235 if( verbose >= 1 )
236 {
237 G4cout << " G4ParameterisationConsPhi no divisions " << fnDiv << " = "
238 << nDiv << G4endl
239 << " Offset " << foffset << " = " << offset << G4endl
240 << " Width " << fwidth << " = " << width << G4endl;
241 }
242#endif
243}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
Definition: G4Cons.hh:75
G4double GetDeltaPhiAngle() const
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const

◆ ~G4ParameterisationConsPhi()

G4ParameterisationConsPhi::~G4ParameterisationConsPhi ( )

Definition at line 246 of file G4ParameterisationCons.cc.

247{
248}

Member Function Documentation

◆ ComputeDimensions()

void G4ParameterisationConsPhi::ComputeDimensions ( G4Cons tubs,
const G4int  copyNo,
const G4VPhysicalVolume physVol 
) const
virtual

Reimplemented from G4VPVParameterisation.

Definition at line 286 of file G4ParameterisationCons.cc.

289{
290 G4Cons* msol = (G4Cons*)(fmotherSolid);
291
292 G4double pRMin1 = msol->GetInnerRadiusMinusZ();
293 G4double pRMax1 = msol->GetOuterRadiusMinusZ();
294 G4double pRMin2 = msol->GetInnerRadiusPlusZ();
295 G4double pRMax2 = msol->GetOuterRadiusPlusZ();
296 G4double pDz = msol->GetZHalfLength();
297
298 //- already rotated double pSPhi = foffset + copyNo*fwidth;
299 G4double pSPhi = foffset + msol->GetStartPhiAngle() + fhgap;
300 G4double pDPhi = fwidth - 2.*fhgap;
301
302 cons.SetInnerRadiusMinusZ( pRMin1 );
303 cons.SetOuterRadiusMinusZ( pRMax1 );
304 cons.SetInnerRadiusPlusZ( pRMin2 );
305 cons.SetOuterRadiusPlusZ( pRMax2 );
306 cons.SetZHalfLength( pDz );
307 cons.SetStartPhiAngle( pSPhi, false );
308 cons.SetDeltaPhiAngle( pDPhi );
309
310#ifdef G4DIVDEBUG
311 if( verbose >= 2 )
312 {
313 G4cout << " G4ParameterisationConsPhi::ComputeDimensions" << G4endl
314 << " pSPhi: " << pSPhi << " - pDPhi: " << pDPhi << G4endl;
315 if( verbose >= 4 ) cons.DumpInfo();
316 }
317#endif
318}
double G4double
Definition: G4Types.hh:64
G4double GetOuterRadiusPlusZ() const
G4double GetStartPhiAngle() const
G4double GetInnerRadiusMinusZ() const
G4double GetInnerRadiusPlusZ() const
G4double GetOuterRadiusMinusZ() const
G4double GetZHalfLength() const

◆ ComputeTransformation()

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

Implements G4VDivisionParameterisation.

Definition at line 259 of file G4ParameterisationCons.cc.

261{
262 //----- translation
263 G4ThreeVector origin(0.,0.,0.);
264 //----- set translation
265 physVol->SetTranslation( origin );
266
267 //----- calculate rotation matrix (so that all volumes point to the centre)
268 G4double posi = foffset + copyNo*fwidth;
269
270#ifdef G4DIVDEBUG
271 if( verbose >= 2 )
272 {
273 G4cout << " G4ParameterisationConsPhi - position: " << posi/deg << G4endl
274 << " Origin: " << origin << " copyNo: " << copyNo
275 << " - foffset: " << foffset/deg
276 << " - fwidth: " << fwidth/deg << G4endl
277 << " - Axis: " << faxis << G4endl;
278 }
279#endif
280
281 ChangeRotMatrix( physVol, -posi );
282}
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.) const
void SetTranslation(const G4ThreeVector &v)

◆ GetMaxParameter()

G4double G4ParameterisationConsPhi::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 251 of file G4ParameterisationCons.cc.

252{
253 G4Cons* msol = (G4Cons*)(fmotherSolid);
254 return msol->GetDeltaPhiAngle();
255}

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