Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PhantomParameterisation Class Reference

#include <G4PhantomParameterisation.hh>

+ Inheritance diagram for G4PhantomParameterisation:

Public Member Functions

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

Protected Attributes

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

Detailed Description

Definition at line 68 of file G4PhantomParameterisation.hh.

Constructor & Destructor Documentation

◆ G4PhantomParameterisation()

G4PhantomParameterisation::G4PhantomParameterisation ( )

Definition at line 42 of file G4PhantomParameterisation.cc.

◆ ~G4PhantomParameterisation()

G4PhantomParameterisation::~G4PhantomParameterisation ( )

Definition at line 49 of file G4PhantomParameterisation.cc.

50{
51}

Member Function Documentation

◆ BuildContainerSolid() [1/2]

void G4PhantomParameterisation::BuildContainerSolid ( G4VPhysicalVolume pPhysicalVol)

◆ BuildContainerSolid() [2/2]

void G4PhantomParameterisation::BuildContainerSolid ( G4VSolid pMotherSolid)

Definition at line 67 of file G4PhantomParameterisation.cc.

69{
70 fContainerSolid = pMotherSolid;
74
75 // CheckVoxelsFillContainer();
76}

◆ CheckVoxelsFillContainer()

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

Definition at line 177 of file G4PhantomParameterisation.cc.

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

◆ 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 std::size_t matIndex = GetMaterialIndex(copyNo);
124
125 return fMaterials[ matIndex ];
126}
std::size_t GetMaterialIndex(std::size_t nx, std::size_t ny, std::size_t nz) const
std::vector< G4Material * > fMaterials

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

◆ 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)

Referenced by G4RegularNavigation::LevelLocate().

◆ GetContainerSolid()

G4VSolid * G4PhantomParameterisation::GetContainerSolid ( ) const
inline

◆ GetMaterial() [1/2]

G4Material * G4PhantomParameterisation::GetMaterial ( std::size_t  copyNo) const

Definition at line 158 of file G4PhantomParameterisation.cc.

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

◆ GetMaterial() [2/2]

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

Definition at line 151 of file G4PhantomParameterisation.cc.

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

Referenced by G4EnergySplitter::SplitEnergyInVolumes().

◆ GetMaterialIndex() [1/2]

std::size_t G4PhantomParameterisation::GetMaterialIndex ( std::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}

◆ GetMaterialIndex() [2/2]

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

Definition at line 141 of file G4PhantomParameterisation.cc.

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

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

◆ GetMaterialIndices()

std::size_t * G4PhantomParameterisation::GetMaterialIndices ( ) const
inline

◆ GetMaterials()

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

◆ GetNoVoxels()

std::size_t G4PhantomParameterisation::GetNoVoxels ( ) const
inline

◆ GetNoVoxelsX()

std::size_t G4PhantomParameterisation::GetNoVoxelsX ( ) const
inline

◆ GetNoVoxelsY()

std::size_t G4PhantomParameterisation::GetNoVoxelsY ( ) const
inline

◆ GetNoVoxelsZ()

std::size_t G4PhantomParameterisation::GetNoVoxelsZ ( ) const
inline

◆ GetReplicaNo()

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

Reimplemented in G4PartialPhantomParameterisation.

Definition at line 224 of file G4PhantomParameterisation.cc.

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

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 std::size_t nx;
98 std::size_t ny;
99 std::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}

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

◆ SetMaterialIndices()

void G4PhantomParameterisation::SetMaterialIndices ( std::size_t *  matInd)
inline

◆ SetMaterials()

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

◆ SetNoVoxels()

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

◆ SetSkipEqualMaterials()

void G4PhantomParameterisation::SetSkipEqualMaterials ( G4bool  skip)

◆ SetVoxelDimensions()

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

◆ SkipEqualMaterials()

G4bool G4PhantomParameterisation::SkipEqualMaterials ( ) const

Member Data Documentation

◆ bSkipEqualMaterials

G4bool G4PhantomParameterisation::bSkipEqualMaterials = true
protected

Definition at line 192 of file G4PhantomParameterisation.hh.

◆ fContainerSolid

G4VSolid* G4PhantomParameterisation::fContainerSolid = nullptr
protected

◆ fContainerWallX

◆ fContainerWallY

◆ fContainerWallZ

◆ fMaterialIndices

std::size_t* G4PhantomParameterisation::fMaterialIndices = nullptr
protected

◆ fMaterials

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

◆ fNoVoxels

std::size_t G4PhantomParameterisation::fNoVoxels = 0
protected

Definition at line 174 of file G4PhantomParameterisation.hh.

◆ fNoVoxelsX

◆ fNoVoxelsXY

std::size_t G4PhantomParameterisation::fNoVoxelsXY = 0
protected

◆ fNoVoxelsY

◆ fNoVoxelsZ

◆ fVoxelHalfX

◆ fVoxelHalfY

◆ fVoxelHalfZ

◆ kCarTolerance

G4double G4PhantomParameterisation::kCarTolerance
protected

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