Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Tet.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 intellectual property of the *
19// * Vanderbilt University Free Electron Laser Center *
20// * Vanderbilt University, Nashville, TN, USA *
21// * Development supported by: *
22// * United States MFEL program under grant FA9550-04-1-0045 *
23// * and NASA under contract number NNG04CT05P *
24// * Written by Marcus H. Mendenhall and Robert A. Weller. *
25// * *
26// * Contributed to the Geant4 Core, January, 2005. *
27// * *
28// ********************************************************************
29//
30// Implementation for G4Tet class
31//
32// 03.09.2004 - Marcus Mendenhall, created
33// 08.01.2020 - Evgueni Tcherniaev, complete revision, speed up
34// --------------------------------------------------------------------
35
36#include "G4Tet.hh"
37
38#if !defined(G4GEOM_USE_UTET)
39
40#include "G4VoxelLimits.hh"
41#include "G4AffineTransform.hh"
42#include "G4BoundingEnvelope.hh"
43
45
46#include "G4QuickRand.hh"
47
48#include "G4VGraphicsScene.hh"
49#include "G4Polyhedron.hh"
50#include "G4VisExtent.hh"
51
52#include "G4AutoLock.hh"
53
54namespace
55{
56 G4Mutex polyhedronMutex = G4MUTEX_INITIALIZER;
57}
58
59using namespace CLHEP;
60
61////////////////////////////////////////////////////////////////////////
62//
63// Constructor - create a tetrahedron
64// A Tet has all of its geometrical information precomputed
65//
67 const G4ThreeVector& p0,
68 const G4ThreeVector& p1,
69 const G4ThreeVector& p2,
70 const G4ThreeVector& p3, G4bool* degeneracyFlag)
71 : G4VSolid(pName)
72{
73 // Check for degeneracy
74 G4bool degenerate = CheckDegeneracy(p0, p1, p2, p3);
75 if (degeneracyFlag != nullptr)
76 {
77 *degeneracyFlag = degenerate;
78 }
79 else if (degenerate)
80 {
81 std::ostringstream message;
82 message << "Degenerate tetrahedron: " << GetName() << " !\n"
83 << " anchor: " << p0 << "\n"
84 << " p1 : " << p1 << "\n"
85 << " p2 : " << p2 << "\n"
86 << " p3 : " << p3 << "\n"
87 << " volume: "
88 << std::abs((p1 - p0).cross(p2 - p0).dot(p3 - p0))/6.;
89 G4Exception("G4Tet::G4Tet()", "GeomSolids0002", FatalException, message);
90 }
91
92 // Define surface thickness
93 halfTolerance = 0.5 * kCarTolerance;
94
95 // Set data members
96 Initialize(p0, p1, p2, p3);
97}
98
99////////////////////////////////////////////////////////////////////////
100//
101// Fake default constructor - sets only member data and allocates memory
102// for usage restricted to object persistency.
103//
104G4Tet::G4Tet( __void__& a )
105 : G4VSolid(a)
106{
107}
108
109////////////////////////////////////////////////////////////////////////
110//
111// Destructor
112//
114{
115 delete fpPolyhedron; fpPolyhedron = nullptr;
116}
117
118////////////////////////////////////////////////////////////////////////
119//
120// Copy constructor
121//
123 : G4VSolid(rhs)
124{
125 halfTolerance = rhs.halfTolerance;
126 fCubicVolume = rhs.fCubicVolume;
127 fSurfaceArea = rhs.fSurfaceArea;
128 for (G4int i = 0; i < 4; ++i) { fVertex[i] = rhs.fVertex[i]; }
129 for (G4int i = 0; i < 4; ++i) { fNormal[i] = rhs.fNormal[i]; }
130 for (G4int i = 0; i < 4; ++i) { fDist[i] = rhs.fDist[i]; }
131 for (G4int i = 0; i < 4; ++i) { fArea[i] = rhs.fArea[i]; }
132 fBmin = rhs.fBmin;
133 fBmax = rhs.fBmax;
134}
135
136////////////////////////////////////////////////////////////////////////
137//
138// Assignment operator
139//
141{
142 // Check assignment to self
143 //
144 if (this == &rhs) { return *this; }
145
146 // Copy base class data
147 //
149
150 // Copy data
151 //
152 halfTolerance = rhs.halfTolerance;
153 fCubicVolume = rhs.fCubicVolume;
154 fSurfaceArea = rhs.fSurfaceArea;
155 for (G4int i = 0; i < 4; ++i) { fVertex[i] = rhs.fVertex[i]; }
156 for (G4int i = 0; i < 4; ++i) { fNormal[i] = rhs.fNormal[i]; }
157 for (G4int i = 0; i < 4; ++i) { fDist[i] = rhs.fDist[i]; }
158 for (G4int i = 0; i < 4; ++i) { fArea[i] = rhs.fArea[i]; }
159 fBmin = rhs.fBmin;
160 fBmax = rhs.fBmax;
161 fRebuildPolyhedron = false;
162 delete fpPolyhedron; fpPolyhedron = nullptr;
163
164 return *this;
165}
166
167////////////////////////////////////////////////////////////////////////
168//
169// Return true if tetrahedron is degenerate
170// Tetrahedron is considered as degenerate in case its minimal
171// height is less than the degeneracy tolerance
172//
174 const G4ThreeVector& p1,
175 const G4ThreeVector& p2,
176 const G4ThreeVector& p3) const
177{
178 G4double hmin = 4. * kCarTolerance; // degeneracy tolerance
179
180 // Calculate volume
181 G4double vol = std::abs((p1 - p0).cross(p2 - p0).dot(p3 - p0));
182
183 // Calculate face areas squared
184 G4double ss[4];
185 ss[0] = ((p1 - p0).cross(p2 - p0)).mag2();
186 ss[1] = ((p2 - p0).cross(p3 - p0)).mag2();
187 ss[2] = ((p3 - p0).cross(p1 - p0)).mag2();
188 ss[3] = ((p2 - p1).cross(p3 - p1)).mag2();
189
190 // Find face with max area
191 G4int k = 0;
192 for (G4int i = 1; i < 4; ++i)
193 {
194 if (ss[i] > ss[k]) { k = i; }
195 }
196
197 // Check: vol^2 / s^2 <= hmin^2
198 return (vol*vol <= ss[k]*hmin*hmin);
199}
200
201////////////////////////////////////////////////////////////////////////
202//
203// Set data members
204//
205void G4Tet::Initialize(const G4ThreeVector& p0,
206 const G4ThreeVector& p1,
207 const G4ThreeVector& p2,
208 const G4ThreeVector& p3)
209{
210 // Set vertices
211 fVertex[0] = p0;
212 fVertex[1] = p1;
213 fVertex[2] = p2;
214 fVertex[3] = p3;
215
216 G4ThreeVector norm[4];
217 norm[0] = (p2 - p0).cross(p1 - p0);
218 norm[1] = (p3 - p0).cross(p2 - p0);
219 norm[2] = (p1 - p0).cross(p3 - p0);
220 norm[3] = (p2 - p1).cross(p3 - p1);
221 G4double volume = norm[0].dot(p3 - p0);
222 if (volume > 0.)
223 {
224 for (auto & i : norm) { i = -i; }
225 }
226
227 // Set normals to face planes
228 for (G4int i = 0; i < 4; ++i) { fNormal[i] = norm[i].unit(); }
229
230 // Set distances to planes
231 for (G4int i = 0; i < 3; ++i) { fDist[i] = fNormal[i].dot(p0); }
232 fDist[3] = fNormal[3].dot(p1);
233
234 // Set face areas
235 for (G4int i = 0; i < 4; ++i) { fArea[i] = 0.5*norm[i].mag(); }
236
237 // Set bounding box
238 for (G4int i = 0; i < 3; ++i)
239 {
240 fBmin[i] = std::min(std::min(std::min(p0[i], p1[i]), p2[i]), p3[i]);
241 fBmax[i] = std::max(std::max(std::max(p0[i], p1[i]), p2[i]), p3[i]);
242 }
243
244 // Set volume and surface area
245 fCubicVolume = std::abs(volume)/6.;
246 fSurfaceArea = fArea[0] + fArea[1] + fArea[2] + fArea[3];
247}
248
249////////////////////////////////////////////////////////////////////////
250//
251// Set vertices
252//
254 const G4ThreeVector& p1,
255 const G4ThreeVector& p2,
256 const G4ThreeVector& p3, G4bool* degeneracyFlag)
257{
258 // Check for degeneracy
259 G4bool degenerate = CheckDegeneracy(p0, p1, p2, p3);
260 if (degeneracyFlag != nullptr)
261 {
262 *degeneracyFlag = degenerate;
263 }
264 else if (degenerate)
265 {
266 std::ostringstream message;
267 message << "Degenerate tetrahedron is not permitted: " << GetName() << " !\n"
268 << " anchor: " << p0 << "\n"
269 << " p1 : " << p1 << "\n"
270 << " p2 : " << p2 << "\n"
271 << " p3 : " << p3 << "\n"
272 << " volume: "
273 << std::abs((p1 - p0).cross(p2 - p0).dot(p3 - p0))/6.;
274 G4Exception("G4Tet::SetVertices()", "GeomSolids0002",
275 FatalException, message);
276 }
277
278 // Set data members
279 Initialize(p0, p1, p2, p3);
280
281 // Set flag to rebuild polyhedron
282 fRebuildPolyhedron = true;
283}
284
285////////////////////////////////////////////////////////////////////////
286//
287// Return four vertices
288//
290 G4ThreeVector& p1,
291 G4ThreeVector& p2,
292 G4ThreeVector& p3) const
293{
294 p0 = fVertex[0];
295 p1 = fVertex[1];
296 p2 = fVertex[2];
297 p3 = fVertex[3];
298}
299
300////////////////////////////////////////////////////////////////////////
301//
302// Return std::vector of vertices
303//
304std::vector<G4ThreeVector> G4Tet::GetVertices() const
305{
306 std::vector<G4ThreeVector> vertices(4);
307 for (G4int i = 0; i < 4; ++i) { vertices[i] = fVertex[i]; }
308 return vertices;
309}
310
311////////////////////////////////////////////////////////////////////////
312//
313// Dispatch to parameterisation for replication mechanism dimension
314// computation & modification.
315//
321
322////////////////////////////////////////////////////////////////////////
323//
324// Set bounding box
325//
327 const G4ThreeVector& pMax)
328{
329 G4int iout[4] = { 0, 0, 0, 0 };
330 for (G4int i = 0; i < 4; ++i)
331 {
332 iout[i] = (G4int)(fVertex[i].x() < pMin.x() ||
333 fVertex[i].y() < pMin.y() ||
334 fVertex[i].z() < pMin.z() ||
335 fVertex[i].x() > pMax.x() ||
336 fVertex[i].y() > pMax.y() ||
337 fVertex[i].z() > pMax.z());
338 }
339 if (iout[0] + iout[1] + iout[2] + iout[3] != 0)
340 {
341 std::ostringstream message;
342 message << "Attempt to set bounding box that does not encapsulate solid: "
343 << GetName() << " !\n"
344 << " Specified bounding box limits:\n"
345 << " pmin: " << pMin << "\n"
346 << " pmax: " << pMax << "\n"
347 << " Tetrahedron vertices:\n"
348 << " anchor " << fVertex[0] << ((iout[0]) != 0 ? " is outside\n" : "\n")
349 << " p1 " << fVertex[1] << ((iout[1]) != 0 ? " is outside\n" : "\n")
350 << " p2 " << fVertex[2] << ((iout[2]) != 0 ? " is outside\n" : "\n")
351 << " p3 " << fVertex[3] << ((iout[3]) != 0 ? " is outside" : "");
352 G4Exception("G4Tet::SetBoundingLimits()", "GeomSolids0002",
353 FatalException, message);
354 }
355 fBmin = pMin;
356 fBmax = pMax;
357}
358
359////////////////////////////////////////////////////////////////////////
360//
361// Return bounding box
362//
364{
365 pMin = fBmin;
366 pMax = fBmax;
367}
368
369////////////////////////////////////////////////////////////////////////
370//
371// Calculate extent under transform and specified limit
372//
374 const G4VoxelLimits& pVoxelLimit,
375 const G4AffineTransform& pTransform,
376 G4double& pMin, G4double& pMax) const
377{
378 G4ThreeVector bmin, bmax;
379
380 // Check bounding box (bbox)
381 //
382 BoundingLimits(bmin,bmax);
383 G4BoundingEnvelope bbox(bmin,bmax);
384
385 // Use simple bounding-box to help in the case of complex 3D meshes
386 //
387 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
388
389#if 0
390 // Precise extent computation (disabled by default for this shape)
391 //
392 G4bool exist;
393 if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
394 {
395 return exist = (pMin < pMax) ? true : false;
396 }
397
398 // Set bounding envelope (benv) and calculate extent
399 //
400 std::vector<G4ThreeVector> vec = GetVertices();
401
402 G4ThreeVectorList anchor(1);
403 anchor[0].set(vec[0].x(),vec[0].y(),vec[0].z());
404
405 G4ThreeVectorList base(3);
406 base[0].set(vec[1].x(),vec[1].y(),vec[1].z());
407 base[1].set(vec[2].x(),vec[2].y(),vec[2].z());
408 base[2].set(vec[3].x(),vec[3].y(),vec[3].z());
409
410 std::vector<const G4ThreeVectorList *> polygons(2);
411 polygons[0] = &anchor;
412 polygons[1] = &base;
413
414 G4BoundingEnvelope benv(bmin,bmax,polygons);
415 return exists = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
416#endif
417}
418
419////////////////////////////////////////////////////////////////////////
420//
421// Return whether point inside/outside/on surface
422//
424{
425 G4double dd[4];
426 for (G4int i = 0; i < 4; ++i) { dd[i] = fNormal[i].dot(p) - fDist[i]; }
427
428 G4double dist = std::max(std::max(std::max(dd[0], dd[1]), dd[2]), dd[3]);
429 return (dist <= -halfTolerance) ?
430 kInside : ((dist <= halfTolerance) ? kSurface : kOutside);
431}
432
433////////////////////////////////////////////////////////////////////////
434//
435// Return unit normal to surface at p
436//
438{
439 G4double k[4];
440 for (G4int i = 0; i < 4; ++i)
441 {
442 k[i] = (G4double)(std::abs(fNormal[i].dot(p) - fDist[i]) <= halfTolerance);
443 }
444 G4double nsurf = k[0] + k[1] + k[2] + k[3];
445 G4ThreeVector norm =
446 k[0]*fNormal[0] + k[1]*fNormal[1] + k[2]*fNormal[2] + k[3]*fNormal[3];
447
448 if (nsurf == 1.)
449 {
450 return norm;
451 }
452 if (nsurf > 1.)
453 {
454 return norm.unit(); // edge or vertex
455 }
456
457#ifdef G4SPECSDEBUG
458 std::ostringstream message;
459 G4long oldprc = message.precision(16);
460 message << "Point p is not on surface (!?) of solid: "
461 << GetName() << "\n";
462 message << "Position:\n";
463 message << " p.x() = " << p.x()/mm << " mm\n";
464 message << " p.y() = " << p.y()/mm << " mm\n";
465 message << " p.z() = " << p.z()/mm << " mm";
466 G4cout.precision(oldprc);
467 G4Exception("G4Tet::SurfaceNormal(p)", "GeomSolids1002",
468 JustWarning, message );
469 DumpInfo();
470#endif
471 return ApproxSurfaceNormal(p);
472}
473
474////////////////////////////////////////////////////////////////////////
475//
476// Find surface nearest to point and return corresponding normal
477// This method normally should not be called
478//
479G4ThreeVector G4Tet::ApproxSurfaceNormal(const G4ThreeVector& p) const
480{
481 G4double dist = -DBL_MAX;
482 G4int iside = 0;
483 for (G4int i = 0; i < 4; ++i)
484 {
485 G4double d = fNormal[i].dot(p) - fDist[i];
486 if (d > dist) { dist = d; iside = i; }
487 }
488 return fNormal[iside];
489}
490
491////////////////////////////////////////////////////////////////////////
492//
493// Calculate distance to surface from outside,
494// return kInfinity if no intersection
495//
497 const G4ThreeVector& v) const
498{
499 G4double tin = -DBL_MAX, tout = DBL_MAX;
500 for (G4int i = 0; i < 4; ++i)
501 {
502 G4double cosa = fNormal[i].dot(v);
503 G4double dist = fNormal[i].dot(p) - fDist[i];
504 if (dist >= -halfTolerance)
505 {
506 if (cosa >= 0.) { return kInfinity; }
507 tin = std::max(tin, -dist/cosa);
508 }
509 else if (cosa > 0.)
510 {
511 tout = std::min(tout, -dist/cosa);
512 }
513 }
514
515 return (tout - tin <= halfTolerance) ?
516 kInfinity : ((tin < halfTolerance) ? 0. : tin);
517}
518
519////////////////////////////////////////////////////////////////////////
520//
521// Estimate safety distance to surface from outside
522//
524{
525 G4double dd[4];
526 for (G4int i = 0; i < 4; ++i) { dd[i] = fNormal[i].dot(p) - fDist[i]; }
527
528 G4double dist = std::max(std::max(std::max(dd[0], dd[1]), dd[2]), dd[3]);
529 return (dist > 0.) ? dist : 0.;
530}
531
532////////////////////////////////////////////////////////////////////////
533//
534// Calcluate distance to surface from inside
535//
537 const G4ThreeVector& v,
538 const G4bool calcNorm,
539 G4bool* validNorm,
540 G4ThreeVector* n) const
541{
542 // Calculate distances and cosines
543 G4double cosa[4], dist[4];
544 G4int ind[4] = {0}, nside = 0;
545 for (G4int i = 0; i < 4; ++i)
546 {
547 G4double tmp = fNormal[i].dot(v);
548 cosa[i] = tmp;
549 ind[nside] = (G4int)(tmp > 0) * i;
550 nside += (G4int)(tmp > 0);
551 dist[i] = fNormal[i].dot(p) - fDist[i];
552 }
553
554 // Find intersection (in most of cases nside == 1)
555 G4double tout = DBL_MAX;
556 G4int iside = 0;
557 for (G4int i = 0; i < nside; ++i)
558 {
559 G4int k = ind[i];
560 // Check: leaving the surface
561 if (dist[k] >= -halfTolerance) { tout = 0.; iside = k; break; }
562 // Compute distance to intersection
563 G4double tmp = -dist[k]/cosa[k];
564 if (tmp < tout) { tout = tmp; iside = k; }
565 }
566
567 // Set normal, if required, and return distance to out
568 if (calcNorm)
569 {
570 *validNorm = true;
571 *n = fNormal[iside];
572 }
573 return tout;
574}
575
576////////////////////////////////////////////////////////////////////////
577//
578// Calculate safety distance to surface from inside
579//
581{
582 G4double dd[4];
583 for (G4int i = 0; i < 4; ++i) { dd[i] = fDist[i] - fNormal[i].dot(p); }
584
585 G4double dist = std::min(std::min(std::min(dd[0], dd[1]), dd[2]), dd[3]);
586 return (dist > 0.) ? dist : 0.;
587}
588
589////////////////////////////////////////////////////////////////////////
590//
591// GetEntityType
592//
594{
595 return {"G4Tet"};
596}
597
598////////////////////////////////////////////////////////////////////////
599//
600// IsFaceted
601//
603{
604 return true;
605}
606
607////////////////////////////////////////////////////////////////////////
608//
609// Make a clone of the object
610//
612{
613 return new G4Tet(*this);
614}
615
616////////////////////////////////////////////////////////////////////////
617//
618// Stream object contents to an output stream
619//
620std::ostream& G4Tet::StreamInfo(std::ostream& os) const
621{
622 G4long oldprc = os.precision(16);
623 os << "-----------------------------------------------------------\n"
624 << " *** Dump for solid - " << GetName() << " ***\n"
625 << " ===================================================\n"
626 << " Solid type: " << GetEntityType() << "\n"
627 << " Parameters: \n"
628 << " anchor: " << fVertex[0]/mm << " mm\n"
629 << " p1 : " << fVertex[1]/mm << " mm\n"
630 << " p2 : " << fVertex[2]/mm << " mm\n"
631 << " p3 : " << fVertex[3]/mm << " mm\n"
632 << "-----------------------------------------------------------\n";
633 os.precision(oldprc);
634 return os;
635}
636
637////////////////////////////////////////////////////////////////////////
638//
639// Return random point on the surface
640//
642{
643 constexpr G4int iface[4][3] = { {0,1,2}, {0,2,3}, {0,3,1}, {1,2,3} };
644
645 // Select face
646 G4double select = fSurfaceArea*G4QuickRand();
647 G4int i = 0;
648 i += (G4int)(select > fArea[0]);
649 i += (G4int)(select > fArea[0] + fArea[1]);
650 i += (G4int)(select > fArea[0] + fArea[1] + fArea[2]);
651
652 // Set selected triangle
653 G4ThreeVector p0 = fVertex[iface[i][0]];
654 G4ThreeVector e1 = fVertex[iface[i][1]] - p0;
655 G4ThreeVector e2 = fVertex[iface[i][2]] - p0;
656
657 // Return random point
658 G4double r1 = G4QuickRand();
659 G4double r2 = G4QuickRand();
660 return (r1 + r2 > 1.) ?
661 p0 + e1*(1. - r1) + e2*(1. - r2) : p0 + e1*r1 + e2*r2;
662}
663
664////////////////////////////////////////////////////////////////////////
665//
666// Return volume of the tetrahedron
667//
669{
670 return fCubicVolume;
671}
672
673////////////////////////////////////////////////////////////////////////
674//
675// Return surface area of the tetrahedron
676//
678{
679 return fSurfaceArea;
680}
681
682////////////////////////////////////////////////////////////////////////
683//
684// Methods for visualisation
685//
687{
688 scene.AddSolid (*this);
689}
690
691////////////////////////////////////////////////////////////////////////
692//
693// Return VisExtent
694//
696{
697 return { fBmin.x(), fBmax.x(),
698 fBmin.y(), fBmax.y(),
699 fBmin.z(), fBmax.z() };
700}
701
702////////////////////////////////////////////////////////////////////////
703//
704// CreatePolyhedron
705//
707{
708 // Check orientation of vertices
709 G4ThreeVector v1 = fVertex[1] - fVertex[0];
710 G4ThreeVector v2 = fVertex[2] - fVertex[0];
711 G4ThreeVector v3 = fVertex[3] - fVertex[0];
712 G4bool invert = v1.cross(v2).dot(v3) < 0.;
713 G4int k2 = (invert) ? 3 : 2;
714 G4int k3 = (invert) ? 2 : 3;
715
716 // Set coordinates of vertices
717 G4double xyz[4][3];
718 for (G4int i = 0; i < 3; ++i)
719 {
720 xyz[0][i] = fVertex[0][i];
721 xyz[1][i] = fVertex[1][i];
722 xyz[2][i] = fVertex[k2][i];
723 xyz[3][i] = fVertex[k3][i];
724 }
725
726 // Create polyhedron
727 G4int faces[4][4] = { {1,3,2,0}, {1,4,3,0}, {1,2,4,0}, {2,3,4,0} };
728 auto ph = new G4Polyhedron;
729 ph->createPolyhedron(4,4,xyz,faces);
730
731 return ph;
732}
733
734////////////////////////////////////////////////////////////////////////
735//
736// GetPolyhedron
737//
739{
740 if (fpPolyhedron == nullptr ||
741 fRebuildPolyhedron ||
742 fpPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() !=
743 fpPolyhedron->GetNumberOfRotationSteps())
744 {
745 G4AutoLock l(&polyhedronMutex);
746 delete fpPolyhedron;
747 fpPolyhedron = CreatePolyhedron();
748 fRebuildPolyhedron = false;
749 l.unlock();
750 }
751 return fpPolyhedron;
752}
753
754#endif
G4TemplateAutoLock< G4Mutex > G4AutoLock
std::vector< G4ThreeVector > G4ThreeVectorList
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4double G4QuickRand(uint32_t seed=0)
#define G4MUTEX_INITIALIZER
std::mutex G4Mutex
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
long G4long
Definition G4Types.hh:87
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4String G4GeometryType
Definition G4VSolid.hh:70
G4GLOB_DLL std::ostream G4cout
double z() const
Hep3Vector unit() const
double x() const
double y() const
Hep3Vector cross(const Hep3Vector &) const
double dot(const Hep3Vector &) const
double mag() 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...
G4bool BoundingBoxVsVoxelLimits(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
G4Polyhedron * CreatePolyhedron() const override
Definition G4Tet.cc:706
G4VSolid * Clone() const override
Definition G4Tet.cc:611
G4Polyhedron * GetPolyhedron() const override
Definition G4Tet.cc:738
G4bool IsFaceted() const override
Definition G4Tet.cc:602
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const override
Definition G4Tet.cc:373
G4GeometryType GetEntityType() const override
Definition G4Tet.cc:593
void DescribeYourselfTo(G4VGraphicsScene &scene) const override
Definition G4Tet.cc:686
G4Tet & operator=(const G4Tet &rhs)
Definition G4Tet.cc:140
G4ThreeVector GetPointOnSurface() const override
Definition G4Tet.cc:641
G4double GetSurfaceArea() override
Definition G4Tet.cc:677
void SetBoundingLimits(const G4ThreeVector &pMin, const G4ThreeVector &pMax)
Definition G4Tet.cc:326
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const override
Definition G4Tet.cc:536
std::ostream & StreamInfo(std::ostream &os) const override
Definition G4Tet.cc:620
G4double GetCubicVolume() override
Definition G4Tet.cc:668
G4Tet(const G4String &pName, const G4ThreeVector &anchor, const G4ThreeVector &p2, const G4ThreeVector &p3, const G4ThreeVector &p4, G4bool *degeneracyFlag=nullptr)
Definition G4Tet.cc:66
G4VisExtent GetExtent() const override
Definition G4Tet.cc:695
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const override
Definition G4Tet.cc:437
G4bool CheckDegeneracy(const G4ThreeVector &p0, const G4ThreeVector &p1, const G4ThreeVector &p2, const G4ThreeVector &p3) const
Definition G4Tet.cc:173
std::vector< G4ThreeVector > GetVertices() const
Definition G4Tet.cc:304
~G4Tet() override
Definition G4Tet.cc:113
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep) override
Definition G4Tet.cc:316
EInside Inside(const G4ThreeVector &p) const override
Definition G4Tet.cc:423
void SetVertices(const G4ThreeVector &anchor, const G4ThreeVector &p1, const G4ThreeVector &p2, const G4ThreeVector &p3, G4bool *degeneracyFlag=nullptr)
Definition G4Tet.cc:253
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override
Definition G4Tet.cc:363
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const override
Definition G4Tet.cc:496
virtual void AddSolid(const G4Box &)=0
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....
G4String GetName() const
G4VSolid(const G4String &name)
Definition G4VSolid.cc:59
void DumpInfo() const
G4double kCarTolerance
Definition G4VSolid.hh:418
G4VSolid & operator=(const G4VSolid &rhs)
Definition G4VSolid.cc:108
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
EInside
Definition geomdefs.hh:67
@ kInside
Definition geomdefs.hh:70
@ kOutside
Definition geomdefs.hh:68
@ kSurface
Definition geomdefs.hh:69
#define DBL_MAX
Definition templates.hh:62