Geant4 9.6.0
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=0)
 
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 (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 (size_t *matInd)
 
void SetVoxelDimensions (G4double halfx, G4double halfy, G4double halfz)
 
void SetNoVoxel (size_t nx, size_t ny, size_t nz)
 
G4double GetVoxelHalfX () const
 
G4double GetVoxelHalfY () const
 
G4double GetVoxelHalfZ () const
 
size_t GetNoVoxelX () const
 
size_t GetNoVoxelY () const
 
size_t GetNoVoxelZ () const
 
size_t GetNoVoxel () const
 
std::vector< G4Material * > GetMaterials () const
 
size_t * GetMaterialIndices () const
 
G4VSolidGetContainerSolid () const
 
G4ThreeVector GetTranslation (const G4int copyNo) const
 
G4bool SkipEqualMaterials () const
 
void SetSkipEqualMaterials (G4bool skip)
 
size_t GetMaterialIndex (size_t nx, size_t ny, size_t nz) const
 
size_t GetMaterialIndex (size_t copyNo) const
 
G4MaterialGetMaterial (size_t nx, size_t ny, size_t nz) const
 
G4MaterialGetMaterial (size_t copyNo) const
 
void CheckVoxelsFillContainer (G4double contX, G4double contY, G4double contZ) const
 
- Public Member Functions inherited from G4VPVParameterisation
 G4VPVParameterisation ()
 
virtual ~G4VPVParameterisation ()
 
virtual void ComputeTransformation (const G4int, G4VPhysicalVolume *) const =0
 
virtual G4VSolidComputeSolid (const G4int, G4VPhysicalVolume *)
 
virtual G4MaterialComputeMaterial (const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=0)
 
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 (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
 
G4double fVoxelHalfY
 
G4double fVoxelHalfZ
 
size_t fNoVoxelX
 
size_t fNoVoxelY
 
size_t fNoVoxelZ
 
size_t fNoVoxelXY
 
size_t fNoVoxel
 
std::vector< G4Material * > fMaterials
 
size_t * fMaterialIndices
 
G4VSolidfContainerSolid
 
G4double fContainerWallX
 
G4double fContainerWallY
 
G4double fContainerWallZ
 
G4double kCarTolerance
 
G4bool bSkipEqualMaterials
 

Detailed Description

Definition at line 72 of file G4PhantomParameterisation.hh.

Constructor & Destructor Documentation

◆ G4PhantomParameterisation()

G4PhantomParameterisation::G4PhantomParameterisation ( )

Definition at line 46 of file G4PhantomParameterisation.cc.

47 : fVoxelHalfX(0.), fVoxelHalfY(0.), fVoxelHalfZ(0.),
52{
54}
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()

◆ ~G4PhantomParameterisation()

G4PhantomParameterisation::~G4PhantomParameterisation ( )

Definition at line 58 of file G4PhantomParameterisation.cc.

59{
60}

Member Function Documentation

◆ BuildContainerSolid() [1/2]

void G4PhantomParameterisation::BuildContainerSolid ( G4VPhysicalVolume pPhysicalVol)

Definition at line 64 of file G4PhantomParameterisation.cc.

66{
67 fContainerSolid = pMotherPhysical->GetLogicalVolume()->GetSolid();
71
72 // CheckVoxelsFillContainer();
73}

◆ BuildContainerSolid() [2/2]

void G4PhantomParameterisation::BuildContainerSolid ( G4VSolid pMotherSolid)

Definition at line 76 of file G4PhantomParameterisation.cc.

78{
79 fContainerSolid = pMotherSolid;
83
84 // CheckVoxelsFillContainer();
85}

◆ CheckVoxelsFillContainer()

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

Definition at line 183 of file G4PhantomParameterisation.cc.

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

◆ ComputeDimensions() [1/12]

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

Reimplemented from G4VPVParameterisation.

Definition at line 88 of file G4PhantomParameterisation.hh.

89 {}

◆ ComputeDimensions() [2/12]

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

Reimplemented from G4VPVParameterisation.

Definition at line 96 of file G4PhantomParameterisation.hh.

97 {}

◆ ComputeDimensions() [3/12]

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

Reimplemented from G4VPVParameterisation.

Definition at line 106 of file G4PhantomParameterisation.hh.

107 {}

◆ ComputeDimensions() [4/12]

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

Reimplemented from G4VPVParameterisation.

Definition at line 98 of file G4PhantomParameterisation.hh.

99 {}

◆ ComputeDimensions() [5/12]

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

Reimplemented from G4VPVParameterisation.

Definition at line 104 of file G4PhantomParameterisation.hh.

105 {}

◆ ComputeDimensions() [6/12]

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

Reimplemented from G4VPVParameterisation.

Definition at line 108 of file G4PhantomParameterisation.hh.

109 {}

◆ ComputeDimensions() [7/12]

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

Reimplemented from G4VPVParameterisation.

Definition at line 110 of file G4PhantomParameterisation.hh.

111 {}

◆ ComputeDimensions() [8/12]

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

Reimplemented from G4VPVParameterisation.

Definition at line 100 of file G4PhantomParameterisation.hh.

101 {}

◆ ComputeDimensions() [9/12]

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

Reimplemented from G4VPVParameterisation.

Definition at line 102 of file G4PhantomParameterisation.hh.

103 {}

◆ ComputeDimensions() [10/12]

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

Reimplemented from G4VPVParameterisation.

Definition at line 94 of file G4PhantomParameterisation.hh.

95 {}

◆ ComputeDimensions() [11/12]

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

Reimplemented from G4VPVParameterisation.

Definition at line 92 of file G4PhantomParameterisation.hh.

93 {}

◆ ComputeDimensions() [12/12]

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

Reimplemented from G4VPVParameterisation.

Definition at line 90 of file G4PhantomParameterisation.hh.

91 {}

◆ ComputeMaterial()

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

Reimplemented from G4VPVParameterisation.

Reimplemented in G4PartialPhantomParameterisation.

Definition at line 128 of file G4PhantomParameterisation.cc.

130{
131 CheckCopyNo( copyNo );
132 size_t matIndex = GetMaterialIndex(copyNo);
133
134 return fMaterials[ matIndex ];
135}
size_t GetMaterialIndex(size_t nx, size_t ny, 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 120 of file G4PhantomParameterisation.cc.

122{
123 return pPhysicalVol->GetLogicalVolume()->GetSolid();
124}
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 89 of file G4PhantomParameterisation.cc.

91{
92 // Voxels cannot be rotated, return translation
93 //
94 G4ThreeVector trans = GetTranslation( copyNo );
95
96 physVol->SetTranslation( trans );
97}
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 ( size_t  copyNo) const

Definition at line 165 of file G4PhantomParameterisation.cc.

166{
167 return fMaterials[GetMaterialIndex(copyNo)];
168}

◆ GetMaterial() [2/2]

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

Definition at line 159 of file G4PhantomParameterisation.cc.

160{
161 return fMaterials[GetMaterialIndex(nx,ny,nz)];
162}

Referenced by G4EnergySplitter::SplitEnergyInVolumes().

◆ GetMaterialIndex() [1/2]

size_t G4PhantomParameterisation::GetMaterialIndex ( size_t  copyNo) const

Definition at line 139 of file G4PhantomParameterisation.cc.

141{
142 CheckCopyNo( copyNo );
143
144 if( !fMaterialIndices ) { return 0; }
145 return *(fMaterialIndices+copyNo);
146}

◆ GetMaterialIndex() [2/2]

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

Definition at line 150 of file G4PhantomParameterisation.cc.

152{
153 size_t copyNo = nx + fNoVoxelX*ny + fNoVoxelXY*nz;
154 return GetMaterialIndex( copyNo );
155}

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

◆ GetMaterialIndices()

size_t * G4PhantomParameterisation::GetMaterialIndices ( ) const
inline

◆ GetMaterials()

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

◆ GetNoVoxel()

size_t G4PhantomParameterisation::GetNoVoxel ( ) const
inline

◆ GetNoVoxelX()

size_t G4PhantomParameterisation::GetNoVoxelX ( ) const
inline

◆ GetNoVoxelY()

size_t G4PhantomParameterisation::GetNoVoxelY ( ) const
inline

◆ GetNoVoxelZ()

size_t G4PhantomParameterisation::GetNoVoxelZ ( ) const
inline

◆ GetReplicaNo()

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

Reimplemented in G4PartialPhantomParameterisation.

Definition at line 230 of file G4PhantomParameterisation.cc.

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

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

◆ GetTranslation()

G4ThreeVector G4PhantomParameterisation::GetTranslation ( const G4int  copyNo) const

Definition at line 101 of file G4PhantomParameterisation.cc.

103{
104 CheckCopyNo( copyNo );
105
106 size_t nx;
107 size_t ny;
108 size_t nz;
109
110 ComputeVoxelIndices( copyNo, nx, ny, nz );
111
112 G4ThreeVector trans( (2*nx+1)*fVoxelHalfX - fContainerWallX,
113 (2*ny+1)*fVoxelHalfY - fContainerWallY,
114 (2*nz+1)*fVoxelHalfZ - fContainerWallZ);
115 return trans;
116}

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 ( size_t *  matInd)
inline

◆ SetMaterials()

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

◆ SetNoVoxel()

void G4PhantomParameterisation::SetNoVoxel ( 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

Member Data Documentation

◆ bSkipEqualMaterials

G4bool G4PhantomParameterisation::bSkipEqualMaterials
protected

Definition at line 194 of file G4PhantomParameterisation.hh.

◆ fContainerSolid

G4VSolid* G4PhantomParameterisation::fContainerSolid
protected

◆ fContainerWallX

◆ fContainerWallY

◆ fContainerWallZ

◆ fMaterialIndices

size_t* G4PhantomParameterisation::fMaterialIndices
protected

◆ fMaterials

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

◆ fNoVoxel

size_t G4PhantomParameterisation::fNoVoxel
protected

Definition at line 176 of file G4PhantomParameterisation.hh.

◆ fNoVoxelX

◆ fNoVoxelXY

size_t G4PhantomParameterisation::fNoVoxelXY
protected

◆ fNoVoxelY

◆ fNoVoxelZ

◆ fVoxelHalfX

◆ fVoxelHalfY

◆ fVoxelHalfZ

◆ kCarTolerance

G4double G4PhantomParameterisation::kCarTolerance
protected

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