Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4UTet.cc
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// Implementation for G4UTet wrapper class
27//
28// 1.11.13 G.Cosmo, CERN
29// --------------------------------------------------------------------
30
31#include "G4Tet.hh"
32#include "G4UTet.hh"
33
34#if ( defined(G4GEOM_USE_USOLIDS) || defined(G4GEOM_USE_PARTIAL_USOLIDS) )
35
36#include "G4AffineTransform.hh"
38#include "G4BoundingEnvelope.hh"
39
40using namespace CLHEP;
41
42////////////////////////////////////////////////////////////////////////
43//
44// Constructor - create a tetrahedron
45// This class is implemented separately from general polyhedra,
46// because the simplex geometry can be computed very quickly,
47// which may become important in situations imported from mesh generators,
48// in which a very large number of G4Tets are created.
49// A Tet has all of its geometrical information precomputed
50//
51G4UTet::G4UTet(const G4String& pName,
52 const G4ThreeVector& anchor,
53 const G4ThreeVector& p1,
54 const G4ThreeVector& p2,
55 const G4ThreeVector& p3, G4bool* degeneracyFlag)
56 : Base_t(pName, U3Vector(anchor.x(),anchor.y(),anchor.z()),
57 U3Vector(p1.x(), p1.y(), p1.z()),
58 U3Vector(p2.x(), p2.y(), p2.z()),
59 U3Vector(p3.x(), p3.y(), p3.z()))
60{
61 // Check for degeneracy
62 G4bool degenerate = CheckDegeneracy(anchor, p1, p2, p3);
63 if(degeneracyFlag != nullptr) *degeneracyFlag = degenerate;
64 else if (degenerate)
65 {
66 G4Exception("G4UTet::G4UTet()", "GeomSolids0002", FatalException,
67 "Degenerate tetrahedron not allowed.");
68 }
69
70 // Set bounding box
71 for (G4int i = 0; i < 3; ++i)
72 {
73 fBmin[i] = std::min(std::min(std::min(anchor[i], p1[i]), p2[i]), p3[i]);
74 fBmax[i] = std::max(std::max(std::max(anchor[i], p1[i]), p2[i]), p3[i]);
75 }
76}
77
78///////////////////////////////////////////////////////////////////////////////
79//
80// Copy constructor
81//
82G4UTet::G4UTet(const G4UTet& rhs)
83 : Base_t(rhs)
84{
85 fBmin = rhs.fBmin;
86 fBmax = rhs.fBmax;
87}
88
89
90///////////////////////////////////////////////////////////////////////////////
91//
92// Assignment operator
93//
94G4UTet& G4UTet::operator = (const G4UTet& rhs)
95{
96 // Check assignment to self
97 if (this == &rhs) { return *this; }
98
99 // Copy base class data
100 Base_t::operator=(rhs);
101
102 // Copy bounding box
103 fBmin = rhs.fBmin;
104 fBmax = rhs.fBmax;
105
106 return *this;
107}
108
109///////////////////////////////////////////////////////////////////////////////
110//
111// Return true if tetrahedron is degenerate
112// Tetrahedron is concidered as degenerate in case if its minimal
113// height is less than the degeneracy tolerance
114//
115G4bool G4UTet::CheckDegeneracy(const G4ThreeVector& p0,
116 const G4ThreeVector& p1,
117 const G4ThreeVector& p2,
118 const G4ThreeVector& p3) const
119{
120 G4double hmin = 4. * kCarTolerance; // degeneracy tolerance
121
122 // Calculate volume
123 G4double vol = std::abs((p1 - p0).cross(p2 - p0).dot(p3 - p0));
124
125 // Calculate face areas squared
126 G4double ss[4];
127 ss[0] = ((p1 - p0).cross(p2 - p0)).mag2();
128 ss[1] = ((p2 - p0).cross(p3 - p0)).mag2();
129 ss[2] = ((p3 - p0).cross(p1 - p0)).mag2();
130 ss[3] = ((p2 - p1).cross(p3 - p1)).mag2();
131
132 // Find face with max area
133 G4int k = 0;
134 for (G4int i = 1; i < 4; ++i) { if (ss[i] > ss[k]) k = i; }
135
136 // Check: vol^2 / s^2 <= hmin^2
137 return (vol*vol <= ss[k]*hmin*hmin);
138}
139
140////////////////////////////////////////////////////////////////////////
141//
142// Dispatch to parameterisation for replication mechanism dimension
143// computation & modification.
144//
145void G4UTet::ComputeDimensions(G4VPVParameterisation*,
146 const G4int,
147 const G4VPhysicalVolume*)
148{
149}
150
151//////////////////////////////////////////////////////////////////////////
152//
153// Make a clone of the object
154//
155G4VSolid* G4UTet::Clone() const
156{
157 return new G4UTet(*this);
158}
159
160///////////////////////////////////////////////////////////////////////////////
161//
162// Modifier
163//
164void G4UTet::SetVertices(const G4ThreeVector& anchor,
165 const G4ThreeVector& p1,
166 const G4ThreeVector& p2,
167 const G4ThreeVector& p3,
168 G4bool* degeneracyFlag)
169{
170 // Check for degeneracy
171 G4bool degenerate = CheckDegeneracy(anchor, p1, p2, p3);
172 if(degeneracyFlag != nullptr) *degeneracyFlag = degenerate;
173 else if (degenerate)
174 {
175 G4Exception("G4UTet::SetVertices()", "GeomSolids0002", FatalException,
176 "Degenerate tetrahedron not allowed.");
177 }
178
179 // Change tetrahedron
180 *this = G4UTet(GetName(), anchor, p1, p2, p3, &degenerate);
181}
182
183///////////////////////////////////////////////////////////////////////////////
184//
185// Accessors
186//
187void G4UTet::GetVertices(G4ThreeVector& anchor,
188 G4ThreeVector& p1,
189 G4ThreeVector& p2,
190 G4ThreeVector& p3) const
191{
192 std::vector<U3Vector> vec(4);
193 Base_t::GetVertices(vec[0], vec[1], vec[2], vec[3]);
194 anchor = G4ThreeVector(vec[0].x(), vec[0].y(), vec[0].z());
195 p1 = G4ThreeVector(vec[1].x(), vec[1].y(), vec[1].z());
196 p2 = G4ThreeVector(vec[2].x(), vec[2].y(), vec[2].z());
197 p3 = G4ThreeVector(vec[3].x(), vec[3].y(), vec[3].z());
198}
199
200std::vector<G4ThreeVector> G4UTet::GetVertices() const
201{
202 std::vector<U3Vector> vec(4);
203 Base_t::GetVertices(vec[0], vec[1], vec[2], vec[3]);
204 std::vector<G4ThreeVector> vertices;
205 for (unsigned int i=0; i<4; ++i)
206 {
207 G4ThreeVector v(vec[i].x(), vec[i].y(), vec[i].z());
208 vertices.push_back(v);
209 }
210 return vertices;
211}
212
213////////////////////////////////////////////////////////////////////////
214//
215// Set bounding box
216//
217void G4UTet::SetBoundingLimits(const G4ThreeVector& pMin,
218 const G4ThreeVector& pMax)
219{
220 G4ThreeVector fVertex[4];
221 GetVertices(fVertex[0], fVertex[1], fVertex[2], fVertex[3]);
222
223 G4int iout[4] = { 0, 0, 0, 0 };
224 for (G4int i = 0; i < 4; ++i)
225 {
226 iout[i] = (G4int)(fVertex[i].x() < pMin.x() ||
227 fVertex[i].y() < pMin.y() ||
228 fVertex[i].z() < pMin.z() ||
229 fVertex[i].x() > pMax.x() ||
230 fVertex[i].y() > pMax.y() ||
231 fVertex[i].z() > pMax.z());
232 }
233 if (iout[0] + iout[1] + iout[2] + iout[3] != 0)
234 {
235 std::ostringstream message;
236 message << "Attempt to set bounding box that does not encapsulate solid: "
237 << GetName() << " !\n"
238 << " Specified bounding box limits:\n"
239 << " pmin: " << pMin << "\n"
240 << " pmax: " << pMax << "\n"
241 << " Tetrahedron vertices:\n"
242 << " anchor " << fVertex[0] << ((iout[0]) != 0 ? " is outside\n" : "\n")
243 << " p1 " << fVertex[1] << ((iout[1]) != 0 ? " is outside\n" : "\n")
244 << " p2 " << fVertex[2] << ((iout[2]) != 0 ? " is outside\n" : "\n")
245 << " p3 " << fVertex[3] << ((iout[3]) != 0 ? " is outside" : "");
246 G4Exception("G4UTet::SetBoundingLimits()", "GeomSolids0002",
247 FatalException, message);
248 }
249 fBmin = pMin;
250 fBmax = pMax;
251}
252
253//////////////////////////////////////////////////////////////////////////
254//
255// Get bounding box
256
257void G4UTet::BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const
258{
259 pMin = fBmin;
260 pMax = fBmax;
261}
262
263//////////////////////////////////////////////////////////////////////////
264//
265// Calculate extent under transform and specified limit
266
267G4bool
268G4UTet::CalculateExtent(const EAxis pAxis,
269 const G4VoxelLimits& pVoxelLimit,
270 const G4AffineTransform& pTransform,
271 G4double& pMin, G4double& pMax) const
272{
273 G4ThreeVector bmin, bmax;
274
275 // Check bounding box (bbox)
276 //
277 BoundingLimits(bmin,bmax);
278 G4BoundingEnvelope bbox(bmin,bmax);
279
280 // Use simple bounding-box to help in the case of complex 3D meshes
281 //
282 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
283
284#if 0
285 // Precise extent computation (disabled by default for this shape)
286 //
287 G4bool exist;
288 if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
289 {
290 return exist = (pMin < pMax) ? true : false;
291 }
292
293 // Set bounding envelope (benv) and calculate extent
294 //
295 std::vector<G4ThreeVector> vec = GetVertices();
296
297 G4ThreeVectorList anchor(1);
298 anchor[0] = vec[0];
299
300 G4ThreeVectorList base(3);
301 base[0] = vec[1];
302 base[1] = vec[2];
303 base[2] = vec[3];
304
305 std::vector<const G4ThreeVectorList *> polygons(2);
306 polygons[0] = &anchor;
307 polygons[1] = &base;
308
309 G4BoundingEnvelope benv(bmin,bmax,polygons);
310 return exists = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
311#endif
312}
313
314////////////////////////////////////////////////////////////////////////
315//
316// CreatePolyhedron
317//
318G4Polyhedron* G4UTet::CreatePolyhedron() const
319{
320 std::vector<U3Vector> vec(4);
321 Base_t::GetVertices(vec[0], vec[1], vec[2], vec[3]);
322
323 G4double xyz[4][3];
324 const G4int faces[4][4] = {{1,3,2,0},{1,4,3,0},{1,2,4,0},{2,3,4,0}};
325 for (unsigned int i=0; i<4; ++i)
326 {
327 xyz[i][0] = vec[i].x();
328 xyz[i][1] = vec[i].y();
329 xyz[i][2] = vec[i].z();
330 }
331
332 auto ph = new G4Polyhedron;
333 ph->createPolyhedron(4,4,xyz,faces);
334 return ph;
335}
336
337#endif // G4GEOM_USE_USOLIDS
const G4double kCarTolerance
std::vector< G4ThreeVector > G4ThreeVectorList
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
double z() const
double x() const
double y() const
G4AffineTransform is a class for geometric affine transformations. It supports efficient arbitrary ro...
G4BoundingEnvelope is a helper class to facilitate calculation of the extent of a solid within the li...
G4VPVParameterisation ia an abstract base class for Parameterisation, able to compute the transformat...
G4VPhysicalVolume is an abstract base class for the representation of a positioned volume....
G4VSolid is an abstract base class for solids, physical shapes that can be tracked through....
Definition G4VSolid.hh:80
G4VoxelLimits represents limitation/restrictions of space, where restrictions are only made perpendic...
G4int createPolyhedron(G4int Nnodes, G4int Nfaces, const G4double xyz[][3], const G4int faces[][4])
EAxis
Definition geomdefs.hh:54
bool exists(std::string const &a_path)
Definition LUPI_file.cc:129