Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PhantomParameterisation.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// class G4PhantomParameterisation
27//
28// Class description:
29//
30// Describes regular parameterisations: a set of boxes of equal dimension
31// in the x, y and z dimensions. The G4PVParameterised volume using this
32// class must be placed inside a volume that is completely filled by these
33// boxes.
34
35// History:
36// - Created: P.Arce, May 2007
37//---------------------------------------------------------------------
38#ifndef G4PhantomParameterisation_HH
39#define G4PhantomParameterisation_HH
40
41#include <vector>
42
43#include "G4Types.hh"
45#include "G4AffineTransform.hh"
46
48class G4VTouchable;
49class G4VSolid;
50class G4Material;
51
52// Dummy forward declarations ...
53
54class G4Box;
55class G4Tubs;
56class G4Trd;
57class G4Trap;
58class G4Cons;
59class G4Orb;
60class G4Sphere;
61class G4Ellipsoid;
62class G4Torus;
63class G4Para;
64class G4Hype;
65class G4Polycone;
66class G4Polyhedra;
67
69{
70 public: // with description
71
74
75 virtual void ComputeTransformation(const G4int, G4VPhysicalVolume *) const;
76
78
79 virtual G4Material* ComputeMaterial(const G4int repNo,
80 G4VPhysicalVolume* currentVol,
81 const G4VTouchable* parentTouch=nullptr);
82 // Dummy declarations ...
83
85 const G4VPhysicalVolume*) const {}
87 const G4VPhysicalVolume*) const {}
89 const G4VPhysicalVolume*) const {}
91 const G4VPhysicalVolume*) const {}
93 const G4VPhysicalVolume*) const {}
95 const G4VPhysicalVolume*) const {}
97 const G4VPhysicalVolume*) const {}
99 const G4VPhysicalVolume*) const {}
101 const G4VPhysicalVolume*) const {}
103 const G4VPhysicalVolume*) const {}
105 const G4VPhysicalVolume*) const {}
107 const G4VPhysicalVolume*) const {}
109 const G4VPhysicalVolume*) const {}
110
111 void BuildContainerSolid( G4VPhysicalVolume* pPhysicalVol );
112 void BuildContainerSolid( G4VSolid* pMotherSolid );
113 // Save as container solid the parent of the voxels. Check that the
114 // voxels fill it completely.
115
116 virtual G4int GetReplicaNo( const G4ThreeVector& localPoint,
117 const G4ThreeVector& localDir );
118 // Get the voxel number corresponding to the point in the container
119 // frame. Use 'localDir' to avoid precision problems at the surfaces.
120
121 // Set and Get methods
122
123 inline void SetMaterials(std::vector<G4Material*>& mates );
124
125 inline void SetMaterialIndices( size_t* matInd );
126
127 void SetVoxelDimensions( G4double halfx, G4double halfy, G4double halfz );
128 void SetNoVoxel( size_t nx, size_t ny, size_t nz );
129
130 inline G4double GetVoxelHalfX() const;
131 inline G4double GetVoxelHalfY() const;
132 inline G4double GetVoxelHalfZ() const;
133 inline size_t GetNoVoxelX() const;
134 inline size_t GetNoVoxelY() const;
135 inline size_t GetNoVoxelZ() const;
136 inline size_t GetNoVoxel() const;
137
138 inline std::vector<G4Material*> GetMaterials() const;
139 inline size_t* GetMaterialIndices() const;
141
142 G4ThreeVector GetTranslation(const G4int copyNo ) const;
143
146
147 size_t GetMaterialIndex( size_t nx, size_t ny, size_t nz) const;
148 size_t GetMaterialIndex( size_t copyNo) const;
149
150 G4Material* GetMaterial( size_t nx, size_t ny, size_t nz) const;
151 G4Material* GetMaterial( size_t copyNo ) const;
152
153 void CheckVoxelsFillContainer( G4double contX, G4double contY,
154 G4double contZ ) const;
155 // Check that the voxels fill it completely.
156
157 private:
158
159 void ComputeVoxelIndices(const G4int copyNo, size_t& nx,
160 size_t& ny, size_t& nz ) const;
161 // Convert the copyNo to voxel numbers in x, y and z.
162
163 void CheckCopyNo( const G4int copyNo ) const;
164 // Check that the copy number is within limits.
165
166 protected:
167
169 // Half dimension of voxels (assume they are boxes).
170 size_t fNoVoxelX = 0, fNoVoxelY = 0, fNoVoxelZ = 0;
171 // Number of voxel in x, y and z dimensions.
172 size_t fNoVoxelXY = 0;
173 // Number of voxels in x times number of voxels in y (for speed-up).
174 size_t fNoVoxel = 0;
175 // Total number of voxels (for speed-up).
176
177 std::vector<G4Material*> fMaterials;
178 // List of materials of the voxels.
179 size_t* fMaterialIndices = nullptr;
180 // Index in fMaterials that correspond to each voxel.
181
183 // Save as container solid the parent of the voxels.
184 // Check that the voxels fill it completely.
185
187 // Save position of container wall for speed-up.
188
190 // Relative surface tolerance.
191
193 // Flag to skip surface when two voxel have same material or not
194};
195
196#include "G4PhantomParameterisation.icc"
197
198#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
Definition: G4Box.hh:56
Definition: G4Cons.hh:78
Definition: G4Hype.hh:69
Definition: G4Orb.hh:56
Definition: G4Para.hh:79
void ComputeDimensions(G4Polycone &, const G4int, const G4VPhysicalVolume *) const
void ComputeDimensions(G4Orb &, const G4int, const G4VPhysicalVolume *) const
void SetVoxelDimensions(G4double halfx, G4double halfy, G4double halfz)
virtual G4int GetReplicaNo(const G4ThreeVector &localPoint, const G4ThreeVector &localDir)
void BuildContainerSolid(G4VPhysicalVolume *pPhysicalVol)
void CheckVoxelsFillContainer(G4double contX, G4double contY, G4double contZ) const
void SetMaterials(std::vector< G4Material * > &mates)
virtual G4Material * ComputeMaterial(const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=nullptr)
G4Material * GetMaterial(size_t nx, size_t ny, size_t nz) const
void ComputeDimensions(G4Sphere &, const G4int, const G4VPhysicalVolume *) const
void ComputeDimensions(G4Para &, const G4int, const G4VPhysicalVolume *) const
void ComputeDimensions(G4Hype &, const G4int, const G4VPhysicalVolume *) const
void ComputeDimensions(G4Cons &, const G4int, const G4VPhysicalVolume *) const
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
size_t GetNoVoxelY() const
G4double GetVoxelHalfZ() const
G4VSolid * GetContainerSolid() const
void ComputeDimensions(G4Torus &, const G4int, const G4VPhysicalVolume *) const
void SetSkipEqualMaterials(G4bool skip)
void ComputeDimensions(G4Trap &, const G4int, const G4VPhysicalVolume *) const
void ComputeDimensions(G4Polyhedra &, const G4int, const G4VPhysicalVolume *) const
G4ThreeVector GetTranslation(const G4int copyNo) const
void SetNoVoxel(size_t nx, size_t ny, size_t nz)
G4double GetVoxelHalfY() const
void ComputeDimensions(G4Ellipsoid &, const G4int, const G4VPhysicalVolume *) const
std::vector< G4Material * > GetMaterials() const
G4bool SkipEqualMaterials() const
G4double GetVoxelHalfX() const
size_t GetMaterialIndex(size_t nx, size_t ny, size_t nz) const
void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
virtual void ComputeTransformation(const G4int, G4VPhysicalVolume *) const
size_t GetNoVoxelZ() const
std::vector< G4Material * > fMaterials
void ComputeDimensions(G4Tubs &, const G4int, const G4VPhysicalVolume *) const
size_t GetNoVoxelX() const
size_t GetNoVoxel() const
void ComputeDimensions(G4Trd &, const G4int, const G4VPhysicalVolume *) const
void SetMaterialIndices(size_t *matInd)
size_t * GetMaterialIndices() const
Definition: G4Trd.hh:63
Definition: G4Tubs.hh:75