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