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

#include <G4PartialPhantomParameterisation.hh>

+ Inheritance diagram for G4PartialPhantomParameterisation:

Public Member Functions

 G4PartialPhantomParameterisation ()
 
 ~G4PartialPhantomParameterisation ()
 
void ComputeTransformation (const G4int, G4VPhysicalVolume *) const
 
G4MaterialComputeMaterial (const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=nullptr)
 
G4int GetReplicaNo (const G4ThreeVector &localPoint, const G4ThreeVector &localDir)
 
G4ThreeVector GetTranslation (const G4int copyNo) const
 
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 SetFilledIDs (std::multimap< G4int, G4int > fid)
 
void SetFilledMins (std::map< G4int, std::map< G4int, G4int > > fmins)
 
void BuildContainerWalls ()
 
- Public Member Functions inherited from G4PhantomParameterisation
 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
 

Additional Inherited Members

- Protected Attributes inherited from G4PhantomParameterisation
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 52 of file G4PartialPhantomParameterisation.hh.

Constructor & Destructor Documentation

◆ G4PartialPhantomParameterisation()

G4PartialPhantomParameterisation::G4PartialPhantomParameterisation ( )

◆ ~G4PartialPhantomParameterisation()

G4PartialPhantomParameterisation::~G4PartialPhantomParameterisation ( )

Definition at line 51 of file G4PartialPhantomParameterisation.cc.

52{
53}

Member Function Documentation

◆ BuildContainerWalls()

◆ ComputeMaterial()

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

Reimplemented from G4PhantomParameterisation.

Definition at line 83 of file G4PartialPhantomParameterisation.cc.

85{
86 CheckCopyNo( copyNo );
87 auto matIndex = GetMaterialIndex(copyNo);
88
89 return fMaterials[ matIndex ];
90}
std::size_t GetMaterialIndex(std::size_t nx, std::size_t ny, std::size_t nz) const
std::vector< G4Material * > fMaterials

◆ ComputeTransformation()

void G4PartialPhantomParameterisation::ComputeTransformation ( const G4int  copyNo,
G4VPhysicalVolume physVol 
) const
virtual

Reimplemented from G4PhantomParameterisation.

Definition at line 56 of file G4PartialPhantomParameterisation.cc.

58{
59 // Voxels cannot be rotated, return translation
60 //
61 G4ThreeVector trans = GetTranslation( copyNo );
62 physVol->SetTranslation( trans );
63}
G4ThreeVector GetTranslation(const G4int copyNo) const
void SetTranslation(const G4ThreeVector &v)

◆ GetMaterial() [1/2]

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

Definition at line 123 of file G4PartialPhantomParameterisation.cc.

125{
126 return fMaterials[GetMaterialIndex(copyNo)];
127}

◆ GetMaterial() [2/2]

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

Definition at line 115 of file G4PartialPhantomParameterisation.cc.

117{
118 return fMaterials[GetMaterialIndex(nx,ny,nz)];
119}

◆ GetMaterialIndex() [1/2]

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

Definition at line 94 of file G4PartialPhantomParameterisation.cc.

96{
97 CheckCopyNo( copyNo );
98
99 if( fMaterialIndices == nullptr ) { return 0; }
100
101 return *(fMaterialIndices+copyNo);
102}

◆ GetMaterialIndex() [2/2]

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

Definition at line 106 of file G4PartialPhantomParameterisation.cc.

108{
109 std::size_t copyNo = nx + fNoVoxelsX*ny + fNoVoxelsXY*nz;
110 return GetMaterialIndex( copyNo );
111}

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

◆ GetReplicaNo()

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

Reimplemented from G4PhantomParameterisation.

Definition at line 159 of file G4PartialPhantomParameterisation.cc.

161{
162 // Check the voxel numbers corresponding to localPoint
163 // When a particle is on a surface, it may be between -kCarTolerance and
164 // +kCartolerance. By a simple distance as:
165 // G4int nx = G4int( (localPoint.x()+)/fVoxelHalfX/2.);
166 // those between -kCartolerance and 0 will be placed on voxel N-1 and those
167 // between 0 and kCarTolerance on voxel N.
168 // To avoid precision problems place the tracks that are on the surface on
169 // voxel N-1 if they have negative direction and on voxel N if they have
170 // positive direction.
171 // Add +kCarTolerance so that they are first placed on voxel N, and then
172 // if the direction is negative substract 1
173
174 G4double fx = (localPoint.x()+fContainerWallX+kCarTolerance)/(fVoxelHalfX*2.);
175 G4int nx = G4int(fx);
176
177 G4double fy = (localPoint.y()+fContainerWallY+kCarTolerance)/(fVoxelHalfY*2.);
178 G4int ny = G4int(fy);
179
180 G4double fz = (localPoint.z()+fContainerWallZ+kCarTolerance)/(fVoxelHalfZ*2.);
181 G4int nz = G4int(fz);
182
183 // If it is on the surface side, check the direction: if direction is
184 // negative place it on the previous voxel (if direction is positive it is
185 // already in the next voxel...).
186 // Correct also cases where n = -1 or n = fNoVoxels. It is always traced to be
187 // due to multiple scattering: track is entering a voxel but multiple
188 // scattering changes the angle towards outside
189 //
190 if( fx - nx < kCarTolerance/fVoxelHalfX )
191 {
192 if( localDir.x() < 0 )
193 {
194 if( nx != 0 )
195 {
196 nx -= 1;
197 }
198 }
199 else
200 {
201 if( nx == G4int(fNoVoxelsX) )
202 {
203 nx -= 1;
204 }
205 }
206 }
207 if( fy - ny < kCarTolerance/fVoxelHalfY )
208 {
209 if( localDir.y() < 0 )
210 {
211 if( ny != 0 )
212 {
213 ny -= 1;
214 }
215 }
216 else
217 {
218 if( ny == G4int(fNoVoxelsY) )
219 {
220 ny -= 1;
221 }
222 }
223 }
224 if( fz - nz < kCarTolerance/fVoxelHalfZ )
225 {
226 if( localDir.z() < 0 )
227 {
228 if( nz != 0 )
229 {
230 nz -= 1;
231 }
232 }
233 else
234 {
235 if( nz == G4int(fNoVoxelsZ) )
236 {
237 nz -= 1;
238 }
239 }
240 }
241
242 // Check if there are still errors
243 //
244 G4bool isOK = true;
245 if( nx < 0 )
246 {
247 nx = 0;
248 isOK = false;
249 }
250 else if( nx >= G4int(fNoVoxelsX) )
251 {
252 nx = G4int(fNoVoxelsX)-1;
253 isOK = false;
254 }
255 if( ny < 0 )
256 {
257 ny = 0;
258 isOK = false;
259 }
260 else if( ny >= G4int(fNoVoxelsY) )
261 {
262 ny = G4int(fNoVoxelsY)-1;
263 isOK = false;
264 }
265 if( nz < 0 )
266 {
267 nz = 0;
268 isOK = false;
269 }
270 else if( nz >= G4int(fNoVoxelsZ) )
271 {
272 nz = G4int(fNoVoxelsZ)-1;
273 isOK = false;
274 }
275 if( !isOK )
276 {
277 std::ostringstream message;
278 message << "Corrected the copy number! It was negative or too big."
279 << G4endl
280 << " LocalPoint: " << localPoint << G4endl
281 << " LocalDir: " << localDir << G4endl
282 << " Voxel container size: " << fContainerWallX
283 << " " << fContainerWallY << " " << fContainerWallZ << G4endl
284 << " LocalPoint - wall: "
285 << localPoint.x()-fContainerWallX << " "
286 << localPoint.y()-fContainerWallY << " "
287 << localPoint.z()-fContainerWallZ;
288 G4Exception("G4PartialPhantomParameterisation::GetReplicaNo()",
289 "GeomNav1002", JustWarning, message);
290 }
291
292 G4int nyz = G4int(nz*fNoVoxelsY+ny);
293 auto ite = fFilledIDs.cbegin();
294/*
295 for( ite = fFilledIDs.cbegin(); ite != fFilledIDs.cend(); ++ite )
296 {
297 G4cout << " G4PartialPhantomParameterisation::GetReplicaNo filled "
298 << (*ite).first << " , " << (*ite).second << std::endl;
299 }
300*/
301
302 advance(ite,nyz);
303 auto iteant = ite; iteant--;
304 G4int copyNo = (*iteant).first + 1 + ( nx - (*ite).second );
305/*
306 G4cout << " G4PartialPhantomParameterisation::GetReplicaNo getting copyNo "
307 << copyNo << " nyz " << nyz << " (*iteant).first "
308 << (*iteant).first << " (*ite).second " << (*ite).second << G4endl;
309
310 G4cout << " G4PartialPhantomParameterisation::GetReplicaNo " << copyNo
311 << " nx " << nx << " ny " << ny << " nz " << nz
312 << " localPoint " << localPoint << " localDir " << localDir << G4endl;
313*/
314 return copyNo;
315}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
double z() const
double x() const
double y() const

◆ GetTranslation()

G4ThreeVector G4PartialPhantomParameterisation::GetTranslation ( const G4int  copyNo) const

Definition at line 67 of file G4PartialPhantomParameterisation.cc.

69{
70 CheckCopyNo( copyNo );
71
72 std::size_t nx, ny, nz;
73 ComputeVoxelIndices( copyNo, nx, ny, nz );
74
77 (2*nz+1)*fVoxelHalfZ - fContainerWallZ);
78 return trans;
79}

Referenced by ComputeTransformation().

◆ SetFilledIDs()

void G4PartialPhantomParameterisation::SetFilledIDs ( std::multimap< G4int, G4int fid)
inline

Definition at line 78 of file G4PartialPhantomParameterisation.hh.

79 {
80 fFilledIDs = fid;
81 }

◆ SetFilledMins()

void G4PartialPhantomParameterisation::SetFilledMins ( std::map< G4int, std::map< G4int, G4int > >  fmins)
inline

Definition at line 83 of file G4PartialPhantomParameterisation.hh.

84 {
85 fFilledMins = fmins;
86 }

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