Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ParameterisationBoxX Class Reference

#include <G4ParameterisationBox.hh>

Inheritance diagram for G4ParameterisationBoxX:
Collaboration diagram for G4ParameterisationBoxX:

Public Member Functions

 G4ParameterisationBoxX (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
 ~G4ParameterisationBoxX ()
 
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
G4RotationMatrix
fRot = 0
 
static const G4int verbose = 5
 

Detailed Description

Definition at line 76 of file G4ParameterisationBox.hh.

Constructor & Destructor Documentation

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

Definition at line 72 of file G4ParameterisationBox.cc.

75  : G4VParameterisationBox( axis, nDiv, width, offset, msolid, divType )
76 {
78  SetType( "DivisionBoxX" );
79 
80  G4Box* mbox = (G4Box*)(fmotherSolid);
81  if( divType == DivWIDTH )
82  {
83  fnDiv = CalculateNDiv( 2*mbox->GetXHalfLength(), width, offset );
84  }
85  else if( divType == DivNDIV )
86  {
87  fwidth = CalculateWidth( 2*mbox->GetXHalfLength(), nDiv, offset );
88  }
89 #ifdef G4DIVDEBUG
90  if( verbose >= 1 )
91  {
92  G4cout << " G4ParameterisationBoxX - no divisions "
93  << fnDiv << " = " << nDiv << G4endl
94  << " Offset " << foffset << " = " << offset << G4endl
95  << " Width " << fwidth << " = " << width << G4endl;
96  }
97 #endif
98 }
G4double GetXHalfLength() const
void SetType(const G4String &type)
Definition: G4Box.hh:64
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) 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

Here is the call graph for this function:

G4ParameterisationBoxX::~G4ParameterisationBoxX ( )

Definition at line 101 of file G4ParameterisationBox.cc.

102 {
103 }

Member Function Documentation

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

Reimplemented from G4VPVParameterisation.

Definition at line 150 of file G4ParameterisationBox.cc.

152 {
153  G4Box* msol = (G4Box*)(fmotherSolid);
154 
155  G4double pDx = fwidth/2. - fhgap;
156  G4double pDy = msol->GetYHalfLength();
157  G4double pDz = msol->GetZHalfLength();
158 
159  box.SetXHalfLength( pDx );
160  box.SetYHalfLength( pDy );
161  box.SetZHalfLength( pDz );
162 
163 #ifdef G4DIVDEBUG
164  if( verbose >= 2 )
165  {
166  G4cout << " G4ParameterisationBoxX::ComputeDimensions()" << G4endl
167  << " pDx: " << pDz << G4endl;
168  box.DumpInfo();
169  }
170 #endif
171 }
void SetZHalfLength(G4double dz)
Definition: G4Box.cc:174
Definition: G4Box.hh:64
G4double GetZHalfLength() const
void DumpInfo() const
G4GLOB_DLL std::ostream G4cout
G4double GetYHalfLength() const
void SetYHalfLength(G4double dy)
Definition: G4Box.cc:154
void SetXHalfLength(G4double dx)
Definition: G4Box.cc:134
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Implements G4VDivisionParameterisation.

Definition at line 115 of file G4ParameterisationBox.cc.

116 {
117  G4Box* msol = (G4Box*)(fmotherSolid );
118  G4double mdx = msol->GetXHalfLength( );
119 
120  //----- translation
121  G4ThreeVector origin(0.,0.,0.);
122  G4double posi = -mdx + foffset+(copyNo+0.5)*fwidth;
123 
124  if( faxis == kXAxis )
125  {
126  origin.setX( posi );
127  }
128  else
129  {
130  std::ostringstream message;
131  message << "Only axes along X are allowed ! Axis: " << faxis;
132  G4Exception("G4ParameterisationBoxX::ComputeTransformation()",
133  "GeomDiv0002", FatalException, message);
134  }
135 #ifdef G4DIVDEBUG
136  if( verbose >= 2 )
137  {
138  G4cout << std::setprecision(8) << " G4ParameterisationBoxX: "
139  << copyNo << G4endl
140  << " Position " << origin << " Axis " << faxis << G4endl;
141  }
142 #endif
143  //----- set translation
144  physVol->SetTranslation( origin );
145 }
G4double GetXHalfLength() const
Definition: G4Box.hh:64
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

Here is the call graph for this function:

G4double G4ParameterisationBoxX::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 106 of file G4ParameterisationBox.cc.

107 {
108  G4Box* msol = (G4Box*)(fmotherSolid);
109  return 2*msol->GetXHalfLength();
110 }
G4double GetXHalfLength() const
Definition: G4Box.hh:64

Here is the call graph for this function:


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