Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4ParameterisationBoxZ Class Reference

#include <G4ParameterisationBox.hh>

Inheritance diagram for G4ParameterisationBoxZ:
G4VParameterisationBox G4VDivisionParameterisation G4VPVParameterisation

Public Member Functions

 G4ParameterisationBoxZ (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
 ~G4ParameterisationBoxZ ()
 
G4double GetMaxParameter () const
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const
 
void ComputeDimensions (G4Box &box, const G4int copyNo, const G4VPhysicalVolume *physVol) const
 
- Public Member Functions inherited from G4VParameterisationBox
 G4VParameterisationBox (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
virtual ~G4VParameterisationBox ()
 
- 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 G4int verbose = 5
 

Detailed Description

Definition at line 165 of file G4ParameterisationBox.hh.

Constructor & Destructor Documentation

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

Definition at line 277 of file G4ParameterisationBox.cc.

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

280  : G4VParameterisationBox( axis, nDiv, width, offset, msolid, divType )
281 {
283  SetType( "DivisionBoxZ" );
284 
285  G4Box* mbox = (G4Box*)(fmotherSolid);
286  if( divType == DivWIDTH )
287  {
288  fnDiv = CalculateNDiv( 2*mbox->GetZHalfLength(), width, offset );
289  }
290  else if ( divType == DivNDIV )
291  {
292  fwidth = CalculateWidth( 2*mbox->GetZHalfLength(), nDiv, offset );
293  }
294 #ifdef G4DIVDEBUG
295  if( verbose >= 1 )
296  {
297  G4cout << " G4ParameterisationBoxZ - no divisions " << fnDiv << " = "
298  << nDiv << ". Offset " << foffset << " = " << offset
299  << ". Width " << fwidth << " = " << width << G4endl;
300  }
301 #endif
302 }
void SetType(const G4String &type)
Definition: G4Box.hh:63
#define width
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4double GetZHalfLength() const
G4GLOB_DLL std::ostream G4cout
G4VParameterisationBox(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const
#define G4endl
Definition: G4ios.hh:61
G4ParameterisationBoxZ::~G4ParameterisationBoxZ ( )

Definition at line 305 of file G4ParameterisationBox.cc.

306 {
307 }

Member Function Documentation

void G4ParameterisationBoxZ::ComputeDimensions ( G4Box box,
const G4int  copyNo,
const G4VPhysicalVolume physVol 
) const
virtual

Reimplemented from G4VPVParameterisation.

Definition at line 354 of file G4ParameterisationBox.cc.

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

356 {
357  G4Box* msol = (G4Box*)(fmotherSolid);
358 
359  G4double pDx = msol->GetXHalfLength();
360  G4double pDy = msol->GetYHalfLength();
361  G4double pDz = fwidth/2. - fhgap;
362 
363  box.SetXHalfLength( pDx );
364  box.SetYHalfLength( pDy );
365  box.SetZHalfLength( pDz );
366 
367 #ifdef G4DIVDEBUG
368  if( verbose >= 2 )
369  {
370  G4cout << " G4ParameterisationBoxZ::ComputeDimensions()" << G4endl
371  << " pDx: " << pDz << G4endl;
372  box.DumpInfo();
373  }
374 #endif
375 }
G4double GetXHalfLength() const
void SetZHalfLength(G4double dz)
Definition: G4Box.cc:172
Definition: G4Box.hh:63
void DumpInfo() const
G4GLOB_DLL std::ostream G4cout
G4double GetYHalfLength() const
void SetYHalfLength(G4double dy)
Definition: G4Box.cc:152
void SetXHalfLength(G4double dx)
Definition: G4Box.cc:132
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void G4ParameterisationBoxZ::ComputeTransformation ( const G4int  copyNo,
G4VPhysicalVolume physVol 
) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 319 of file G4ParameterisationBox.cc.

References FatalException, G4VDivisionParameterisation::faxis, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4Exception(), G4Box::GetZHalfLength(), kZAxis, G4VDivisionParameterisation::OffsetZ(), G4VPhysicalVolume::SetTranslation(), CLHEP::Hep3Vector::setZ(), and G4VDivisionParameterisation::verbose.

320 {
321  G4Box* msol = (G4Box*)(fmotherSolid );
322  G4double mdz = msol->GetZHalfLength();
323 
324  //----- translation
325  G4ThreeVector origin(0.,0.,0.);
326  G4double posi = -mdz + OffsetZ() + (copyNo+0.5)*fwidth;
327 
328  if( faxis == kZAxis )
329  {
330  origin.setZ( posi );
331  }
332  else
333  {
334  std::ostringstream message;
335  message << "Only axes along Z are allowed ! Axis: " << faxis;
336  G4Exception("G4ParameterisationBoxZ::ComputeTransformation()",
337  "GeomDiv0002", FatalException, message);
338  }
339 #ifdef G4DIVDEBUG
340  if( verbose >= 2 )
341  {
342  G4cout << std::setprecision(8) << " G4ParameterisationBoxZ: "
343  << copyNo << G4endl
344  << " Position " << origin << " Axis " << faxis << G4endl;
345  }
346 #endif
347  //----- set translation
348  physVol->SetTranslation( origin );
349 }
Definition: G4Box.hh:63
G4double GetZHalfLength() const
G4GLOB_DLL std::ostream G4cout
void SetTranslation(const G4ThreeVector &v)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double G4ParameterisationBoxZ::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 310 of file G4ParameterisationBox.cc.

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

311 {
312  G4Box* msol = (G4Box*)(fmotherSolid);
313  return 2*msol->GetZHalfLength();
314 }
Definition: G4Box.hh:63
G4double GetZHalfLength() const

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