ROOT  6.06/09
Reference Guide
TGeoPolygon.cxx
Go to the documentation of this file.
1 // @(#)root/geom:$Id$
2 // Author: Mihaela Gheata 5/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 //____________________________________________________________________________
13 // TGeoPolygon - Arbitrary polygon class
14 //____________________________________________________________________________
15 //
16 // A polygon is a 2D shape defined by vertices in the XY plane. It is used by
17 // TGeoXtru class for computing Contains() and Safety(). Only the pointers to
18 // the actual lists of XY values are used - these are not owned by the class.
19 //
20 // To check if a point in XY plane is contained by a polygon, this is split
21 // into an outscribed convex polygon and the remaining polygons of its subtracton
22 // from the outscribed one. A point is INSIDE if it is
23 // contained by the outscribed polygon but NOT by the remaining ones. Since these
24 // can also be arbitrary polygons at their turn, a tree structure is formed:
25 //
26 // P = Pconvex - (Pconvex-P) where (-) means 'subtraction'
27 // Pconvex-P = P1 + P2 + ... where (+) means 'union'
28 //
29 // *Note that P1, P2, ... do not intersect each other and they are defined
30 // by subsets of the list of vertices of P. They can be split in the same
31 // way as P*
32 //
33 // Therefore, if C(P) represents the Boolean : 'does P contains a given point?',
34 // then:
35 //
36 // C(P) = C(Pconvex) .and. not(C(P1) | C(P2) | ...)
37 //
38 // For creating a polygon without TGeoXtru class, one has to call the constructor
39 // TGeoPolygon(nvert) and then SetXY(Double_t *x, Double_t *y) providing the
40 // arrays of X and Y vertex positions (defined clockwise) that have to 'live' longer
41 // than the polygon they will describe. This complication is due to efficiency reasons.
42 // At the end one has to call the FinishPolygon() method.
43 
44 #include "TObjArray.h"
45 #include "TGeoPolygon.h"
46 #include "TMath.h"
47 #include "TGeoShape.h"
48 #include "TGeoManager.h"
49 #include "TVirtualGeoPainter.h"
50 
52 
53 ////////////////////////////////////////////////////////////////////////////////
54 /// Dummy constructor.
55 
57 {
58  fNvert = 0;
59  fNconvex = 0;
60  fInd = 0;
61  fIndc = 0;
62  fX = 0;
63  fY = 0;
64  fDaughters = 0;
65  SetConvex(kFALSE);
66  TObject::SetBit(kGeoFinishPolygon, kFALSE);
67 }
68 
69 ////////////////////////////////////////////////////////////////////////////////
70 /// Default constructor.
71 
73 {
74  if (nvert<3) {
75  Fatal("Ctor", "Invalid number of vertices %i", nvert);
76  return;
77  }
78  fNvert = nvert;
79  fNconvex = 0;
80  fInd = new Int_t[nvert];
81  fIndc = 0;
82  fX = 0;
83  fY = 0;
84  fDaughters = 0;
87  SetNextIndex();
88 }
89 
90 ////////////////////////////////////////////////////////////////////////////////
91 /// Destructor
92 
94 {
95  if (fInd) delete [] fInd;
96  if (fIndc) delete [] fIndc;
97  if (fDaughters) {
98  fDaughters->Delete();
99  delete fDaughters;
100  }
101 }
102 ////////////////////////////////////////////////////////////////////////////////
103 /// Computes area of the polygon in [length^2].
104 
106 {
107  Int_t ic,i,j;
108  Double_t area = 0;
109  // Compute area of the convex part
110  for (ic=0; ic<fNvert; ic++) {
111  i = fInd[ic];
112  j = fInd[(ic+1)%fNvert];
113  area += 0.5*(fX[i]*fY[j]-fX[j]*fY[i]);
114  }
115  return TMath::Abs(area);
116 }
117 
118 ////////////////////////////////////////////////////////////////////////////////
119 /// Check if a point given by X = point[0], Y = point[1] is inside the polygon.
120 
122 {
123  Int_t i;
124  TGeoPolygon *poly;
125  for (i=0; i<fNconvex; i++)
126  if (!IsRightSided(point, fIndc[i], fIndc[(i+1)%fNconvex])) return kFALSE;
127  if (!fDaughters) return kTRUE;
129  for (i=0; i<nd; i++) {
130  poly = (TGeoPolygon*)fDaughters->UncheckedAt(i);
131  if (poly->Contains(point)) return kFALSE;
132  }
133  return kTRUE;
134 }
135 
136 ////////////////////////////////////////////////////////////////////////////////
137 /// Check polygon convexity.
138 
140 {
141  if (fNvert==3) {
142  SetConvex();
143  return;
144  }
145  Int_t j,k;
146  Double_t point[2];
147  for (Int_t i=0; i<fNvert; i++) {
148  j = (i+1)%fNvert;
149  k = (i+2)%fNvert;
150  point[0] = fX[fInd[k]];
151  point[1] = fY[fInd[k]];
152  if (!IsRightSided(point, fInd[i], fInd[j])) return;
153  }
154  SetConvex();
155 }
156 
157 ////////////////////////////////////////////////////////////////////////////////
158 /// Draw the polygon.
159 
161 {
162  if (!gGeoManager) return;
164 }
165 
166 ////////////////////////////////////////////////////////////////////////////////
167 /// Decompose polygon in a convex outscribed part and a list of daughter
168 /// polygons that have to be substracted to get the actual one.
169 
171 {
173  // check convexity
174  ConvexCheck();
175  // find outscribed convex polygon indices
177  if (IsConvex()) {
178 // printf(" -> polygon convex -> same indices\n");
179  memcpy(fIndc, fInd, fNvert*sizeof(Int_t));
180  return;
181  }
182 // printf(" -> polygon NOT convex\n");
183  // make daughters if necessary
184  if (IsConvex()) return;
185  // ... algorithm here
186  if (!fDaughters) fDaughters = new TObjArray();
187  TGeoPolygon *poly = 0;
188  Int_t indconv = 0;
189  Int_t indnext, indback;
190  Int_t nskip;
191  while (indconv < fNconvex) {
192  indnext = (indconv+1)%fNconvex;
193  nskip = fIndc[indnext]-fIndc[indconv];
194  if (nskip<0) nskip+=fNvert;
195  if (nskip==1) {
196  indconv++;
197  continue;
198  }
199  // gap -> make polygon
200  poly = new TGeoPolygon(nskip+1);
201  poly->SetXY(fX,fY);
202  poly->SetNextIndex(fInd[fIndc[indconv]]);
203  poly->SetNextIndex(fInd[fIndc[indnext]]);
204  indback = fIndc[indnext]-1;
205  if (indback < 0) indback+=fNvert;
206  while (indback != fIndc[indconv]) {
207  poly->SetNextIndex(fInd[indback]);
208  indback--;
209  if (indback < 0) indback+=fNvert;
210  }
211  poly->FinishPolygon();
212  fDaughters->Add(poly);
213  indconv++;
214  }
215  for (indconv=0; indconv<fNconvex; indconv++) fIndc[indconv] = fInd[fIndc[indconv]];
216 }
217 
218 ////////////////////////////////////////////////////////////////////////////////
219 /// Fill list of vertices into provided arrays.
220 
222 {
223  memcpy(x, fX, fNvert*sizeof(Double_t));
224  memcpy(y, fY, fNvert*sizeof(Double_t));
225 }
226 
227 ////////////////////////////////////////////////////////////////////////////////
228 /// Fill list of vertices of the convex outscribed polygon into provided arrays.
229 
231 {
232  for (Int_t ic=0; ic<fNconvex; ic++) {
233  x[ic] = fX[fIndc[ic]];
234  y[ic] = fY[fIndc[ic]];
235  }
236 }
237 
238 ////////////////////////////////////////////////////////////////////////////////
239 /// Check if POINT is right-sided with respect to the segment defined by IND1 and IND2.
240 
241 Bool_t TGeoPolygon::IsRightSided(const Double_t *point, Int_t ind1, Int_t ind2) const
242 {
243  Double_t dot = (point[0]-fX[ind1])*(fY[ind2]-fY[ind1]) -
244  (point[1]-fY[ind1])*(fX[ind2]-fX[ind1]);
245  if (!IsClockwise()) dot = -dot;
246  if (dot<-1.E-10) return kFALSE;
247  return kTRUE;
248 }
249 
250 ////////////////////////////////////////////////////////////////////////////////
251 /// Check if a segment [0..fNvert-1] belongs to the outscribed convex pgon.
252 
254 {
255  if (i2<0) i2=(i1+1)%fNvert;
256  Double_t point[2];
257  for (Int_t i=0; i<fNvert; i++) {
258  if (i==i1 || i==i2) continue;
259  point[0] = fX[fInd[i]];
260  point[1] = fY[fInd[i]];
261  if (!IsRightSided(point, fInd[i1], fInd[i2])) return kFALSE;
262  }
263  return kTRUE;
264 }
265 
266 ////////////////////////////////////////////////////////////////////////////////
267 /// Check for illegal crossings between non-consecutive segments
268 
270 {
271  if (fNvert<4) return kFALSE;
272  Bool_t is_illegal = kFALSE;
273  Double_t x1,y1,x2,y2,x3,y3,x4,y4;
274  for (Int_t i=0; i<fNvert-2; i++) {
275  // Check segment i
276  for (Int_t j=i+2; j<fNvert; j++) {
277  // Versus segment j
278  if (i==0 && j==(fNvert-1)) continue;
279  x1 = fX[i];
280  y1 = fY[i];
281  x2 = fX[i+1];
282  y2 = fY[i+1];
283  x3 = fX[j];
284  y3 = fY[j];
285  x4 = fX[(j+1)%fNvert];
286  y4 = fY[(j+1)%fNvert];
287  if (TGeoShape::IsSegCrossing(x1,y1,x2,y2,x3,y3,x4,y4)) {
288  Error("IsIllegalCheck", "Illegal crossing of segment %d vs. segment %d", i,j);
289  is_illegal = kTRUE;
290  }
291  }
292  }
293  return is_illegal;
294 }
295 
296 ////////////////////////////////////////////////////////////////////////////////
297 /// Compute indices for the outscribed convex polygon.
298 
300 {
301  fNconvex = 0;
302  Int_t iseg = 0;
303  Int_t ivnew;
304  Bool_t conv;
305  Int_t *indconv = new Int_t[fNvert];
306  memset(indconv, 0, fNvert*sizeof(Int_t));
307  while (iseg<fNvert) {
308  if (!IsSegConvex(iseg)) {
309  if (iseg+2 > fNvert) break;
310  ivnew = (iseg+2)%fNvert;
311  conv = kFALSE;
312  // check iseg with next vertices
313  while (ivnew) {
314  if (IsSegConvex(iseg, ivnew)) {
315  conv = kTRUE;
316  break;
317  }
318  ivnew = (ivnew+1)%fNvert;
319  }
320  if (!conv) {
321 // Error("OutscribedConvex","NO convex line connection to vertex %d\n", iseg);
322  iseg++;
323  continue;
324  }
325  } else {
326  ivnew = (iseg+1)%fNvert;
327  }
328  // segment belonging to convex outscribed poligon
329  if (!fNconvex) indconv[fNconvex++] = iseg;
330  else if (indconv[fNconvex-1] != iseg) indconv[fNconvex++] = iseg;
331  if (iseg<fNvert-1) indconv[fNconvex++] = ivnew;
332  if (ivnew<iseg) break;
333  iseg = ivnew;
334  }
335  if (!fNconvex) {
336  delete [] indconv;
337  Fatal("OutscribedConvex","cannot build outscribed convex");
338  return;
339  }
340  fIndc = new Int_t[fNvert];
341  memcpy(fIndc, indconv, fNconvex*sizeof(Int_t)); // does not contain real indices yet
342  delete [] indconv;
343 }
344 
345 ////////////////////////////////////////////////////////////////////////////////
346 /// Compute minimum distance from POINT to any segment. Returns segment index.
347 
348 Double_t TGeoPolygon::Safety(const Double_t *point, Int_t &isegment) const
349 {
350  Int_t i1, i2;
351  Double_t p1[2], p2[3];
352  Double_t lsq, ssq, dx, dy, dpx, dpy, u;
353  Double_t safe=1E30;
354  Int_t isegmin=0;
355  for (i1=0; i1<fNvert; i1++) {
356  if (TGeoShape::IsSameWithinTolerance(safe,0)) {
357  isegment = isegmin;
358  return 0.;
359  }
360  i2 = (i1+1)%fNvert;
361  p1[0] = fX[i1];
362  p1[1] = fY[i1];
363  p2[0] = fX[i2];
364  p2[1] = fY[i2];
365 
366  dx = p2[0] - p1[0];
367  dy = p2[1] - p1[1];
368  dpx = point[0] - p1[0];
369  dpy = point[1] - p1[1];
370 
371  lsq = dx*dx + dy*dy;
373  ssq = dpx*dpx + dpy*dpy;
374  if (ssq < safe) {
375  safe = ssq;
376  isegmin = i1;
377  }
378  continue;
379  }
380  u = (dpx*dx + dpy*dy)/lsq;
381  if (u>1) {
382  dpx = point[0]-p2[0];
383  dpy = point[1]-p2[1];
384  } else {
385  if (u>=0) {
386  dpx -= u*dx;
387  dpy -= u*dy;
388  }
389  }
390  ssq = dpx*dpx + dpy*dpy;
391  if (ssq < safe) {
392  safe = ssq;
393  isegmin = i1;
394  }
395  }
396  isegment = isegmin;
397  safe = TMath::Sqrt(safe);
398 // printf("== segment %d: (%f, %f) - (%f, %f) safe=%f\n", isegment, fX[isegment],fY[isegment],fX[(isegment+1)%fNvert],fY[(isegment+1)%fNvert],safe);
399  return safe;
400 }
401 
402 ////////////////////////////////////////////////////////////////////////////////
403 /// Sets the next polygone index. If index<0 sets all indices consecutive
404 /// in increasing order.
405 
407 {
408  Int_t i;
409  if (index <0) {
410  for (i=0; i<fNvert; i++) fInd[i] = i;
411  return;
412  }
413  if (fNconvex >= fNvert) {
414  Error("SetNextIndex", "all indices already set");
415  return;
416  }
417  fInd[fNconvex++] = index;
418  if (fNconvex == fNvert) {
419  if (!fX || !fY) return;
420  Double_t area = 0.0;
421  for (i=0; i<fNvert; i++) area += fX[fInd[i]]*fY[fInd[(i+1)%fNvert]]-fX[fInd[(i+1)%fNvert]]*fY[fInd[i]];
422  if (area<0) TObject::SetBit(kGeoACW, kFALSE);
424  }
425 }
426 
427 ////////////////////////////////////////////////////////////////////////////////
428 /// Set X/Y array pointer for the polygon and daughters.
429 
431 {
432  Int_t i;
433  fX = x;
434  fY = y;
435  Double_t area = 0.0;
436  for (i=0; i<fNvert; i++) area += fX[fInd[i]]*fY[fInd[(i+1)%fNvert]]-fX[fInd[(i+1)%fNvert]]*fY[fInd[i]];
437  if (area<0) TObject::SetBit(kGeoACW, kFALSE);
439 
440  if (!fDaughters) return;
441  TGeoPolygon *poly;
443  for (i=0; i<nd; i++) {
444  poly = (TGeoPolygon*)fDaughters->At(i);
445  if (poly) poly->SetXY(x,y);
446  }
447 }
Double_t Area() const
Computes area of the polygon in [length^2].
Bool_t IsRightSided(const Double_t *point, Int_t ind1, Int_t ind2) const
Check if POINT is right-sided with respect to the segment defined by IND1 and IND2.
virtual void Draw(Option_t *option="")
Draw the polygon.
An array of TObjects.
Definition: TObjArray.h:39
ClassImp(TGeoPolygon) TGeoPolygon
Dummy constructor.
Definition: TGeoPolygon.cxx:51
void SetNextIndex(Int_t index=-1)
Sets the next polygone index.
const char Option_t
Definition: RtypesCore.h:62
virtual void Delete(Option_t *option="")
Remove all objects from the array AND delete all heap based objects.
Definition: TObjArray.cxx:328
void ConvexCheck()
Check polygon convexity.
void GetVertices(Double_t *x, Double_t *y) const
Fill list of vertices into provided arrays.
Int_t * fIndc
Definition: TGeoPolygon.h:44
Bool_t IsSegConvex(Int_t i1, Int_t i2=-1) const
Check if a segment [0..fNvert-1] belongs to the outscribed convex pgon.
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
Int_t GetEntriesFast() const
Definition: TObjArray.h:66
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:325
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:732
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition: TObject.cxx:946
Bool_t IsIllegalCheck() const
Check for illegal crossings between non-consecutive segments.
static const double x2[5]
Double_t dot(const TVector2 &v1, const TVector2 &v2)
Definition: CsgOps.cxx:333
Double_t x[n]
Definition: legend1.C:17
static double p2(double t, double a, double b, double c)
static const double x4[22]
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:918
void GetConvexVertices(Double_t *x, Double_t *y) const
Fill list of vertices of the convex outscribed polygon into provided arrays.
TObject * UncheckedAt(Int_t i) const
Definition: TObjArray.h:91
static Bool_t IsSegCrossing(Double_t x1, Double_t y1, Double_t x2, Double_t y2, Double_t x3, Double_t y3, Double_t x4, Double_t y4)
Check if segments (A,B) and (C,D) are crossing, where: A(x1,y1), B(x2,y2), C(x3,y3), D(x4,y4)
Definition: TGeoShape.cxx:335
void SetXY(Double_t *x, Double_t *y)
Set X/Y array pointer for the polygon and daughters.
Double_t E()
Definition: TMath.h:54
static double p1(double t, double a, double b)
Double_t * fX
Definition: TGeoPolygon.h:45
void SetConvex(Bool_t flag=kTRUE)
Definition: TGeoPolygon.h:75
void OutscribedConvex()
Compute indices for the outscribed convex polygon.
static const double x1[5]
Bool_t IsConvex() const
Definition: TGeoPolygon.h:71
R__EXTERN TGeoManager * gGeoManager
Definition: TGeoManager.h:556
Double_t Safety(const Double_t *point, Int_t &isegment) const
Compute minimum distance from POINT to any segment. Returns segment index.
double Double_t
Definition: RtypesCore.h:55
Double_t * fY
pointer to list of current X coordinates of vertices
Definition: TGeoPolygon.h:46
void FinishPolygon()
Decompose polygon in a convex outscribed part and a list of daughter polygons that have to be substra...
Int_t * fInd
Definition: TGeoPolygon.h:43
Double_t y[n]
Definition: legend1.C:17
int nskip
Definition: mixmax.cxx:39
Bool_t IsClockwise() const
Definition: TGeoPolygon.h:70
virtual void DrawPolygon(const TGeoPolygon *poly)=0
TObjArray * fDaughters
pointer to list of current Y coordinates of vertices
Definition: TGeoPolygon.h:47
Int_t fNconvex
Definition: TGeoPolygon.h:42
Bool_t Contains(const Double_t *point) const
Check if a point given by X = point[0], Y = point[1] is inside the polygon.
void Add(TObject *obj)
Definition: TObjArray.h:75
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
TObject * At(Int_t idx) const
Definition: TObjArray.h:167
Int_t fNvert
Definition: TGeoPolygon.h:41
const Bool_t kTRUE
Definition: Rtypes.h:91
virtual ~TGeoPolygon()
Destructor.
Definition: TGeoPolygon.cxx:93
TVirtualGeoPainter * GetGeomPainter()
Make a default painter if none present. Returns pointer to it.
static const double x3[11]