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

#include <G4ParameterisationPolyhedra.hh>

Inheritance diagram for G4ParameterisationPolyhedraRho:
Collaboration diagram for G4ParameterisationPolyhedraRho:

Public Member Functions

 G4ParameterisationPolyhedraRho (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
 
 ~G4ParameterisationPolyhedraRho ()
 
void CheckParametersValidity ()
 
G4double GetMaxParameter () const
 
void ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const
 
void ComputeDimensions (G4Polyhedra &phedra, const G4int copyNo, const G4VPhysicalVolume *physVol) const
 
- Public Member Functions inherited from G4VParameterisationPolyhedra
 G4VParameterisationPolyhedra (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
 
virtual ~G4VParameterisationPolyhedra ()
 
- 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
 
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 86 of file G4ParameterisationPolyhedra.hh.

Constructor & Destructor Documentation

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

Definition at line 137 of file G4ParameterisationPolyhedra.cc.

140  : G4VParameterisationPolyhedra( axis, nDiv, width, offset, msolid, divType )
141 {
142 
144  SetType( "DivisionPolyhedraRho" );
145 
147  G4PolyhedraHistorical* original_pars = msol->GetOriginalParameters();
148 
149  if( divType == DivWIDTH )
150  {
151  fnDiv = CalculateNDiv( original_pars->Rmax[0]
152  - original_pars->Rmin[0], width, offset );
153  }
154  else if( divType == DivNDIV )
155  {
156  fwidth = CalculateWidth( original_pars->Rmax[0]
157  - original_pars->Rmin[0], nDiv, offset );
158  }
159 
160 #ifdef G4DIVDEBUG
161  if( verbose >= 1 )
162  {
163  G4cout << " G4ParameterisationPolyhedraRho - # divisions " << fnDiv
164  << " = " << nDiv << G4endl
165  << " Offset " << foffset << " = " << offset << G4endl
166  << " Width " << fwidth << " = " << width << G4endl;
167  }
168 #endif
169 }
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4VParameterisationPolyhedra(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
G4GLOB_DLL std::ostream G4cout
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const
G4PolyhedraHistorical * GetOriginalParameters() const
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

G4ParameterisationPolyhedraRho::~G4ParameterisationPolyhedraRho ( )

Definition at line 172 of file G4ParameterisationPolyhedra.cc.

173 {
174 }

Member Function Documentation

void G4ParameterisationPolyhedraRho::CheckParametersValidity ( )
virtual

Reimplemented from G4VDivisionParameterisation.

Definition at line 177 of file G4ParameterisationPolyhedra.cc.

178 {
180 
182 
184  {
185  std::ostringstream message;
186  message << "In solid " << msol->GetName() << G4endl
187  << "Division along R will be done with a width "
188  << "different for each solid section." << G4endl
189  << "WIDTH will not be used !";
190  G4Exception("G4ParameterisationPolyhedraRho::CheckParametersValidity()",
191  "GeomDiv1001", JustWarning, message);
192  }
193  if( foffset != 0. )
194  {
195  std::ostringstream message;
196  message << "In solid " << msol->GetName() << G4endl
197  << "Division along R will be done with a width "
198  << "different for each solid section." << G4endl
199  << "OFFSET will not be used !";
200  G4Exception("G4ParameterisationPolyhedraRho::CheckParametersValidity()",
201  "GeomDiv1001", JustWarning, message);
202  }
203 }
G4String GetName() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

Here is the caller graph for this function:

void G4ParameterisationPolyhedraRho::ComputeDimensions ( G4Polyhedra phedra,
const G4int  copyNo,
const G4VPhysicalVolume physVol 
) const
virtual

Reimplemented from G4VPVParameterisation.

Definition at line 252 of file G4ParameterisationPolyhedra.cc.

254 {
256 
257  G4PolyhedraHistorical* origparamMother = msol->GetOriginalParameters();
258  G4PolyhedraHistorical origparam( *origparamMother );
259  G4int nZplanes = origparamMother->Num_z_planes;
260 
261  G4double width = 0.;
262  for( G4int ii = 0; ii < nZplanes; ii++ )
263  {
264  width = CalculateWidth( origparamMother->Rmax[ii]
265  - origparamMother->Rmin[ii], fnDiv, foffset );
266  origparam.Rmin[ii] = origparamMother->Rmin[ii]+foffset+width*copyNo;
267  origparam.Rmax[ii] = origparamMother->Rmin[ii]+foffset+width*(copyNo+1);
268  }
269 
270  phedra.SetOriginalParameters(&origparam); // copy values & transfer pointers
271  phedra.Reset(); // reset to new solid parameters
272 
273 #ifdef G4DIVDEBUG
274  if( verbose >= -2 )
275  {
276  G4cout << "G4ParameterisationPolyhedraRho::ComputeDimensions()" << G4endl
277  << "-- Parametrised phedra copy-number: " << copyNo << G4endl;
278  phedra.DumpInfo();
279  }
280 #endif
281 }
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4bool Reset()
Definition: G4Polyhedra.cc:486
int G4int
Definition: G4Types.hh:78
void DumpInfo() const
G4GLOB_DLL std::ostream G4cout
G4PolyhedraHistorical * GetOriginalParameters() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void SetOriginalParameters(G4PolyhedraHistorical *pars)

Here is the call graph for this function:

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

Implements G4VDivisionParameterisation.

Definition at line 216 of file G4ParameterisationPolyhedra.cc.

217 {
218  //----- translation
219  G4ThreeVector origin(0.,0.,0.);
220 
221  //----- set translation
222  physVol->SetTranslation( origin );
223 
224  //----- calculate rotation matrix: unit
225 
226 #ifdef G4DIVDEBUG
227  if( verbose >= 2 )
228  {
229  G4cout << " G4ParameterisationPolyhedraRho " << G4endl
230  << " foffset: " << foffset/deg
231  << " - fwidth: " << fwidth/deg << G4endl;
232  }
233 #endif
234 
235  ChangeRotMatrix( physVol );
236 
237 #ifdef G4DIVDEBUG
238  if( verbose >= 2 )
239  {
240  G4cout << std::setprecision(8) << " G4ParameterisationPolyhedraRho "
241  << G4endl
242  << " Position: " << origin
243  << " - Width: " << fwidth
244  << " - Axis: " << faxis << G4endl;
245  }
246 #endif
247 }
G4GLOB_DLL std::ostream G4cout
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.) const
void SetTranslation(const G4ThreeVector &v)
#define G4endl
Definition: G4ios.hh:61
static constexpr double deg
Definition: G4SIunits.hh:152

Here is the call graph for this function:

G4double G4ParameterisationPolyhedraRho::GetMaxParameter ( ) const
virtual

Implements G4VDivisionParameterisation.

Definition at line 206 of file G4ParameterisationPolyhedra.cc.

207 {
209  G4PolyhedraHistorical* original_pars = msol->GetOriginalParameters();
210  return original_pars->Rmax[0] - original_pars->Rmin[0];
211 }
G4PolyhedraHistorical * GetOriginalParameters() const

Here is the call graph for this function:


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