Geant4-11
Public Member Functions | Protected Attributes | Private Member Functions
G4PhantomParameterisation Class Reference

#include <G4PhantomParameterisation.hh>

Inheritance diagram for G4PhantomParameterisation:
G4VPVParameterisation G4PartialPhantomParameterisation

Public Member Functions

void BuildContainerSolid (G4VPhysicalVolume *pPhysicalVol)
 
void BuildContainerSolid (G4VSolid *pMotherSolid)
 
void CheckVoxelsFillContainer (G4double contX, G4double contY, G4double contZ) const
 
void ComputeDimensions (G4Box &, const G4int, const G4VPhysicalVolume *) const
 
void ComputeDimensions (G4Cons &, const G4int, const G4VPhysicalVolume *) const
 
void ComputeDimensions (G4Ellipsoid &, const G4int, const G4VPhysicalVolume *) const
 
void ComputeDimensions (G4Hype &, const G4int, const G4VPhysicalVolume *) const
 
void ComputeDimensions (G4Orb &, const G4int, const G4VPhysicalVolume *) const
 
void ComputeDimensions (G4Para &, const G4int, const G4VPhysicalVolume *) const
 
void ComputeDimensions (G4Polycone &, const G4int, const G4VPhysicalVolume *) const
 
void ComputeDimensions (G4Polyhedra &, const G4int, const G4VPhysicalVolume *) const
 
void ComputeDimensions (G4Sphere &, const G4int, const G4VPhysicalVolume *) const
 
void ComputeDimensions (G4Torus &, const G4int, const G4VPhysicalVolume *) const
 
void ComputeDimensions (G4Trap &, const G4int, const G4VPhysicalVolume *) const
 
void ComputeDimensions (G4Trd &, const G4int, const G4VPhysicalVolume *) const
 
void ComputeDimensions (G4Tubs &, const G4int, const G4VPhysicalVolume *) const
 
virtual G4MaterialComputeMaterial (const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=nullptr)
 
virtual G4VSolidComputeSolid (const G4int, G4VPhysicalVolume *)
 
virtual void ComputeTransformation (const G4int, G4VPhysicalVolume *) const
 
 G4PhantomParameterisation ()
 
G4VSolidGetContainerSolid () const
 
G4MaterialGetMaterial (size_t copyNo) const
 
G4MaterialGetMaterial (size_t nx, size_t ny, size_t nz) const
 
size_t GetMaterialIndex (size_t copyNo) const
 
size_t GetMaterialIndex (size_t nx, size_t ny, size_t nz) const
 
size_t * GetMaterialIndices () const
 
std::vector< G4Material * > GetMaterials () const
 
virtual G4VVolumeMaterialScannerGetMaterialScanner ()
 
size_t GetNoVoxels () const
 
size_t GetNoVoxelsX () const
 
size_t GetNoVoxelsY () const
 
size_t GetNoVoxelsZ () const
 
virtual G4int GetReplicaNo (const G4ThreeVector &localPoint, const G4ThreeVector &localDir)
 
G4ThreeVector GetTranslation (const G4int copyNo) const
 
G4double GetVoxelHalfX () const
 
G4double GetVoxelHalfY () const
 
G4double GetVoxelHalfZ () const
 
virtual G4bool IsNested () const
 
void SetMaterialIndices (size_t *matInd)
 
void SetMaterials (std::vector< G4Material * > &mates)
 
void SetNoVoxels (size_t nx, size_t ny, size_t nz)
 
void SetSkipEqualMaterials (G4bool skip)
 
void SetVoxelDimensions (G4double halfx, G4double halfy, G4double halfz)
 
G4bool SkipEqualMaterials () const
 
 ~G4PhantomParameterisation ()
 

Protected Attributes

G4bool bSkipEqualMaterials = true
 
G4VSolidfContainerSolid = nullptr
 
G4double fContainerWallX =0.0
 
G4double fContainerWallY =0.0
 
G4double fContainerWallZ =0.0
 
size_t * fMaterialIndices = nullptr
 
std::vector< G4Material * > fMaterials
 
size_t fNoVoxels = 0
 
size_t fNoVoxelsX = 0
 
size_t fNoVoxelsXY = 0
 
size_t fNoVoxelsY = 0
 
size_t fNoVoxelsZ = 0
 
G4double fVoxelHalfX = 0.0
 
G4double fVoxelHalfY = 0.0
 
G4double fVoxelHalfZ = 0.0
 
G4double kCarTolerance
 

Private Member Functions

void CheckCopyNo (const G4int copyNo) const
 
void ComputeVoxelIndices (const G4int copyNo, size_t &nx, size_t &ny, size_t &nz) const
 

Detailed Description

Definition at line 68 of file G4PhantomParameterisation.hh.

Constructor & Destructor Documentation

◆ G4PhantomParameterisation()

G4PhantomParameterisation::G4PhantomParameterisation ( )

◆ ~G4PhantomParameterisation()

G4PhantomParameterisation::~G4PhantomParameterisation ( )

Definition at line 49 of file G4PhantomParameterisation.cc.

50{
51}

Member Function Documentation

◆ BuildContainerSolid() [1/2]

void G4PhantomParameterisation::BuildContainerSolid ( G4VPhysicalVolume pPhysicalVol)

Definition at line 55 of file G4PhantomParameterisation.cc.

References fContainerSolid, fContainerWallX, fContainerWallY, fContainerWallZ, fNoVoxelsX, fNoVoxelsY, fNoVoxelsZ, fVoxelHalfX, fVoxelHalfY, fVoxelHalfZ, G4VPhysicalVolume::GetLogicalVolume(), and G4LogicalVolume::GetSolid().

◆ BuildContainerSolid() [2/2]

void G4PhantomParameterisation::BuildContainerSolid ( G4VSolid pMotherSolid)

◆ CheckCopyNo()

void G4PhantomParameterisation::CheckCopyNo ( const G4int  copyNo) const
private

Definition at line 390 of file G4PhantomParameterisation.cc.

391{
392 if( copyNo < 0 || copyNo >= G4int(fNoVoxels) )
393 {
394 std::ostringstream message;
395 message << "Copy number is negative or too big!" << G4endl
396 << " Copy number: " << copyNo << G4endl
397 << " Total number of voxels: " << fNoVoxels;
398 G4Exception("G4PhantomParameterisation::CheckCopyNo()",
399 "GeomNav0002", FatalErrorInArgument, message);
400 }
401}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57

References FatalErrorInArgument, fNoVoxels, G4endl, and G4Exception().

Referenced by ComputeMaterial(), ComputeVoxelIndices(), GetMaterialIndex(), and GetTranslation().

◆ CheckVoxelsFillContainer()

void G4PhantomParameterisation::CheckVoxelsFillContainer ( G4double  contX,
G4double  contY,
G4double  contZ 
) const

Definition at line 175 of file G4PhantomParameterisation.cc.

177{
178 G4double toleranceForWarning = 0.25*kCarTolerance;
179
180 // Any bigger value than 0.25*kCarTolerance will give a warning in
181 // G4NormalNavigation::ComputeStep(), because the Inverse of a container
182 // translation that is Z+epsilon gives -Z+epsilon (and the maximum tolerance
183 // in G4Box::Inside is 0.5*kCarTolerance
184 //
185 G4double toleranceForError = 1.*kCarTolerance;
186
187 // Any bigger value than kCarTolerance will give an error in GetReplicaNo()
188 //
189 if( std::fabs(contX-fNoVoxelsX*fVoxelHalfX) >= toleranceForError
190 || std::fabs(contY-fNoVoxelsY*fVoxelHalfY) >= toleranceForError
191 || std::fabs(contZ-fNoVoxelsZ*fVoxelHalfZ) >= toleranceForError )
192 {
193 std::ostringstream message;
194 message << "Voxels do not fully fill the container: "
196 << " DiffX= " << contX-fNoVoxelsX*fVoxelHalfX << G4endl
197 << " DiffY= " << contY-fNoVoxelsY*fVoxelHalfY << G4endl
198 << " DiffZ= " << contZ-fNoVoxelsZ*fVoxelHalfZ << G4endl
199 << " Maximum difference is: " << toleranceForError;
200 G4Exception("G4PhantomParameterisation::CheckVoxelsFillContainer()",
201 "GeomNav0002", FatalException, message);
202
203 }
204 else if( std::fabs(contX-fNoVoxelsX*fVoxelHalfX) >= toleranceForWarning
205 || std::fabs(contY-fNoVoxelsY*fVoxelHalfY) >= toleranceForWarning
206 || std::fabs(contZ-fNoVoxelsZ*fVoxelHalfZ) >= toleranceForWarning )
207 {
208 std::ostringstream message;
209 message << "Voxels do not fully fill the container: "
211 << " DiffX= " << contX-fNoVoxelsX*fVoxelHalfX << G4endl
212 << " DiffY= " << contY-fNoVoxelsY*fVoxelHalfY << G4endl
213 << " DiffZ= " << contZ-fNoVoxelsZ*fVoxelHalfZ << G4endl
214 << " Maximum difference is: " << toleranceForWarning;
215 G4Exception("G4PhantomParameterisation::CheckVoxelsFillContainer()",
216 "GeomNav1002", JustWarning, message);
217 }
218}
@ JustWarning
@ FatalException
double G4double
Definition: G4Types.hh:83
G4String GetName() const

References FatalException, fContainerSolid, fNoVoxelsX, fNoVoxelsY, fNoVoxelsZ, fVoxelHalfX, fVoxelHalfY, fVoxelHalfZ, G4endl, G4Exception(), G4VSolid::GetName(), JustWarning, and kCarTolerance.

◆ ComputeDimensions() [1/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Box ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 84 of file G4PhantomParameterisation.hh.

85 {}

◆ ComputeDimensions() [2/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Cons ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 92 of file G4PhantomParameterisation.hh.

93 {}

◆ ComputeDimensions() [3/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Ellipsoid ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 98 of file G4PhantomParameterisation.hh.

99 {}

◆ ComputeDimensions() [4/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Hype ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 104 of file G4PhantomParameterisation.hh.

105 {}

◆ ComputeDimensions() [5/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Orb ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 94 of file G4PhantomParameterisation.hh.

95 {}

◆ ComputeDimensions() [6/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Para ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 102 of file G4PhantomParameterisation.hh.

103 {}

◆ ComputeDimensions() [7/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Polycone ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 106 of file G4PhantomParameterisation.hh.

107 {}

◆ ComputeDimensions() [8/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Polyhedra ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 108 of file G4PhantomParameterisation.hh.

109 {}

◆ ComputeDimensions() [9/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Sphere ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 96 of file G4PhantomParameterisation.hh.

97 {}

◆ ComputeDimensions() [10/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Torus ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 100 of file G4PhantomParameterisation.hh.

101 {}

◆ ComputeDimensions() [11/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Trap ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 90 of file G4PhantomParameterisation.hh.

91 {}

◆ ComputeDimensions() [12/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Trd ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 88 of file G4PhantomParameterisation.hh.

89 {}

◆ ComputeDimensions() [13/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Tubs ,
const  G4int,
const G4VPhysicalVolume  
) const
inlinevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 86 of file G4PhantomParameterisation.hh.

87 {}

◆ ComputeMaterial()

G4Material * G4PhantomParameterisation::ComputeMaterial ( const G4int  repNo,
G4VPhysicalVolume currentVol,
const G4VTouchable parentTouch = nullptr 
)
virtual

Reimplemented from G4VPVParameterisation.

Reimplemented in G4PartialPhantomParameterisation.

Definition at line 119 of file G4PhantomParameterisation.cc.

121{
122 CheckCopyNo( copyNo );
123 size_t matIndex = GetMaterialIndex(copyNo);
124
125 return fMaterials[ matIndex ];
126}
size_t GetMaterialIndex(size_t nx, size_t ny, size_t nz) const
std::vector< G4Material * > fMaterials
void CheckCopyNo(const G4int copyNo) const

References CheckCopyNo(), fMaterials, and GetMaterialIndex().

Referenced by G4GMocrenFileSceneHandler::AddSolid(), G4RegularNavigation::ComputeStepSkippingEqualMaterials(), and G4RegularNavigation::LevelLocate().

◆ ComputeSolid()

G4VSolid * G4PhantomParameterisation::ComputeSolid ( const  G4int,
G4VPhysicalVolume pPhysicalVol 
)
virtual

Reimplemented from G4VPVParameterisation.

Definition at line 111 of file G4PhantomParameterisation.cc.

113{
114 return pPhysicalVol->GetLogicalVolume()->GetSolid();
115}
G4VSolid * GetSolid() const
G4LogicalVolume * GetLogicalVolume() const

References G4VPhysicalVolume::GetLogicalVolume(), and G4LogicalVolume::GetSolid().

◆ ComputeTransformation()

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

Implements G4VPVParameterisation.

Reimplemented in G4PartialPhantomParameterisation.

Definition at line 80 of file G4PhantomParameterisation.cc.

82{
83 // Voxels cannot be rotated, return translation
84 //
85 G4ThreeVector trans = GetTranslation( copyNo );
86
87 physVol->SetTranslation( trans );
88}
G4ThreeVector GetTranslation(const G4int copyNo) const
void SetTranslation(const G4ThreeVector &v)

References GetTranslation(), and G4VPhysicalVolume::SetTranslation().

Referenced by G4RegularNavigation::LevelLocate().

◆ ComputeVoxelIndices()

void G4PhantomParameterisation::ComputeVoxelIndices ( const G4int  copyNo,
size_t &  nx,
size_t &  ny,
size_t &  nz 
) const
private

Definition at line 163 of file G4PhantomParameterisation.cc.

166{
167 CheckCopyNo( copyNo );
168 nx = size_t(copyNo%fNoVoxelsX);
169 ny = size_t( (copyNo/fNoVoxelsX)%fNoVoxelsY );
170 nz = size_t(copyNo/fNoVoxelsXY);
171}

References CheckCopyNo(), fNoVoxelsX, fNoVoxelsXY, and fNoVoxelsY.

Referenced by GetTranslation().

◆ GetContainerSolid()

G4VSolid * G4PhantomParameterisation::GetContainerSolid ( ) const
inline

◆ GetMaterial() [1/2]

G4Material * G4PhantomParameterisation::GetMaterial ( size_t  copyNo) const

Definition at line 157 of file G4PhantomParameterisation.cc.

158{
159 return fMaterials[GetMaterialIndex(copyNo)];
160}

References fMaterials, and GetMaterialIndex().

◆ GetMaterial() [2/2]

G4Material * G4PhantomParameterisation::GetMaterial ( size_t  nx,
size_t  ny,
size_t  nz 
) const

Definition at line 151 of file G4PhantomParameterisation.cc.

152{
153 return fMaterials[GetMaterialIndex(nx,ny,nz)];
154}

References fMaterials, and GetMaterialIndex().

Referenced by G4EnergySplitter::SplitEnergyInVolumes().

◆ GetMaterialIndex() [1/2]

size_t G4PhantomParameterisation::GetMaterialIndex ( size_t  copyNo) const

Definition at line 130 of file G4PhantomParameterisation.cc.

132{
133 CheckCopyNo( copyNo );
134
135 if( fMaterialIndices == nullptr ) { return 0; }
136 return *(fMaterialIndices+copyNo);
137}

References CheckCopyNo(), and fMaterialIndices.

◆ GetMaterialIndex() [2/2]

size_t G4PhantomParameterisation::GetMaterialIndex ( size_t  nx,
size_t  ny,
size_t  nz 
) const

Definition at line 141 of file G4PhantomParameterisation.cc.

143{
144 size_t copyNo = nx + fNoVoxelsX*ny + fNoVoxelsXY*nz;
145 return GetMaterialIndex( copyNo );
146}

References fNoVoxelsX, fNoVoxelsXY, and GetMaterialIndex().

Referenced by ComputeMaterial(), GetMaterial(), and GetMaterialIndex().

◆ GetMaterialIndices()

size_t * G4PhantomParameterisation::GetMaterialIndices ( ) const
inline

◆ GetMaterials()

std::vector< G4Material * > G4PhantomParameterisation::GetMaterials ( ) const
inline

◆ GetMaterialScanner()

G4VVolumeMaterialScanner * G4VPVParameterisation::GetMaterialScanner ( )
virtualinherited

Reimplemented in G4VNestedParameterisation.

Definition at line 62 of file G4VPVParameterisation.cc.

63{
64 return nullptr;
65}

Referenced by G4Region::ScanVolumeTree().

◆ GetNoVoxels()

size_t G4PhantomParameterisation::GetNoVoxels ( ) const
inline

◆ GetNoVoxelsX()

size_t G4PhantomParameterisation::GetNoVoxelsX ( ) const
inline

◆ GetNoVoxelsY()

size_t G4PhantomParameterisation::GetNoVoxelsY ( ) const
inline

◆ GetNoVoxelsZ()

size_t G4PhantomParameterisation::GetNoVoxelsZ ( ) const
inline

◆ GetReplicaNo()

G4int G4PhantomParameterisation::GetReplicaNo ( const G4ThreeVector localPoint,
const G4ThreeVector localDir 
)
virtual

Reimplemented in G4PartialPhantomParameterisation.

Definition at line 222 of file G4PhantomParameterisation.cc.

224{
225
226 // Check first that point is really inside voxels
227 //
228 if( fContainerSolid->Inside( localPoint ) == kOutside )
229 {
230 if( std::fabs(localPoint.x()) - fContainerWallX > kCarTolerance
231 && std::fabs(localPoint.y()) - fContainerWallY > kCarTolerance
232 && std::fabs(localPoint.z()) - fContainerWallZ > kCarTolerance )
233 {
234 std::ostringstream message;
235 message << "Point outside voxels!" << G4endl
236 << " localPoint - " << localPoint
237 << " - is outside container solid: "
239 << "DIFFERENCE WITH PHANTOM WALLS X: "
240 << std::fabs(localPoint.x()) - fContainerWallX
241 << " Y: " << std::fabs(localPoint.y()) - fContainerWallY
242 << " Z: " << std::fabs(localPoint.z()) - fContainerWallZ;
243 G4Exception("G4PhantomParameterisation::GetReplicaNo()", "GeomNav0003",
244 FatalErrorInArgument, message);
245 }
246 }
247
248 // Check the voxel numbers corresponding to localPoint
249 // When a particle is on a surface, it may be between -kCarTolerance and
250 // +kCartolerance. By a simple distance as:
251 // G4int nx = G4int( (localPoint.x()+)/fVoxelHalfX/2.);
252 // those between -kCartolerance and 0 will be placed on voxel N-1 and those
253 // between 0 and kCarTolerance on voxel N.
254 // To avoid precision problems place the tracks that are on the surface on
255 // voxel N-1 if they have negative direction and on voxel N if they have
256 // positive direction.
257 // Add +kCarTolerance so that they are first placed on voxel N, and then
258 // if the direction is negative substract 1
259
260 G4double fx = (localPoint.x()+fContainerWallX+kCarTolerance)/(fVoxelHalfX*2.);
261 G4int nx = G4int(fx);
262
263 G4double fy = (localPoint.y()+fContainerWallY+kCarTolerance)/(fVoxelHalfY*2.);
264 G4int ny = G4int(fy);
265
266 G4double fz = (localPoint.z()+fContainerWallZ+kCarTolerance)/(fVoxelHalfZ*2.);
267 G4int nz = G4int(fz);
268
269 // If it is on the surface side, check the direction: if direction is
270 // negative place it in the previous voxel (if direction is positive it is
271 // already in the next voxel).
272 // Correct also cases where n = -1 or n = fNoVoxels. It is always traced to be
273 // due to multiple scattering: track is entering a voxel but multiple
274 // scattering changes the angle towards outside
275 //
276 if( fx - nx < kCarTolerance*fVoxelHalfX )
277 {
278 if( localDir.x() < 0 )
279 {
280 if( nx != 0 )
281 {
282 nx -= 1;
283 }
284 }
285 else
286 {
287 if( nx == G4int(fNoVoxelsX) )
288 {
289 nx -= 1;
290 }
291 }
292 }
293 if( fy - ny < kCarTolerance*fVoxelHalfY )
294 {
295 if( localDir.y() < 0 )
296 {
297 if( ny != 0 )
298 {
299 ny -= 1;
300 }
301 }
302 else
303 {
304 if( ny == G4int(fNoVoxelsY) )
305 {
306 ny -= 1;
307 }
308 }
309 }
310 if( fz - nz < kCarTolerance*fVoxelHalfZ )
311 {
312 if( localDir.z() < 0 )
313 {
314 if( nz != 0 )
315 {
316 nz -= 1;
317 }
318 }
319 else
320 {
321 if( nz == G4int(fNoVoxelsZ) )
322 {
323 nz -= 1;
324 }
325 }
326 }
327
328 G4int copyNo = nx + fNoVoxelsX*ny + fNoVoxelsXY*nz;
329
330 // Check if there are still errors
331 //
332 G4bool isOK = true;
333 if( nx < 0 )
334 {
335 nx = 0;
336 isOK = false;
337 }
338 else if( nx >= G4int(fNoVoxelsX) )
339 {
340 nx = fNoVoxelsX-1;
341 isOK = false;
342 }
343 if( ny < 0 )
344 {
345 ny = 0;
346 isOK = false;
347 }
348 else if( ny >= G4int(fNoVoxelsY) )
349 {
350 ny = fNoVoxelsY-1;
351 isOK = false;
352 }
353 if( nz < 0 )
354 {
355 nz = 0;
356 isOK = false;
357 }
358 else if( nz >= G4int(fNoVoxelsZ) )
359 {
360 nz = fNoVoxelsZ-1;
361 isOK = false;
362 }
363 if( !isOK )
364 {
365 if( std::fabs(localPoint.x()-fContainerWallX) > kCarTolerance &&
366 std::fabs(localPoint.y()-fContainerWallY) > kCarTolerance &&
367 std::fabs(localPoint.z()-fContainerWallZ) > kCarTolerance ){ // only if difference is big
368 std::ostringstream message;
369 message << "Corrected the copy number! It was negative or too big" << G4endl
370 << " LocalPoint: " << localPoint << G4endl
371 << " LocalDir: " << localDir << G4endl
372 << " Voxel container size: " << fContainerWallX
373 << " " << fContainerWallY << " " << fContainerWallZ << G4endl
374 << " LocalPoint - wall: "
375 << localPoint.x()-fContainerWallX << " "
376 << localPoint.y()-fContainerWallY << " "
377 << localPoint.z()-fContainerWallZ;
378 G4Exception("G4PhantomParameterisation::GetReplicaNo()",
379 "GeomNav1002", JustWarning, message);
380 }
381
382 copyNo = nx + fNoVoxelsX*ny + fNoVoxelsXY*nz;
383 }
384
385 return copyNo;
386}
bool G4bool
Definition: G4Types.hh:86
double z() const
double x() const
double y() const
virtual EInside Inside(const G4ThreeVector &p) const =0
@ kOutside
Definition: geomdefs.hh:68

References FatalErrorInArgument, fContainerSolid, fContainerWallX, fContainerWallY, fContainerWallZ, fNoVoxelsX, fNoVoxelsXY, fNoVoxelsY, fNoVoxelsZ, fVoxelHalfX, fVoxelHalfY, fVoxelHalfZ, G4endl, G4Exception(), G4VSolid::GetName(), G4VSolid::Inside(), JustWarning, kCarTolerance, kOutside, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by G4RegularNavigation::ComputeStep(), G4RegularNavigation::ComputeStepSkippingEqualMaterials(), and G4RegularNavigation::LevelLocate().

◆ GetTranslation()

G4ThreeVector G4PhantomParameterisation::GetTranslation ( const G4int  copyNo) const

Definition at line 92 of file G4PhantomParameterisation.cc.

94{
95 CheckCopyNo( copyNo );
96
97 size_t nx;
98 size_t ny;
99 size_t nz;
100
101 ComputeVoxelIndices( copyNo, nx, ny, nz );
102
103 G4ThreeVector trans( (2*nx+1)*fVoxelHalfX - fContainerWallX,
104 (2*ny+1)*fVoxelHalfY - fContainerWallY,
105 (2*nz+1)*fVoxelHalfZ - fContainerWallZ);
106 return trans;
107}
void ComputeVoxelIndices(const G4int copyNo, size_t &nx, size_t &ny, size_t &nz) const

References CheckCopyNo(), ComputeVoxelIndices(), fContainerWallX, fContainerWallY, fContainerWallZ, fVoxelHalfX, fVoxelHalfY, and fVoxelHalfZ.

Referenced by G4RegularNavigation::ComputeStep(), G4RegularNavigation::ComputeStepSkippingEqualMaterials(), and ComputeTransformation().

◆ GetVoxelHalfX()

G4double G4PhantomParameterisation::GetVoxelHalfX ( ) const
inline

◆ GetVoxelHalfY()

G4double G4PhantomParameterisation::GetVoxelHalfY ( ) const
inline

◆ GetVoxelHalfZ()

G4double G4PhantomParameterisation::GetVoxelHalfZ ( ) const
inline

◆ IsNested()

G4bool G4VPVParameterisation::IsNested ( ) const
virtualinherited

◆ SetMaterialIndices()

void G4PhantomParameterisation::SetMaterialIndices ( size_t *  matInd)
inline

◆ SetMaterials()

void G4PhantomParameterisation::SetMaterials ( std::vector< G4Material * > &  mates)
inline

◆ SetNoVoxels()

void G4PhantomParameterisation::SetNoVoxels ( size_t  nx,
size_t  ny,
size_t  nz 
)

◆ SetSkipEqualMaterials()

void G4PhantomParameterisation::SetSkipEqualMaterials ( G4bool  skip)

◆ SetVoxelDimensions()

void G4PhantomParameterisation::SetVoxelDimensions ( G4double  halfx,
G4double  halfy,
G4double  halfz 
)

◆ SkipEqualMaterials()

G4bool G4PhantomParameterisation::SkipEqualMaterials ( ) const

Field Documentation

◆ bSkipEqualMaterials

G4bool G4PhantomParameterisation::bSkipEqualMaterials = true
protected

Definition at line 192 of file G4PhantomParameterisation.hh.

◆ fContainerSolid

G4VSolid* G4PhantomParameterisation::fContainerSolid = nullptr
protected

◆ fContainerWallX

G4double G4PhantomParameterisation::fContainerWallX =0.0
protected

◆ fContainerWallY

G4double G4PhantomParameterisation::fContainerWallY =0.0
protected

◆ fContainerWallZ

G4double G4PhantomParameterisation::fContainerWallZ =0.0
protected

◆ fMaterialIndices

size_t* G4PhantomParameterisation::fMaterialIndices = nullptr
protected

◆ fMaterials

std::vector<G4Material*> G4PhantomParameterisation::fMaterials
protected

◆ fNoVoxels

size_t G4PhantomParameterisation::fNoVoxels = 0
protected

◆ fNoVoxelsX

size_t G4PhantomParameterisation::fNoVoxelsX = 0
protected

◆ fNoVoxelsXY

size_t G4PhantomParameterisation::fNoVoxelsXY = 0
protected

◆ fNoVoxelsY

size_t G4PhantomParameterisation::fNoVoxelsY = 0
protected

◆ fNoVoxelsZ

size_t G4PhantomParameterisation::fNoVoxelsZ = 0
protected

◆ fVoxelHalfX

G4double G4PhantomParameterisation::fVoxelHalfX = 0.0
protected

◆ fVoxelHalfY

G4double G4PhantomParameterisation::fVoxelHalfY = 0.0
protected

◆ fVoxelHalfZ

G4double G4PhantomParameterisation::fVoxelHalfZ = 0.0
protected

◆ kCarTolerance

G4double G4PhantomParameterisation::kCarTolerance
protected

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