Geant4 10.7.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=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 (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=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
 
size_t fNoVoxelX = 0
 
size_t fNoVoxelY = 0
 
size_t fNoVoxelZ = 0
 
size_t fNoVoxelXY = 0
 
size_t fNoVoxel = 0
 
std::vector< G4Material * > fMaterials
 
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]

◆ 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 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-fNoVoxelX*fVoxelHalfX) >= toleranceForError
190 || std::fabs(contY-fNoVoxelY*fVoxelHalfY) >= toleranceForError
191 || std::fabs(contZ-fNoVoxelZ*fVoxelHalfZ) >= toleranceForError )
192 {
193 std::ostringstream message;
194 message << "Voxels do not fully fill the container: "
196 << " DiffX= " << contX-fNoVoxelX*fVoxelHalfX << G4endl
197 << " DiffY= " << contY-fNoVoxelY*fVoxelHalfY << G4endl
198 << " DiffZ= " << contZ-fNoVoxelZ*fVoxelHalfZ << G4endl
199 << " Maximum difference is: " << toleranceForError;
200 G4Exception("G4PhantomParameterisation::CheckVoxelsFillContainer()",
201 "GeomNav0002", FatalException, message);
202
203 }
204 else if( std::fabs(contX-fNoVoxelX*fVoxelHalfX) >= toleranceForWarning
205 || std::fabs(contY-fNoVoxelY*fVoxelHalfY) >= toleranceForWarning
206 || std::fabs(contZ-fNoVoxelZ*fVoxelHalfZ) >= toleranceForWarning )
207 {
208 std::ostringstream message;
209 message << "Voxels do not fully fill the container: "
211 << " DiffX= " << contX-fNoVoxelX*fVoxelHalfX << G4endl
212 << " DiffY= " << contY-fNoVoxelY*fVoxelHalfY << G4endl
213 << " DiffZ= " << contZ-fNoVoxelZ*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)
Definition: G4Exception.cc:35
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 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

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 ( 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 ( 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}

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}

◆ 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 + fNoVoxelX*ny + fNoVoxelXY*nz;
145 return GetMaterialIndex( copyNo );
146}

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

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 = true
protected

Definition at line 192 of file G4PhantomParameterisation.hh.

◆ fContainerSolid

G4VSolid* G4PhantomParameterisation::fContainerSolid = nullptr
protected

◆ fContainerWallX

◆ fContainerWallY

◆ fContainerWallZ

◆ fMaterialIndices

size_t* G4PhantomParameterisation::fMaterialIndices = nullptr
protected

◆ fMaterials

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

◆ fNoVoxel

size_t G4PhantomParameterisation::fNoVoxel = 0
protected

Definition at line 174 of file G4PhantomParameterisation.hh.

◆ fNoVoxelX

◆ fNoVoxelXY

size_t G4PhantomParameterisation::fNoVoxelXY = 0
protected

◆ fNoVoxelY

◆ fNoVoxelZ

◆ fVoxelHalfX

◆ fVoxelHalfY

◆ fVoxelHalfZ

◆ kCarTolerance

G4double G4PhantomParameterisation::kCarTolerance
protected

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