Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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$
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"
50#include "G4AffineTransform.hh"
51
53class G4VTouchable;
54class G4VSolid;
55class G4Material;
56
57// Dummy forward declarations ...
58
59class G4Box;
60class G4Tubs;
61class G4Trd;
62class G4Trap;
63class G4Cons;
64class G4Orb;
65class G4Sphere;
66class G4Torus;
67class G4Para;
68class G4Hype;
69class G4Polycone;
70class G4Polyhedra;
71
73{
74 public: // with description
75
78
79 virtual void ComputeTransformation(const G4int, G4VPhysicalVolume *) const;
80
82
83 virtual G4Material* ComputeMaterial(const G4int repNo,
84 G4VPhysicalVolume *currentVol,
85 const G4VTouchable *parentTouch=0);
86 // Dummy declarations ...
87
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 {}
111 const G4VPhysicalVolume*) const {}
112
113 void BuildContainerSolid( G4VPhysicalVolume *pPhysicalVol );
114 void BuildContainerSolid( G4VSolid *pMotherSolid );
115 // Save as container solid the parent of the voxels. Check that the
116 // voxels fill it completely.
117
118 virtual G4int GetReplicaNo( const G4ThreeVector& localPoint,
119 const G4ThreeVector& localDir );
120 // Get the voxel number corresponding to the point in the container
121 // frame. Use 'localDir' to avoid precision problems at the surfaces.
122
123 // Set and Get methods
124
125 inline void SetMaterials(std::vector<G4Material*>& mates );
126
127 inline void SetMaterialIndices( size_t* matInd );
128
129 void SetVoxelDimensions( G4double halfx, G4double halfy, G4double halfz );
130 void SetNoVoxel( size_t nx, size_t ny, size_t nz );
131
132 inline G4double GetVoxelHalfX() const;
133 inline G4double GetVoxelHalfY() const;
134 inline G4double GetVoxelHalfZ() const;
135 inline size_t GetNoVoxelX() const;
136 inline size_t GetNoVoxelY() const;
137 inline size_t GetNoVoxelZ() const;
138 inline size_t GetNoVoxel() const;
139
140 inline std::vector<G4Material*> GetMaterials() const;
141 inline size_t* GetMaterialIndices() const;
143
144 G4ThreeVector GetTranslation(const G4int copyNo ) const;
145
148
149 size_t GetMaterialIndex( size_t nx, size_t ny, size_t nz) const;
150 size_t GetMaterialIndex( size_t copyNo) const;
151
152 G4Material* GetMaterial( size_t nx, size_t ny, size_t nz) const;
153 G4Material* GetMaterial( size_t copyNo ) const;
154
155 void CheckVoxelsFillContainer( G4double contX, G4double contY,
156 G4double contZ ) const;
157 // Check that the voxels fill it completely.
158
159 private:
160
161 void ComputeVoxelIndices(const G4int copyNo, size_t& nx,
162 size_t& ny, size_t& nz ) const;
163 // Convert the copyNo to voxel numbers in x, y and z.
164
165 void CheckCopyNo( const G4int copyNo ) const;
166 // Check that the copy number is within limits.
167
168 protected:
169
171 // Half dimension of voxels (assume they are boxes).
173 // Number of voxel in x, y and z dimensions.
175 // Number of voxels in x times number of voxels in y (for speed-up).
176 size_t fNoVoxel;
177 // Total number of voxels (for speed-up).
178
179 std::vector<G4Material*> fMaterials;
180 // List of materials of the voxels.
182 // Index in fMaterials that correspond to each voxel.
183
185 // Save as container solid the parent of the voxels.
186 // Check that the voxels fill it completely.
187
189 // Save position of container wall for speed-up.
190
192 // Relative surface tolerance.
193
195 // Flag to skip surface when two voxel have same material or not
196};
197
198#include "G4PhantomParameterisation.icc"
199
200#endif
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
Definition: G4Box.hh:55
Definition: G4Cons.hh:75
Definition: G4Hype.hh:67
Definition: G4Orb.hh:52
Definition: G4Para.hh:77
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)
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
std::vector< G4Material * > GetMaterials() const
virtual G4Material * ComputeMaterial(const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=0)
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:77