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
5/*************************************************************************
6 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
7 * All rights reserved. *
8 * *
9 * For the licensing terms see $ROOTSYS/LICENSE. *
10 * For the list of contributors see $ROOTSYS/README/CREDITS. *
11 *************************************************************************/
12
13/** \class TGeoTessellated
14\ingroup Geometry_classes
15
16Tessellated solid class. It is composed by a set of planar faces having triangular or
17quadrilateral shape. The class does not provide navigation functionality, it just wraps the data
18for the composing faces.
19*/
20
21#include <iostream>
22#include <sstream>
23
24#include "TGeoManager.h"
25#include "TGeoMatrix.h"
26#include "TGeoVolume.h"
27#include "TVirtualGeoPainter.h"
28#include "TGeoTessellated.h"
29#include "TBuffer3D.h"
30#include "TBuffer3DTypes.h"
31#include "TMath.h"
32
33#include <array>
34#include <vector>
35
37
39
40////////////////////////////////////////////////////////////////////////////////
41/// Compact consecutive equal vertices
42
43int TGeoFacet::CompactFacet(Vertex_t *vert, int nvertices)
44{
45 // Compact the common vertices and return new facet
46 if (nvertices < 2)
47 return nvertices;
48 int nvert = nvertices;
49 int i = 0;
50 while (i < nvert) {
51 if (vert[(i + 1) % nvert] == vert[i]) {
52 // shift last vertices left by one element
53 for (int j = i + 2; j < nvert; ++j)
54 vert[j - 1] = vert[j];
55 nvert--;
56 }
57 i++;
58 }
59 return nvert;
60}
61
62////////////////////////////////////////////////////////////////////////////////
63/// Check if a connected neighbour facet has compatible normal
64
65bool TGeoFacet::IsNeighbour(const TGeoFacet &other, bool &flip) const
66{
67
68 // Find a connecting segment
69 bool neighbour = false;
70 int line1[2], line2[2];
71 int npoints = 0;
72 for (int i = 0; i < fNvert; ++i) {
73 auto ivert = fIvert[i];
74 // Check if the other facet has the same vertex
75 for (int j = 0; j < other.GetNvert(); ++j) {
76 if (ivert == other[j]) {
77 line1[npoints] = i;
78 line2[npoints] = j;
79 if (++npoints == 2) {
80 neighbour = true;
81 bool order1 = line1[1] == line1[0] + 1;
82 bool order2 = line2[1] == (line2[0] + 1) % other.GetNvert();
83 flip = (order1 == order2);
84 return neighbour;
85 }
86 }
87 }
88 }
89 return neighbour;
90}
91
92////////////////////////////////////////////////////////////////////////////////
93/// Constructor. In case nfacets is zero, it is user's responsibility to
94/// call CloseShape once all faces are defined.
95
96TGeoTessellated::TGeoTessellated(const char *name, int nfacets) : TGeoBBox(name, 0, 0, 0)
97{
98 fNfacets = nfacets;
99 if (nfacets)
100 fFacets.reserve(nfacets);
101}
102
103////////////////////////////////////////////////////////////////////////////////
104/// Constructor providing directly the array of vertices. Facets have to be added
105/// providing vertex indices rather than coordinates.
106
107TGeoTessellated::TGeoTessellated(const char *name, const std::vector<Vertex_t> &vertices) : TGeoBBox(name, 0, 0, 0)
108{
109 fVertices = vertices;
110 fNvert = fVertices.size();
111}
112
113////////////////////////////////////////////////////////////////////////////////
114/// Add a vertex checking for duplicates, returning the vertex index
115
117{
118 constexpr double tolerance = 1.e-10;
119 auto vertexHash = [&](Vertex_t const &vertex) {
120 // Compute hash for the vertex
121 long hash = 0;
122 // helper function to generate hash from integer numbers
123 auto hash_combine = [](long seed, const long value) {
124 return seed ^ (std::hash<long>{}(value) + 0x9e3779b9 + (seed << 6) + (seed >> 2));
125 };
126 for (int i = 0; i < 3; i++) {
127 // use tolerance to generate int with the desired precision from a real number for hashing
128 hash = hash_combine(hash, std::roundl(vertex[i] / tolerance));
129 }
130 return hash;
131 };
132
133 auto hash = vertexHash(vert);
134 bool isAdded = false;
135 int ivert = -1;
136 // Get the compatible vertices
137 auto range = fVerticesMap.equal_range(hash);
138 for (auto it = range.first; it != range.second; ++it) {
139 ivert = it->second;
140 if (fVertices[ivert] == vert) {
141 isAdded = true;
142 break;
143 }
144 }
145 if (!isAdded) {
146 ivert = fVertices.size();
147 fVertices.push_back(vert);
148 fVerticesMap.insert(std::make_pair(hash, ivert));
149 }
150 return ivert;
151}
152
153////////////////////////////////////////////////////////////////////////////////
154/// Adding a triangular facet from vertex positions in absolute coordinates
155
156bool TGeoTessellated::AddFacet(const Vertex_t &pt0, const Vertex_t &pt1, const Vertex_t &pt2)
157{
158 if (fDefined) {
159 Error("AddFacet", "Shape %s already fully defined. Not adding", GetName());
160 return false;
161 }
162
163 Vertex_t vert[3];
164 vert[0] = pt0;
165 vert[1] = pt1;
166 vert[2] = pt2;
167 int nvert = TGeoFacet::CompactFacet(vert, 3);
168 if (nvert < 3) {
169 Error("AddFacet", "Triangular facet at index %d degenerated. Not adding.", GetNfacets());
170 return false;
171 }
172 int ind[3];
173 for (auto i = 0; i < 3; ++i)
174 ind[i] = AddVertex(vert[i]);
175 fNseg += 3;
176 fFacets.emplace_back(ind[0], ind[1], ind[2]);
177
178 if (fNfacets > 0 && GetNfacets() == fNfacets)
179 CloseShape();
180 return true;
181}
182
183////////////////////////////////////////////////////////////////////////////////
184/// Adding a triangular facet from indices of vertices
185
186bool TGeoTessellated::AddFacet(int i0, int i1, int i2)
187{
188 if (fDefined) {
189 Error("AddFacet", "Shape %s already fully defined. Not adding", GetName());
190 return false;
191 }
192 if (fVertices.empty()) {
193 Error("AddFacet", "Shape %s Cannot add facets by indices without vertices. Not adding", GetName());
194 return false;
195 }
196
197 fNseg += 3;
198 fFacets.emplace_back(i0, i1, i2);
199 return true;
200}
201
202////////////////////////////////////////////////////////////////////////////////
203/// Adding a quadrilateral facet from vertex positions in absolute coordinates
204
205bool TGeoTessellated::AddFacet(const Vertex_t &pt0, const Vertex_t &pt1, const Vertex_t &pt2, const Vertex_t &pt3)
206{
207 if (fDefined) {
208 Error("AddFacet", "Shape %s already fully defined. Not adding", GetName());
209 return false;
210 }
211 Vertex_t vert[4];
212 vert[0] = pt0;
213 vert[1] = pt1;
214 vert[2] = pt2;
215 vert[3] = pt3;
216 int nvert = TGeoFacet::CompactFacet(vert, 4);
217 if (nvert < 3) {
218 Error("AddFacet", "Quadrilateral facet at index %d degenerated. Not adding.", GetNfacets());
219 return false;
220 }
221
222 int ind[4];
223 for (auto i = 0; i < nvert; ++i)
224 ind[i] = AddVertex(vert[i]);
225 fNseg += nvert;
226 if (nvert == 3)
227 fFacets.emplace_back(ind[0], ind[1], ind[2]);
228 else
229 fFacets.emplace_back(ind[0], ind[1], ind[2], ind[3]);
230
231 if (fNfacets > 0 && GetNfacets() == fNfacets)
232 CloseShape(false);
233 return true;
234}
235
236////////////////////////////////////////////////////////////////////////////////
237/// Adding a quadrilateral facet from indices of vertices
238
239bool TGeoTessellated::AddFacet(int i0, int i1, int i2, int i3)
240{
241 if (fDefined) {
242 Error("AddFacet", "Shape %s already fully defined. Not adding", GetName());
243 return false;
244 }
245 if (fVertices.empty()) {
246 Error("AddFacet", "Shape %s Cannot add facets by indices without vertices. Not adding", GetName());
247 return false;
248 }
249
250 fNseg += 4;
251 fFacets.emplace_back(i0, i1, i2, i3);
252 return true;
253}
254
255////////////////////////////////////////////////////////////////////////////////
256/// Compute normal for a given facet
257
258Vertex_t TGeoTessellated::FacetComputeNormal(int ifacet, bool &degenerated) const
259{
260 // Compute normal using non-zero segments
261 constexpr double kTolerance = 1.e-20;
262 auto const &facet = fFacets[ifacet];
263 int nvert = facet.GetNvert();
264 degenerated = true;
265 Vertex_t normal;
266 for (int i = 0; i < nvert - 1; ++i) {
267 Vertex_t e1 = fVertices[facet[i + 1]] - fVertices[facet[i]];
268 if (e1.Mag2() < kTolerance)
269 continue;
270 for (int j = i + 1; j < nvert; ++j) {
271 Vertex_t e2 = fVertices[facet[(j + 1) % nvert]] - fVertices[facet[j]];
272 if (e2.Mag2() < kTolerance)
273 continue;
274 normal = Vertex_t::Cross(e1, e2);
275 // e1 and e2 may be colinear
276 if (normal.Mag2() < kTolerance)
277 continue;
278 normal.Normalize();
279 degenerated = false;
280 break;
281 }
282 if (!degenerated)
283 break;
284 }
285 return normal;
286}
287
288////////////////////////////////////////////////////////////////////////////////
289/// Check validity of facet
290
291bool TGeoTessellated::FacetCheck(int ifacet) const
292{
293 constexpr double kTolerance = 1.e-10;
294 auto const &facet = fFacets[ifacet];
295 int nvert = facet.GetNvert();
296 bool degenerated = true;
297 FacetComputeNormal(ifacet, degenerated);
298 if (degenerated) {
299 std::cout << "Facet: " << ifacet << " is degenerated\n";
300 return false;
301 }
302
303 // Compute surface area
304 double surfaceArea = 0.;
305 for (int i = 1; i < nvert - 1; ++i) {
306 Vertex_t e1 = fVertices[facet[i]] - fVertices[facet[0]];
307 Vertex_t e2 = fVertices[facet[i + 1]] - fVertices[facet[0]];
308 surfaceArea += 0.5 * Vertex_t::Cross(e1, e2).Mag();
309 }
310 if (surfaceArea < kTolerance) {
311 std::cout << "Facet: " << ifacet << " has zero surface area\n";
312 return false;
313 }
314
315 return true;
316}
317
318////////////////////////////////////////////////////////////////////////////////
319/// Close the shape: calculate bounding box and compact vertices
320
321void TGeoTessellated::CloseShape(bool check, bool fixFlipped, bool verbose)
322{
323 // Compute bounding box
324 fDefined = true;
325 fNvert = fVertices.size();
326 fNfacets = fFacets.size();
327 ComputeBBox();
328 // Cleanup the vertex map
329 std::multimap<long, int>().swap(fVerticesMap);
330
331 if (fVertices.size() > 0) {
332 if (!check)
333 return;
334
335 // Check facets
336 for (auto i = 0; i < fNfacets; ++i)
337 FacetCheck(i);
338
339 fClosedBody = CheckClosure(fixFlipped, verbose);
340 }
341}
342
343////////////////////////////////////////////////////////////////////////////////
344/// Check closure of the solid and check/fix flipped normals
345
346bool TGeoTessellated::CheckClosure(bool fixFlipped, bool verbose)
347{
348 int *nn = new int[fNfacets];
349 bool *flipped = new bool[fNfacets];
350 bool hasorphans = false;
351 bool hasflipped = false;
352 for (int i = 0; i < fNfacets; ++i) {
353 nn[i] = 0;
354 flipped[i] = false;
355 }
356
357 for (int icrt = 0; icrt < fNfacets; ++icrt) {
358 // all neighbours checked?
359 if (nn[icrt] >= fFacets[icrt].GetNvert())
360 continue;
361 for (int i = icrt + 1; i < fNfacets; ++i) {
362 bool isneighbour = fFacets[icrt].IsNeighbour(fFacets[i], flipped[i]);
363 if (isneighbour) {
364 if (flipped[icrt])
365 flipped[i] = !flipped[i];
366 if (flipped[i])
367 hasflipped = true;
368 nn[icrt]++;
369 nn[i]++;
370 if (nn[icrt] == fFacets[icrt].GetNvert())
371 break;
372 }
373 }
374 if (nn[icrt] < fFacets[icrt].GetNvert())
375 hasorphans = true;
376 }
377
378 if (hasorphans && verbose) {
379 Error("Check", "Tessellated solid %s has following not fully connected facets:", GetName());
380 for (int icrt = 0; icrt < fNfacets; ++icrt) {
381 if (nn[icrt] < fFacets[icrt].GetNvert())
382 std::cout << icrt << " (" << fFacets[icrt].GetNvert() << " edges, " << nn[icrt] << " neighbours)\n";
383 }
384 }
385 fClosedBody = !hasorphans;
386 int nfixed = 0;
387 if (hasflipped) {
388 if (verbose)
389 Warning("Check", "Tessellated solid %s has following facets with flipped normals:", GetName());
390 for (int icrt = 0; icrt < fNfacets; ++icrt) {
391 if (flipped[icrt]) {
392 if (verbose)
393 std::cout << icrt << "\n";
394 if (fixFlipped) {
395 fFacets[icrt].Flip();
396 nfixed++;
397 }
398 }
399 }
400 if (nfixed && verbose)
401 Info("Check", "Automatically flipped %d facets to match first defined facet", nfixed);
402 }
403 delete[] nn;
404 delete[] flipped;
405
406 return !hasorphans;
407}
408
409////////////////////////////////////////////////////////////////////////////////
410/// Compute bounding box
411
413{
414 const double kBig = TGeoShape::Big();
415 double vmin[3] = {kBig, kBig, kBig};
416 double vmax[3] = {-kBig, -kBig, -kBig};
417 for (const auto &facet : fFacets) {
418 for (int i = 0; i < facet.GetNvert(); ++i) {
419 for (int j = 0; j < 3; ++j) {
420 vmin[j] = TMath::Min(vmin[j], fVertices[facet[i]].operator[](j));
421 vmax[j] = TMath::Max(vmax[j], fVertices[facet[i]].operator[](j));
422 }
423 }
424 }
425 fDX = 0.5 * (vmax[0] - vmin[0]);
426 fDY = 0.5 * (vmax[1] - vmin[1]);
427 fDZ = 0.5 * (vmax[2] - vmin[2]);
428 for (int i = 0; i < 3; ++i)
429 fOrigin[i] = 0.5 * (vmax[i] + vmin[i]);
430}
431
432////////////////////////////////////////////////////////////////////////////////
433/// Returns numbers of vertices, segments and polygons composing the shape mesh.
434
435void TGeoTessellated::GetMeshNumbers(int &nvert, int &nsegs, int &npols) const
436{
437 nvert = fNvert;
438 nsegs = fNseg;
439 npols = GetNfacets();
440}
441
442////////////////////////////////////////////////////////////////////////////////
443/// Creates a TBuffer3D describing *this* shape.
444/// Coordinates are in local reference frame.
445
447{
448 const int nvert = fNvert;
449 const int nsegs = fNseg;
450 const int npols = GetNfacets();
451 auto buff = new TBuffer3D(TBuffer3DTypes::kGeneric, nvert, 3 * nvert, nsegs, 3 * nsegs, npols, 6 * npols);
452 if (buff) {
453 SetPoints(buff->fPnts);
454 SetSegsAndPols(*buff);
455 }
456 return buff;
457}
458
459////////////////////////////////////////////////////////////////////////////////
460/// Prints basic info
461
463{
464 std::cout << "=== Tessellated shape " << GetName() << " having " << GetNvertices() << " vertices and "
465 << GetNfacets() << " facets\n";
466}
467
468////////////////////////////////////////////////////////////////////////////////
469/// Fills TBuffer3D structure for segments and polygons.
470
472{
473 const int c = GetBasicColor();
474 int *segs = buff.fSegs;
475 int *pols = buff.fPols;
476
477 int indseg = 0; // segment internal data index
478 int indpol = 0; // polygon internal data index
479 int sind = 0; // segment index
480 for (const auto &facet : fFacets) {
481 auto nvert = facet.GetNvert();
482 pols[indpol++] = c;
483 pols[indpol++] = nvert;
484 for (auto j = 0; j < nvert; ++j) {
485 int k = (j + 1) % nvert;
486 // segment made by next consecutive points
487 segs[indseg++] = c;
488 segs[indseg++] = facet[j];
489 segs[indseg++] = facet[k];
490 // add segment to current polygon and increment segment index
491 pols[indpol + nvert - j - 1] = sind++;
492 }
493 indpol += nvert;
494 }
495}
496
497////////////////////////////////////////////////////////////////////////////////
498/// Fill tessellated points to an array.
499
501{
502 int ind = 0;
503 for (const auto &vertex : fVertices) {
504 vertex.CopyTo(&points[ind]);
505 ind += 3;
506 }
507}
508
509////////////////////////////////////////////////////////////////////////////////
510/// Fill tessellated points in float.
511
513{
514 int ind = 0;
515 for (const auto &vertex : fVertices) {
516 points[ind++] = vertex.x();
517 points[ind++] = vertex.y();
518 points[ind++] = vertex.z();
519 }
520}
521
522////////////////////////////////////////////////////////////////////////////////
523/// Resize the shape by scaling vertices within maxsize and center to origin
524
526{
527 using Vector3_t = Vertex_t;
528
529 if (!fDefined) {
530 Error("ResizeCenter", "Not all faces are defined");
531 return;
532 }
533 Vector3_t origin(fOrigin[0], fOrigin[1], fOrigin[2]);
534 double maxedge = TMath::Max(TMath::Max(fDX, fDY), fDZ);
535 double scale = maxsize / maxedge;
536 for (size_t i = 0; i < fVertices.size(); ++i) {
537 fVertices[i] = scale * (fVertices[i] - origin);
538 }
539 fOrigin[0] = fOrigin[1] = fOrigin[2] = 0;
540 fDX *= scale;
541 fDY *= scale;
542 fDZ *= scale;
543}
544
545////////////////////////////////////////////////////////////////////////////////
546/// Fills a static 3D buffer and returns a reference.
547
548const TBuffer3D &TGeoTessellated::GetBuffer3D(int reqSections, Bool_t localFrame) const
549{
550 static TBuffer3D buffer(TBuffer3DTypes::kGeneric);
551
552 FillBuffer3D(buffer, reqSections, localFrame);
553
554 const int nvert = fNvert;
555 const int nsegs = fNseg;
556 const int npols = GetNfacets();
557
558 if (reqSections & TBuffer3D::kRawSizes) {
559 if (buffer.SetRawSizes(nvert, 3 * nvert, nsegs, 3 * nsegs, npols, 6 * npols)) {
561 }
562 }
563 if ((reqSections & TBuffer3D::kRaw) && buffer.SectionsValid(TBuffer3D::kRawSizes)) {
564 SetPoints(buffer.fPnts);
565 if (!buffer.fLocalFrame) {
566 TransformPoints(buffer.fPnts, buffer.NbPnts());
567 }
568
569 SetSegsAndPols(buffer);
571 }
572
573 return buffer;
574}
575
576////////////////////////////////////////////////////////////////////////////////
577/// Reads a single tessellated solid from an .obj file.
578
579TGeoTessellated *TGeoTessellated::ImportFromObjFormat(const char *objfile, bool check, bool verbose)
580{
581 using std::vector, std::string, std::ifstream, std::stringstream, std::endl;
582
583 vector<Vertex_t> vertices;
584 vector<string> sfacets;
585
586 struct FacetInd_t {
587 int i0 = -1;
588 int i1 = -1;
589 int i2 = -1;
590 int i3 = -1;
591 int nvert = 0;
592 FacetInd_t(int a, int b, int c)
593 {
594 i0 = a;
595 i1 = b;
596 i2 = c;
597 nvert = 3;
598 };
599 FacetInd_t(int a, int b, int c, int d)
600 {
601 i0 = a;
602 i1 = b;
603 i2 = c;
604 i3 = d;
605 nvert = 4;
606 };
607 };
608
609 vector<FacetInd_t> facets;
610 // List of geometric vertices, with (x, y, z [,w]) coordinates, w is optional and defaults to 1.0.
611 // struct vtx_t { double x = 0; double y = 0; double z = 0; double w = 1; };
612
613 // Texture coordinates in u, [,v ,w]) coordinates, these will vary between 0 and 1. v, w are optional and default to
614 // 0.
615 // struct tex_t { double u; double v; double w; };
616
617 // List of vertex normals in (x,y,z) form; normals might not be unit vectors.
618 // struct vn_t { double x; double y; double z; };
619
620 // Parameter space vertices in ( u [,v] [,w] ) form; free form geometry statement
621 // struct vp_t { double u; double v; double w; };
622
623 // Faces are defined using lists of vertex, texture and normal indices which start at 1.
624 // Polygons such as quadrilaterals can be defined by using more than three vertex/texture/normal indices.
625 // f v1//vn1 v2//vn2 v3//vn3 ...
626
627 // Records starting with the letter "l" specify the order of the vertices which build a polyline.
628 // l v1 v2 v3 v4 v5 v6 ...
629
630 string line;
631 int ind[4] = {0};
632 ifstream file(objfile);
633 if (!file.is_open()) {
634 ::Error("TGeoTessellated::ImportFromObjFormat", "Unable to open %s", objfile);
635 return nullptr;
636 }
637
638 while (getline(file, line)) {
639 stringstream ss(line);
640 string tag;
641
642 // We ignore everything which is not a vertex or a face
643 if (line.rfind('v', 0) == 0 && line.rfind("vt", 0) != 0 && line.rfind("vn", 0) != 0 && line.rfind("vn", 0) != 0) {
644 // Decode the vertex
645 double pos[4] = {0, 0, 0, 1};
646 ss >> tag >> pos[0] >> pos[1] >> pos[2] >> pos[3];
647 vertices.emplace_back(pos[0] * pos[3], pos[1] * pos[3], pos[2] * pos[3]);
648 }
649
650 else if (line.rfind('f', 0) == 0) {
651 // Decode the face
652 ss >> tag;
653 string word;
654 sfacets.clear();
655 while (ss >> word)
656 sfacets.push_back(word);
657 if (sfacets.size() > 4 || sfacets.size() < 3) {
658 ::Error("TGeoTessellated::ImportFromObjFormat", "Detected face having unsupported %zu vertices",
659 sfacets.size());
660 return nullptr;
661 }
662 int nvert = 0;
663 for (auto &sword : sfacets) {
664 stringstream ssword(sword);
665 string token;
666 getline(ssword, token, '/'); // just need the vertex index, which is the first token
667 // Convert string token to integer
668
669 ind[nvert++] = stoi(token) - 1;
670 if (ind[nvert - 1] < 0) {
671 ::Error("TGeoTessellated::ImportFromObjFormat", "Unsupported relative vertex index definition in %s",
672 objfile);
673 return nullptr;
674 }
675 }
676 if (nvert == 3)
677 facets.emplace_back(ind[0], ind[1], ind[2]);
678 else
679 facets.emplace_back(ind[0], ind[1], ind[2], ind[3]);
680 }
681 }
682
683 int nvertices = (int)vertices.size();
684 int nfacets = (int)facets.size();
685 if (nfacets < 3) {
686 ::Error("TGeoTessellated::ImportFromObjFormat", "Not enough faces detected in %s", objfile);
687 return nullptr;
688 }
689
690 string sobjfile(objfile);
691 if (verbose)
692 std::cout << "Read " << nvertices << " vertices and " << nfacets << " facets from " << sobjfile << endl;
693
694 auto tsl = new TGeoTessellated(sobjfile.erase(sobjfile.find_last_of('.')).c_str(), vertices);
695
696 for (int i = 0; i < nfacets; ++i) {
697 auto facet = facets[i];
698 if (facet.nvert == 3)
699 tsl->AddFacet(facet.i0, facet.i1, facet.i2);
700 else
701 tsl->AddFacet(facet.i0, facet.i1, facet.i2, facet.i3);
702 }
703 tsl->CloseShape(check, true, verbose);
704 tsl->Print();
705 return tsl;
706}
#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
float Float_t
Definition RtypesCore.h:57
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:382
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:110
segment * segs
Definition X3DBuffer.c:23
Generic 3D primitive description class.
Definition TBuffer3D.h:18
Int_t * fPols
Definition TBuffer3D.h:115
UInt_t NbPnts() const
Definition TBuffer3D.h:80
Bool_t SectionsValid(UInt_t mask) const
Definition TBuffer3D.h:67
void SetSectionsValid(UInt_t mask)
Definition TBuffer3D.h:65
Int_t * fSegs
Definition TBuffer3D.h:114
Bool_t fLocalFrame
Definition TBuffer3D.h:90
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:113
void FillBuffer3D(TBuffer3D &buffer, Int_t reqSections, Bool_t localFrame) const override
Fill the supplied buffer, with sections in desired frame See TBuffer3D.h for explanation of sections,...
Double_t fDX
Definition TGeoBBox.h:20
Double_t fOrigin[3]
Definition TGeoBBox.h:23
Double_t fDY
Definition TGeoBBox.h:21
Double_t fDZ
Definition TGeoBBox.h:22
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.
int GetNvert() const
static Double_t Big()
Definition TGeoShape.h:87
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.
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 FacetCheck(int ifacet) const
Check validity of facet.
void Print(Option_t *option="") const override
Prints basic info.
void SetSegsAndPols(TBuffer3D &buff) const override
Fills TBuffer3D structure for segments and polygons.
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
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
bool AddFacet(const Vertex_t &pt0, const Vertex_t &pt1, const Vertex_t &pt2)
Adding a triangular facet from vertex positions in absolute coordinates.
std::vector< Vertex_t > fVertices
int GetNfacets() const
bool fClosedBody
Shape fully defined.
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:979
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:993
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition TObject.cxx:967
TLine * line
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:250
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:198
ROOT::Geom::Vertex_t Vertex_t
double Mag() const
double Mag2() const
static Vertex_t Cross(Vertex_t const &left, Vertex_t const &right)
The cross (vector) product of two Vector3D<T> objects.
void Normalize()
Normalizes the vector by dividing each entry by the length.