Geant4 11.3.0
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 ()=default
 
 ~G4PartialPhantomParameterisation () override=default
 
void ComputeTransformation (const G4int, G4VPhysicalVolume *) const override
 
G4MaterialComputeMaterial (const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=nullptr) override
 
G4int GetReplicaNo (const G4ThreeVector &localPoint, const G4ThreeVector &localDir) override
 
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 () 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)
 
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 ()
 

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 ( )
default

◆ ~G4PartialPhantomParameterisation()

G4PartialPhantomParameterisation::~G4PartialPhantomParameterisation ( )
overridedefault

Member Function Documentation

◆ BuildContainerWalls()

◆ ComputeMaterial()

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

Reimplemented from G4VPVParameterisation.

Definition at line 71 of file G4PartialPhantomParameterisation.cc.

73{
74 CheckCopyNo( copyNo );
75 auto matIndex = GetMaterialIndex(copyNo);
76
77 return fMaterials[ matIndex ];
78}
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
overridevirtual

Implements G4VPVParameterisation.

Definition at line 44 of file G4PartialPhantomParameterisation.cc.

46{
47 // Voxels cannot be rotated, return translation
48 //
49 G4ThreeVector trans = GetTranslation( copyNo );
50 physVol->SetTranslation( trans );
51}
CLHEP::Hep3Vector G4ThreeVector
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 111 of file G4PartialPhantomParameterisation.cc.

113{
114 return fMaterials[GetMaterialIndex(copyNo)];
115}

◆ GetMaterial() [2/2]

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

Definition at line 103 of file G4PartialPhantomParameterisation.cc.

105{
106 return fMaterials[GetMaterialIndex(nx,ny,nz)];
107}

◆ GetMaterialIndex() [1/2]

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

Definition at line 82 of file G4PartialPhantomParameterisation.cc.

84{
85 CheckCopyNo( copyNo );
86
87 if( fMaterialIndices == nullptr ) { return 0; }
88
89 return *(fMaterialIndices+copyNo);
90}

◆ GetMaterialIndex() [2/2]

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

Definition at line 94 of file G4PartialPhantomParameterisation.cc.

96{
97 std::size_t copyNo = nx + fNoVoxelsX*ny + fNoVoxelsXY*nz;
98 return GetMaterialIndex( copyNo );
99}

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

◆ GetReplicaNo()

G4int G4PartialPhantomParameterisation::GetReplicaNo ( const G4ThreeVector & localPoint,
const G4ThreeVector & localDir )
overridevirtual

Reimplemented from G4PhantomParameterisation.

Definition at line 147 of file G4PartialPhantomParameterisation.cc.

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

◆ GetTranslation()

G4ThreeVector G4PartialPhantomParameterisation::GetTranslation ( const G4int copyNo) const

Definition at line 55 of file G4PartialPhantomParameterisation.cc.

57{
58 CheckCopyNo( copyNo );
59
60 std::size_t nx, ny, nz;
61 ComputeVoxelIndices( copyNo, nx, ny, nz );
62
65 (2*nz+1)*fVoxelHalfZ - fContainerWallZ);
66 return trans;
67}

Referenced by ComputeTransformation().

◆ SetFilledIDs()

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

Definition at line 78 of file G4PartialPhantomParameterisation.hh.

79 {
80 fFilledIDs = std::move(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 = std::move(fmins);
86 }

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