Geant4  10.02.p03
G4ReplicatedSlice Class Reference

#include <G4ReplicatedSlice.hh>

Inheritance diagram for G4ReplicatedSlice:
Collaboration diagram for G4ReplicatedSlice:

Public Member Functions

 G4ReplicatedSlice (const G4String &pName, G4LogicalVolume *pLogical, G4LogicalVolume *pMotherLogical, const EAxis pAxis, const G4int nReplicas, const G4double width, const G4double half_gap, const G4double offset)
 
 G4ReplicatedSlice (const G4String &pName, G4LogicalVolume *pLogical, G4LogicalVolume *pMotherLogical, const EAxis pAxis, const G4int nReplicas, const G4double half_gap, const G4double offset)
 
 G4ReplicatedSlice (const G4String &pName, G4LogicalVolume *pLogical, G4LogicalVolume *pMotherLogical, const EAxis pAxis, const G4double width, const G4double half_gap, const G4double offset)
 
 G4ReplicatedSlice (const G4String &pName, G4LogicalVolume *pLogical, G4VPhysicalVolume *pMotherPhysical, const EAxis pAxis, const G4int nReplicas, const G4double width, const G4double half_gap, const G4double offset)
 
 G4ReplicatedSlice (const G4String &pName, G4LogicalVolume *pLogical, G4VPhysicalVolume *pMotherPhysical, const EAxis pAxis, const G4int nReplicas, const G4double half_gap, const G4double offset)
 
 G4ReplicatedSlice (const G4String &pName, G4LogicalVolume *pLogical, G4VPhysicalVolume *pMotherPhysical, const EAxis pAxis, const G4double width, const G4double half_gap, const G4double offset)
 
virtual ~G4ReplicatedSlice ()
 
virtual G4bool IsMany () const
 
virtual G4int GetCopyNo () const
 
virtual void SetCopyNo (G4int CopyNo)
 
virtual G4bool IsReplicated () const
 
virtual G4VPVParameterisationGetParameterisation () const
 
virtual void GetReplicationData (EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const
 
EAxis GetDivisionAxis () const
 
G4bool IsParameterised () const
 
G4bool IsRegularStructure () const
 
G4int GetRegularStructureId () const
 
- Public Member Functions inherited from G4VPhysicalVolume
 G4VPhysicalVolume (G4RotationMatrix *pRot, const G4ThreeVector &tlate, const G4String &pName, G4LogicalVolume *pLogical, G4VPhysicalVolume *pMother)
 
virtual ~G4VPhysicalVolume ()
 
G4bool operator== (const G4VPhysicalVolume &p) const
 
G4RotationMatrixGetObjectRotation () const
 
G4RotationMatrix GetObjectRotationValue () const
 
G4ThreeVector GetObjectTranslation () const
 
const G4RotationMatrixGetFrameRotation () const
 
G4ThreeVector GetFrameTranslation () const
 
const G4ThreeVectorGetTranslation () const
 
const G4RotationMatrixGetRotation () const
 
void SetTranslation (const G4ThreeVector &v)
 
G4RotationMatrixGetRotation ()
 
void SetRotation (G4RotationMatrix *)
 
G4LogicalVolumeGetLogicalVolume () const
 
void SetLogicalVolume (G4LogicalVolume *pLogical)
 
G4LogicalVolumeGetMotherLogical () const
 
void SetMotherLogical (G4LogicalVolume *pMother)
 
const G4StringGetName () const
 
void SetName (const G4String &pName)
 
EVolume VolumeType () const
 
virtual G4int GetMultiplicity () const
 
virtual G4bool CheckOverlaps (G4int res=1000, G4double tol=0., G4bool verbose=true, G4int errMax=1)
 
 G4VPhysicalVolume (__void__ &)
 
G4int GetInstanceID () const
 

Protected Attributes

EAxis faxis
 
EAxis fdivAxis
 
G4int fnReplicas
 
G4double fwidth
 
G4double foffset
 
G4int fcopyNo
 
G4VDivisionParameterisationfparam
 
- Protected Attributes inherited from G4VPhysicalVolume
G4int instanceID
 

Private Member Functions

void CheckAndSetParameters (const EAxis pAxis, const G4int nDivs, const G4double width, const G4double half_gap, const G4double offset, DivisionType divType, G4LogicalVolume *pMotherLogical, const G4LogicalVolume *pLogical)
 
void SetParameterisation (G4LogicalVolume *motherLogical, const EAxis pAxis, const G4int nReplicas, const G4double width, const G4double half_gap, const G4double offset, DivisionType divType)
 
void ErrorInAxis (EAxis axis, G4VSolid *solid)
 
 G4ReplicatedSlice (const G4ReplicatedSlice &)
 
const G4ReplicatedSliceoperator= (const G4ReplicatedSlice &)
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VPhysicalVolume
static const G4PVManagerGetSubInstanceManager ()
 
- Protected Member Functions inherited from G4VPhysicalVolume
void InitialiseWorker (G4VPhysicalVolume *pMasterObject, G4RotationMatrix *pRot, const G4ThreeVector &tlate)
 
void TerminateWorker (G4VPhysicalVolume *pMasterObject)
 
- Static Protected Attributes inherited from G4VPhysicalVolume
static G4GEOM_DLL G4PVManager subInstanceManager
 

Detailed Description

Definition at line 71 of file G4ReplicatedSlice.hh.

Constructor & Destructor Documentation

◆ G4ReplicatedSlice() [1/7]

G4ReplicatedSlice::G4ReplicatedSlice ( const G4String pName,
G4LogicalVolume pLogical,
G4LogicalVolume pMotherLogical,
const EAxis  pAxis,
const G4int  nReplicas,
const G4double  width,
const G4double  half_gap,
const G4double  offset 
)

Definition at line 44 of file G4ReplicatedSlice.cc.

52  : G4VPhysicalVolume(0,G4ThreeVector(),pName,pLogical,0), fcopyNo(-1)
53 {
54  CheckAndSetParameters(pAxis, nDivs, width, half_gap, offset,
55  DivNDIVandWIDTH, pMotherLogical, pLogical);
56 }
G4VPhysicalVolume(G4RotationMatrix *pRot, const G4ThreeVector &tlate, const G4String &pName, G4LogicalVolume *pLogical, G4VPhysicalVolume *pMother)
CLHEP::Hep3Vector G4ThreeVector
#define width
void CheckAndSetParameters(const EAxis pAxis, const G4int nDivs, const G4double width, const G4double half_gap, const G4double offset, DivisionType divType, G4LogicalVolume *pMotherLogical, const G4LogicalVolume *pLogical)
Here is the call graph for this function:

◆ G4ReplicatedSlice() [2/7]

G4ReplicatedSlice::G4ReplicatedSlice ( const G4String pName,
G4LogicalVolume pLogical,
G4LogicalVolume pMotherLogical,
const EAxis  pAxis,
const G4int  nReplicas,
const G4double  half_gap,
const G4double  offset 
)

Definition at line 59 of file G4ReplicatedSlice.cc.

66  : G4VPhysicalVolume(0,G4ThreeVector(),pName,pLogical,0), fcopyNo(-1)
67 {
68  CheckAndSetParameters(pAxis, nDivs, 0., half_gap, offset,
69  DivNDIV, pMotherLogical, pLogical);
70 }
G4VPhysicalVolume(G4RotationMatrix *pRot, const G4ThreeVector &tlate, const G4String &pName, G4LogicalVolume *pLogical, G4VPhysicalVolume *pMother)
CLHEP::Hep3Vector G4ThreeVector
void CheckAndSetParameters(const EAxis pAxis, const G4int nDivs, const G4double width, const G4double half_gap, const G4double offset, DivisionType divType, G4LogicalVolume *pMotherLogical, const G4LogicalVolume *pLogical)
Here is the call graph for this function:

◆ G4ReplicatedSlice() [3/7]

G4ReplicatedSlice::G4ReplicatedSlice ( const G4String pName,
G4LogicalVolume pLogical,
G4LogicalVolume pMotherLogical,
const EAxis  pAxis,
const G4double  width,
const G4double  half_gap,
const G4double  offset 
)

Definition at line 73 of file G4ReplicatedSlice.cc.

80  : G4VPhysicalVolume(0,G4ThreeVector(),pName,pLogical,0), fcopyNo(-1)
81 {
82  CheckAndSetParameters(pAxis, 0, width, half_gap, offset,
83  DivWIDTH, pMotherLogical, pLogical);
84 }
G4VPhysicalVolume(G4RotationMatrix *pRot, const G4ThreeVector &tlate, const G4String &pName, G4LogicalVolume *pLogical, G4VPhysicalVolume *pMother)
CLHEP::Hep3Vector G4ThreeVector
#define width
void CheckAndSetParameters(const EAxis pAxis, const G4int nDivs, const G4double width, const G4double half_gap, const G4double offset, DivisionType divType, G4LogicalVolume *pMotherLogical, const G4LogicalVolume *pLogical)
Here is the call graph for this function:

◆ G4ReplicatedSlice() [4/7]

G4ReplicatedSlice::G4ReplicatedSlice ( const G4String pName,
G4LogicalVolume pLogical,
G4VPhysicalVolume pMotherPhysical,
const EAxis  pAxis,
const G4int  nReplicas,
const G4double  width,
const G4double  half_gap,
const G4double  offset 
)

Definition at line 87 of file G4ReplicatedSlice.cc.

95  : G4VPhysicalVolume(0,G4ThreeVector(),pName,pLogical,0), fcopyNo(-1)
96 {
97  CheckAndSetParameters(pAxis, nDivs, width, half_gap, offset,
98  DivNDIVandWIDTH, pMotherPhysical->GetLogicalVolume(), pLogical);
99 }
G4VPhysicalVolume(G4RotationMatrix *pRot, const G4ThreeVector &tlate, const G4String &pName, G4LogicalVolume *pLogical, G4VPhysicalVolume *pMother)
CLHEP::Hep3Vector G4ThreeVector
#define width
void CheckAndSetParameters(const EAxis pAxis, const G4int nDivs, const G4double width, const G4double half_gap, const G4double offset, DivisionType divType, G4LogicalVolume *pMotherLogical, const G4LogicalVolume *pLogical)
G4LogicalVolume * GetLogicalVolume() const
Here is the call graph for this function:

◆ G4ReplicatedSlice() [5/7]

G4ReplicatedSlice::G4ReplicatedSlice ( const G4String pName,
G4LogicalVolume pLogical,
G4VPhysicalVolume pMotherPhysical,
const EAxis  pAxis,
const G4int  nReplicas,
const G4double  half_gap,
const G4double  offset 
)

Definition at line 102 of file G4ReplicatedSlice.cc.

109  : G4VPhysicalVolume(0,G4ThreeVector(),pName,pLogical,0), fcopyNo(-1)
110 {
111  CheckAndSetParameters(pAxis, nDivs, 0., half_gap, offset,
112  DivNDIV, pMotherPhysical->GetLogicalVolume(), pLogical);
113 }
G4VPhysicalVolume(G4RotationMatrix *pRot, const G4ThreeVector &tlate, const G4String &pName, G4LogicalVolume *pLogical, G4VPhysicalVolume *pMother)
CLHEP::Hep3Vector G4ThreeVector
void CheckAndSetParameters(const EAxis pAxis, const G4int nDivs, const G4double width, const G4double half_gap, const G4double offset, DivisionType divType, G4LogicalVolume *pMotherLogical, const G4LogicalVolume *pLogical)
G4LogicalVolume * GetLogicalVolume() const
Here is the call graph for this function:

◆ G4ReplicatedSlice() [6/7]

G4ReplicatedSlice::G4ReplicatedSlice ( const G4String pName,
G4LogicalVolume pLogical,
G4VPhysicalVolume pMotherPhysical,
const EAxis  pAxis,
const G4double  width,
const G4double  half_gap,
const G4double  offset 
)

Definition at line 116 of file G4ReplicatedSlice.cc.

123  : G4VPhysicalVolume(0,G4ThreeVector(),pName,pLogical,0), fcopyNo(-1)
124 {
125  CheckAndSetParameters(pAxis, 0, width, half_gap, offset,
126  DivWIDTH, pMotherPhysical->GetLogicalVolume(), pLogical);
127 }
G4VPhysicalVolume(G4RotationMatrix *pRot, const G4ThreeVector &tlate, const G4String &pName, G4LogicalVolume *pLogical, G4VPhysicalVolume *pMother)
CLHEP::Hep3Vector G4ThreeVector
#define width
void CheckAndSetParameters(const EAxis pAxis, const G4int nDivs, const G4double width, const G4double half_gap, const G4double offset, DivisionType divType, G4LogicalVolume *pMotherLogical, const G4LogicalVolume *pLogical)
G4LogicalVolume * GetLogicalVolume() const
Here is the call graph for this function:

◆ ~G4ReplicatedSlice()

G4ReplicatedSlice::~G4ReplicatedSlice ( )
virtual

Definition at line 249 of file G4ReplicatedSlice.cc.

250 {
251  delete GetRotation();
252 }
const G4RotationMatrix * GetRotation() const
Here is the call graph for this function:

◆ G4ReplicatedSlice() [7/7]

G4ReplicatedSlice::G4ReplicatedSlice ( const G4ReplicatedSlice )
private

Member Function Documentation

◆ CheckAndSetParameters()

void G4ReplicatedSlice::CheckAndSetParameters ( const EAxis  pAxis,
const G4int  nDivs,
const G4double  width,
const G4double  half_gap,
const G4double  offset,
DivisionType  divType,
G4LogicalVolume pMotherLogical,
const G4LogicalVolume pLogical 
)
private

!!!! axis has to be x/y/z in G4VoxelLimits::GetMinExtent

Definition at line 131 of file G4ReplicatedSlice.cc.

139 {
140  if(!pMotherLogical)
141  {
142  std::ostringstream message;
143  message << "Invalid setup." << G4endl
144  << "NULL pointer specified as mother! Volume: " << GetName();
145  G4Exception("G4ReplicatedSlice::CheckAndSetParameters()", "GeomDiv0002",
146  FatalException, message);
147  }
148  if(pLogical == pMotherLogical)
149  {
150  std::ostringstream message;
151  message << "Invalid setup." << G4endl
152  << "Cannot place a volume inside itself! Volume: " << GetName();
153  G4Exception("G4ReplicatedSlice::CheckAndSetParameters()", "GeomDiv0002",
154  FatalException, message);
155  }
156 
157  //----- Check that mother solid is of the same type as
158  // daughter solid (otherwise, the corresponding
159  // Parameterisation::ComputeDimension() will not be called)
160  //
161  G4String msolType = pMotherLogical->GetSolid()->GetEntityType();
162  G4String dsolType = pLogical->GetSolid()->GetEntityType();
163  if( msolType != dsolType && ( msolType != "G4Trd" || dsolType != "G4Trap" ) )
164  {
165  std::ostringstream message;
166  message << "Invalid setup." << G4endl
167  << "Incorrect solid type for division of volume: "
168  << GetName() << G4endl
169  << " It is: " << msolType
170  << ", while it should be: " << dsolType;
171  G4Exception("G4ReplicatedSlice::CheckAndSetParameters()",
172  "GeomDiv0002", FatalException, message);
173  }
174 
175  pMotherLogical->AddDaughter(this);
176  SetMotherLogical(pMotherLogical);
177  SetParameterisation(pMotherLogical, pAxis, nDivs,
178  width, half_gap, offset, divType);
179 
180  if( divType == DivWIDTH )
181  {
183  }
184  else
185  {
186  fnReplicas = nDivs;
187  }
188  if (fnReplicas < 1 )
189  {
190  G4Exception("G4ReplicatedSlice::CheckAndSetParameters()", "GeomDiv0002",
191  FatalException, "Illegal number of replicas!");
192  }
193  if( divType != DivNDIV)
194  {
195  fwidth = fparam->GetWidth();
196  }
197  else
198  {
199  fwidth = width;
200  }
201  if( fwidth < 0 )
202  {
203  G4Exception("G4ReplicatedSlice::CheckAndSetParameters()", "GeomDiv0002",
204  FatalException, "Width must be positive!");
205  }
206  if( fwidth < 2.*half_gap )
207  {
208  G4Exception("G4ReplicatedSlice::CheckAndSetParameters()", "GeomDiv0002",
209  FatalException, "Half_gap is too large!");
210  }
211 
212  foffset = offset;
213  fdivAxis = pAxis;
214 
216  //
217  if( pAxis == kRho || pAxis == kRadial3D || pAxis == kPhi )
218  {
219  faxis = kZAxis;
220  }
221  else
222  {
223  faxis = pAxis;
224  }
225 
226  switch (faxis)
227  {
228  case kPhi:
229  case kRho:
230  case kXAxis:
231  case kYAxis:
232  case kZAxis:
233  break;
234  default:
235  G4Exception("G4ReplicatedSlice::CheckAndSetParameters()", "GeomDiv0002",
236  FatalException, "Unknown axis of replication.");
237  break;
238  }
239 
240  // Create rotation matrix: for phi axis it will be changed
241  // in G4VPVParameterisation::ComputeTransformation, for others
242  // it will stay the unity
243  //
244  G4RotationMatrix *pRMat = new G4RotationMatrix();
245  SetRotation(pRMat);
246 }
Definition: geomdefs.hh:54
CLHEP::HepRotation G4RotationMatrix
G4double GetWidth() const
#define width
virtual G4GeometryType GetEntityType() const =0
void SetRotation(G4RotationMatrix *)
G4VDivisionParameterisation * fparam
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
const G4String & GetName() const
void SetParameterisation(G4LogicalVolume *motherLogical, const EAxis pAxis, const G4int nReplicas, const G4double width, const G4double half_gap, const G4double offset, DivisionType divType)
#define G4endl
Definition: G4ios.hh:61
void SetMotherLogical(G4LogicalVolume *pMother)
Definition: geomdefs.hh:54
void AddDaughter(G4VPhysicalVolume *p)
G4VSolid * GetSolid() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ErrorInAxis()

void G4ReplicatedSlice::ErrorInAxis ( EAxis  axis,
G4VSolid solid 
)
private

Definition at line 498 of file G4ReplicatedSlice.cc.

499 {
500  G4String error = "Trying to divide solid " + solid->GetName()
501  + " of type " + solid->GetEntityType() + " along axis ";
502  switch( axis )
503  {
504  case kXAxis:
505  error += "X.";
506  break;
507  case kYAxis:
508  error += "Y.";
509  break;
510  case kZAxis:
511  error += "Z.";
512  break;
513  case kRho:
514  error += "Rho.";
515  break;
516  case kRadial3D:
517  error += "Radial3D.";
518  break;
519  case kPhi:
520  error += "Phi.";
521  break;
522  default:
523  break;
524  }
525  G4Exception("G4ReplicatedSlice::ErrorInAxis()", "GeomDiv0002",
526  FatalException, error);
527 }
Definition: geomdefs.hh:54
virtual G4GeometryType GetEntityType() const =0
G4String GetName() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static PROLOG_HANDLER error
Definition: xmlrole.cc:112
Definition: geomdefs.hh:54
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetCopyNo()

G4int G4ReplicatedSlice::GetCopyNo ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 273 of file G4ReplicatedSlice.cc.

274 {
275  return fcopyNo;
276 }

◆ GetDivisionAxis()

EAxis G4ReplicatedSlice::GetDivisionAxis ( ) const

Definition at line 255 of file G4ReplicatedSlice.cc.

256 {
257  return fdivAxis;
258 }

◆ GetParameterisation()

G4VPVParameterisation * G4ReplicatedSlice::GetParameterisation ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 291 of file G4ReplicatedSlice.cc.

292 {
293  return fparam;
294 }
G4VDivisionParameterisation * fparam

◆ GetRegularStructureId()

G4int G4ReplicatedSlice::GetRegularStructureId ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 543 of file G4ReplicatedSlice.cc.

544 {
545  return 0;
546 }

◆ GetReplicationData()

void G4ReplicatedSlice::GetReplicationData ( EAxis axis,
G4int nReplicas,
G4double width,
G4double offset,
G4bool consuming 
) const
virtual

Implements G4VPhysicalVolume.

Definition at line 297 of file G4ReplicatedSlice.cc.

302 {
303  axis=faxis;
304  nDivs=fnReplicas;
305  width=fwidth;
306  offset=foffset;
307  consuming=false;
308 }
#define width

◆ IsMany()

G4bool G4ReplicatedSlice::IsMany ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 267 of file G4ReplicatedSlice.cc.

268 {
269  return false;
270 }

◆ IsParameterised()

G4bool G4ReplicatedSlice::IsParameterised ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 261 of file G4ReplicatedSlice.cc.

262 {
263  return true;
264 }

◆ IsRegularStructure()

G4bool G4ReplicatedSlice::IsRegularStructure ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 535 of file G4ReplicatedSlice.cc.

536 {
537  return false;
538 }

◆ IsReplicated()

G4bool G4ReplicatedSlice::IsReplicated ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 285 of file G4ReplicatedSlice.cc.

286 {
287  return true;
288 }

◆ operator=()

const G4ReplicatedSlice& G4ReplicatedSlice::operator= ( const G4ReplicatedSlice )
private

◆ SetCopyNo()

void G4ReplicatedSlice::SetCopyNo ( G4int  CopyNo)
virtual

Implements G4VPhysicalVolume.

Definition at line 279 of file G4ReplicatedSlice.cc.

280 {
281  fcopyNo= newCopyNo;
282 }

◆ SetParameterisation()

void G4ReplicatedSlice::SetParameterisation ( G4LogicalVolume motherLogical,
const EAxis  pAxis,
const G4int  nReplicas,
const G4double  width,
const G4double  half_gap,
const G4double  offset,
DivisionType  divType 
)
private

Definition at line 312 of file G4ReplicatedSlice.cc.

319 {
320  G4VSolid* mSolid = motherLogical->GetSolid();
321  G4String mSolidType = mSolid->GetEntityType();
322  fparam = 0;
323 
324  // If the solid is a reflected one, update type to its
325  // real constituent solid.
326  //
327  if (mSolidType == "G4ReflectedSolid")
328  {
329  mSolidType = ((G4ReflectedSolid*)mSolid)->GetConstituentMovedSolid()
330  ->GetEntityType();
331  }
332 
333  // Parameterisation type depend of mother solid type and axis of division
334  //
335  if( mSolidType == "G4Box" )
336  {
337  switch( axis )
338  {
339  case kXAxis:
340  fparam = new G4ParameterisationBoxX( axis, nDivs, width,
341  offset, mSolid, divType );
342  break;
343  case kYAxis:
344  fparam = new G4ParameterisationBoxY( axis, nDivs, width,
345  offset, mSolid, divType );
346  break;
347  case kZAxis:
348  fparam = new G4ParameterisationBoxZ( axis, nDivs, width,
349  offset, mSolid, divType );
350  break;
351  default:
352  ErrorInAxis( axis, mSolid );
353  break;
354  }
355  }
356  else if( mSolidType == "G4Tubs" )
357  {
358  switch( axis )
359  {
360  case kRho:
361  fparam = new G4ParameterisationTubsRho( axis, nDivs, width,
362  offset, mSolid, divType );
363  break;
364  case kPhi:
365  fparam = new G4ParameterisationTubsPhi( axis, nDivs, width,
366  offset, mSolid, divType );
367  break;
368  case kZAxis:
369  fparam = new G4ParameterisationTubsZ( axis, nDivs, width,
370  offset, mSolid, divType );
371  break;
372  default:
373  ErrorInAxis( axis, mSolid );
374  break;
375  }
376  }
377  else if( mSolidType == "G4Cons" )
378  {
379  switch( axis )
380  {
381  case kRho:
382  fparam = new G4ParameterisationConsRho( axis, nDivs, width,
383  offset, mSolid, divType );
384  break;
385  case kPhi:
386  fparam = new G4ParameterisationConsPhi( axis, nDivs, width,
387  offset, mSolid, divType );
388  break;
389  case kZAxis:
390  fparam = new G4ParameterisationConsZ( axis, nDivs, width,
391  offset, mSolid, divType );
392  break;
393  default:
394  ErrorInAxis( axis, mSolid );
395  break;
396  }
397  }
398  else if( mSolidType == "G4Trd" )
399  {
400  switch( axis )
401  {
402  case kXAxis:
403  fparam = new G4ParameterisationTrdX( axis, nDivs, width,
404  offset, mSolid, divType );
405  break;
406  case kYAxis:
407  fparam = new G4ParameterisationTrdY( axis, nDivs, width,
408  offset, mSolid, divType );
409  break;
410  case kZAxis:
411  fparam = new G4ParameterisationTrdZ( axis, nDivs, width,
412  offset, mSolid, divType );
413  break;
414  default:
415  ErrorInAxis( axis, mSolid );
416  break;
417  }
418  }
419  else if( mSolidType == "G4Para" )
420  {
421  switch( axis )
422  {
423  case kXAxis:
424  fparam = new G4ParameterisationParaX( axis, nDivs, width,
425  offset, mSolid, divType );
426  break;
427  case kYAxis:
428  fparam = new G4ParameterisationParaY( axis, nDivs, width,
429  offset, mSolid, divType );
430  break;
431  case kZAxis:
432  fparam = new G4ParameterisationParaZ( axis, nDivs, width,
433  offset, mSolid, divType );
434  break;
435  default:
436  ErrorInAxis( axis, mSolid );
437  break;
438  }
439  }
440 // else if( mSolidType == "G4Trap" )
441 // {
442 // }
443 // else if( mSolidType == "G4Polycone" )
444 // {
445 // switch( axis )
446 // {
447 // case kRho:
448 // fparam = new G4ParameterisationPolyconeRho( axis, nDivs, width,
449 // offset, mSolid, divType );
450 // break;
451 // case kPhi:
452 // fparam = new G4ParameterisationPolyconePhi( axis, nDivs, width,
453 // offset, mSolid, divType );
454 // break;
455 // case kZAxis:
456 // fparam = new G4ParameterisationPolyconeZ( axis, nDivs, width,
457 // offset, mSolid, divType );
458 // break;
459 // default:
460 // ErrorInAxis( axis, mSolid );
461 // break;
462 // }
463 // }
464 // else if( mSolidType == "G4Polyhedra" )
465 // {
466 // switch( axis )
467 // {
468 // case kRho:
469 // fparam = new G4ParameterisationPolyhedraRho( axis, nDivs, width,
470 // offset, mSolid, divType );
471 // break;
472 // case kPhi:
473 // fparam = new G4ParameterisationPolyhedraPhi( axis, nDivs, width,
474 // offset, mSolid, divType );
475 // break;
476 // case kZAxis:
477 // fparam = new G4ParameterisationPolyhedraZ( axis, nDivs, width,
478 // offset, mSolid, divType );
479 // break;
480 // default:
481 // ErrorInAxis( axis, mSolid );
482 // break;
483 // }
484 // }
485  else
486  {
487  std::ostringstream message;
488  message << "Solid type not supported: " << mSolidType << "." << G4endl
489  << "Divisions for " << mSolidType << " not implemented.";
490  G4Exception("G4ReplicatedSlice::SetParameterisation()", "GeomDiv0001",
491  FatalException, message);
492  }
493 
494  fparam->SetHalfGap(half_gap);
495 }
Definition: geomdefs.hh:54
#define width
void SetHalfGap(G4double hg)
virtual G4GeometryType GetEntityType() const =0
G4VDivisionParameterisation * fparam
void ErrorInAxis(EAxis axis, G4VSolid *solid)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
Definition: geomdefs.hh:54
G4VSolid * GetSolid() const
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ faxis

EAxis G4ReplicatedSlice::faxis
protected

Definition at line 186 of file G4ReplicatedSlice.hh.

◆ fcopyNo

G4int G4ReplicatedSlice::fcopyNo
protected

Definition at line 190 of file G4ReplicatedSlice.hh.

◆ fdivAxis

EAxis G4ReplicatedSlice::fdivAxis
protected

Definition at line 187 of file G4ReplicatedSlice.hh.

◆ fnReplicas

G4int G4ReplicatedSlice::fnReplicas
protected

Definition at line 188 of file G4ReplicatedSlice.hh.

◆ foffset

G4double G4ReplicatedSlice::foffset
protected

Definition at line 189 of file G4ReplicatedSlice.hh.

◆ fparam

G4VDivisionParameterisation* G4ReplicatedSlice::fparam
protected

Definition at line 191 of file G4ReplicatedSlice.hh.

◆ fwidth

G4double G4ReplicatedSlice::fwidth
protected

Definition at line 189 of file G4ReplicatedSlice.hh.


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