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

G4PhantomParameterisation describes regular parameterisations: a set of boxes of equal dimension in the x, y and z dimensions. The G4PVParameterised volume using this class must be placed inside a volume that is completely filled by these boxes. More...

#include <G4PhantomParameterisation.hh>

Inheritance diagram for G4PhantomParameterisation:

Public Member Functions

 G4PhantomParameterisation ()
 ~G4PhantomParameterisation () override
void ComputeTransformation (const G4int, G4VPhysicalVolume *) const override
G4VSolidComputeSolid (const G4int, G4VPhysicalVolume *) override
G4MaterialComputeMaterial (const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=nullptr) override
void ComputeDimensions (G4Box &, const G4int, const G4VPhysicalVolume *) const override
void ComputeDimensions (G4Tubs &, const G4int, const G4VPhysicalVolume *) const override
void ComputeDimensions (G4Trd &, const G4int, const G4VPhysicalVolume *) const override
void ComputeDimensions (G4Trap &, const G4int, const G4VPhysicalVolume *) const override
void ComputeDimensions (G4Cons &, const G4int, const G4VPhysicalVolume *) const override
void ComputeDimensions (G4Orb &, const G4int, const G4VPhysicalVolume *) const override
void ComputeDimensions (G4Sphere &, const G4int, const G4VPhysicalVolume *) const override
void ComputeDimensions (G4Ellipsoid &, const G4int, const G4VPhysicalVolume *) const override
void ComputeDimensions (G4Torus &, const G4int, const G4VPhysicalVolume *) const override
void ComputeDimensions (G4Para &, const G4int, const G4VPhysicalVolume *) const override
void ComputeDimensions (G4Hype &, const G4int, const G4VPhysicalVolume *) const override
void ComputeDimensions (G4Polycone &, const G4int, const G4VPhysicalVolume *) const override
void ComputeDimensions (G4Polyhedra &, const G4int, const G4VPhysicalVolume *) const override
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 G4bool IsNested () const
virtual G4VVolumeMaterialScannerGetMaterialScanner ()

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

G4PhantomParameterisation describes regular parameterisations: a set of boxes of equal dimension in the x, y and z dimensions. The G4PVParameterised volume using this class must be placed inside a volume that is completely filled by these boxes.

Definition at line 74 of file G4PhantomParameterisation.hh.

Constructor & Destructor Documentation

◆ G4PhantomParameterisation()

G4PhantomParameterisation::G4PhantomParameterisation ( )

Constructor and Destructor.

Definition at line 41 of file G4PhantomParameterisation.cc.

◆ ~G4PhantomParameterisation()

G4PhantomParameterisation::~G4PhantomParameterisation ( )
overridedefault

Member Function Documentation

◆ BuildContainerSolid() [1/2]

void G4PhantomParameterisation::BuildContainerSolid ( G4VPhysicalVolume * pPhysicalVol)

Saves as container solid the parent of the voxels. Checks that the voxels fill it completely.

Definition at line 52 of file G4PhantomParameterisation.cc.

◆ BuildContainerSolid() [2/2]

void G4PhantomParameterisation::BuildContainerSolid ( G4VSolid * pMotherSolid)

Definition at line 65 of file G4PhantomParameterisation.cc.

67{
68 fContainerSolid = pMotherSolid;
72
73 // CheckVoxelsFillContainer();
74}

◆ CheckVoxelsFillContainer()

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

Checks that the voxels fill it completely.

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: "
195 << fContainerSolid->GetName() << G4endl
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: "
210 << fContainerSolid->GetName() << G4endl
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
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
double G4double
Definition G4Types.hh:83
#define G4endl
Definition G4ios.hh:67

◆ ComputeDimensions() [1/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Box & ,
const G4int ,
const G4VPhysicalVolume *  ) const
inlineoverridevirtual

Dummy declarations ...

Reimplemented from G4VPVParameterisation.

Definition at line 94 of file G4PhantomParameterisation.hh.

95 {}

◆ ComputeDimensions() [2/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Cons & ,
const G4int ,
const G4VPhysicalVolume *  ) const
inlineoverridevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 102 of file G4PhantomParameterisation.hh.

103 {}

◆ ComputeDimensions() [3/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Ellipsoid & ,
const G4int ,
const G4VPhysicalVolume *  ) const
inlineoverridevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 108 of file G4PhantomParameterisation.hh.

109 {}

◆ ComputeDimensions() [4/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Hype & ,
const G4int ,
const G4VPhysicalVolume *  ) const
inlineoverridevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 114 of file G4PhantomParameterisation.hh.

115 {}

◆ ComputeDimensions() [5/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Orb & ,
const G4int ,
const G4VPhysicalVolume *  ) const
inlineoverridevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 104 of file G4PhantomParameterisation.hh.

105 {}

◆ ComputeDimensions() [6/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Para & ,
const G4int ,
const G4VPhysicalVolume *  ) const
inlineoverridevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 112 of file G4PhantomParameterisation.hh.

113 {}

◆ ComputeDimensions() [7/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Polycone & ,
const G4int ,
const G4VPhysicalVolume *  ) const
inlineoverridevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 116 of file G4PhantomParameterisation.hh.

117 {}

◆ ComputeDimensions() [8/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Polyhedra & ,
const G4int ,
const G4VPhysicalVolume *  ) const
inlineoverridevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 118 of file G4PhantomParameterisation.hh.

119 {}

◆ ComputeDimensions() [9/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Sphere & ,
const G4int ,
const G4VPhysicalVolume *  ) const
inlineoverridevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 106 of file G4PhantomParameterisation.hh.

107 {}

◆ ComputeDimensions() [10/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Torus & ,
const G4int ,
const G4VPhysicalVolume *  ) const
inlineoverridevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 110 of file G4PhantomParameterisation.hh.

111 {}

◆ ComputeDimensions() [11/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Trap & ,
const G4int ,
const G4VPhysicalVolume *  ) const
inlineoverridevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 100 of file G4PhantomParameterisation.hh.

101 {}

◆ ComputeDimensions() [12/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Trd & ,
const G4int ,
const G4VPhysicalVolume *  ) const
inlineoverridevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 98 of file G4PhantomParameterisation.hh.

99 {}

◆ ComputeDimensions() [13/13]

void G4PhantomParameterisation::ComputeDimensions ( G4Tubs & ,
const G4int ,
const G4VPhysicalVolume *  ) const
inlineoverridevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 96 of file G4PhantomParameterisation.hh.

97 {}

◆ ComputeMaterial()

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

Computes the material for the 'currentVol' and replica number 'repNo'. Must cope with 'parentTouch' for navigator's SetupHierarchy() when used for nested parameterisations.

Parameters
[in]currentVolPointer to the current physical volume.
[in]repNoThe copy number index.
[in]parentTouchPointer to the touchable of the parent volume.
Returns
A pointer to the associated material.

Reimplemented from G4VPVParameterisation.

Definition at line 117 of file G4PhantomParameterisation.cc.

119{
120 CheckCopyNo( copyNo );
121 std::size_t matIndex = GetMaterialIndex(copyNo);
122
123 return fMaterials[ matIndex ];
124}
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(), and G4RegularNavigation::LevelLocate().

◆ ComputeSolid()

G4VSolid * G4PhantomParameterisation::ComputeSolid ( const G4int no,
G4VPhysicalVolume * pv )
overridevirtual

Computes the solid for the 'pv' volume and replica number 'no'. To be optionally defined in derived classes, for parameterisation of the solid type.

Parameters
[in]noThe copy number index.
[in]pvPointer to the current physical volume.

Reimplemented from G4VPVParameterisation.

Definition at line 109 of file G4PhantomParameterisation.cc.

111{
112 return pPhysicalVol->GetLogicalVolume()->GetSolid();
113}

◆ ComputeTransformation()

void G4PhantomParameterisation::ComputeTransformation ( const G4int no,
G4VPhysicalVolume * pv ) const
overridevirtual

Computes the transformation for the 'pv' volume and replica number 'no'. It is a required method, as it is the reason for this class.

Parameters
[in]pvPointer to the current physical volume.
[in]noThe copy number index.

Implements G4VPVParameterisation.

Definition at line 78 of file G4PhantomParameterisation.cc.

80{
81 // Voxels cannot be rotated, return translation
82 //
83 G4ThreeVector trans = GetTranslation( copyNo );
84
85 physVol->SetTranslation( trans );
86}
CLHEP::Hep3Vector G4ThreeVector
G4ThreeVector GetTranslation(const G4int copyNo) const

Referenced by G4RegularNavigation::LevelLocate().

◆ GetContainerSolid()

G4VSolid * G4PhantomParameterisation::GetContainerSolid ( ) const
inline

◆ GetMaterial() [1/2]

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

Definition at line 156 of file G4PhantomParameterisation.cc.

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

◆ GetMaterial() [2/2]

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

Definition at line 149 of file G4PhantomParameterisation.cc.

150{
151 return fMaterials[GetMaterialIndex(nx,ny,nz)];
152}

◆ GetMaterialIndex() [1/2]

std::size_t G4PhantomParameterisation::GetMaterialIndex ( std::size_t copyNo) const

Definition at line 128 of file G4PhantomParameterisation.cc.

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

◆ GetMaterialIndex() [2/2]

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

Definition at line 139 of file G4PhantomParameterisation.cc.

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

Referenced by ComputeMaterial(), GetMaterial(), 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

Gets the voxel number corresponding to the point in the container frame. Uses 'localDir' to avoid precision problems at the surfaces.

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: "
238 << fContainerSolid->GetName() << G4endl
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 auto nx = G4int(fx);
262
263 G4double fy = (localPoint.y()+fContainerWallY+kCarTolerance)/(fVoxelHalfY*2.);
264 auto ny = G4int(fy);
265
266 G4double fz = (localPoint.z()+fContainerWallZ+kCarTolerance)/(fVoxelHalfZ*2.);
267 auto 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 auto copyNo = G4int(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 = G4int(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 = G4int(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 = G4int(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 = G4int(nx + fNoVoxelsX*ny + fNoVoxelsXY*nz);
383 }
384
385 return copyNo;
386}
@ FatalErrorInArgument
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
double z() const
double x() const
double y() const
@ kOutside
Definition geomdefs.hh:68

Referenced by G4RegularNavigation::LevelLocate().

◆ GetTranslation()

G4ThreeVector G4PhantomParameterisation::GetTranslation ( const G4int copyNo) const

Definition at line 90 of file G4PhantomParameterisation.cc.

92{
93 CheckCopyNo( copyNo );
94
95 std::size_t nx;
96 std::size_t ny;
97 std::size_t nz;
98
99 ComputeVoxelIndices( copyNo, nx, ny, nz );
100
101 G4ThreeVector trans( (2*nx+1)*fVoxelHalfX - fContainerWallX,
102 (2*ny+1)*fVoxelHalfY - fContainerWallY,
103 (2*nz+1)*fVoxelHalfZ - fContainerWallZ);
104 return trans;
105}

Referenced by 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

Set and Get methods

◆ 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

Flag to skip surface when two voxel have same material or not.

Definition at line 217 of file G4PhantomParameterisation.hh.

◆ fContainerSolid

G4VSolid* G4PhantomParameterisation::fContainerSolid = nullptr
protected

Saves as container solid the parent of the voxels.

Definition at line 208 of file G4PhantomParameterisation.hh.

Referenced by BuildContainerSolid(), BuildContainerSolid(), CheckVoxelsFillContainer(), and GetReplicaNo().

◆ fContainerWallX

◆ fContainerWallY

◆ fContainerWallZ

◆ fMaterialIndices

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

Index in fMaterials that correspond to each voxel.

Definition at line 205 of file G4PhantomParameterisation.hh.

Referenced by G4PartialPhantomParameterisation::GetMaterialIndex(), and GetMaterialIndex().

◆ fMaterials

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

◆ fNoVoxels

std::size_t G4PhantomParameterisation::fNoVoxels = 0
protected

Total number of voxels (for speed-up).

Definition at line 199 of file G4PhantomParameterisation.hh.

◆ fNoVoxelsX

◆ fNoVoxelsXY

std::size_t G4PhantomParameterisation::fNoVoxelsXY = 0
protected

Number of voxels in x times number of voxels in y (for speed-up).

Definition at line 196 of file G4PhantomParameterisation.hh.

Referenced by G4PartialPhantomParameterisation::GetMaterialIndex(), GetMaterialIndex(), and GetReplicaNo().

◆ fNoVoxelsY

◆ fNoVoxelsZ

◆ fVoxelHalfX

◆ fVoxelHalfY

◆ fVoxelHalfZ

◆ kCarTolerance

G4double G4PhantomParameterisation::kCarTolerance
protected

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