G4ParameterisationConsZ Class Reference

#include <G4ParameterisationCons.hh>

Inheritance diagram for G4ParameterisationConsZ:

G4VParameterisationCons G4VDivisionParameterisation G4VPVParameterisation

Public Member Functions

 G4ParameterisationConsZ (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
 ~G4ParameterisationConsZ ()
G4double GetMaxParameter () const
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const
void ComputeDimensions (G4Cons &tubs, const G4int copyNo, const G4VPhysicalVolume *physVol) const

Detailed Description

Definition at line 160 of file G4ParameterisationCons.hh.


Constructor & Destructor Documentation

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

Definition at line 322 of file G4ParameterisationCons.cc.

References G4VDivisionParameterisation::CalculateNDiv(), G4VDivisionParameterisation::CalculateWidth(), G4VDivisionParameterisation::CheckParametersValidity(), DivNDIV, DivWIDTH, G4VDivisionParameterisation::faxis, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::fnDiv, G4VDivisionParameterisation::foffset, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4Cons::GetZHalfLength(), G4VDivisionParameterisation::SetType(), and G4VDivisionParameterisation::verbose.

00325   :  G4VParameterisationCons( axis, nDiv, width, offset, msolid, divType )
00326 { 
00327   CheckParametersValidity();
00328   SetType( "DivisionConsZ" );
00329 
00330   G4Cons* msol = (G4Cons*)(fmotherSolid);
00331   if( divType == DivWIDTH )
00332   {
00333     fnDiv = CalculateNDiv( 2*msol->GetZHalfLength(), width, offset );
00334   }
00335   else if( divType == DivNDIV )
00336   {
00337     fwidth = CalculateWidth( 2*msol->GetZHalfLength(), nDiv, offset );
00338   }
00339 
00340 #ifdef G4DIVDEBUG
00341   if( verbose >= 1 )
00342   {
00343     G4cout << " G4ParameterisationConsZ: # divisions " << fnDiv << " = "
00344            << nDiv << G4endl
00345            << " Offset " << foffset << " = " << offset << G4endl
00346            << " Width " << fwidth << " = " << width << G4endl
00347            << " - Axis: " << faxis << G4endl;
00348   }
00349 #endif
00350 }

G4ParameterisationConsZ::~G4ParameterisationConsZ (  ) 

Definition at line 353 of file G4ParameterisationCons.cc.

00354 {
00355 }


Member Function Documentation

void G4ParameterisationConsZ::ComputeDimensions ( G4Cons tubs,
const G4int  copyNo,
const G4VPhysicalVolume physVol 
) const [virtual]

Reimplemented from G4VPVParameterisation.

Definition at line 395 of file G4ParameterisationCons.cc.

References G4VSolid::DumpInfo(), G4VDivisionParameterisation::fhgap, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4Cons::GetDeltaPhiAngle(), G4Cons::GetInnerRadiusMinusZ(), G4Cons::GetInnerRadiusPlusZ(), G4Cons::GetOuterRadiusMinusZ(), G4Cons::GetOuterRadiusPlusZ(), G4Cons::GetStartPhiAngle(), G4Cons::GetZHalfLength(), G4VDivisionParameterisation::OffsetZ(), G4Cons::SetDeltaPhiAngle(), G4Cons::SetInnerRadiusMinusZ(), G4Cons::SetInnerRadiusPlusZ(), G4Cons::SetOuterRadiusMinusZ(), G4Cons::SetOuterRadiusPlusZ(), G4Cons::SetStartPhiAngle(), G4Cons::SetZHalfLength(), and G4VDivisionParameterisation::verbose.

00397 {
00398   G4Cons* msol = (G4Cons*)(fmotherSolid);
00399 
00400   G4double mHalfLength = msol->GetZHalfLength() - fhgap;
00401   G4double aRInner = (msol->GetInnerRadiusPlusZ()
00402                    - msol->GetInnerRadiusMinusZ()) / (2*mHalfLength);
00403   G4double bRInner = (msol->GetInnerRadiusPlusZ()
00404                    + msol->GetInnerRadiusMinusZ()) / 2;
00405   G4double aROuter = (msol->GetOuterRadiusPlusZ()
00406                    - msol->GetOuterRadiusMinusZ()) / (2*mHalfLength);
00407   G4double bROuter = (msol->GetOuterRadiusPlusZ()
00408                    + msol->GetOuterRadiusMinusZ()) / 2;
00409   G4double xMinusZ = -mHalfLength + OffsetZ() + fwidth*copyNo + fhgap;
00410   G4double xPlusZ  = -mHalfLength + OffsetZ() + fwidth*(copyNo+1) - fhgap;
00411   cons.SetInnerRadiusMinusZ( aRInner * xMinusZ + bRInner );
00412   cons.SetOuterRadiusMinusZ( aROuter * xMinusZ + bROuter );
00413   cons.SetInnerRadiusPlusZ( aRInner * xPlusZ + bRInner );
00414   cons.SetOuterRadiusPlusZ( aROuter * xPlusZ + bROuter );
00415  
00416   G4double pDz = fwidth / 2. - fhgap;
00417   G4double pSPhi = msol->GetStartPhiAngle();
00418   G4double pDPhi = msol->GetDeltaPhiAngle();
00419 
00420   cons.SetZHalfLength( pDz );
00421   cons.SetStartPhiAngle( pSPhi, false );
00422   cons.SetDeltaPhiAngle( pDPhi );
00423 
00424 #ifdef G4DIVDEBUG
00425   if( verbose >= 2 )
00426   {
00427     G4cout << " G4ParameterisationConsZ::ComputeDimensions()" << G4endl
00428            << " pDz: " << pDz << G4endl;
00429     if( verbose >= 4 ) cons.DumpInfo();
00430   }
00431 #endif
00432 
00433 }

void G4ParameterisationConsZ::ComputeTransformation ( const G4int  copyNo,
G4VPhysicalVolume physVol 
) const [virtual]

Implements G4VDivisionParameterisation.

Definition at line 367 of file G4ParameterisationCons.cc.

References G4VDivisionParameterisation::ChangeRotMatrix(), G4VDivisionParameterisation::foffset, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4VDivisionParameterisation::GetMotherSolid(), G4Cons::GetZHalfLength(), G4VDivisionParameterisation::OffsetZ(), G4VPhysicalVolume::SetTranslation(), and G4VDivisionParameterisation::verbose.

00368 {
00369   //----- set translation: along Z axis
00370   G4Cons* motherCons = (G4Cons*)(GetMotherSolid());
00371   G4double posi = - motherCons->GetZHalfLength() + OffsetZ()
00372                   + fwidth/2 + copyNo*fwidth;
00373   G4ThreeVector origin(0.,0.,posi); 
00374   physVol->SetTranslation( origin );
00375 
00376   //----- calculate rotation matrix: unit
00377 
00378 #ifdef G4DIVDEBUG
00379   if( verbose >= 2 )
00380   {
00381     G4cout << " G4ParameterisationConsZ::ComputeTransformation()" << G4endl
00382            << " Origin: " << origin << " - copyNo: " << copyNo << G4endl
00383            << " foffset: " << foffset << " - fwidth: " << fwidth
00384            << G4endl;
00385   }
00386 #endif
00387 
00388   ChangeRotMatrix( physVol );
00389 }

G4double G4ParameterisationConsZ::GetMaxParameter (  )  const [virtual]

Implements G4VDivisionParameterisation.

Definition at line 358 of file G4ParameterisationCons.cc.

References G4VDivisionParameterisation::fmotherSolid, and G4Cons::GetZHalfLength().

00359 {
00360   G4Cons* msol = (G4Cons*)(fmotherSolid);
00361   return 2*msol->GetZHalfLength();
00362 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:52:49 2013 for Geant4 by  doxygen 1.4.7