Geant4 11.2.2
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 () 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

Definition at line 68 of file G4PhantomParameterisation.hh.

Constructor & Destructor Documentation

◆ G4PhantomParameterisation()

G4PhantomParameterisation::G4PhantomParameterisation ( )

◆ ~G4PhantomParameterisation()

G4PhantomParameterisation::~G4PhantomParameterisation ( )
overridedefault

Member Function Documentation

◆ BuildContainerSolid() [1/2]

◆ BuildContainerSolid() [2/2]

void G4PhantomParameterisation::BuildContainerSolid ( G4VSolid * pMotherSolid)

Definition at line 66 of file G4PhantomParameterisation.cc.

68{
69 fContainerSolid = pMotherSolid;
73
74 // CheckVoxelsFillContainer();
75}

◆ CheckVoxelsFillContainer()

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

Definition at line 176 of file G4PhantomParameterisation.cc.

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

◆ ComputeDimensions() [1/13]

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

Reimplemented from G4VPVParameterisation.

Definition at line 84 of file G4PhantomParameterisation.hh.

85 {}

◆ ComputeDimensions() [2/13]

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

Reimplemented from G4VPVParameterisation.

Definition at line 92 of file G4PhantomParameterisation.hh.

93 {}

◆ ComputeDimensions() [3/13]

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

Reimplemented from G4VPVParameterisation.

Definition at line 98 of file G4PhantomParameterisation.hh.

99 {}

◆ ComputeDimensions() [4/13]

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

Reimplemented from G4VPVParameterisation.

Definition at line 104 of file G4PhantomParameterisation.hh.

105 {}

◆ ComputeDimensions() [5/13]

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

Reimplemented from G4VPVParameterisation.

Definition at line 94 of file G4PhantomParameterisation.hh.

95 {}

◆ ComputeDimensions() [6/13]

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

Reimplemented from G4VPVParameterisation.

Definition at line 102 of file G4PhantomParameterisation.hh.

103 {}

◆ ComputeDimensions() [7/13]

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

Reimplemented from G4VPVParameterisation.

Definition at line 106 of file G4PhantomParameterisation.hh.

107 {}

◆ ComputeDimensions() [8/13]

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

Reimplemented from G4VPVParameterisation.

Definition at line 108 of file G4PhantomParameterisation.hh.

109 {}

◆ ComputeDimensions() [9/13]

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

Reimplemented from G4VPVParameterisation.

Definition at line 96 of file G4PhantomParameterisation.hh.

97 {}

◆ ComputeDimensions() [10/13]

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

Reimplemented from G4VPVParameterisation.

Definition at line 100 of file G4PhantomParameterisation.hh.

101 {}

◆ ComputeDimensions() [11/13]

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

Reimplemented from G4VPVParameterisation.

Definition at line 90 of file G4PhantomParameterisation.hh.

91 {}

◆ ComputeDimensions() [12/13]

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

Reimplemented from G4VPVParameterisation.

Definition at line 88 of file G4PhantomParameterisation.hh.

89 {}

◆ ComputeDimensions() [13/13]

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

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

Reimplemented from G4VPVParameterisation.

Definition at line 118 of file G4PhantomParameterisation.cc.

120{
121 CheckCopyNo( copyNo );
122 std::size_t matIndex = GetMaterialIndex(copyNo);
123
124 return fMaterials[ matIndex ];
125}
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 ,
G4VPhysicalVolume * pPhysicalVol )
overridevirtual

Reimplemented from G4VPVParameterisation.

Definition at line 110 of file G4PhantomParameterisation.cc.

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

◆ ComputeTransformation()

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

Implements G4VPVParameterisation.

Definition at line 79 of file G4PhantomParameterisation.cc.

81{
82 // Voxels cannot be rotated, return translation
83 //
84 G4ThreeVector trans = GetTranslation( copyNo );
85
86 physVol->SetTranslation( trans );
87}
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 157 of file G4PhantomParameterisation.cc.

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

◆ GetMaterial() [2/2]

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

Definition at line 150 of file G4PhantomParameterisation.cc.

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

Referenced by G4EnergySplitter::SplitEnergyInVolumes().

◆ GetMaterialIndex() [1/2]

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

Definition at line 129 of file G4PhantomParameterisation.cc.

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

◆ GetMaterialIndex() [2/2]

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

Definition at line 140 of file G4PhantomParameterisation.cc.

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

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

Reimplemented in G4PartialPhantomParameterisation.

Definition at line 223 of file G4PhantomParameterisation.cc.

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

◆ GetTranslation()

G4ThreeVector G4PhantomParameterisation::GetTranslation ( const G4int copyNo) const

Definition at line 91 of file G4PhantomParameterisation.cc.

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

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

◆ 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

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