Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeoTessellated.cxx
Go to the documentation of this file.
1// @(#)root/geom:$Id$// Author: Andrei Gheata 24/10/01
2
3// Contains() and DistFromOutside/Out() implemented by Mihaela Gheata
4// 2026-01: Revision to use BVH for navigation functions by Sandro Wenzel
5
6/*************************************************************************
7 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
8 * All rights reserved. *
9 * *
10 * For the licensing terms see $ROOTSYS/LICENSE. *
11 * For the list of contributors see $ROOTSYS/README/CREDITS. *
12 *************************************************************************/
13
14/** \class TGeoTessellated
15\ingroup Geometry_classes
16
17Tessellated solid class. It is composed by a set of planar faces having triangular or
18quadrilateral shape. The class does not provide navigation functionality, it just wraps the data
19for the composing faces.
20*/
21
22#include <iostream>
23#include <sstream>
24
25#include "TGeoManager.h"
26#include "TGeoMatrix.h"
27#include "TGeoVolume.h"
28#include "TVirtualGeoPainter.h"
29#include "TGeoTessellated.h"
30#include "TBuffer3D.h"
31#include "TBuffer3DTypes.h"
32#include "TMath.h"
33#include "TBuffer.h"
34
35#include <array>
36#include <vector>
37
38// include the Third-party BVH headers
39#include <bvh2_third_party.h>
40// some kernels on top of BVH
41#include <bvh2_extra_kernels.h>
42
43#include <cmath>
44#include <limits>
45
47
49
50////////////////////////////////////////////////////////////////////////////////
51/// Compact consecutive equal vertices
52
54{
55 // Compact the common vertices and return new facet
56 if (nvertices < 2)
57 return nvertices;
58 int nvert = nvertices;
59 int i = 0;
60 while (i < nvert) {
61 if (vert[(i + 1) % nvert] == vert[i]) {
62 // shift last vertices left by one element
63 for (int j = i + 2; j < nvert; ++j)
64 vert[j - 1] = vert[j];
65 nvert--;
66 }
67 i++;
68 }
69 return nvert;
70}
71
72////////////////////////////////////////////////////////////////////////////////
73/// Check if a connected neighbour facet has compatible normal
74
75bool TGeoFacet::IsNeighbour(const TGeoFacet &other, bool &flip) const
76{
77
78 // Find a connecting segment
79 bool neighbour = false;
80 int line1[2], line2[2];
81 int npoints = 0;
82 for (int i = 0; i < fNvert; ++i) {
83 auto ivert = fIvert[i];
84 // Check if the other facet has the same vertex
85 for (int j = 0; j < other.GetNvert(); ++j) {
86 if (ivert == other[j]) {
87 line1[npoints] = i;
88 line2[npoints] = j;
89 if (++npoints == 2) {
90 neighbour = true;
91 bool order1 = line1[1] == line1[0] + 1;
92 bool order2 = line2[1] == (line2[0] + 1) % other.GetNvert();
93 flip = (order1 == order2);
94 return neighbour;
95 }
96 }
97 }
98 }
99 return neighbour;
100}
101
102////////////////////////////////////////////////////////////////////////////////
103/// Constructor. In case nfacets is zero, it is user's responsibility to
104/// call CloseShape once all faces are defined.
105
107{
109 if (nfacets)
110 fFacets.reserve(nfacets);
111}
112
113////////////////////////////////////////////////////////////////////////////////
114/// Constructor providing directly the array of vertices. Facets have to be added
115/// providing vertex indices rather than coordinates.
116
117TGeoTessellated::TGeoTessellated(const char *name, const std::vector<Vertex_t> &vertices) : TGeoBBox(name, 0, 0, 0)
118{
119 fVertices = vertices;
120 fNvert = fVertices.size();
121}
122
123////////////////////////////////////////////////////////////////////////////////
124/// Add a vertex checking for duplicates, returning the vertex index
125
127{
128 constexpr double tolerance = 1.e-10;
129 auto vertexHash = [&](Vertex_t const &vertex) {
130 // Compute hash for the vertex
131 long hash = 0;
132 // helper function to generate hash from integer numbers
133 auto hash_combine = [](long seed, const long value) {
134 return seed ^ (std::hash<long>{}(value) + 0x9e3779b9 + (seed << 6) + (seed >> 2));
135 };
136 for (int i = 0; i < 3; i++) {
137 // use tolerance to generate int with the desired precision from a real number for hashing
138 hash = hash_combine(hash, std::roundl(vertex[i] / tolerance));
139 }
140 return hash;
141 };
142
143 auto hash = vertexHash(vert);
144 bool isAdded = false;
145 int ivert = -1;
146 // Get the compatible vertices
147 auto range = fVerticesMap.equal_range(hash);
148 for (auto it = range.first; it != range.second; ++it) {
149 ivert = it->second;
150 if (fVertices[ivert] == vert) {
151 isAdded = true;
152 break;
153 }
154 }
155 if (!isAdded) {
156 ivert = fVertices.size();
157 fVertices.push_back(vert);
158 fVerticesMap.insert(std::make_pair(hash, ivert));
159 }
160 return ivert;
161}
162
163////////////////////////////////////////////////////////////////////////////////
164/// Adding a triangular facet from vertex positions in absolute coordinates
165
167{
168 if (fDefined) {
169 Error("AddFacet", "Shape %s already fully defined. Not adding", GetName());
170 return false;
171 }
172
173 Vertex_t vert[3];
174 vert[0] = pt0;
175 vert[1] = pt1;
176 vert[2] = pt2;
178 if (nvert < 3) {
179 Error("AddFacet", "Triangular facet at index %d degenerated. Not adding.", GetNfacets());
180 return false;
181 }
182 int ind[3];
183 for (auto i = 0; i < 3; ++i)
184 ind[i] = AddVertex(vert[i]);
185 fNseg += 3;
186 fFacets.emplace_back(ind[0], ind[1], ind[2]);
187
188 return true;
189}
190
191////////////////////////////////////////////////////////////////////////////////
192/// Adding a triangular facet from indices of vertices
193
195{
196 if (fDefined) {
197 Error("AddFacet", "Shape %s already fully defined. Not adding", GetName());
198 return false;
199 }
200 if (fVertices.empty()) {
201 Error("AddFacet", "Shape %s Cannot add facets by indices without vertices. Not adding", GetName());
202 return false;
203 }
204
205 fNseg += 3;
206 fFacets.emplace_back(i0, i1, i2);
207 return true;
208}
209
210////////////////////////////////////////////////////////////////////////////////
211/// Adding a quadrilateral facet from vertex positions in absolute coordinates
212
214{
215 if (fDefined) {
216 Error("AddFacet", "Shape %s already fully defined. Not adding", GetName());
217 return false;
218 }
219 Vertex_t vert[4];
220 vert[0] = pt0;
221 vert[1] = pt1;
222 vert[2] = pt2;
223 vert[3] = pt3;
225 if (nvert < 3) {
226 Error("AddFacet", "Quadrilateral facet at index %d degenerated. Not adding.", GetNfacets());
227 return false;
228 }
229
230 int ind[4];
231 for (auto i = 0; i < nvert; ++i)
232 ind[i] = AddVertex(vert[i]);
233 fNseg += nvert;
234 if (nvert == 3)
235 fFacets.emplace_back(ind[0], ind[1], ind[2]);
236 else
237 fFacets.emplace_back(ind[0], ind[1], ind[2], ind[3]);
238
239 if (fNfacets > 0 && GetNfacets() == fNfacets)
240 CloseShape(false);
241 return true;
242}
243
244////////////////////////////////////////////////////////////////////////////////
245/// Adding a quadrilateral facet from indices of vertices
246
247bool TGeoTessellated::AddFacet(int i0, int i1, int i2, int i3)
248{
249 if (fDefined) {
250 Error("AddFacet", "Shape %s already fully defined. Not adding", GetName());
251 return false;
252 }
253 if (fVertices.empty()) {
254 Error("AddFacet", "Shape %s Cannot add facets by indices without vertices. Not adding", GetName());
255 return false;
256 }
257
258 fNseg += 4;
259 fFacets.emplace_back(i0, i1, i2, i3);
260 return true;
261}
262
263////////////////////////////////////////////////////////////////////////////////
264/// Compute normal for a given facet
265
267{
268 // Compute normal using non-zero segments
269 constexpr double kTolerance = 1.e-20;
270 auto const &facet = fFacets[ifacet];
271 int nvert = facet.GetNvert();
272 degenerated = true;
274 for (int i = 0; i < nvert - 1; ++i) {
275 Vertex_t e1 = fVertices[facet[i + 1]] - fVertices[facet[i]];
276 if (e1.Mag2() < kTolerance)
277 continue;
278 for (int j = i + 1; j < nvert; ++j) {
279 Vertex_t e2 = fVertices[facet[(j + 1) % nvert]] - fVertices[facet[j]];
280 if (e2.Mag2() < kTolerance)
281 continue;
283 // e1 and e2 may be colinear
284 if (normal.Mag2() < kTolerance)
285 continue;
286 normal.Normalize();
287 degenerated = false;
288 break;
289 }
290 if (!degenerated)
291 break;
292 }
293 return normal;
294}
295
296////////////////////////////////////////////////////////////////////////////////
297/// Check validity of facet
298
300{
301 constexpr double kTolerance = 1.e-10;
302 auto const &facet = fFacets[ifacet];
303 int nvert = facet.GetNvert();
304 bool degenerated = true;
306 if (degenerated) {
307 std::cout << "Facet: " << ifacet << " is degenerated\n";
308 return false;
309 }
310
311 // Compute surface area
312 double surfaceArea = 0.;
313 for (int i = 1; i < nvert - 1; ++i) {
315 Vertex_t e2 = fVertices[facet[i + 1]] - fVertices[facet[0]];
316 surfaceArea += 0.5 * Vertex_t::Cross(e1, e2).Mag();
317 }
318 if (surfaceArea < kTolerance) {
319 std::cout << "Facet: " << ifacet << " has zero surface area\n";
320 return false;
321 }
322
323 return true;
324}
325
326////////////////////////////////////////////////////////////////////////////////
327/// Close the shape: calculate bounding box and compact vertices
328
329void TGeoTessellated::CloseShape(bool check, bool fixFlipped, bool verbose)
330{
331 if (fIsClosed && fBVH) {
332 return;
333 }
334 // Compute bounding box
335 fDefined = true;
336 fNvert = fVertices.size();
337 fNfacets = fFacets.size();
338 ComputeBBox();
339
340 BuildBVH();
341 if (fOutwardNormals.size() == 0) {
343 } else {
344 // short check if the normal container is of correct size
345 if (fOutwardNormals.size() != fFacets.size()) {
346 std::cerr << "Inconsistency in normal container";
347 }
348 }
349 fIsClosed = true;
350
351 // Cleanup the vertex map
352 std::multimap<long, int>().swap(fVerticesMap);
353
354 if (fVertices.size() > 0) {
355 if (!check)
356 return;
357
358 // Check facets
359 for (auto i = 0; i < fNfacets; ++i)
360 FacetCheck(i);
361
363 }
364}
365
366////////////////////////////////////////////////////////////////////////////////
367/// Check closure of the solid and check/fix flipped normals
368
370{
371 int *nn = new int[fNfacets];
372 bool *flipped = new bool[fNfacets];
373 bool hasorphans = false;
374 bool hasflipped = false;
375 for (int i = 0; i < fNfacets; ++i) {
376 nn[i] = 0;
377 flipped[i] = false;
378 }
379
380 for (int icrt = 0; icrt < fNfacets; ++icrt) {
381 // all neighbours checked?
382 if (nn[icrt] >= fFacets[icrt].GetNvert())
383 continue;
384 for (int i = icrt + 1; i < fNfacets; ++i) {
385 bool isneighbour = fFacets[icrt].IsNeighbour(fFacets[i], flipped[i]);
386 if (isneighbour) {
387 if (flipped[icrt])
388 flipped[i] = !flipped[i];
389 if (flipped[i])
390 hasflipped = true;
391 nn[icrt]++;
392 nn[i]++;
393 if (nn[icrt] == fFacets[icrt].GetNvert())
394 break;
395 }
396 }
397 if (nn[icrt] < fFacets[icrt].GetNvert())
398 hasorphans = true;
399 }
400
401 if (hasorphans && verbose) {
402 Error("Check", "Tessellated solid %s has following not fully connected facets:", GetName());
403 for (int icrt = 0; icrt < fNfacets; ++icrt) {
404 if (nn[icrt] < fFacets[icrt].GetNvert())
405 std::cout << icrt << " (" << fFacets[icrt].GetNvert() << " edges, " << nn[icrt] << " neighbours)\n";
406 }
407 }
409 int nfixed = 0;
410 if (hasflipped) {
411 if (verbose)
412 Warning("Check", "Tessellated solid %s has following facets with flipped normals:", GetName());
413 for (int icrt = 0; icrt < fNfacets; ++icrt) {
414 if (flipped[icrt]) {
415 if (verbose)
416 std::cout << icrt << "\n";
417 if (fixFlipped) {
418 fFacets[icrt].Flip();
419 nfixed++;
420 }
421 }
422 }
423 if (nfixed && verbose)
424 Info("Check", "Automatically flipped %d facets to match first defined facet", nfixed);
425 }
426 delete[] nn;
427 delete[] flipped;
428
429 return !hasorphans;
430}
431
432////////////////////////////////////////////////////////////////////////////////
433/// Compute bounding box
434
436{
437 const double kBig = TGeoShape::Big();
438 double vmin[3] = {kBig, kBig, kBig};
439 double vmax[3] = {-kBig, -kBig, -kBig};
440 for (const auto &facet : fFacets) {
441 for (int i = 0; i < facet.GetNvert(); ++i) {
442 for (int j = 0; j < 3; ++j) {
443 vmin[j] = TMath::Min(vmin[j], fVertices[facet[i]].operator[](j));
444 vmax[j] = TMath::Max(vmax[j], fVertices[facet[i]].operator[](j));
445 }
446 }
447 }
448 fDX = 0.5 * (vmax[0] - vmin[0]);
449 fDY = 0.5 * (vmax[1] - vmin[1]);
450 fDZ = 0.5 * (vmax[2] - vmin[2]);
451 for (int i = 0; i < 3; ++i)
452 fOrigin[i] = 0.5 * (vmax[i] + vmin[i]);
453}
454
455////////////////////////////////////////////////////////////////////////////////
456/// Returns numbers of vertices, segments and polygons composing the shape mesh.
457
459{
460 nvert = fNvert;
461 nsegs = fNseg;
462 npols = GetNfacets();
463}
464
465////////////////////////////////////////////////////////////////////////////////
466/// Creates a TBuffer3D describing *this* shape.
467/// Coordinates are in local reference frame.
468
470{
471 const int nvert = fNvert;
472 const int nsegs = fNseg;
473 const int npols = GetNfacets();
475 if (buff) {
476 SetPoints(buff->fPnts);
478 }
479 return buff;
480}
481
482////////////////////////////////////////////////////////////////////////////////
483/// Prints basic info
484
486{
487 std::cout << "=== Tessellated shape " << GetName() << " having " << GetNvertices() << " vertices and "
488 << GetNfacets() << " facets\n";
489}
490
491////////////////////////////////////////////////////////////////////////////////
492/// Fills TBuffer3D structure for segments and polygons.
493
495{
496 const int c = GetBasicColor();
497 int *segs = buff.fSegs;
498 int *pols = buff.fPols;
499
500 int indseg = 0; // segment internal data index
501 int indpol = 0; // polygon internal data index
502 int sind = 0; // segment index
503 for (const auto &facet : fFacets) {
504 auto nvert = facet.GetNvert();
505 pols[indpol++] = c;
506 pols[indpol++] = nvert;
507 for (auto j = 0; j < nvert; ++j) {
508 int k = (j + 1) % nvert;
509 // segment made by next consecutive points
510 segs[indseg++] = c;
511 segs[indseg++] = facet[j];
512 segs[indseg++] = facet[k];
513 // add segment to current polygon and increment segment index
514 pols[indpol + nvert - j - 1] = sind++;
515 }
516 indpol += nvert;
517 }
518}
519
520////////////////////////////////////////////////////////////////////////////////
521/// Fill tessellated points to an array.
522
524{
525 int ind = 0;
526 for (const auto &vertex : fVertices) {
527 vertex.CopyTo(&points[ind]);
528 ind += 3;
529 }
530}
531
532////////////////////////////////////////////////////////////////////////////////
533/// Fill tessellated points in float.
534
536{
537 int ind = 0;
538 for (const auto &vertex : fVertices) {
539 points[ind++] = vertex.x();
540 points[ind++] = vertex.y();
541 points[ind++] = vertex.z();
542 }
543}
544
545////////////////////////////////////////////////////////////////////////////////
546/// Resize the shape by scaling vertices within maxsize and center to origin
547
549{
550 using Vector3_t = Vertex_t;
551
552 if (!fDefined) {
553 Error("ResizeCenter", "Not all faces are defined");
554 return;
555 }
557 double maxedge = TMath::Max(TMath::Max(fDX, fDY), fDZ);
558 double scale = maxsize / maxedge;
559 constexpr double kTol = 1e-12;
560 const bool modified = (std::abs(scale - 1.0) > kTol) || (std::abs(origin[0]) > kTol) ||
561 (std::abs(origin[1]) > kTol) || (std::abs(origin[2]) > kTol);
562 for (size_t i = 0; i < fVertices.size(); ++i) {
563 fVertices[i] = scale * (fVertices[i] - origin);
564 }
565 fOrigin[0] = fOrigin[1] = fOrigin[2] = 0;
566 fDX *= scale;
567 fDY *= scale;
568 fDZ *= scale;
569 if (modified) {
570 BuildBVH();
572 }
573}
574
575////////////////////////////////////////////////////////////////////////////////
576/// Fills a static 3D buffer and returns a reference.
577
579{
580 static TBuffer3D buffer(TBuffer3DTypes::kGeneric);
581
583
584 const int nvert = fNvert;
585 const int nsegs = fNseg;
586 const int npols = GetNfacets();
587
589 if (buffer.SetRawSizes(nvert, 3 * nvert, nsegs, 3 * nsegs, npols, 6 * npols)) {
591 }
592 }
594 SetPoints(buffer.fPnts);
595 if (!buffer.fLocalFrame) {
596 TransformPoints(buffer.fPnts, buffer.NbPnts());
597 }
598
599 SetSegsAndPols(buffer);
601 }
602
603 return buffer;
604}
605
606////////////////////////////////////////////////////////////////////////////////
607/// Reads a single tessellated solid from an .obj file.
608
610{
611 using std::vector, std::string, std::ifstream, std::stringstream, std::endl;
612
613 vector<Vertex_t> vertices;
615
616 struct FacetInd_t {
617 int i0 = -1;
618 int i1 = -1;
619 int i2 = -1;
620 int i3 = -1;
621 int nvert = 0;
622 FacetInd_t(int a, int b, int c)
623 {
624 i0 = a;
625 i1 = b;
626 i2 = c;
627 nvert = 3;
628 };
629 FacetInd_t(int a, int b, int c, int d)
630 {
631 i0 = a;
632 i1 = b;
633 i2 = c;
634 i3 = d;
635 nvert = 4;
636 };
637 };
638
640 // List of geometric vertices, with (x, y, z [,w]) coordinates, w is optional and defaults to 1.0.
641 // struct vtx_t { double x = 0; double y = 0; double z = 0; double w = 1; };
642
643 // Texture coordinates in u, [,v ,w]) coordinates, these will vary between 0 and 1. v, w are optional and default to
644 // 0.
645 // struct tex_t { double u; double v; double w; };
646
647 // List of vertex normals in (x,y,z) form; normals might not be unit vectors.
648 // struct vn_t { double x; double y; double z; };
649
650 // Parameter space vertices in ( u [,v] [,w] ) form; free form geometry statement
651 // struct vp_t { double u; double v; double w; };
652
653 // Faces are defined using lists of vertex, texture and normal indices which start at 1.
654 // Polygons such as quadrilaterals can be defined by using more than three vertex/texture/normal indices.
655 // f v1//vn1 v2//vn2 v3//vn3 ...
656
657 // Records starting with the letter "l" specify the order of the vertices which build a polyline.
658 // l v1 v2 v3 v4 v5 v6 ...
659
660 string line;
661 int ind[4] = {0};
662 ifstream file(objfile);
663 if (!file.is_open()) {
664 ::Error("TGeoTessellated::ImportFromObjFormat", "Unable to open %s", objfile);
665 return nullptr;
666 }
667
668 while (getline(file, line)) {
669 stringstream ss(line);
670 string tag;
671
672 // We ignore everything which is not a vertex or a face
673 if (line.rfind('v', 0) == 0 && line.rfind("vt", 0) != 0 && line.rfind("vn", 0) != 0 && line.rfind("vn", 0) != 0) {
674 // Decode the vertex
675 double pos[4] = {0, 0, 0, 1};
676 ss >> tag >> pos[0] >> pos[1] >> pos[2] >> pos[3];
677 vertices.emplace_back(pos[0] * pos[3], pos[1] * pos[3], pos[2] * pos[3]);
678 }
679
680 else if (line.rfind('f', 0) == 0) {
681 // Decode the face
682 ss >> tag;
683 string word;
684 sfacets.clear();
685 while (ss >> word)
686 sfacets.push_back(word);
687 if (sfacets.size() > 4 || sfacets.size() < 3) {
688 ::Error("TGeoTessellated::ImportFromObjFormat", "Detected face having unsupported %zu vertices",
689 sfacets.size());
690 return nullptr;
691 }
692 int nvert = 0;
693 for (auto &sword : sfacets) {
694 stringstream ssword(sword);
695 string token;
696 getline(ssword, token, '/'); // just need the vertex index, which is the first token
697 // Convert string token to integer
698
699 ind[nvert++] = stoi(token) - 1;
700 if (ind[nvert - 1] < 0) {
701 ::Error("TGeoTessellated::ImportFromObjFormat", "Unsupported relative vertex index definition in %s",
702 objfile);
703 return nullptr;
704 }
705 }
706 if (nvert == 3)
707 facets.emplace_back(ind[0], ind[1], ind[2]);
708 else
709 facets.emplace_back(ind[0], ind[1], ind[2], ind[3]);
710 }
711 }
712
713 int nvertices = (int)vertices.size();
714 int nfacets = (int)facets.size();
715 if (nfacets < 3) {
716 ::Error("TGeoTessellated::ImportFromObjFormat", "Not enough faces detected in %s", objfile);
717 return nullptr;
718 }
719
720 string sobjfile(objfile);
721 if (verbose)
722 std::cout << "Read " << nvertices << " vertices and " << nfacets << " facets from " << sobjfile << endl;
723
724 auto tsl = new TGeoTessellated(sobjfile.erase(sobjfile.find_last_of('.')).c_str(), vertices);
725
726 for (int i = 0; i < nfacets; ++i) {
727 auto facet = facets[i];
728 if (facet.nvert == 3)
729 tsl->AddFacet(facet.i0, facet.i1, facet.i2);
730 else
731 tsl->AddFacet(facet.i0, facet.i1, facet.i2, facet.i3);
732 }
733 tsl->CloseShape(check, true, verbose);
734 tsl->Print();
735 return tsl;
736}
737
738// implementation of some geometry helper functions in anonymous namespace
739namespace {
740
742// The classic Moeller-Trumbore ray triangle-intersection kernel:
743// - Compute triangle edges e1, e2
744// - Compute determinant det
745// - Reject parallel rays
746// - Compute barycentric coordinates u, v
747// - Compute ray parameter t
748double rayTriangle(const Vertex_t &orig, const Vertex_t &dir, const Vertex_t &v0, const Vertex_t &v1,
749 const Vertex_t &v2, double rayEPS = 1e-8)
750{
751 constexpr double EPS = 1e-8;
752 const double INF = std::numeric_limits<double>::infinity();
753 Vertex_t e1{v1[0] - v0[0], v1[1] - v0[1], v1[2] - v0[2]};
754 Vertex_t e2{v2[0] - v0[0], v2[1] - v0[1], v2[2] - v0[2]};
755 auto p = Vertex_t::Cross(dir, e2);
756 auto det = e1.Dot(p);
757 if (std::abs(det) <= EPS) {
758 return INF;
759 }
760
761 Vertex_t tvec{orig[0] - v0[0], orig[1] - v0[1], orig[2] - v0[2]};
762 auto invDet = 1.0 / det;
763 auto u = tvec.Dot(p) * invDet;
764 if (u < 0.0 || u > 1.0) {
765 return INF;
766 }
767 auto q = Vertex_t::Cross(tvec, e1);
768 auto v = dir.Dot(q) * invDet;
769 if (v < 0.0 || u + v > 1.0) {
770 return INF;
771 }
772 auto t = e2.Dot(q) * invDet;
773 return (t > rayEPS) ? t : INF;
774}
775
776inline double rayFacet(const Vertex_t &orig, const Vertex_t &dir, const TGeoFacet &facet,
777 const std::vector<Vertex_t> &vertices, double rayEPS = 1e-8)
778{
779 // Keep the stored facet topology intact and triangulate quads only for geometric queries.
780 const auto &v0 = vertices[facet[0]];
781 const auto &v1 = vertices[facet[1]];
782 const auto &v2 = vertices[facet[2]];
783 auto t = rayTriangle(orig, dir, v0, v1, v2, rayEPS);
784 if (facet.GetNvert() == 3)
785 return t;
786 const auto &v3 = vertices[facet[3]];
787 auto t2 = rayTriangle(orig, dir, v0, v2, v3, rayEPS);
788 return std::min(t, t2);
789}
790
791inline bool rayFacetHit(const Vertex_t &orig, const Vertex_t &dir, const TGeoFacet &facet,
792 const std::vector<Vertex_t> &vertices, double rayEPS = 1e-8)
793{
794 // Contains/parity checks only need a boolean hit, so keep the finite-distance test in one place.
795 return rayFacet(orig, dir, facet, vertices, rayEPS) != std::numeric_limits<double>::infinity();
796}
797
798template <typename T = float>
799struct Vec3f {
800 T x, y, z;
801 Vec3f(T x_, T y_, T z_) : x(x_), y(y_), z(z_){};
802};
803
804template <typename T>
805inline Vec3f<T> operator-(const Vec3f<T> &a, const Vec3f<T> &b)
806{
807 return {a.x - b.x, a.y - b.y, a.z - b.z};
808}
809
810template <typename T>
811inline Vec3f<T> cross(const Vec3f<T> &a, const Vec3f<T> &b)
812{
813 return {a.y * b.z - a.z * b.y, a.z * b.x - a.x * b.z, a.x * b.y - a.y * b.x};
814}
815
816template <typename T>
817inline T dot(const Vec3f<T> &a, const Vec3f<T> &b)
818{
819 return a.x * b.x + a.y * b.y + a.z * b.z;
820}
821
822// Kernel to get closest/shortest distance between a point and a triangl (a,b,c).
823// Performed by default in float since Safety can be approximate.
824// Project point onto triangle plane
825// If projection lies inside → distance to plane
826// Otherwise compute min distance to the three edges
827// Return squared distance
828template <typename T = float>
829T pointTriangleDistSq(const Vec3f<T> &p, const Vec3f<T> &a, const Vec3f<T> &b, const Vec3f<T> &c)
830{
831 // Edges
832 Vec3f<T> ab = b - a;
833 Vec3f<T> ac = c - a;
834 Vec3f<T> ap = p - a;
835
836 auto d1 = dot(ab, ap);
837 auto d2 = dot(ac, ap);
838 if (d1 <= T(0.0) && d2 <= T(0.0)) {
839 return dot(ap, ap); // barycentric (1,0,0)
840 }
841
842 Vec3f<T> bp = p - b;
843 auto d3 = dot(ab, bp);
844 auto d4 = dot(ac, bp);
845 if (d3 >= T(0.0) && d4 <= d3) {
846 return dot(bp, bp); // (0,1,0)
847 }
848
849 T vc = d1 * d4 - d3 * d2;
850 if (vc <= 0.0f && d1 >= 0.0f && d3 <= 0.0f) {
851 T v = d1 / (d1 - d3);
852 Vec3f<T> proj = {a.x + v * ab.x, a.y + v * ab.y, a.z + v * ab.z};
853 Vec3f<T> d = p - proj;
854 return dot(d, d); // edge AB
855 }
856
857 Vec3f<T> cp = p - c;
858 T d5 = dot(ab, cp);
859 T d6 = dot(ac, cp);
860 if (d6 >= T(0.0f) && d5 <= d6) {
861 return dot(cp, cp); // (0,0,1)
862 }
863
864 T vb = d5 * d2 - d1 * d6;
865 if (vb <= 0.0f && d2 >= 0.0f && d6 <= 0.0f) {
866 T w = d2 / (d2 - d6);
867 Vec3f<T> proj = {a.x + w * ac.x, a.y + w * ac.y, a.z + w * ac.z};
868 Vec3f<T> d = p - proj;
869 return dot(d, d); // edge AC
870 }
871
872 T va = d3 * d6 - d5 * d4;
873 if (va <= 0.0f && (d4 - d3) >= 0.0f && (d5 - d6) >= 0.0f) {
874 T w = (d4 - d3) / ((d4 - d3) + (d5 - d6));
875 Vec3f<T> proj = {b.x + w * (c.x - b.x), b.y + w * (c.y - b.y), b.z + w * (c.z - b.z)};
876 Vec3f<T> d = p - proj;
877 return dot(d, d); // edge BC
878 }
879
880 // Inside face region
881 T denom = T(1.0f) / (va + vb + vc);
882 T v = vb * denom;
883 T w = vc * denom;
884
885 Vec3f<T> proj = {a.x + ab.x * v + ac.x * w, a.y + ab.y * v + ac.y * w, a.z + ab.z * v + ac.z * w};
886
887 Vec3f<T> d = p - proj;
888 return dot(d, d);
889}
890
891template <typename T = float>
892T pointFacetDistSq(const Vec3f<T> &p, const TGeoFacet &facet, const std::vector<Vertex_t> &vertices)
893{
894 // Safety uses the same on-the-fly split as ray queries so triangles and quads stay consistent.
895 const auto &v0 = vertices[facet[0]];
896 const auto &v1 = vertices[facet[1]];
897 const auto &v2 = vertices[facet[2]];
898 auto d = pointTriangleDistSq(p, Vec3f<T>(v0[0], v0[1], v0[2]), Vec3f<T>(v1[0], v1[1], v1[2]),
899 Vec3f<T>(v2[0], v2[1], v2[2]));
900 if (facet.GetNvert() == 3)
901 return d;
902 const auto &v3 = vertices[facet[3]];
903 auto d2 = pointTriangleDistSq(p, Vec3f<T>(v0[0], v0[1], v0[2]), Vec3f<T>(v2[0], v2[1], v2[2]),
904 Vec3f<T>(v3[0], v3[1], v3[2]));
905 return std::min(d, d2);
906}
907
908} // end anonymous namespace
909
910////////////////////////////////////////////////////////////////////////////////
911/// DistFromOutside
912
914 Double_t * /*safe*/) const
915{
916 // use the BVH intersector in combination with leaf ray-triangle testing
917 double local_step = Big(); // we need this otherwise the lambda get's confused
918
919 using Scalar = float;
921 using Node = bvh::v2::Node<Scalar, 3>;
922 using Bvh = bvh::v2::Bvh<Node>;
923 using Ray = bvh::v2::Ray<Scalar, 3>;
924
925 // let's fetch the bvh
926 auto mybvh = (Bvh *)fBVH;
927 if (!mybvh) {
928 assert(false);
929 return -1.;
930 }
931
932 auto truncate_roundup = [](double orig) {
933 float epsilon = std::numeric_limits<float>::epsilon() * std::fabs(orig);
934 // Add the bias to x before assigning it to y
935 return static_cast<float>(orig + epsilon);
936 };
937
938 // let's do very quick checks against the top node
939 const auto topnode_bbox = mybvh->get_root().get_bbox();
940 if ((-point[0] + topnode_bbox.min[0]) > stepmax) {
941 return Big();
942 }
943 if ((-point[1] + topnode_bbox.min[1]) > stepmax) {
944 return Big();
945 }
946 if ((-point[2] + topnode_bbox.min[2]) > stepmax) {
947 return Big();
948 }
949 if ((point[0] - topnode_bbox.max[0]) > stepmax) {
950 return Big();
951 }
952 if ((point[1] - topnode_bbox.max[1]) > stepmax) {
953 return Big();
954 }
955 if ((point[2] - topnode_bbox.max[2]) > stepmax) {
956 return Big();
957 }
958
959 // the ray used for bvh interaction
960 Ray ray(Vec3(point[0], point[1], point[2]), // origin
961 Vec3(dir[0], dir[1], dir[2]), // direction
962 0.0f, // minimum distance (could give stepmax ?)
964
965 static constexpr bool use_robust_traversal = true;
966
967 Vertex_t dir_v{dir[0], dir[1], dir[2]};
968 // Traverse the BVH and apply concrete object intersection in BVH leafs
970 mybvh->intersect<false, use_robust_traversal>(ray, mybvh->get_root().index, stack, [&](size_t begin, size_t end) {
971 for (size_t prim_id = begin; prim_id < end; ++prim_id) {
972 auto objectid = mybvh->prim_ids[prim_id];
973 const auto &facet = fFacets[objectid];
974 const auto &n = fOutwardNormals[objectid];
975
976 // quick normal test. Coming from outside, the dot product must be negative
977 if (n.Dot(dir_v) > 0.) {
978 continue;
979 }
980
981 auto thisdist = rayFacet(Vertex_t(point[0], point[1], point[2]), dir_v, facet, fVertices, 0.);
982
983 if (thisdist < local_step) {
985 }
986 }
987 return false; // go on after this
988 });
989
990 return local_step;
991}
992
993////////////////////////////////////////////////////////////////////////////////
994/// DistFromOutside
995
997 Double_t /*stepmax*/, Double_t * /*safe*/) const
998{
999 // use the BVH intersector in combination with leaf ray-triangle testing
1000 double local_step = Big(); // we need this otherwise the lambda get's confused
1001
1002 using Scalar = float;
1004 using Node = bvh::v2::Node<Scalar, 3>;
1005 using Bvh = bvh::v2::Bvh<Node>;
1006 using Ray = bvh::v2::Ray<Scalar, 3>;
1007
1008 // let's fetch the bvh
1009 auto mybvh = (Bvh *)fBVH;
1010 if (!mybvh) {
1011 assert(false);
1012 return -1.;
1013 }
1014
1015 auto truncate_roundup = [](double orig) {
1016 float epsilon = std::numeric_limits<float>::epsilon() * std::fabs(orig);
1017 // Add the bias to x before assigning it to y
1018 return static_cast<float>(orig + epsilon);
1019 };
1020
1021 // the ray used for bvh interaction
1022 Ray ray(Vec3(point[0], point[1], point[2]), // origin
1023 Vec3(dir[0], dir[1], dir[2]), // direction
1024 0., // minimum distance (could give stepmax ?)
1026
1027 static constexpr bool use_robust_traversal = true;
1028
1029 Vertex_t dir_v{dir[0], dir[1], dir[2]};
1030 // Traverse the BVH and apply concrete object intersection in BVH leafs
1032 mybvh->intersect<false, use_robust_traversal>(ray, mybvh->get_root().index, stack, [&](size_t begin, size_t end) {
1033 for (size_t prim_id = begin; prim_id < end; ++prim_id) {
1034 auto objectid = mybvh->prim_ids[prim_id];
1035 auto facet = fFacets[objectid];
1036 const auto &n = fOutwardNormals[objectid];
1037
1038 // Only exiting surfaces are relevant (from inside--> dot product must be positive)
1039 if (n.Dot(dir_v) <= 0.) {
1040 continue;
1041 }
1042
1043 const double t = rayFacet(Vertex_t{point[0], point[1], point[2]}, dir_v, facet, fVertices, 0.);
1044 if (t < local_step) {
1045 local_step = t;
1046 }
1047 }
1048 return false; // go on after this
1049 });
1050
1051 return local_step;
1052}
1053
1054////////////////////////////////////////////////////////////////////////////////
1055/// Capacity
1056
1058{
1059 // For explanation of the following algorithm see:
1060 // https://en.wikipedia.org/wiki/Polyhedron#Volume
1061 // http://wwwf.imperial.ac.uk/~rn/centroid.pdf
1062
1063 double vol = 0.0;
1064 for (size_t i = 0; i < fFacets.size(); ++i) {
1065 auto &facet = fFacets[i];
1066 auto a = fVertices[facet[0]];
1067 auto b = fVertices[facet[1]];
1068 auto c = fVertices[facet[2]];
1069 vol +=
1070 a[0] * (b[1] * c[2] - b[2] * c[1]) + b[0] * (c[1] * a[2] - c[2] * a[1]) + c[0] * (a[1] * b[2] - a[2] * b[1]);
1071 }
1072 return vol / 6.0;
1073}
1074
1075////////////////////////////////////////////////////////////////////////////////
1076/// BuildBVH
1077
1079{
1080 using Scalar = float;
1081 using BBox = bvh::v2::BBox<Scalar, 3>;
1083 using Node = bvh::v2::Node<Scalar, 3>;
1084 using Bvh = bvh::v2::Bvh<Node>;
1085
1086 // helper determining axis aligned bounding box from a facet;
1087 auto GetBoundingBox = [this](TGeoFacet const &facet) {
1088 const auto nvertices = facet.GetNvert();
1089 if (nvertices != 3 && nvertices != 4)
1090 Fatal("BuildBVH", "only facets with 3 or 4 vertices supported");
1091 const auto &v1 = fVertices[facet[0]];
1092 const auto &v2 = fVertices[facet[1]];
1093 const auto &v3 = fVertices[facet[2]];
1094 const auto &v4 = (nvertices == 4) ? fVertices[facet[3]] : v3;
1095 BBox bbox;
1096 // A quad needs all four vertices in the BVH bounds even though traversal tests two triangles later on.
1097 bbox.min[0] = std::min(std::min(std::min(v1[0], v2[0]), v3[0]), v4[0]) - 0.001f;
1098 bbox.min[1] = std::min(std::min(std::min(v1[1], v2[1]), v3[1]), v4[1]) - 0.001f;
1099 bbox.min[2] = std::min(std::min(std::min(v1[2], v2[2]), v3[2]), v4[2]) - 0.001f;
1100 bbox.max[0] = std::max(std::max(std::max(v1[0], v2[0]), v3[0]), v4[0]) + 0.001f;
1101 bbox.max[1] = std::max(std::max(std::max(v1[1], v2[1]), v3[1]), v4[1]) + 0.001f;
1102 bbox.max[2] = std::max(std::max(std::max(v1[2], v2[2]), v3[2]), v4[2]) + 0.001f;
1103 return bbox;
1104 };
1105
1106 // we need bounding boxes enclosing the primitives and centers of primitives
1107 // (replaced here by centers of bounding boxes) to build the bvh
1108 std::vector<BBox> bboxes;
1109 std::vector<Vec3> centers;
1110
1111 int nd = fFacets.size();
1112 bboxes.reserve(nd);
1113 centers.reserve(nd);
1114
1115 for (int i = 0; i < nd; ++i) {
1116 auto &facet = fFacets[i];
1117
1118 // fetch the bounding box of this node and add to the vector of bounding boxes
1119 (bboxes).push_back(GetBoundingBox(facet));
1120 centers.emplace_back((bboxes).back().get_center());
1121 }
1122
1123 // check if some previous object is registered and delete if necessary
1124 if (fBVH) {
1125 delete (Bvh *)fBVH;
1126 fBVH = nullptr;
1127 }
1128
1129 // create the bvh
1132 auto bvh = bvh::v2::DefaultBuilder<Node>::build(bboxes, centers, config);
1133 auto bvhptr = new Bvh;
1134 *bvhptr = std::move(bvh); // copy structure
1135 fBVH = (void *)(bvhptr);
1136}
1137
1138////////////////////////////////////////////////////////////////////////////////
1139/// Contains
1140
1141bool TGeoTessellated::Contains(Double_t const *point) const
1142{
1143 // we do the parity test
1144 using Scalar = float;
1146 using Node = bvh::v2::Node<Scalar, 3>;
1147 using Bvh = bvh::v2::Bvh<Node>;
1148 using Ray = bvh::v2::Ray<Scalar, 3>;
1149
1150 // let's fetch the bvh
1151 auto mybvh = (Bvh *)fBVH;
1152 if (!mybvh) {
1153 assert(false);
1154 return false;
1155 }
1156
1157 auto truncate_roundup = [](double orig) {
1158 float epsilon = std::numeric_limits<float>::epsilon() * std::fabs(orig);
1159 // Add the bias to x before assigning it to y
1160 return static_cast<float>(orig + epsilon);
1161 };
1162
1163 // let's do very quick checks against the top node
1164 if (!TGeoBBox::Contains(point)) {
1165 return false;
1166 }
1167
1168 // An arbitrary test direction.
1169 // Doesn't need to be normalized and probes all normals. Also ensuring to be skewed somewhat
1170 // without evident symmetries.
1171 Vertex_t test_dir{1.0, 1.41421356237, 1.73205080757};
1172
1173 double local_step = Big();
1174 // the ray used for bvh interaction
1175 Ray ray(Vec3(point[0], point[1], point[2]), // origin
1176 Vec3(test_dir[0], test_dir[1], test_dir[2]), // direction
1177 0.0f, // minimum distance (could give stepmax ?)
1179
1180 static constexpr bool use_robust_traversal = true;
1181
1182 // Traverse the BVH and apply concrete object intersection in BVH leafs
1184 size_t crossings = 0;
1185 mybvh->intersect<false, use_robust_traversal>(ray, mybvh->get_root().index, stack, [&](size_t begin, size_t end) {
1186 for (size_t prim_id = begin; prim_id < end; ++prim_id) {
1187 auto objectid = mybvh->prim_ids[prim_id];
1188 auto &facet = fFacets[objectid];
1189
1190 // for the parity test, we probe all crossing surfaces
1191 if (rayFacetHit(Vertex_t(point[0], point[1], point[2]), test_dir, facet, fVertices, 0.)) {
1192 ++crossings;
1193 }
1194 }
1195 return false;
1196 });
1197
1198 return crossings & 1;
1199}
1200
1201namespace {
1202
1203// Helper classes/structs used for priority queue - BVH traversal
1204// structure keeping cost (value) for a BVH index
1205struct BVHPrioElement {
1206 size_t bvh_node_id;
1207 float value;
1208};
1209
1210// A priority queue for BVHPrioElement with an additional clear method
1211// for quick reset. We intentionally derive from std::priority_queue here to expose a
1212// clear() convenience method via access to the protected container `c`.
1213// This is internal, non-polymorphic code and relies on standard-library
1214// implementation details that are stable across supported platforms.
1215template <typename Comparator>
1216class BVHPrioQueue : public std::priority_queue<BVHPrioElement, std::vector<BVHPrioElement>, Comparator> {
1217public:
1218 using std::priority_queue<BVHPrioElement, std::vector<BVHPrioElement>,
1219 Comparator>::priority_queue; // constructor inclusion
1220
1221 // convenience method to quickly clear/reset the queue (instead of having to pop one by one)
1222 void clear() { this->c.clear(); }
1223};
1224
1225} // namespace
1226
1227/// a reusable safety kernel, which optionally returns the closest face
1228template <bool returnFace>
1230{
1231 // This is the classic traversal/pruning of a BVH based on priority queue search
1232
1234
1235 using Scalar = float;
1237 using Node = bvh::v2::Node<Scalar, 3>;
1238 using Bvh = bvh::v2::Bvh<Node>;
1239
1240 // let's fetch the bvh
1241 auto mybvh = (Bvh *)fBVH;
1242
1243 // testpoint object in float for quick BVH interaction
1244 Vec3 testpoint(point[0], point[1], point[2]);
1245
1246 auto currnode = mybvh->nodes[0]; // we start from the top BVH node
1247 // we do a quick check on the top node (in case we are outside shape)
1248 bool outside_top = false;
1249 if (!in) {
1251 if (outside_top) {
1253 // we simply return safety to the outer bounding box as an estimate
1254 return std::sqrt(safety_sq_to_top);
1255 }
1256 }
1257
1258 // comparator bringing out "smallest" value on top
1259 auto cmp = [](BVHPrioElement a, BVHPrioElement b) { return a.value > b.value; };
1260 static thread_local BVHPrioQueue<decltype(cmp)> queue(cmp);
1261 queue.clear();
1262
1263 // algorithm is based on standard iterative tree traversal with priority queues
1264 float current_safety_to_node_sq = 0.f;
1265
1266 if (returnFace) {
1267 *closest_facet_id = -1;
1268 }
1269
1270 do {
1271 if (currnode.is_leaf()) {
1272 // we are in a leaf node and actually talk to a face primitive
1273 const auto begin_prim_id = currnode.index.first_id();
1274 const auto end_prim_id = begin_prim_id + currnode.index.prim_count();
1275
1276 for (auto p_id = begin_prim_id; p_id < end_prim_id; p_id++) {
1277 const auto object_id = mybvh->prim_ids[p_id];
1278
1279 const auto &facet = fFacets[object_id];
1280 auto thissafetySQ =
1281 pointFacetDistSq(Vec3f<float>(point[0], point[1], point[2]), facet, fVertices);
1282
1285 if (returnFace) {
1287 }
1288 }
1289 }
1290 } else {
1291 // not a leave node ... for further traversal,
1292 // we inject the children into priority queue based on distance to it's bounding box
1293 const auto leftchild_id = currnode.index.first_id();
1294 const auto rightchild_id = leftchild_id + 1;
1295
1296 for (size_t childid : {leftchild_id, rightchild_id}) {
1297 if (childid >= mybvh->nodes.size()) {
1298 continue;
1299 }
1300
1301 const auto &node = mybvh->nodes[childid];
1302 const auto inside = bvh::v2::extra::contains(node.get_bbox(), testpoint);
1303
1304 if (inside) {
1305 // this must be further considered because we are inside the bounding box
1306 queue.push(BVHPrioElement{childid, -1.});
1307 } else {
1310 // this should be further considered
1311 queue.push(BVHPrioElement{childid, safety_to_node_square});
1312 }
1313 }
1314 }
1315 }
1316
1317 if (queue.size() > 0) {
1318 auto currElement = queue.top();
1319 currnode = mybvh->nodes[currElement.bvh_node_id];
1321 queue.pop();
1322 } else {
1323 break;
1324 }
1326
1327 return std::nextafter(std::sqrt(smallest_safety_sq), 0.0f);
1328}
1329
1330////////////////////////////////////////////////////////////////////////////////
1331/// Safety
1332
1334{
1335 // we could use some caching here (in future) since queries to the solid will likely
1336 // be made with some locality
1337
1338 // fall-back to precise safety kernel
1339 return SafetyKernel<false>(point, in);
1340}
1341
1342////////////////////////////////////////////////////////////////////////////////
1343/// ComputeNormal interface
1344
1345void TGeoTessellated::ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm) const
1346{
1347 // We take the approach to identify closest facet to the point via safety
1348 // and returning the normal from this face.
1349
1350 // TODO: Before doing that we could check for cached points from other queries
1351
1352 // use safety kernel
1353 int closest_face_id = -1;
1354 SafetyKernel<true>(point, true, &closest_face_id);
1355
1356 if (closest_face_id < 0) {
1357 norm[0] = 1.;
1358 norm[1] = 0.;
1359 norm[2] = 0.;
1360 return;
1361 }
1362
1363 const auto &n = fOutwardNormals[closest_face_id];
1364 norm[0] = n[0];
1365 norm[1] = n[1];
1366 norm[2] = n[2];
1367
1368 // change sign depending on dir
1369 if (norm[0] * dir[0] + norm[1] * dir[1] + norm[2] * dir[2] < 0) {
1370 norm[0] = -norm[0];
1371 norm[1] = -norm[1];
1372 norm[2] = -norm[2];
1373 }
1374 return;
1375}
1376
1377////////////////////////////////////////////////////////////////////////////////
1378/// Custom streamer which performs Closing on read.
1379/// Recalculation of BVH and normals is fast
1380
1382{
1383 if (b.IsReading()) {
1384 b.ReadClassBuffer(TGeoTessellated::Class(), this);
1385 CloseShape(false); // close shape but do not re-perform checks
1386 } else {
1387 b.WriteClassBuffer(TGeoTessellated::Class(), this);
1388 }
1389}
1390
1391////////////////////////////////////////////////////////////////////////////////
1392/// Calculate the normals
1393
1395{
1396 fOutwardNormals.clear();
1397 fOutwardNormals.reserve(fFacets.size());
1398 for (size_t i = 0; i < fFacets.size(); ++i) {
1399 // Reuse the generic facet normal code so triangles and quads follow the same path.
1400 bool degenerated = false;
1401 auto norm = FacetComputeNormal(static_cast<int>(i), degenerated);
1402 fOutwardNormals.emplace_back(norm);
1403 }
1404}
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
#define e(i)
Definition RSha256.hxx:103
float Float_t
Float 4 bytes (float)
Definition RtypesCore.h:71
const char Option_t
Option string (const char)
Definition RtypesCore.h:80
#define ClassImp(name)
Definition Rtypes.h:376
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t points
char name[80]
Definition TGX11.cxx:157
float * q
TTime operator-(const TTime &t1, const TTime &t2)
Definition TTime.h:83
segment * segs
Definition X3DBuffer.c:23
Generic 3D primitive description class.
Definition TBuffer3D.h:18
UInt_t NbPnts() const
Definition TBuffer3D.h:89
Bool_t SectionsValid(UInt_t mask) const
Definition TBuffer3D.h:76
void SetSectionsValid(UInt_t mask)
Definition TBuffer3D.h:74
Bool_t fLocalFrame
Definition TBuffer3D.h:99
Bool_t SetRawSizes(UInt_t reqPnts, UInt_t reqPntsCapacity, UInt_t reqSegs, UInt_t reqSegsCapacity, UInt_t reqPols, UInt_t reqPolsCapacity)
Set kRaw tessellation section of buffer with supplied sizes.
Double_t * fPnts
Definition TBuffer3D.h:122
Buffer base class used for serializing objects.
Definition TBuffer.h:43
Box class.
Definition TGeoBBox.h:18
void FillBuffer3D(TBuffer3D &buffer, Int_t reqSections, Bool_t localFrame) const override
Fills the supplied buffer, with sections in desired frame See TBuffer3D.h for explanation of sections...
Double_t fDX
Definition TGeoBBox.h:21
Double_t fOrigin[3]
Definition TGeoBBox.h:24
Bool_t Contains(const Double_t *point) const override
Test if point is inside this shape.
Definition TGeoBBox.cxx:322
Double_t fDY
Definition TGeoBBox.h:22
Double_t fDZ
Definition TGeoBBox.h:23
bool IsNeighbour(const TGeoFacet &other, bool &flip) const
Check if a connected neighbour facet has compatible normal.
static int CompactFacet(Vertex_t *vert, int nvertices)
Compact consecutive equal vertices.
static Double_t Big()
Definition TGeoShape.h:95
Int_t GetBasicColor() const
Get the basic color (0-7).
void TransformPoints(Double_t *points, UInt_t NbPoints) const
Tranform a set of points (LocalToMaster)
const char * GetName() const override
Get the shape name.
Tessellated solid class.
Double_t Safety(const Double_t *point, Bool_t in=kTRUE) const override
Safety.
void ResizeCenter(double maxsize)
Resize and center the shape in a box of size maxsize.
int AddVertex(const Vertex_t &vert)
Add a vertex checking for duplicates, returning the vertex index.
bool Contains(const Double_t *point) const override
Contains.
void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm) const override
ComputeNormal interface.
bool FacetCheck(int ifacet) const
Check validity of facet.
void Streamer(TBuffer &) override
Custom streamer which performs Closing on read.
Double_t SafetyKernel(const Double_t *point, bool in, int *closest_facet_id=nullptr) const
a reusable safety kernel, which optionally returns the closest face
void Print(Option_t *option="") const override
Prints basic info.
Double_t Capacity() const override
Capacity.
Double_t DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=nullptr) const override
DistFromOutside.
void * fBVH
to know if shape still needs closure/initialization
Double_t DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=nullptr) const override
DistFromOutside.
void SetSegsAndPols(TBuffer3D &buff) const override
Fills TBuffer3D structure for segments and polygons.
void BuildBVH()
BuildBVH.
const TBuffer3D & GetBuffer3D(int reqSections, Bool_t localFrame) const override
Fills a static 3D buffer and returns a reference.
void SetPoints(double *points) const override
Fill tessellated points to an array.
bool CheckClosure(bool fixFlipped=true, bool verbose=true)
Check closure of the solid and check/fix flipped normals.
int GetNvertices() const
bool fDefined
! Shape fully defined
Vertex_t FacetComputeNormal(int ifacet, bool &degenerated) const
Compute normal for a given facet.
void CloseShape(bool check=true, bool fixFlipped=true, bool verbose=true)
Close the shape: calculate bounding box and compact vertices.
Tessellated::Vertex_t Vertex_t
void GetMeshNumbers(int &nvert, int &nsegs, int &npols) const override
Returns numbers of vertices, segments and polygons composing the shape mesh.
std::vector< TGeoFacet > fFacets
static TGeoTessellated * ImportFromObjFormat(const char *objfile, bool check=false, bool verbose=false)
Reader from .obj format.
TBuffer3D * MakeBuffer3D() const override
Creates a TBuffer3D describing this shape.
void ComputeBBox() override
Compute bounding box.
std::multimap< long, int > fVerticesMap
! Temporary map used to deduplicate vertices
std::vector< Vertex_t > fOutwardNormals
bool AddFacet(const Vertex_t &pt0, const Vertex_t &pt1, const Vertex_t &pt2)
Adding a triangular facet from vertex positions in absolute coordinates.
void CalculateNormals()
Calculate the normals.
static TClass * Class()
std::vector< Vertex_t > fVertices
int GetNfacets() const
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:1081
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1095
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition TObject.cxx:1123
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition TObject.cxx:1069
This builder is only a wrapper around all the other builders, which selects the best builder dependin...
static BVH_ALWAYS_INLINE Bvh< Node > build(ThreadPool &thread_pool, std::span< const BBox > bboxes, std::span< const Vec > centers, const Config &config={})
Build a BVH in parallel using the given thread pool.
TLine * line
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:249
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:197
ROOT::Geom::Vertex_t Vertex_t
bool contains(bvh::v2::BBox< T, 3 > const &box, bvh::v2::Vec< T, 3 > const &p)
auto SafetySqToNode(bvh::v2::BBox< T, 3 > const &box, bvh::v2::Vec< T, 3 > const &p)
BVH_ALWAYS_INLINE Vec< T, 3 > cross(const Vec< T, 3 > &a, const Vec< T, 3 > &b)
Definition vec.h:104
BVH_ALWAYS_INLINE T dot(const Vec< T, N > &a, const Vec< T, N > &b)
Definition vec.h:98
Definition bbox.h:9
static Vertex_t Cross(Vertex_t const &left, Vertex_t const &right)
The cross (vector) product of two Vector3D<T> objects.
static double Dot(Vertex_t const &left, Vertex_t const &right)
The dot product of two vector objects.
Definition TGeoVector3.h:97
Quality quality
The quality of the BVH produced by the builder.
Growing stack that can be used for BVH traversal.
Definition stack.h:34
Binary BVH node, containing its bounds and an index into its children or the primitives it contains.
Definition node.h:23