G4ParameterisationBoxY Class Reference

#include <G4ParameterisationBox.hh>

Inheritance diagram for G4ParameterisationBoxY:

G4VParameterisationBox G4VDivisionParameterisation G4VPVParameterisation

Public Member Functions

 G4ParameterisationBoxY (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 ~G4ParameterisationBoxY ()
G4double GetMaxParameter () const
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const
void ComputeDimensions (G4Box &box, const G4int copyNo, const G4VPhysicalVolume *physVol) const

Detailed Description

Definition at line 117 of file G4ParameterisationBox.hh.


Constructor & Destructor Documentation

G4ParameterisationBoxY::G4ParameterisationBoxY ( EAxis  axis,
G4int  nCopies,
G4double  offset,
G4double  step,
G4VSolid msolid,
DivisionType  divType 
)

Definition at line 175 of file G4ParameterisationBox.cc.

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

00178   :  G4VParameterisationBox( axis, nDiv, width, offset, msolid, divType )
00179 {
00180   CheckParametersValidity();
00181   SetType( "DivisionBoxY" );
00182 
00183   G4Box* mbox = (G4Box*)(fmotherSolid);
00184   if( divType == DivWIDTH )
00185   {
00186     fnDiv = CalculateNDiv( 2*mbox->GetYHalfLength(), width, offset );
00187   }
00188   else if( divType == DivNDIV )
00189   {
00190     fwidth = CalculateWidth( 2*mbox->GetYHalfLength(), nDiv, offset );
00191   }
00192 
00193 #ifdef G4DIVDEBUG
00194   if( verbose >= 1 )
00195   {
00196     G4cout << " G4ParameterisationBoxY - no divisions " << fnDiv << " = "
00197            << nDiv << ". Offset " << foffset << " = " << offset
00198            << ". Width " << fwidth << " = " << width << G4endl;
00199   }
00200 #endif
00201 }

G4ParameterisationBoxY::~G4ParameterisationBoxY (  ) 

Definition at line 204 of file G4ParameterisationBox.cc.

00205 {
00206 }


Member Function Documentation

void G4ParameterisationBoxY::ComputeDimensions ( G4Box box,
const G4int  copyNo,
const G4VPhysicalVolume physVol 
) const [virtual]

Reimplemented from G4VPVParameterisation.

Definition at line 252 of file G4ParameterisationBox.cc.

References G4VSolid::DumpInfo(), G4VDivisionParameterisation::fhgap, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4Box::GetXHalfLength(), G4Box::GetZHalfLength(), G4Box::SetXHalfLength(), G4Box::SetYHalfLength(), G4Box::SetZHalfLength(), and G4VDivisionParameterisation::verbose.

00254 {
00255   G4Box* msol = (G4Box*)(fmotherSolid);
00256 
00257   G4double pDx = msol->GetXHalfLength();
00258   G4double pDy = fwidth/2. - fhgap;
00259   G4double pDz = msol->GetZHalfLength();
00260 
00261   box.SetXHalfLength( pDx );
00262   box.SetYHalfLength( pDy );
00263   box.SetZHalfLength( pDz );
00264 
00265 #ifdef G4DIVDEBUG
00266   if( verbose >= 2 )
00267   {
00268     G4cout << " G4ParameterisationBoxY::ComputeDimensions()" << G4endl
00269            << " pDx: " << pDz << G4endl;
00270     box.DumpInfo();
00271   }
00272 #endif
00273 }

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

Implements G4VDivisionParameterisation.

Definition at line 218 of file G4ParameterisationBox.cc.

References FatalException, G4VDivisionParameterisation::faxis, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::foffset, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4Exception(), G4Box::GetYHalfLength(), kYAxis, G4VPhysicalVolume::SetTranslation(), and G4VDivisionParameterisation::verbose.

00219 {
00220   G4Box* msol = (G4Box*)(fmotherSolid);
00221   G4double mdy = msol->GetYHalfLength();
00222 
00223   //----- translation 
00224   G4ThreeVector origin(0.,0.,0.); 
00225   G4double posi = -mdy + foffset + (copyNo+0.5)*fwidth;
00226   if( faxis == kYAxis )
00227   {
00228     origin.setY( posi ); 
00229   }
00230   else
00231   {
00232     std::ostringstream message;
00233     message << "Only axes along Y are allowed !  Axis: " << faxis;
00234     G4Exception("G4ParameterisationBoxY::ComputeTransformation()",
00235                 "GeomDiv0002", FatalException, message);
00236   }
00237 #ifdef G4DIVDEBUG
00238   if( verbose >= 2 )
00239   {
00240     G4cout << std::setprecision(8) << " G4ParameterisationBoxY: "
00241            << copyNo << G4endl
00242            << " Position " << origin << " Axis " << faxis << G4endl;
00243   }
00244 #endif
00245    //----- set translation 
00246   physVol->SetTranslation( origin );
00247 }

G4double G4ParameterisationBoxY::GetMaxParameter (  )  const [virtual]

Implements G4VDivisionParameterisation.

Definition at line 209 of file G4ParameterisationBox.cc.

References G4VDivisionParameterisation::fmotherSolid, and G4Box::GetYHalfLength().

00210 {
00211   G4Box* msol = (G4Box*)(fmotherSolid);
00212   return 2*msol->GetYHalfLength();
00213 }


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