Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
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
 
- 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 120 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(), G4VDivisionParameterisation::verbose, and width.

178  : G4VParameterisationBox( axis, nDiv, width, offset, msolid, divType )
179 {
181  SetType( "DivisionBoxY" );
182 
183  G4Box* mbox = (G4Box*)(fmotherSolid);
184  if( divType == DivWIDTH )
185  {
186  fnDiv = CalculateNDiv( 2*mbox->GetYHalfLength(), width, offset );
187  }
188  else if( divType == DivNDIV )
189  {
190  fwidth = CalculateWidth( 2*mbox->GetYHalfLength(), nDiv, offset );
191  }
192 
193 #ifdef G4DIVDEBUG
194  if( verbose >= 1 )
195  {
196  G4cout << " G4ParameterisationBoxY - no divisions " << fnDiv << " = "
197  << nDiv << ". Offset " << foffset << " = " << offset
198  << ". Width " << fwidth << " = " << width << G4endl;
199  }
200 #endif
201 }
void SetType(const G4String &type)
Definition: G4Box.hh:63
#define width
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4GLOB_DLL std::ostream G4cout
G4double GetYHalfLength() const
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
G4ParameterisationBoxY::~G4ParameterisationBoxY ( )

Definition at line 204 of file G4ParameterisationBox.cc.

205 {
206 }

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.

254 {
255  G4Box* msol = (G4Box*)(fmotherSolid);
256 
257  G4double pDx = msol->GetXHalfLength();
258  G4double pDy = fwidth/2. - fhgap;
259  G4double pDz = msol->GetZHalfLength();
260 
261  box.SetXHalfLength( pDx );
262  box.SetYHalfLength( pDy );
263  box.SetZHalfLength( pDz );
264 
265 #ifdef G4DIVDEBUG
266  if( verbose >= 2 )
267  {
268  G4cout << " G4ParameterisationBoxY::ComputeDimensions()" << G4endl
269  << " pDx: " << pDz << G4endl;
270  box.DumpInfo();
271  }
272 #endif
273 }
G4double GetXHalfLength() const
void SetZHalfLength(G4double dz)
Definition: G4Box.cc:172
Definition: G4Box.hh:63
G4double GetZHalfLength() const
void DumpInfo() const
G4GLOB_DLL std::ostream G4cout
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 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(), CLHEP::Hep3Vector::setY(), and G4VDivisionParameterisation::verbose.

219 {
220  G4Box* msol = (G4Box*)(fmotherSolid);
221  G4double mdy = msol->GetYHalfLength();
222 
223  //----- translation
224  G4ThreeVector origin(0.,0.,0.);
225  G4double posi = -mdy + foffset + (copyNo+0.5)*fwidth;
226  if( faxis == kYAxis )
227  {
228  origin.setY( posi );
229  }
230  else
231  {
232  std::ostringstream message;
233  message << "Only axes along Y are allowed ! Axis: " << faxis;
234  G4Exception("G4ParameterisationBoxY::ComputeTransformation()",
235  "GeomDiv0002", FatalException, message);
236  }
237 #ifdef G4DIVDEBUG
238  if( verbose >= 2 )
239  {
240  G4cout << std::setprecision(8) << " G4ParameterisationBoxY: "
241  << copyNo << G4endl
242  << " Position " << origin << " Axis " << faxis << G4endl;
243  }
244 #endif
245  //----- set translation
246  physVol->SetTranslation( origin );
247 }
Definition: G4Box.hh:63
G4GLOB_DLL std::ostream G4cout
void SetTranslation(const G4ThreeVector &v)
G4double GetYHalfLength() const
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 G4ParameterisationBoxY::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 209 of file G4ParameterisationBox.cc.

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

210 {
211  G4Box* msol = (G4Box*)(fmotherSolid);
212  return 2*msol->GetYHalfLength();
213 }
Definition: G4Box.hh:63
G4double GetYHalfLength() const

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