Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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) *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// Fake default constructor - sets only member data and allocates memory
81// for usage restricted to object persistency.
82//
83G4UTet::G4UTet( __void__& a )
84 : Base_t(a)
85{
86}
87
88//////////////////////////////////////////////////////////////////////////
89//
90// Destructor
91//
92G4UTet::~G4UTet()
93{
94}
95
96///////////////////////////////////////////////////////////////////////////////
97//
98// Copy constructor
99//
100G4UTet::G4UTet(const G4UTet& rhs)
101 : Base_t(rhs)
102{
103 fBmin = rhs.fBmin;
104 fBmax = rhs.fBmax;
105}
106
107
108///////////////////////////////////////////////////////////////////////////////
109//
110// Assignment operator
111//
112G4UTet& G4UTet::operator = (const G4UTet& rhs)
113{
114 // Check assignment to self
115 if (this == &rhs) { return *this; }
116
117 // Copy base class data
118 Base_t::operator=(rhs);
119
120 // Copy bounding box
121 fBmin = rhs.fBmin;
122 fBmax = rhs.fBmax;
123
124 return *this;
125}
126
127///////////////////////////////////////////////////////////////////////////////
128//
129// Return true if tetrahedron is degenerate
130// Tetrahedron is concidered as degenerate in case if its minimal
131// height is less than the degeneracy tolerance
132//
133G4bool G4UTet::CheckDegeneracy(const G4ThreeVector& p0,
134 const G4ThreeVector& p1,
135 const G4ThreeVector& p2,
136 const G4ThreeVector& p3) const
137{
138 G4double hmin = 4. * kCarTolerance; // degeneracy tolerance
139
140 // Calculate volume
141 G4double vol = std::abs((p1 - p0).cross(p2 - p0).dot(p3 - p0));
142
143 // Calculate face areas squared
144 G4double ss[4];
145 ss[0] = ((p1 - p0).cross(p2 - p0)).mag2();
146 ss[1] = ((p2 - p0).cross(p3 - p0)).mag2();
147 ss[2] = ((p3 - p0).cross(p1 - p0)).mag2();
148 ss[3] = ((p2 - p1).cross(p3 - p1)).mag2();
149
150 // Find face with max area
151 G4int k = 0;
152 for (G4int i = 1; i < 4; ++i) { if (ss[i] > ss[k]) k = i; }
153
154 // Check: vol^2 / s^2 <= hmin^2
155 return (vol*vol <= ss[k]*hmin*hmin);
156}
157
158////////////////////////////////////////////////////////////////////////
159//
160// Dispatch to parameterisation for replication mechanism dimension
161// computation & modification.
162//
163void G4UTet::ComputeDimensions(G4VPVParameterisation*,
164 const G4int,
165 const G4VPhysicalVolume*)
166{
167}
168
169//////////////////////////////////////////////////////////////////////////
170//
171// Make a clone of the object
172//
173G4VSolid* G4UTet::Clone() const
174{
175 return new G4UTet(*this);
176}
177
178///////////////////////////////////////////////////////////////////////////////
179//
180// Modifier
181//
182void G4UTet::SetVertices(const G4ThreeVector& anchor,
183 const G4ThreeVector& p1,
184 const G4ThreeVector& p2,
185 const G4ThreeVector& p3,
186 G4bool* degeneracyFlag)
187{
188 // Check for degeneracy
189 G4bool degenerate = CheckDegeneracy(anchor, p1, p2, p3);
190 if(degeneracyFlag) *degeneracyFlag = degenerate;
191 else if (degenerate)
192 {
193 G4Exception("G4UTet::SetVertices()", "GeomSolids0002", FatalException,
194 "Degenerate tetrahedron not allowed.");
195 }
196
197 // Change tetrahedron
198 *this = G4UTet(GetName(), anchor, p1, p2, p3, &degenerate);
199}
200
201///////////////////////////////////////////////////////////////////////////////
202//
203// Accessors
204//
205void G4UTet::GetVertices(G4ThreeVector& anchor,
206 G4ThreeVector& p1,
207 G4ThreeVector& p2,
208 G4ThreeVector& p3) const
209{
210 std::vector<U3Vector> vec(4);
211 Base_t::GetVertices(vec[0], vec[1], vec[2], vec[3]);
212 anchor = G4ThreeVector(vec[0].x(), vec[0].y(), vec[0].z());
213 p1 = G4ThreeVector(vec[1].x(), vec[1].y(), vec[1].z());
214 p2 = G4ThreeVector(vec[2].x(), vec[2].y(), vec[2].z());
215 p3 = G4ThreeVector(vec[3].x(), vec[3].y(), vec[3].z());
216}
217
218std::vector<G4ThreeVector> G4UTet::GetVertices() const
219{
220 std::vector<U3Vector> vec(4);
221 Base_t::GetVertices(vec[0], vec[1], vec[2], vec[3]);
222 std::vector<G4ThreeVector> vertices;
223 for (unsigned int i=0; i<4; ++i)
224 {
225 G4ThreeVector v(vec[i].x(), vec[i].y(), vec[i].z());
226 vertices.push_back(v);
227 }
228 return vertices;
229}
230
231////////////////////////////////////////////////////////////////////////
232//
233// Set bounding box
234//
235void G4UTet::SetBoundingLimits(const G4ThreeVector& pMin,
236 const G4ThreeVector& pMax)
237{
238 G4ThreeVector fVertex[4];
239 GetVertices(fVertex[0], fVertex[1], fVertex[2], fVertex[3]);
240
241 G4int iout[4] = { 0, 0, 0, 0 };
242 for (G4int i = 0; i < 4; ++i)
243 {
244 iout[i] = (fVertex[i].x() < pMin.x() ||
245 fVertex[i].y() < pMin.y() ||
246 fVertex[i].z() < pMin.z() ||
247 fVertex[i].x() > pMax.x() ||
248 fVertex[i].y() > pMax.y() ||
249 fVertex[i].z() > pMax.z());
250 }
251 if (iout[0] + iout[1] + iout[2] + iout[3] != 0)
252 {
253 std::ostringstream message;
254 message << "Attempt to set bounding box that does not encapsulate solid: "
255 << GetName() << " !\n"
256 << " Specified bounding box limits:\n"
257 << " pmin: " << pMin << "\n"
258 << " pmax: " << pMax << "\n"
259 << " Tetrahedron vertices:\n"
260 << " anchor " << fVertex[0] << ((iout[0]) ? " is outside\n" : "\n")
261 << " p1 " << fVertex[1] << ((iout[1]) ? " is outside\n" : "\n")
262 << " p2 " << fVertex[2] << ((iout[2]) ? " is outside\n" : "\n")
263 << " p3 " << fVertex[3] << ((iout[3]) ? " is outside" : "");
264 G4Exception("G4UTet::SetBoundingLimits()", "GeomSolids0002",
265 FatalException, message);
266 }
267 fBmin = pMin;
268 fBmax = pMax;
269}
270
271//////////////////////////////////////////////////////////////////////////
272//
273// Get bounding box
274
275void G4UTet::BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const
276{
277 pMin = fBmin;
278 pMax = fBmax;
279}
280
281//////////////////////////////////////////////////////////////////////////
282//
283// Calculate extent under transform and specified limit
284
285G4bool
286G4UTet::CalculateExtent(const EAxis pAxis,
287 const G4VoxelLimits& pVoxelLimit,
288 const G4AffineTransform& pTransform,
289 G4double& pMin, G4double& pMax) const
290{
291 G4ThreeVector bmin, bmax;
292
293 // Check bounding box (bbox)
294 //
295 BoundingLimits(bmin,bmax);
296 G4BoundingEnvelope bbox(bmin,bmax);
297
298 // Use simple bounding-box to help in the case of complex 3D meshes
299 //
300 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
301
302#if 0
303 // Precise extent computation (disabled by default for this shape)
304 //
305 G4bool exist;
306 if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
307 {
308 return exist = (pMin < pMax) ? true : false;
309 }
310
311 // Set bounding envelope (benv) and calculate extent
312 //
313 std::vector<G4ThreeVector> vec = GetVertices();
314
315 G4ThreeVectorList anchor(1);
316 anchor[0] = vec[0];
317
318 G4ThreeVectorList base(3);
319 base[0] = vec[1];
320 base[1] = vec[2];
321 base[2] = vec[3];
322
323 std::vector<const G4ThreeVectorList *> polygons(2);
324 polygons[0] = &anchor;
325 polygons[1] = &base;
326
327 G4BoundingEnvelope benv(bmin,bmax,polygons);
328 return exists = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
329#endif
330}
331
332////////////////////////////////////////////////////////////////////////
333//
334// CreatePolyhedron
335//
336G4Polyhedron* G4UTet::CreatePolyhedron() const
337{
338 std::vector<U3Vector> vec(4);
339 Base_t::GetVertices(vec[0], vec[1], vec[2], vec[3]);
340
341 G4double xyz[4][3];
342 const G4int faces[4][4] = {{1,3,2,0},{1,4,3,0},{1,2,4,0},{2,3,4,0}};
343 for (unsigned int i=0; i<4; ++i)
344 {
345 xyz[i][0] = vec[i].x();
346 xyz[i][1] = vec[i].y();
347 xyz[i][2] = vec[i].z();
348 }
349
350 G4Polyhedron* ph = new G4Polyhedron;
351 ph->createPolyhedron(4,4,xyz,faces);
352 return ph;
353}
354
355#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)
Definition: G4Exception.cc:59
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
G4int createPolyhedron(G4int Nnodes, G4int Nfaces, const G4double xyz[][3], const G4int faces[][4])
EAxis
Definition: geomdefs.hh:54
Definition: DoubConv.h:17