Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 //
27 // $Id: G4PhantomParameterisation.hh 73435 2013-08-27 11:06:54Z gcosmo $
28 //
29 //
30 // class G4PhantomParameterisation
31 //
32 // Class description:
33 //
34 // Describes regular parameterisations: a set of boxes of equal dimension
35 // in the x, y and z dimensions. The G4PVParameterised volume using this
36 // class must be placed inside a volume that is completely filled by these
37 // boxes.
38 
39 // History:
40 // - Created. P. Arce, May 2007
41 // *********************************************************************
42 
43 #ifndef G4PhantomParameterisation_HH
44 #define G4PhantomParameterisation_HH
45 
46 #include <vector>
47 
48 #include "G4Types.hh"
49 #include "G4VPVParameterisation.hh"
50 #include "G4AffineTransform.hh"
51 
52 class G4VPhysicalVolume;
53 class G4VTouchable;
54 class G4VSolid;
55 class G4Material;
56 
57 // Dummy forward declarations ...
58 
59 class G4Box;
60 class G4Tubs;
61 class G4Trd;
62 class G4Trap;
63 class G4Cons;
64 class G4Orb;
65 class G4Sphere;
66 class G4Ellipsoid;
67 class G4Torus;
68 class G4Para;
69 class G4Hype;
70 class G4Polycone;
71 class G4Polyhedra;
72 
74 {
75  public: // with description
76 
79 
80  virtual void ComputeTransformation(const G4int, G4VPhysicalVolume *) const;
81 
82  virtual G4VSolid* ComputeSolid(const G4int, G4VPhysicalVolume *);
83 
84  virtual G4Material* ComputeMaterial(const G4int repNo,
85  G4VPhysicalVolume *currentVol,
86  const G4VTouchable *parentTouch=0);
87  // Dummy declarations ...
88 
89  void ComputeDimensions (G4Box &, const G4int,
90  const G4VPhysicalVolume*) const {}
92  const G4VPhysicalVolume*) const {}
93  void ComputeDimensions (G4Trd&, const G4int,
94  const G4VPhysicalVolume*) const {}
96  const G4VPhysicalVolume*) const {}
98  const G4VPhysicalVolume*) const {}
99  void ComputeDimensions (G4Orb&, const G4int,
100  const G4VPhysicalVolume*) const {}
102  const G4VPhysicalVolume*) const {}
104  const G4VPhysicalVolume*) const {}
106  const G4VPhysicalVolume*) const {}
108  const G4VPhysicalVolume*) const {}
110  const G4VPhysicalVolume*) const {}
112  const G4VPhysicalVolume*) const {}
114  const G4VPhysicalVolume*) const {}
115 
116  void BuildContainerSolid( G4VPhysicalVolume *pPhysicalVol );
117  void BuildContainerSolid( G4VSolid *pMotherSolid );
118  // Save as container solid the parent of the voxels. Check that the
119  // voxels fill it completely.
120 
121  virtual G4int GetReplicaNo( const G4ThreeVector& localPoint,
122  const G4ThreeVector& localDir );
123  // Get the voxel number corresponding to the point in the container
124  // frame. Use 'localDir' to avoid precision problems at the surfaces.
125 
126  // Set and Get methods
127 
128  inline void SetMaterials(std::vector<G4Material*>& mates );
129 
130  inline void SetMaterialIndices( size_t* matInd );
131 
132  void SetVoxelDimensions( G4double halfx, G4double halfy, G4double halfz );
133  void SetNoVoxel( size_t nx, size_t ny, size_t nz );
134 
135  inline G4double GetVoxelHalfX() const;
136  inline G4double GetVoxelHalfY() const;
137  inline G4double GetVoxelHalfZ() const;
138  inline size_t GetNoVoxelX() const;
139  inline size_t GetNoVoxelY() const;
140  inline size_t GetNoVoxelZ() const;
141  inline size_t GetNoVoxel() const;
142 
143  inline std::vector<G4Material*> GetMaterials() const;
144  inline size_t* GetMaterialIndices() const;
145  inline G4VSolid* GetContainerSolid() const;
146 
147  G4ThreeVector GetTranslation(const G4int copyNo ) const;
148 
149  G4bool SkipEqualMaterials() const;
150  void SetSkipEqualMaterials( G4bool skip );
151 
152  size_t GetMaterialIndex( size_t nx, size_t ny, size_t nz) const;
153  size_t GetMaterialIndex( size_t copyNo) const;
154 
155  G4Material* GetMaterial( size_t nx, size_t ny, size_t nz) const;
156  G4Material* GetMaterial( size_t copyNo ) const;
157 
158  void CheckVoxelsFillContainer( G4double contX, G4double contY,
159  G4double contZ ) const;
160  // Check that the voxels fill it completely.
161 
162  private:
163 
164  void ComputeVoxelIndices(const G4int copyNo, size_t& nx,
165  size_t& ny, size_t& nz ) const;
166  // Convert the copyNo to voxel numbers in x, y and z.
167 
168  void CheckCopyNo( const G4int copyNo ) const;
169  // Check that the copy number is within limits.
170 
171  protected:
172 
174  // Half dimension of voxels (assume they are boxes).
176  // Number of voxel in x, y and z dimensions.
177  size_t fNoVoxelXY;
178  // Number of voxels in x times number of voxels in y (for speed-up).
179  size_t fNoVoxel;
180  // Total number of voxels (for speed-up).
181 
182  std::vector<G4Material*> fMaterials;
183  // List of materials of the voxels.
185  // Index in fMaterials that correspond to each voxel.
186 
188  // Save as container solid the parent of the voxels.
189  // Check that the voxels fill it completely.
190 
192  // Save position of container wall for speed-up.
193 
195  // Relative surface tolerance.
196 
198  // Flag to skip surface when two voxel have same material or not
199 };
200 
201 #include "G4PhantomParameterisation.icc"
202 
203 #endif
size_t * GetMaterialIndices() const
void ComputeDimensions(G4Tubs &, const G4int, const G4VPhysicalVolume *) const
Definition: G4Para.hh:76
void ComputeDimensions(G4Polycone &, const G4int, const G4VPhysicalVolume *) const
Definition: G4Box.hh:63
Definition: G4Tubs.hh:84
virtual void ComputeTransformation(const G4int, G4VPhysicalVolume *) const
Definition: G4Trd.hh:71
void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
virtual G4Material * ComputeMaterial(const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=0)
size_t GetMaterialIndex(size_t nx, size_t ny, size_t nz) const
int G4int
Definition: G4Types.hh:78
virtual G4int GetReplicaNo(const G4ThreeVector &localPoint, const G4ThreeVector &localDir)
void ComputeDimensions(G4Orb &, const G4int, const G4VPhysicalVolume *) const
void BuildContainerSolid(G4VPhysicalVolume *pPhysicalVol)
size_t GetNoVoxelZ() const
G4bool SkipEqualMaterials() const
void SetVoxelDimensions(G4double halfx, G4double halfy, G4double halfz)
Definition: G4Hype.hh:66
void ComputeDimensions(G4Torus &, const G4int, const G4VPhysicalVolume *) const
bool G4bool
Definition: G4Types.hh:79
void ComputeDimensions(G4Sphere &, const G4int, const G4VPhysicalVolume *) const
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
Definition: G4Cons.hh:82
void ComputeDimensions(G4Ellipsoid &, const G4int, const G4VPhysicalVolume *) const
std::vector< G4Material * > GetMaterials() const
void ComputeDimensions(G4Trd &, const G4int, const G4VPhysicalVolume *) const
G4double GetVoxelHalfY() const
Definition: G4Orb.hh:60
void CheckVoxelsFillContainer(G4double contX, G4double contY, G4double contZ) const
G4Material * GetMaterial(size_t nx, size_t ny, size_t nz) const
size_t GetNoVoxelY() const
G4double GetVoxelHalfX() const
size_t GetNoVoxel() const
void SetNoVoxel(size_t nx, size_t ny, size_t nz)
void ComputeDimensions(G4Polyhedra &, const G4int, const G4VPhysicalVolume *) const
G4VSolid * GetContainerSolid() const
G4double GetVoxelHalfZ() const
void SetMaterialIndices(size_t *matInd)
size_t GetNoVoxelX() const
void ComputeDimensions(G4Cons &, const G4int, const G4VPhysicalVolume *) const
void ComputeDimensions(G4Para &, const G4int, const G4VPhysicalVolume *) const
double G4double
Definition: G4Types.hh:76
void SetMaterials(std::vector< G4Material * > &mates)
void ComputeDimensions(G4Hype &, const G4int, const G4VPhysicalVolume *) const
std::vector< G4Material * > fMaterials
void SetSkipEqualMaterials(G4bool skip)
void ComputeDimensions(G4Trap &, const G4int, const G4VPhysicalVolume *) const
G4ThreeVector GetTranslation(const G4int copyNo) const