Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeoXtru.cxx
Go to the documentation of this file.
1// @(#)root/geom:$Id$
2// Author: Mihaela Gheata 24/01/04
3
4/*************************************************************************
5 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12/** \class TGeoXtru
13\ingroup Shapes_classes
14A TGeoXtru shape is represented by the extrusion of an arbitrary
15polygon with fixed outline between several Z sections. Each Z section is
16a scaled version of the same "blueprint" polygon. Different global XY
17translations are allowed from section to section. Corresponding polygon
18vertices from consecutive sections are connected.
19
20An extruded polygon can be created using the constructor:
21
22~~~ {.cpp}
23TGeoXtru::TGeoXtru(Int_t nplanes);
24~~~
25
26 - `nplanes: `number of Z sections (minimum 2)
27
28Begin_Macro
29{
30 new TGeoManager("xtru", "poza12");
31 TGeoMaterial *mat = new TGeoMaterial("Al", 26.98,13,2.7);
32 TGeoMedium *med = new TGeoMedium("MED",1,mat);
33 TGeoVolume *top = gGeoManager->MakeBox("TOP",med,100,100,100);
34 gGeoManager->SetTopVolume(top);
35 TGeoVolume *vol = gGeoManager->MakeXtru("XTRU",med,4);
36 TGeoXtru *xtru = (TGeoXtru*)vol->GetShape();
37 Double_t x[8] = {-30,-30,30,30,15,15,-15,-15};
38 Double_t y[8] = {-30,30,30,-30,-30,15,15,-30};
39 xtru->DefinePolygon(8,x,y);
40 xtru->DefineSection(0,-40, -20., 10., 1.5);
41 xtru->DefineSection(1, 10, 0., 0., 0.5);
42 xtru->DefineSection(2, 10, 0., 0., 0.7);
43 xtru->DefineSection(3, 40, 10., 20., 0.9);
44 top->AddNode(vol,1);
45 gGeoManager->CloseGeometry();
46 gGeoManager->SetNsegments(80);
47 top->Draw();
48 if (gPad) {
49 TView *view = gPad->GetView();
50 if (view) view->ShowAxis();
51 }
52}
53End_Macro
54
55The lists of X and Y positions for all vertices have to be provided for
56the "blueprint" polygon:
57
58~~~{.cpp}
59TGeoXtru::DefinePolygon (Int_t nvertices, Double_t *xv,
60Double_t *yv);
61~~~
62
63 - `nvertices: `number of vertices of the polygon
64 - `xv,yv: `arrays of X and Y coordinates for polygon vertices
65
66The method creates an object of the class **`TGeoPolygon`** for which
67the convexity is automatically determined . The polygon is decomposed
68into convex polygons if needed.
69
70Next step is to define the Z positions for each section plane as well as
71the XY offset and scaling for the corresponding polygons.
72
73~~~{.cpp}
74TGeoXtru::DefineSection(Int_t snum,Double_t zsection,Double_t x0,
75Double_t y0, Double_t scale);
76~~~
77
78 - `snum: `Z section index (0, nplanes-1). The section with
79 `snum = nplanes-1` must be defined last and triggers the computation
80 of the bounding box for the whole shape
81 - `zsection: `Z position of section `snum`. Sections must be defined
82 in increasing order of Z (e.g. `snum=0` correspond to the minimum Z
83 and `snum=nplanes-1` to the maximum one).
84 - `x0,y0: `offset of section `snum` with respect to the local shape
85 reference frame `T`
86 - `scale: `factor that multiplies the X/Y coordinates of each vertex
87 of the polygon at section `snum`:
88 - `x[ivert] = x0 + scale*xv[ivert]`
89 - `y[ivert] = y0 + scale*yv[ivert]`
90*/
91
92#include "TGeoXtru.h"
93
94#include <iostream>
95
96#include "TBuffer3D.h"
97#include "TBuffer3DTypes.h"
98#include "TMath.h"
99
100#include "TVirtualGeoPainter.h"
101#include "TGeoManager.h"
102#include "TGeoVolume.h"
103#include "TGeoPolygon.h"
104
106
107////////////////////////////////////////////////////////////////////////////////
108/// Constructor.
109
110TGeoXtru::ThreadData_t::ThreadData_t() : fSeg(0), fIz(0), fXc(nullptr), fYc(nullptr), fPoly(nullptr) {}
111
112////////////////////////////////////////////////////////////////////////////////
113/// Destructor.
114
116{
117 delete[] fXc;
118 delete[] fYc;
119 delete fPoly;
120}
121
122////////////////////////////////////////////////////////////////////////////////
123
125{
126 if (!fThreadSize)
127 ((TGeoXtru *)this)->CreateThreadData(1);
129 return *fThreadData[tid];
130}
131
132////////////////////////////////////////////////////////////////////////////////
133
135{
136 std::lock_guard<std::mutex> guard(fMutex);
137 std::vector<ThreadData_t *>::iterator i = fThreadData.begin();
138 while (i != fThreadData.end()) {
139 delete *i;
140 ++i;
141 }
142 fThreadData.clear();
143 fThreadSize = 0;
144}
145
146////////////////////////////////////////////////////////////////////////////////
147/// Create thread data for n threads max.
148
150{
151 std::lock_guard<std::mutex> guard(fMutex);
152 fThreadData.resize(nthreads);
154 for (Int_t tid = 0; tid < nthreads; tid++) {
155 if (fThreadData[tid] == nullptr) {
158 td.fXc = new Double_t[fNvert];
159 td.fYc = new Double_t[fNvert];
160 memcpy(td.fXc, fX, fNvert * sizeof(Double_t));
161 memcpy(td.fYc, fY, fNvert * sizeof(Double_t));
162 td.fPoly = new TGeoPolygon(fNvert);
163 td.fPoly->SetXY(td.fXc, td.fYc); // initialize with current coordinates
164 td.fPoly->FinishPolygon();
165 if (tid == 0 && td.fPoly->IsIllegalCheck()) {
166 Error("DefinePolygon", "Shape %s of type XTRU has an illegal polygon.", GetName());
167 }
168 }
169 }
170}
171
172////////////////////////////////////////////////////////////////////////////////
173/// Set current z-plane.
174
176{
177 GetThreadData().fIz = iz;
178}
179////////////////////////////////////////////////////////////////////////////////
180/// Set current segment.
181
186
187////////////////////////////////////////////////////////////////////////////////
188/// dummy ctor
189
191 : TGeoBBox(),
192 fNvert(0),
193 fNz(0),
194 fZcurrent(0.),
195 fX(nullptr),
196 fY(nullptr),
197 fZ(nullptr),
198 fScale(nullptr),
199 fX0(nullptr),
200 fY0(nullptr),
201 fThreadData(0),
202 fThreadSize(0)
203{
205}
206
207////////////////////////////////////////////////////////////////////////////////
208/// Default constructor
209
211 : TGeoBBox(0, 0, 0),
212 fNvert(0),
213 fNz(nz),
214 fZcurrent(0.),
215 fX(nullptr),
216 fY(nullptr),
217 fZ(new Double_t[nz]),
218 fScale(new Double_t[nz]),
219 fX0(new Double_t[nz]),
220 fY0(new Double_t[nz]),
221 fThreadData(0),
222 fThreadSize(0)
223{
225 if (nz < 2) {
226 Error("ctor", "Cannot create TGeoXtru %s with less than 2 Z planes", GetName());
228 return;
229 }
230}
231
232////////////////////////////////////////////////////////////////////////////////
233/// Default constructor in GEANT3 style
234/// - param[0] = nz // number of z planes
235///
236/// - param[1] = z1 // Z position of first plane
237/// - param[2] = x1 // X position of first plane
238/// - param[3] = y1 // Y position of first plane
239/// - param[4] = scale1 // scale factor for first plane
240/// ...
241/// - param[4*(nz-1]+1] = zn
242/// - param[4*(nz-1)+2] = xn
243/// - param[4*(nz-1)+3] = yn
244/// - param[4*(nz-1)+4] = scalen
245
247 : TGeoBBox(0, 0, 0),
248 fNvert(0),
249 fNz(0),
250 fZcurrent(0.),
251 fX(nullptr),
252 fY(nullptr),
253 fZ(nullptr),
254 fScale(nullptr),
255 fX0(nullptr),
256 fY0(nullptr),
257 fThreadData(0),
258 fThreadSize(0)
259{
261 SetDimensions(param);
262}
263
264////////////////////////////////////////////////////////////////////////////////
265/// destructor
266
268{
269 if (fX) {
270 delete[] fX;
271 fX = nullptr;
272 }
273 if (fY) {
274 delete[] fY;
275 fY = nullptr;
276 }
277 if (fZ) {
278 delete[] fZ;
279 fZ = nullptr;
280 }
281 if (fScale) {
282 delete[] fScale;
283 fScale = nullptr;
284 }
285 if (fX0) {
286 delete[] fX0;
287 fX0 = nullptr;
288 }
289 if (fY0) {
290 delete[] fY0;
291 fY0 = nullptr;
292 }
294}
295
296////////////////////////////////////////////////////////////////////////////////
297/// Compute capacity [length^3] of this shape.
298
300{
302 Int_t iz;
303 Double_t capacity = 0;
305 TGeoXtru *xtru = (TGeoXtru *)this;
306 xtru->SetCurrentVertices(0., 0., 1.);
307 area = td.fPoly->Area();
308 for (iz = 0; iz < fNz - 1; iz++) {
309 dz = fZ[iz + 1] - fZ[iz];
311 continue;
312 sc1 = fScale[iz];
313 sc2 = fScale[iz + 1];
314 capacity += (area * dz / 3.) * (sc1 * sc1 + sc1 * sc2 + sc2 * sc2);
315 }
316 return capacity;
317}
318
319////////////////////////////////////////////////////////////////////////////////
320/// compute bounding box of the pcon
321
323{
325 if (!fX || !fZ || !fNvert) {
326 Error("ComputeBBox", "In shape %s polygon not defined", GetName());
328 return;
329 }
330 Double_t zmin = fZ[0];
331 Double_t zmax = fZ[fNz - 1];
336 for (Int_t i = 0; i < fNz; i++) {
337 SetCurrentVertices(fX0[i], fY0[i], fScale[i]);
338 for (Int_t j = 0; j < fNvert; j++) {
339 if (td.fXc[j] < xmin)
340 xmin = td.fXc[j];
341 if (td.fXc[j] > xmax)
342 xmax = td.fXc[j];
343 if (td.fYc[j] < ymin)
344 ymin = td.fYc[j];
345 if (td.fYc[j] > ymax)
346 ymax = td.fYc[j];
347 }
348 }
349 fOrigin[0] = 0.5 * (xmin + xmax);
350 fOrigin[1] = 0.5 * (ymin + ymax);
351 fOrigin[2] = 0.5 * (zmin + zmax);
352 fDX = 0.5 * (xmax - xmin);
353 fDY = 0.5 * (ymax - ymin);
354 fDZ = 0.5 * (zmax - zmin);
355}
356
357////////////////////////////////////////////////////////////////////////////////
358/// Compute normal to closest surface from POINT.
359
360void TGeoXtru::ComputeNormal(const Double_t * /*point*/, const Double_t *dir, Double_t *norm) const
361{
363 if (td.fIz < 0) {
364 memset(norm, 0, 3 * sizeof(Double_t));
365 norm[2] = (dir[2] > 0) ? 1 : -1;
366 return;
367 }
368 Double_t vert[12];
369 GetPlaneVertices(td.fIz, td.fSeg, vert);
371 Double_t ndotd = norm[0] * dir[0] + norm[1] * dir[1] + norm[2] * dir[2];
372 if (ndotd < 0) {
373 norm[0] = -norm[0];
374 norm[1] = -norm[1];
375 norm[2] = -norm[2];
376 }
377}
378
379////////////////////////////////////////////////////////////////////////////////
380/// test if point is inside this shape
381
383{
385 // Check Z range
386 TGeoXtru *xtru = (TGeoXtru *)this;
387 if (point[2] < fZ[0])
388 return kFALSE;
389 if (point[2] > fZ[fNz - 1])
390 return kFALSE;
391 Int_t iz = TMath::BinarySearch(fNz, fZ, point[2]);
392 if (iz < 0 || iz == fNz - 1)
393 return kFALSE;
394 if (TGeoShape::IsSameWithinTolerance(point[2], fZ[iz])) {
395 xtru->SetIz(-1);
396 xtru->SetCurrentVertices(fX0[iz], fY0[iz], fScale[iz]);
397 if (td.fPoly->Contains(point))
398 return kTRUE;
399 if (iz > 1 && TGeoShape::IsSameWithinTolerance(fZ[iz], fZ[iz - 1])) {
400 xtru->SetCurrentVertices(fX0[iz - 1], fY0[iz - 1], fScale[iz - 1]);
401 return td.fPoly->Contains(point);
402 } else if (iz < fNz - 2 && TGeoShape::IsSameWithinTolerance(fZ[iz], fZ[iz + 1])) {
403 xtru->SetCurrentVertices(fX0[iz + 1], fY0[iz + 1], fScale[iz + 1]);
404 return td.fPoly->Contains(point);
405 }
406 }
407 xtru->SetCurrentZ(point[2], iz);
408 if (TMath::Abs(point[2] - fZ[iz]) < TGeoShape::Tolerance() ||
409 TMath::Abs(fZ[iz + 1] - point[2]) < TGeoShape::Tolerance())
410 xtru->SetIz(-1);
411 // Now td.fXc,fYc represent the vertices of the section at point[2]
412 return td.fPoly->Contains(point);
413}
414
415////////////////////////////////////////////////////////////////////////////////
416/// compute closest distance from point px,py to each corner
417
419{
420 const Int_t numPoints = fNvert * fNz;
421 return ShapeDistancetoPrimitive(numPoints, px, py);
422}
423
424////////////////////////////////////////////////////////////////////////////////
425/// Draw the section polygon.
426
428{
430 if (td.fPoly)
431 td.fPoly->Draw(option);
432}
433
434////////////////////////////////////////////////////////////////////////////////
435/// Compute distance to a Xtru lateral surface.
436
438 Bool_t in) const
439{
442 Double_t vert[12];
443 Double_t norm[3];
445 Double_t pt[3];
447 if (TGeoShape::IsSameWithinTolerance(fZ[iz], fZ[iz + 1]) && !in) {
448 TGeoXtru *xtru = (TGeoXtru *)this;
449 snext = (fZ[iz] - point[2]) / dir[2];
450 if (snext < 0)
451 return TGeoShape::Big();
452 pt[0] = point[0] + snext * dir[0];
453 pt[1] = point[1] + snext * dir[1];
454 pt[2] = point[2] + snext * dir[2];
455 if (dir[2] < 0.)
456 xtru->SetCurrentVertices(fX0[iz], fY0[iz], fScale[iz]);
457 else
458 xtru->SetCurrentVertices(fX0[iz + 1], fY0[iz + 1], fScale[iz + 1]);
459 if (!td.fPoly->Contains(pt))
460 return TGeoShape::Big();
461 return snext;
462 }
465 Double_t ndotd = norm[0] * dir[0] + norm[1] * dir[1] + norm[2] * dir[2];
466 if (in) {
467 if (ndotd <= 0)
468 return TGeoShape::Big();
469 safe = (vert[0] - point[0]) * norm[0] + (vert[1] - point[1]) * norm[1] + (vert[2] - point[2]) * norm[2];
470 if (safe < -1.E-8)
471 return TGeoShape::Big(); // direction outwards plane
472 } else {
473 ndotd = -ndotd;
474 if (ndotd <= 0)
475 return TGeoShape::Big();
476 safe = (point[0] - vert[0]) * norm[0] + (point[1] - vert[1]) * norm[1] + (point[2] - vert[2]) * norm[2];
477 if (safe < -1.E-8)
478 return TGeoShape::Big(); // direction outwards plane
479 }
480 snext = safe / ndotd;
481 if (snext > stepmax)
482 return TGeoShape::Big();
483 if (fZ[iz] < fZ[iz + 1]) {
484 znew = point[2] + snext * dir[2];
485 if (znew < fZ[iz])
486 return TGeoShape::Big();
487 if (znew > fZ[iz + 1])
488 return TGeoShape::Big();
489 }
490 pt[0] = point[0] + snext * dir[0];
491 pt[1] = point[1] + snext * dir[1];
492 pt[2] = point[2] + snext * dir[2];
494 return TGeoShape::Big();
495 return TMath::Max(snext, 0.);
496}
497
498////////////////////////////////////////////////////////////////////////////////
499/// compute distance from inside point to surface of the polycone
500/// locate Z segment
501
504{
506 if (iact < 3 && safe) {
507 *safe = Safety(point, kTRUE);
508 if (iact == 0)
509 return TGeoShape::Big();
510 if (iact == 1 && step < *safe)
511 return TGeoShape::Big();
512 }
513 TGeoXtru *xtru = (TGeoXtru *)this;
514 Int_t iz = TMath::BinarySearch(fNz, fZ, point[2]);
515 if (iz < 0) {
516 if (dir[2] <= 0) {
517 xtru->SetIz(-1);
518 return 0.;
519 }
520 iz = 0;
521 }
522 if (iz == fNz - 1) {
523 if (dir[2] >= 0) {
524 xtru->SetIz(-1);
525 return 0.;
526 }
527 iz--;
528 } else {
529 if (iz > 0) {
530 if (TGeoShape::IsSameWithinTolerance(point[2], fZ[iz])) {
531 if (TGeoShape::IsSameWithinTolerance(fZ[iz], fZ[iz + 1]) && dir[2] < 0)
532 iz++;
533 else if (TGeoShape::IsSameWithinTolerance(fZ[iz], fZ[iz - 1]) && dir[2] > 0)
534 iz--;
535 }
536 }
537 }
538 Bool_t convex = td.fPoly->IsConvex();
539 // Double_t stepmax = step;
540 // if (stepmax>TGeoShape::Big()) stepmax = TGeoShape::Big();
542 Double_t dist, sz;
543 Double_t pt[3];
544 Int_t iv, ipl, inext;
545 // we treat the special case when dir[2]=0
546 if (TGeoShape::IsSameWithinTolerance(dir[2], 0)) {
547 for (iv = 0; iv < fNvert; iv++) {
548 xtru->SetIz(-1);
549 dist = DistToPlane(point, dir, iz, iv, TGeoShape::Big(), kTRUE);
550 if (dist < snext) {
551 snext = dist;
552 xtru->SetSeg(iv);
553 if (convex)
554 return snext;
555 }
556 }
557 if (snext < 1.E10)
558 return snext;
559 return TGeoShape::Tolerance();
560 }
561
562 // normal case
563 Int_t incseg = (dir[2] > 0) ? 1 : -1;
564 Int_t iznext = iz;
566 while (iz >= 0 && iz < fNz - 1) {
567 // find the distance to current segment end Z surface
568 ipl = iz + ((incseg + 1) >> 1); // next plane
569 inext = ipl + incseg; // next next plane
570 sz = (fZ[ipl] - point[2]) / dir[2];
571 if (sz < snext) {
572 iznext += incseg;
573 // we cross the next Z section before stepmax
574 pt[0] = point[0] + sz * dir[0];
575 pt[1] = point[1] + sz * dir[1];
576 xtru->SetCurrentVertices(fX0[ipl], fY0[ipl], fScale[ipl]);
577 if (td.fPoly->Contains(pt)) {
578 // ray gets through next polygon - is it the last one?
579 if (ipl == 0 || ipl == fNz - 1) {
580 xtru->SetIz(-1);
581 if (convex)
582 return sz;
583 zexit = kTRUE;
584 snext = sz;
585 }
586 // maybe a Z discontinuity - check this
588 xtru->SetCurrentVertices(fX0[inext], fY0[inext], fScale[inext]);
589 // if we do not cross the next polygone, we are out
590 if (!td.fPoly->Contains(pt)) {
591 xtru->SetIz(-1);
592 if (convex)
593 return sz;
594 zexit = kTRUE;
595 snext = sz;
596 } else {
597 iznext = inext;
598 }
599 }
600 }
601 } else {
602 iznext = fNz - 1; // stop
603 }
604 // ray may cross the lateral surfaces of section iz
605 for (iv = 0; iv < fNvert; iv++) {
606 dist = DistToPlane(point, dir, iz, iv, TGeoShape::Big(), kTRUE);
607 if (dist < snext) {
608 xtru->SetIz(iz);
609 xtru->SetSeg(iv);
610 snext = dist;
611 if (convex)
612 return snext;
613 zexit = kTRUE;
614 }
615 }
616 if (zexit)
617 return snext;
618 iz = iznext;
619 }
620 return TGeoShape::Tolerance();
621}
622
623////////////////////////////////////////////////////////////////////////////////
624/// compute distance from outside point to surface of the tube
625/// Warning("DistFromOutside", "not implemented");
626
629{
631 if (iact < 3 && safe) {
632 *safe = Safety(point, kTRUE);
633 if (iact == 0)
634 return TGeoShape::Big();
635 if (iact == 1 && step < *safe)
636 return TGeoShape::Big();
637 }
638 // Check if the bounding box is crossed within the requested distance
639 Double_t sdist = TGeoBBox::DistFromOutside(point, dir, fDX, fDY, fDZ, fOrigin, step);
640 if (sdist >= step)
641 return TGeoShape::Big();
642 Double_t stepmax = step;
643 if (stepmax > TGeoShape::Big())
645 Double_t snext = 0.;
646 Int_t i, iv;
647 Double_t pt[3];
648 memcpy(pt, point, 3 * sizeof(Double_t));
649 TGeoXtru *xtru = (TGeoXtru *)this;
650 // We might get out easy with Z checks
651 Int_t iz = TMath::BinarySearch(fNz, fZ, point[2]);
652 if (iz < 0) {
653 if (dir[2] <= 0)
654 return TGeoShape::Big();
655 // propagate to first Z plane
656 snext = (fZ[0] - point[2]) / dir[2];
657 if (snext > stepmax)
658 return TGeoShape::Big();
659 for (i = 0; i < 3; i++)
660 pt[i] = point[i] + snext * dir[i];
661 xtru->SetCurrentVertices(fX0[0], fY0[0], fScale[0]);
662 if (td.fPoly->Contains(pt)) {
663 xtru->SetIz(-1);
664 return snext;
665 }
666 iz = 0; // valid starting value = first segment
667 stepmax -= snext;
668 } else {
669 if (iz == fNz - 1) {
670 if (dir[2] >= 0)
671 return TGeoShape::Big();
672 // propagate to last Z plane
673 snext = (fZ[fNz - 1] - point[2]) / dir[2];
674 if (snext > stepmax)
675 return TGeoShape::Big();
676 for (i = 0; i < 3; i++)
677 pt[i] = point[i] + snext * dir[i];
678 xtru->SetCurrentVertices(fX0[fNz - 1], fY0[fNz - 1], fScale[fNz - 1]);
679 if (td.fPoly->Contains(pt)) {
680 xtru->SetIz(-1);
681 return snext;
682 }
683 iz = fNz - 2; // valid value = last segment
684 stepmax -= snext;
685 }
686 }
687 // Check if the bounding box is missed by the track
688 if (!TGeoBBox::Contains(pt)) {
689 Double_t dist = TGeoBBox::DistFromOutside(pt, dir, 3);
690 if (dist > stepmax)
691 return TGeoShape::Big();
692 if (dist > 1E-6)
693 dist -= 1E-6; // decrease snext to make sure we do not cross the xtru
694 else
695 dist = 0;
696 for (i = 0; i < 3; i++)
697 pt[i] += dist * dir[i]; // we are now closer
698 iz = TMath::BinarySearch(fNz, fZ, pt[2]);
699 if (iz < 0)
700 iz = 0;
701 else if (iz == fNz - 1)
702 iz = fNz - 2;
703 snext += dist;
704 stepmax -= dist;
705 }
706 // not the case - we have to do some work...
707 // Start tracking from current iz
708 // - first solve particular case dir[2]=0
709 Bool_t convex = td.fPoly->IsConvex();
710 Bool_t hit = kFALSE;
711 if (TGeoShape::IsSameWithinTolerance(dir[2], 0)) {
712 // loop lateral planes to see if we cross something
713 xtru->SetIz(iz);
714 for (iv = 0; iv < fNvert; iv++) {
715 Double_t dist = DistToPlane(pt, dir, iz, iv, stepmax, kFALSE);
716 if (dist < stepmax) {
717 xtru->SetSeg(iv);
718 if (convex)
719 return (snext + dist);
720 stepmax = dist;
721 hit = kTRUE;
722 }
723 }
724 if (hit)
725 return (snext + stepmax);
726 return TGeoShape::Big();
727 }
728 // general case
729 Int_t incseg = (dir[2] > 0) ? 1 : -1;
730 while (iz >= 0 && iz < fNz - 1) {
731 // compute distance to lateral planes
732 xtru->SetIz(iz);
733 if (TGeoShape::IsSameWithinTolerance(fZ[iz], fZ[iz + 1]))
734 xtru->SetIz(-1);
735 for (iv = 0; iv < fNvert; iv++) {
736 Double_t dist = DistToPlane(pt, dir, iz, iv, stepmax, kFALSE);
737 if (dist < stepmax) {
738 // HIT
739 xtru->SetSeg(iv);
740 if (convex)
741 return (snext + dist);
742 stepmax = dist;
743 hit = kTRUE;
744 }
745 }
746 if (hit)
747 return (snext + stepmax);
748 iz += incseg;
749 }
750 return TGeoShape::Big();
751}
752
753////////////////////////////////////////////////////////////////////////////////
754/// Creates the polygon representing the blueprint of any Xtru section.
755/// - nvert = number of vertices >2
756/// - xv[nvert] = array of X vertex positions
757/// - yv[nvert] = array of Y vertex positions
758///
759/// *NOTE* should be called before DefineSection or ctor with 'param'
760
762{
763 if (nvert < 3) {
764 Error("DefinePolygon", "In shape %s cannot create polygon with less than 3 vertices", GetName());
766 return kFALSE;
767 }
768 for (Int_t i = 0; i < nvert - 1; i++) {
769 for (Int_t j = i + 1; j < nvert; j++) {
771 Error("DefinePolygon", "In shape %s 2 vertices cannot be identical", GetName());
773 // return kFALSE;
774 }
775 }
776 }
777 fNvert = nvert;
778 if (fX)
779 delete[] fX;
780 fX = new Double_t[nvert];
781 if (fY)
782 delete[] fY;
783 fY = new Double_t[nvert];
784 memcpy(fX, xv, nvert * sizeof(Double_t));
785 memcpy(fY, yv, nvert * sizeof(Double_t));
786
788
789 return kTRUE;
790}
791
792////////////////////////////////////////////////////////////////////////////////
793/// defines z position of a section plane, rmin and rmax at this z.
794
796{
797 if ((snum < 0) || (snum >= fNz))
798 return;
799 fZ[snum] = z;
800 fX0[snum] = x0;
801 fY0[snum] = y0;
802 fScale[snum] = scale;
803 if (snum) {
804 if (fZ[snum] < fZ[snum - 1]) {
805 Warning("DefineSection",
806 "In shape: %s, Z position of section "
807 "%i, z=%e, not in increasing order, %i, z=%e",
808 GetName(), snum, fZ[snum], snum - 1, fZ[snum - 1]);
809 return;
810 }
811 }
812 if (snum == (fNz - 1)) {
813 ComputeBBox();
815 InspectShape();
816 }
817}
818
819////////////////////////////////////////////////////////////////////////////////
820/// Return the Z coordinate for segment ipl.
821
823{
824 if (ipl < 0 || ipl > (fNz - 1)) {
825 Error("GetZ", "In shape %s, ipl=%i out of range (0,%i)", GetName(), ipl, fNz - 1);
826 return 0.;
827 }
828 return fZ[ipl];
829}
830////////////////////////////////////////////////////////////////////////////////
831/// Returns normal vector to the planar quadrilateral defined by vector VERT.
832/// The normal points outwards the xtru.
833
835{
836 Double_t cross = 0.;
837 Double_t v1[3], v2[3];
838 v1[0] = vert[9] - vert[0];
839 v1[1] = vert[10] - vert[1];
840 v1[2] = vert[11] - vert[2];
841 v2[0] = vert[3] - vert[0];
842 v2[1] = vert[4] - vert[1];
843 v2[2] = vert[5] - vert[2];
844 norm[0] = v1[1] * v2[2] - v1[2] * v2[1];
845 cross += norm[0] * norm[0];
846 norm[1] = v1[2] * v2[0] - v1[0] * v2[2];
847 cross += norm[1] * norm[1];
848 norm[2] = v1[0] * v2[1] - v1[1] * v2[0];
849 cross += norm[2] * norm[2];
850 if (cross < TGeoShape::Tolerance())
851 return;
852 cross = 1. / TMath::Sqrt(cross);
853 for (Int_t i = 0; i < 3; i++)
854 norm[i] *= cross;
855}
856
857////////////////////////////////////////////////////////////////////////////////
858/// Returns (x,y,z) of 3 vertices of the surface defined by Z sections (iz, iz+1)
859/// and polygon vertices (ivert, ivert+1). No range check.
860
862{
864 Double_t x, y, z1, z2;
865 Int_t iv1 = (ivert + 1) % fNvert;
866 Int_t icrt = 0;
867 z1 = fZ[iz];
868 z2 = fZ[iz + 1];
869 if (td.fPoly->IsClockwise()) {
870 x = fX[ivert] * fScale[iz] + fX0[iz];
871 y = fY[ivert] * fScale[iz] + fY0[iz];
872 vert[icrt++] = x;
873 vert[icrt++] = y;
874 vert[icrt++] = z1;
875 x = fX[iv1] * fScale[iz] + fX0[iz];
876 y = fY[iv1] * fScale[iz] + fY0[iz];
877 vert[icrt++] = x;
878 vert[icrt++] = y;
879 vert[icrt++] = z1;
880 x = fX[iv1] * fScale[iz + 1] + fX0[iz + 1];
881 y = fY[iv1] * fScale[iz + 1] + fY0[iz + 1];
882 vert[icrt++] = x;
883 vert[icrt++] = y;
884 vert[icrt++] = z2;
885 x = fX[ivert] * fScale[iz + 1] + fX0[iz + 1];
886 y = fY[ivert] * fScale[iz + 1] + fY0[iz + 1];
887 vert[icrt++] = x;
888 vert[icrt++] = y;
889 vert[icrt++] = z2;
890 } else {
891 x = fX[iv1] * fScale[iz] + fX0[iz];
892 y = fY[iv1] * fScale[iz] + fY0[iz];
893 vert[icrt++] = x;
894 vert[icrt++] = y;
895 vert[icrt++] = z1;
896 x = fX[ivert] * fScale[iz] + fX0[iz];
897 y = fY[ivert] * fScale[iz] + fY0[iz];
898 vert[icrt++] = x;
899 vert[icrt++] = y;
900 vert[icrt++] = z1;
901 x = fX[ivert] * fScale[iz + 1] + fX0[iz + 1];
902 y = fY[ivert] * fScale[iz + 1] + fY0[iz + 1];
903 vert[icrt++] = x;
904 vert[icrt++] = y;
905 vert[icrt++] = z2;
906 x = fX[iv1] * fScale[iz + 1] + fX0[iz + 1];
907 y = fY[iv1] * fScale[iz + 1] + fY0[iz + 1];
908 vert[icrt++] = x;
909 vert[icrt++] = y;
910 vert[icrt++] = z2;
911 }
912}
913////////////////////////////////////////////////////////////////////////////////
914/// Check if the quadrilateral defined by VERT contains a coplanar POINT.
915
917{
918 Double_t v1[3], v2[3];
919 Double_t cross;
920 Int_t j, k;
921 for (Int_t i = 0; i < 4; i++) { // loop vertices
922 j = 3 * i;
923 k = 3 * ((i + 1) % 4);
924 v1[0] = point[0] - vert[j];
925 v1[1] = point[1] - vert[j + 1];
926 v1[2] = point[2] - vert[j + 2];
927 v2[0] = vert[k] - vert[j];
928 v2[1] = vert[k + 1] - vert[j + 1];
929 v2[2] = vert[k + 2] - vert[j + 2];
930 cross = (v1[1] * v2[2] - v1[2] * v2[1]) * norm[0] + (v1[2] * v2[0] - v1[0] * v2[2]) * norm[1] +
931 (v1[0] * v2[1] - v1[1] * v2[0]) * norm[2];
932 if (cross < 0)
933 return kFALSE;
934 }
935 return kTRUE;
936}
937
938////////////////////////////////////////////////////////////////////////////////
939/// Print actual Xtru parameters.
940
942{
943 printf("*** Shape %s: TGeoXtru ***\n", GetName());
944 printf(" Nz = %i\n", fNz);
945 printf(" List of (x,y) of polygon vertices:\n");
946 for (Int_t ivert = 0; ivert < fNvert; ivert++)
947 printf(" x = %11.5f y = %11.5f\n", fX[ivert], fY[ivert]);
948 for (Int_t ipl = 0; ipl < fNz; ipl++)
949 printf(" plane %i: z=%11.5f x0=%11.5f y0=%11.5f scale=%11.5f\n", ipl, fZ[ipl], fX0[ipl], fY0[ipl],
950 fScale[ipl]);
951 printf(" Bounding box:\n");
953}
954
955////////////////////////////////////////////////////////////////////////////////
956/// Creates a TBuffer3D describing *this* shape.
957/// Coordinates are in local reference frame.
958
960{
961 Int_t nz = GetNz();
962 Int_t nvert = GetNvert();
963 Int_t nbPnts = nz * nvert;
964 Int_t nbSegs = nvert * (2 * nz - 1);
965 Int_t nbPols = nvert * (nz - 1) + 2;
966
968 6 * (nbPols - 2) + 2 * (2 + nvert));
969 if (buff) {
970 SetPoints(buff->fPnts);
972 }
973
974 return buff;
975}
976
977////////////////////////////////////////////////////////////////////////////////
978/// Fill TBuffer3D structure for segments and polygons.
979
981{
982 Int_t nz = GetNz();
983 Int_t nvert = GetNvert();
985
986 Int_t i, j;
987 Int_t indx = 0, indx2, k;
988 for (i = 0; i < nz; i++) {
989 // loop Z planes
990 indx2 = i * nvert;
991 // loop polygon segments
992 for (j = 0; j < nvert; j++) {
993 k = (j + 1) % nvert;
994 buff.fSegs[indx++] = c;
995 buff.fSegs[indx++] = indx2 + j;
996 buff.fSegs[indx++] = indx2 + k;
997 }
998 } // total: nz*nvert polygon segments
999 for (i = 0; i < nz - 1; i++) {
1000 // loop Z planes
1001 indx2 = i * nvert;
1002 // loop polygon segments
1003 for (j = 0; j < nvert; j++) {
1004 k = j + nvert;
1005 buff.fSegs[indx++] = c;
1006 buff.fSegs[indx++] = indx2 + j;
1007 buff.fSegs[indx++] = indx2 + k;
1008 }
1009 } // total (nz-1)*nvert lateral segments
1010
1011 indx = 0;
1012
1013 // fill lateral polygons
1014 for (i = 0; i < nz - 1; i++) {
1015 indx2 = i * nvert;
1016 for (j = 0; j < nvert; j++) {
1017 k = (j + 1) % nvert;
1018 buff.fPols[indx++] = c + j % 3;
1019 buff.fPols[indx++] = 4;
1020 buff.fPols[indx++] = indx2 + j;
1021 buff.fPols[indx++] = nz * nvert + indx2 + k;
1022 buff.fPols[indx++] = indx2 + nvert + j;
1023 buff.fPols[indx++] = nz * nvert + indx2 + j;
1024 }
1025 } // total (nz-1)*nvert polys
1026 buff.fPols[indx++] = c + 2;
1027 buff.fPols[indx++] = nvert;
1028 indx2 = 0;
1029 for (j = nvert - 1; j >= 0; --j) {
1030 buff.fPols[indx++] = indx2 + j;
1031 }
1032
1033 buff.fPols[indx++] = c;
1034 buff.fPols[indx++] = nvert;
1035 indx2 = (nz - 1) * nvert;
1036
1037 for (j = 0; j < nvert; j++) {
1038 buff.fPols[indx++] = indx2 + j;
1039 }
1040}
1041
1042////////////////////////////////////////////////////////////////////////////////
1043/// Compute safety to sector iz, returning also the closest segment index.
1044
1046{
1049 Bool_t in1, in2;
1050 Int_t iseg;
1051 // segment-break case
1052 if (TGeoShape::IsSameWithinTolerance(fZ[iz], fZ[iz + 1])) {
1053 safz = TMath::Abs(point[2] - fZ[iz]);
1054 if (safz > safmin)
1055 return TGeoShape::Big();
1056 SetCurrentVertices(fX0[iz], fY0[iz], fScale[iz]);
1057 saf1 = td.fPoly->Safety(point, iseg);
1058 in1 = td.fPoly->Contains(point);
1059 // if (!in1 && saf1>safmin) return TGeoShape::Big();
1060 SetCurrentVertices(fX0[iz + 1], fY0[iz + 1], fScale[iz + 1]);
1061 saf2 = td.fPoly->Safety(point, iseg);
1062 in2 = td.fPoly->Contains(point);
1063 if ((in1 & !in2) | (in2 & !in1)) {
1064 safe = safz;
1065 } else {
1068 }
1069 if (safe > safmin)
1070 return TGeoShape::Big();
1071 return safe;
1072 }
1073 // normal case
1074 safz = fZ[iz] - point[2];
1075 if (safz > safmin)
1076 return TGeoShape::Big();
1077 if (safz < 0) {
1078 saf1 = point[2] - fZ[iz + 1];
1079 if (saf1 > safmin)
1080 return TGeoShape::Big();
1081 if (saf1 < 0) {
1082 safz = TMath::Max(safz, saf1); // we are in between the 2 Z segments - we ignore safz
1083 } else {
1084 safz = saf1;
1085 }
1086 }
1087
1088 // loop segments
1089 Bool_t found = kFALSE;
1090 Double_t vert[12];
1091 Double_t norm[3];
1092 // printf("plane %d: safz=%f in=%d\n", iz, safz, in);
1093 for (iseg = 0; iseg < fNvert; iseg++) {
1096 saf1 = (point[0] - vert[0]) * norm[0] + (point[1] - vert[1]) * norm[1] + (point[2] - vert[2]) * norm[2];
1097 if (in)
1098 saf1 = -saf1;
1099 // printf("segment %d: (%f,%f)-(%f,%f) norm=(%f,%f,%f): saf1=%f\n", iseg,
1100 // vert[0],vert[1],vert[3],vert[4],norm[0],norm[1],norm[2],saf1);
1101 if (saf1 < -1.E-8)
1102 continue;
1104 safe = TMath::Abs(safe);
1105 if (safe > safmin)
1106 continue;
1107 safmin = safe;
1108 found = kTRUE;
1109 }
1110 if (found)
1111 return safmin;
1112 return TGeoShape::Big();
1113}
1114
1115////////////////////////////////////////////////////////////////////////////////
1116/// computes the closest distance from given point to this shape, according
1117/// to option. The matching point on the shape is stored in spoint.
1118///> localize the Z segment
1119
1121{
1123 Double_t safe;
1125 TGeoXtru *xtru = (TGeoXtru *)this;
1126 Int_t iz;
1127 if (in) {
1128 safmin = TMath::Min(point[2] - fZ[0], fZ[fNz - 1] - point[2]);
1129 for (iz = 0; iz < fNz - 1; iz++) {
1130 safe = xtru->SafetyToSector(point, iz, safmin, in);
1131 if (safe < safmin)
1132 safmin = safe;
1133 }
1134 return safmin;
1135 }
1136 // Accurate safety is expensive, use the bounding box
1137 if (!TGeoBBox::Contains(point))
1138 return TGeoBBox::Safety(point, in);
1139 iz = TMath::BinarySearch(fNz, fZ, point[2]);
1140 if (iz < 0) {
1141 iz = 0;
1142 safz = fZ[0] - point[2];
1143 } else {
1144 if (iz == fNz - 1) {
1145 iz = fNz - 2;
1146 safz = point[2] - fZ[fNz - 1];
1147 }
1148 }
1149 // loop segments from iz up
1150 Int_t i;
1151 for (i = iz; i < fNz - 1; i++) {
1152 safe = xtru->SafetyToSector(point, i, safmin, in);
1153 if (safe < safmin)
1154 safmin = safe;
1155 }
1156 // loop segments from iz-1 down
1157 for (i = iz - 1; i >= 0; i--) {
1158 safe = xtru->SafetyToSector(point, i, safmin, in);
1159 if (safe < safmin)
1160 safmin = safe;
1161 }
1163 return safe;
1164}
1165
1166////////////////////////////////////////////////////////////////////////////////
1167/// Save a primitive as a C++ statement(s) on output stream "out".
1168
1169void TGeoXtru::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
1170{
1172 return;
1173 out << " // Shape: " << GetName() << " type: " << ClassName() << std::endl;
1174 out << " auto " << GetPointerName() << " = new TGeoXtru(" << fNz << ");" << std::endl;
1175 out << " " << GetPointerName() << "->SetName(\"" << GetName() << "\");" << std::endl;
1176 for (Int_t i = 0; i < fNvert; i++) {
1177 out << " xvert[" << i << "] = " << fX[i] << "; yvert[" << i << "] = " << fY[i] << ";" << std::endl;
1178 }
1179 out << " " << GetPointerName() << "->DefinePolygon(" << fNvert << ", xvert, yvert);" << std::endl;
1180 for (Int_t i = 0; i < fNz; i++)
1181 out << " " << GetPointerName() << "->DefineSection(" << i << ", " << fZ[i] << ", " << fX0[i] << ", " << fY0[i]
1182 << ", " << fScale[i] << ");" << std::endl;
1184}
1185
1186////////////////////////////////////////////////////////////////////////////////
1187/// Recompute current section vertices for a given Z position within range of section iz.
1188
1190{
1191 Double_t x0, y0, scale, a, b;
1192 Int_t ind1, ind2;
1193 ind1 = iz;
1194 ind2 = iz + 1;
1195 Double_t invdz = 1. / (fZ[ind2] - fZ[ind1]);
1196 a = (fX0[ind1] * fZ[ind2] - fX0[ind2] * fZ[ind1]) * invdz;
1197 b = (fX0[ind2] - fX0[ind1]) * invdz;
1198 x0 = a + b * z;
1199 a = (fY0[ind1] * fZ[ind2] - fY0[ind2] * fZ[ind1]) * invdz;
1200 b = (fY0[ind2] - fY0[ind1]) * invdz;
1201 y0 = a + b * z;
1202 a = (fScale[ind1] * fZ[ind2] - fScale[ind2] * fZ[ind1]) * invdz;
1203 b = (fScale[ind2] - fScale[ind1]) * invdz;
1204 scale = a + b * z;
1206}
1207
1208////////////////////////////////////////////////////////////////////////////////
1209/// Set current vertex coordinates according X0, Y0 and SCALE.
1210
1212{
1214 for (Int_t i = 0; i < fNvert; i++) {
1215 td.fXc[i] = scale * fX[i] + x0;
1216 td.fYc[i] = scale * fY[i] + y0;
1217 }
1218}
1219
1220////////////////////////////////////////////////////////////////////////////////
1221/// - param[0] = nz // number of z planes
1222///
1223/// - param[1] = z1 // Z position of first plane
1224/// - param[2] = x1 // X position of first plane
1225/// - param[3] = y1 // Y position of first plane
1226/// - param[4] = scale1 // scale factor for first plane
1227/// ...
1228/// - param[4*(nz-1]+1] = zn
1229/// - param[4*(nz-1)+2] = xn
1230/// - param[4*(nz-1)+3] = yn
1231/// - param[4*(nz-1)+4] = scalen
1232
1234{
1235 fNz = (Int_t)param[0];
1236 if (fNz < 2) {
1237 Error("SetDimensions", "Cannot create TGeoXtru %s with less than 2 Z planes", GetName());
1239 return;
1240 }
1241 if (fZ)
1242 delete[] fZ;
1243 if (fScale)
1244 delete[] fScale;
1245 if (fX0)
1246 delete[] fX0;
1247 if (fY0)
1248 delete[] fY0;
1249 fZ = new Double_t[fNz];
1250 fScale = new Double_t[fNz];
1251 fX0 = new Double_t[fNz];
1252 fY0 = new Double_t[fNz];
1253
1254 for (Int_t i = 0; i < fNz; i++)
1255 DefineSection(i, param[1 + 4 * i], param[2 + 4 * i], param[3 + 4 * i], param[4 + 4 * i]);
1256}
1257
1258////////////////////////////////////////////////////////////////////////////////
1259/// create polycone mesh points
1260
1262{
1264 Int_t i, j;
1265 Int_t indx = 0;
1266 TGeoXtru *xtru = (TGeoXtru *)this;
1267 if (points) {
1268 for (i = 0; i < fNz; i++) {
1269 xtru->SetCurrentVertices(fX0[i], fY0[i], fScale[i]);
1270 if (td.fPoly->IsClockwise()) {
1271 for (j = 0; j < fNvert; j++) {
1272 points[indx++] = td.fXc[j];
1273 points[indx++] = td.fYc[j];
1274 points[indx++] = fZ[i];
1275 }
1276 } else {
1277 for (j = 0; j < fNvert; j++) {
1278 points[indx++] = td.fXc[fNvert - 1 - j];
1279 points[indx++] = td.fYc[fNvert - 1 - j];
1280 points[indx++] = fZ[i];
1281 }
1282 }
1283 }
1284 }
1285}
1286
1287////////////////////////////////////////////////////////////////////////////////
1288/// create polycone mesh points
1289
1291{
1293 Int_t i, j;
1294 Int_t indx = 0;
1295 TGeoXtru *xtru = (TGeoXtru *)this;
1296 if (points) {
1297 for (i = 0; i < fNz; i++) {
1298 xtru->SetCurrentVertices(fX0[i], fY0[i], fScale[i]);
1299 if (td.fPoly->IsClockwise()) {
1300 for (j = 0; j < fNvert; j++) {
1301 points[indx++] = td.fXc[j];
1302 points[indx++] = td.fYc[j];
1303 points[indx++] = fZ[i];
1304 }
1305 } else {
1306 for (j = 0; j < fNvert; j++) {
1307 points[indx++] = td.fXc[fNvert - 1 - j];
1308 points[indx++] = td.fYc[fNvert - 1 - j];
1309 points[indx++] = fZ[i];
1310 }
1311 }
1312 }
1313 }
1314}
1315
1316////////////////////////////////////////////////////////////////////////////////
1317/// Returns numbers of vertices, segments and polygons composing the shape mesh.
1318
1320{
1321 Int_t nz = GetNz();
1322 Int_t nv = GetNvert();
1323 nvert = nz * nv;
1324 nsegs = nv * (2 * nz - 1);
1325 npols = nv * (nz - 1) + 2;
1326}
1327
1328////////////////////////////////////////////////////////////////////////////////
1329/// Return number of vertices of the mesh representation
1330
1332{
1333 Int_t numPoints = fNz * fNvert;
1334 return numPoints;
1335}
1336
1337////////////////////////////////////////////////////////////////////////////////
1338/// fill size of this 3-D object
1339
1340void TGeoXtru::Sizeof3D() const {}
1341
1342////////////////////////////////////////////////////////////////////////////////
1343/// Fills a static 3D buffer and returns a reference.
1344
1346{
1347 static TBuffer3D buffer(TBuffer3DTypes::kGeneric);
1348
1350
1352 Int_t nz = GetNz();
1353 Int_t nvert = GetNvert();
1354 Int_t nbPnts = nz * nvert;
1355 Int_t nbSegs = nvert * (2 * nz - 1);
1356 Int_t nbPols = nvert * (nz - 1) + 2;
1357 if (buffer.SetRawSizes(nbPnts, 3 * nbPnts, nbSegs, 3 * nbSegs, nbPols, 6 * (nbPols - 2) + 2 * (2 + nvert))) {
1359 }
1360 }
1361 // TODO: Push down to TGeoShape?
1363 SetPoints(buffer.fPnts);
1364 if (!buffer.fLocalFrame) {
1365 TransformPoints(buffer.fPnts, buffer.NbPnts());
1366 }
1367
1368 SetSegsAndPols(buffer);
1370 }
1371
1372 return buffer;
1373}
1374
1375////////////////////////////////////////////////////////////////////////////////
1376/// Check the inside status for each of the points in the array.
1377/// Input: Array of point coordinates + vector size
1378/// Output: Array of Booleans for the inside of each point
1379
1381{
1382 for (Int_t i = 0; i < vecsize; i++)
1383 inside[i] = Contains(&points[3 * i]);
1384}
1385
1386////////////////////////////////////////////////////////////////////////////////
1387/// Compute the normal for an array o points so that norm.dot.dir is positive
1388/// Input: Arrays of point coordinates and directions + vector size
1389/// Output: Array of normal directions
1390
1392{
1393 for (Int_t i = 0; i < vecsize; i++)
1394 ComputeNormal(&points[3 * i], &dirs[3 * i], &norms[3 * i]);
1395}
1396
1397////////////////////////////////////////////////////////////////////////////////
1398/// Compute distance from array of input points having directions specified by dirs. Store output in dists
1399
1401 Double_t *step) const
1402{
1403 for (Int_t i = 0; i < vecsize; i++)
1404 dists[i] = DistFromInside(&points[3 * i], &dirs[3 * i], 3, step[i]);
1405}
1406
1407////////////////////////////////////////////////////////////////////////////////
1408/// Compute distance from array of input points having directions specified by dirs. Store output in dists
1409
1411 Double_t *step) const
1412{
1413 for (Int_t i = 0; i < vecsize; i++)
1414 dists[i] = DistFromOutside(&points[3 * i], &dirs[3 * i], 3, step[i]);
1415}
1416
1417////////////////////////////////////////////////////////////////////////////////
1418/// Compute safe distance from each of the points in the input array.
1419/// Input: Array of point coordinates, array of statuses for these points, size of the arrays
1420/// Output: Safety values
1421
1423{
1424 for (Int_t i = 0; i < vecsize; i++)
1425 safe[i] = Safety(&points[3 * i], inside[i]);
1426}
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
int Int_t
Definition RtypesCore.h:45
float Float_t
Definition RtypesCore.h:57
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
double Double_t
Definition RtypesCore.h:59
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:374
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t points
float xmin
float ymin
float xmax
float ymax
Generic 3D primitive description class.
Definition TBuffer3D.h:18
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
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
Box class.
Definition TGeoBBox.h:17
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:20
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
Compute distance from outside point to surface of the box.
Definition TGeoBBox.cxx:433
Double_t fOrigin[3]
Definition TGeoBBox.h:23
void InspectShape() const override
Prints shape parameters.
Definition TGeoBBox.cxx:811
Bool_t Contains(const Double_t *point) const override
Test if point is inside this shape.
Definition TGeoBBox.cxx:324
Double_t Safety(const Double_t *point, Bool_t in=kTRUE) const override
Computes the closest distance from given point to this shape.
Definition TGeoBBox.cxx:920
Double_t fDY
Definition TGeoBBox.h:21
Double_t fDZ
Definition TGeoBBox.h:22
static Int_t ThreadId()
Translates the current thread id to an ordinal number.
An arbitrary polygon defined by vertices.
Definition TGeoPolygon.h:19
static Double_t Big()
Definition TGeoShape.h:94
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)
void SetShapeBit(UInt_t f, Bool_t set)
Equivalent of TObject::SetBit.
static Bool_t IsSameWithinTolerance(Double_t a, Double_t b)
Check if two numbers differ with less than a tolerance.
const char * GetPointerName() const
Provide a pointer name containing uid.
Int_t ShapeDistancetoPrimitive(Int_t numpoints, Int_t px, Int_t py) const
Returns distance to shape primitive mesh.
const char * GetName() const override
Get the shape name.
@ kGeoSavePrimitive
Definition TGeoShape.h:64
static Double_t Tolerance()
Definition TGeoShape.h:97
Bool_t TestShapeBit(UInt_t f) const
Definition TGeoShape.h:175
A TGeoXtru shape is represented by the extrusion of an arbitrary polygon with fixed outline between s...
Definition TGeoXtru.h:22
Double_t * fScale
Definition TGeoXtru.h:46
Int_t GetNmeshVertices() const override
Return number of vertices of the mesh representation.
Double_t * fX0
Definition TGeoXtru.h:47
void DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const override
Compute distance from array of input points having directions specified by dirs. Store output in dist...
TGeoXtru()
dummy ctor
Definition TGeoXtru.cxx:190
ThreadData_t & GetThreadData() const
Definition TGeoXtru.cxx:124
void ClearThreadData() const override
Definition TGeoXtru.cxx:134
Double_t Capacity() const override
Compute capacity [length^3] of this shape.
Definition TGeoXtru.cxx:299
Int_t fNz
Definition TGeoXtru.h:41
Double_t * fZ
Definition TGeoXtru.h:45
void ComputeBBox() override
compute bounding box of the pcon
Definition TGeoXtru.cxx:322
void Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const override
Check the inside status for each of the points in the array.
Bool_t Contains(const Double_t *point) const override
test if point is inside this shape
Definition TGeoXtru.cxx:382
Int_t fThreadSize
Navigation data per thread.
Definition TGeoXtru.h:51
Int_t fNvert
Definition TGeoXtru.h:40
void SetCurrentVertices(Double_t x0, Double_t y0, Double_t scale)
Set current vertex coordinates according X0, Y0 and SCALE.
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save a primitive as a C++ statement(s) on output stream "out".
void GetPlaneVertices(Int_t iz, Int_t ivert, Double_t *vert) const
Returns (x,y,z) of 3 vertices of the surface defined by Z sections (iz, iz+1) and polygon vertices (i...
Definition TGeoXtru.cxx:861
void InspectShape() const override
Print actual Xtru parameters.
Definition TGeoXtru.cxx:941
Int_t DistancetoPrimitive(Int_t px, Int_t py) override
compute closest distance from point px,py to each corner
Definition TGeoXtru.cxx:418
Double_t * GetZ() const
Definition TGeoXtru.h:102
Bool_t IsPointInsidePlane(const Double_t *point, Double_t *vert, Double_t *norm) const
Check if the quadrilateral defined by VERT contains a coplanar POINT.
Definition TGeoXtru.cxx:916
void Sizeof3D() const override
fill size of this 3-D object
~TGeoXtru() override
destructor
Definition TGeoXtru.cxx:267
const TBuffer3D & GetBuffer3D(Int_t reqSections, Bool_t localFrame) const override
Fills a static 3D buffer and returns a reference.
Int_t GetNvert() const
Definition TGeoXtru.h:96
Double_t * fY0
Definition TGeoXtru.h:48
Bool_t DefinePolygon(Int_t nvert, const Double_t *xv, const Double_t *yv)
Creates the polygon representing the blueprint of any Xtru section.
Definition TGeoXtru.cxx:761
void SetCurrentZ(Double_t z, Int_t iz)
Recompute current section vertices for a given Z position within range of section iz.
virtual void DefineSection(Int_t snum, Double_t z, Double_t x0=0., Double_t y0=0., Double_t scale=1.)
defines z position of a section plane, rmin and rmax at this z.
Definition TGeoXtru.cxx:795
Double_t Safety(const Double_t *point, Bool_t in=kTRUE) const override
computes the closest distance from given point to this shape, according to option.
void ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize) override
Compute the normal for an array o points so that norm.dot.dir is positive Input: Arrays of point coor...
void SetDimensions(Double_t *param) override
void DrawPolygon(Option_t *option="")
Draw the section polygon.
Definition TGeoXtru.cxx:427
void SetPoints(Double_t *points) const override
create polycone mesh points
Double_t SafetyToSector(const Double_t *point, Int_t iz, Double_t safmin, Bool_t in)
Compute safety to sector iz, returning also the closest segment index.
Double_t * fY
Definition TGeoXtru.h:44
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
compute distance from inside point to surface of the polycone locate Z segment
Definition TGeoXtru.cxx:503
Double_t * fX
Definition TGeoXtru.h:43
void DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const override
Compute distance from array of input points having directions specified by dirs. Store output in dist...
void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm) const override
Compute normal to closest surface from POINT.
Definition TGeoXtru.cxx:360
TBuffer3D * MakeBuffer3D() const override
Creates a TBuffer3D describing this shape.
Definition TGeoXtru.cxx:959
std::vector< ThreadData_t * > fThreadData
Definition TGeoXtru.h:50
Double_t DistToPlane(const Double_t *point, const Double_t *dir, Int_t iz, Int_t ivert, Double_t stepmax, Bool_t in) const
Compute distance to a Xtru lateral surface.
Definition TGeoXtru.cxx:437
void SetSegsAndPols(TBuffer3D &buff) const override
Fill TBuffer3D structure for segments and polygons.
Definition TGeoXtru.cxx:980
void SetSeg(Int_t iseg)
Set current segment.
Definition TGeoXtru.cxx:182
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
compute distance from outside point to surface of the tube Warning("DistFromOutside",...
Definition TGeoXtru.cxx:628
std::mutex fMutex
size of thread-specific array
Definition TGeoXtru.h:52
void GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const override
Returns numbers of vertices, segments and polygons composing the shape mesh.
void GetPlaneNormal(const Double_t *vert, Double_t *norm) const
Returns normal vector to the planar quadrilateral defined by vector VERT.
Definition TGeoXtru.cxx:834
void SetIz(Int_t iz)
Set current z-plane.
Definition TGeoXtru.cxx:175
Int_t GetNz() const
Definition TGeoXtru.h:95
void Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const override
Compute safe distance from each of the points in the input array.
void CreateThreadData(Int_t nthreads) override
Create thread data for n threads max.
Definition TGeoXtru.cxx:149
Double_t fZcurrent
Definition TGeoXtru.h:42
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:201
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:225
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:1006
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:813
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1020
TPaveText * pt
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:250
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:666
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:198
Long64_t BinarySearch(Long64_t n, const T *array, T value)
Binary search in an array of n values to locate value.
Definition TMathBase.h:347
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123
~ThreadData_t()
Destructor.
Definition TGeoXtru.cxx:115
ThreadData_t()
Constructor.
Definition TGeoXtru.cxx:110