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