Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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/** \class TGeoPolygon
13\ingroup Shapes_classes
14
15An arbitrary polygon defined by vertices. The vertices
16have to be defined CLOCKWISE in the XY plane, making either a convex
17or concave polygon. No test for malformed polygons is performed.
18
19A polygon is a 2D shape defined by vertices in the XY plane. It is used by
20TGeoXtru class for computing Contains() and Safety(). Only the pointers to
21the actual lists of XY values are used - these are not owned by the class.
22
23To check if a point in XY plane is contained by a polygon, this is split
24into an outscribed convex polygon and the remaining polygons of its subtraction
25from the outscribed one. A point is INSIDE if it is
26contained by the outscribed polygon but NOT by the remaining ones. Since these
27can also be arbitrary polygons at their turn, a tree structure is formed:
28
29~~~ {.cpp}
30 P = Pconvex - (Pconvex-P) where (-) means 'subtraction'
31 Pconvex-P = P1 + P2 + ... where (+) means 'union'
32~~~
33
34*Note that P1, P2, ... do not intersect each other and they are defined
35by subsets of the list of vertices of P. They can be split in the same
36way as P*
37
38Therefore, if C(P) represents the Boolean : 'does P contains a given point?',
39then:
40
41~~~ {.cpp}
42C(P) = C(Pconvex) .and. not(C(P1) | C(P2) | ...)
43~~~
44
45For creating a polygon without TGeoXtru class, one has to call the constructor
46TGeoPolygon(nvert) and then SetXY(Double_t *x, Double_t *y) providing the
47arrays of X and Y vertex positions (defined clockwise) that have to 'live' longer
48than the polygon they will describe. This complication is due to efficiency reasons.
49At the end one has to call the FinishPolygon() method.
50*/
51
52#include "TGeoPolygon.h"
53
54#include "TObjArray.h"
55#include "TMath.h"
56#include "TGeoShape.h"
57#include "TGeoManager.h"
58#include "TVirtualGeoPainter.h"
59
60
61////////////////////////////////////////////////////////////////////////////////
62/// Dummy constructor.
63
65{
66 fNvert = 0;
67 fNconvex = 0;
68 fInd = nullptr;
69 fIndc = nullptr;
70 fX = nullptr;
71 fY = nullptr;
72 fDaughters = nullptr;
75}
76
77////////////////////////////////////////////////////////////////////////////////
78/// Default constructor.
79
81{
82 if (nvert < 3) {
83 Fatal("Ctor", "Invalid number of vertices %i", nvert);
84 return;
85 }
86 fNvert = nvert;
87 fNconvex = 0;
88 fInd = new Int_t[nvert];
89 fIndc = nullptr;
90 fX = nullptr;
91 fY = nullptr;
92 fDaughters = nullptr;
96}
97
98////////////////////////////////////////////////////////////////////////////////
99/// Destructor
100
102{
103 if (fInd)
104 delete[] fInd;
105 if (fIndc)
106 delete[] fIndc;
107 if (fDaughters) {
109 delete fDaughters;
110 }
111}
112////////////////////////////////////////////////////////////////////////////////
113/// Computes area of the polygon in [length^2].
114
116{
117 Int_t ic, i, j;
118 Double_t area = 0;
119 // Compute area of the convex part
120 for (ic = 0; ic < fNvert; ic++) {
121 i = fInd[ic];
122 j = fInd[(ic + 1) % fNvert];
123 area += 0.5 * (fX[i] * fY[j] - fX[j] * fY[i]);
124 }
125 return TMath::Abs(area);
126}
127
128////////////////////////////////////////////////////////////////////////////////
129/// Check if a point given by X = point[0], Y = point[1] is inside the polygon.
130
132{
133 Int_t i;
135 for (i = 0; i < fNconvex; i++)
136 if (!IsRightSided(point, fIndc[i], fIndc[(i + 1) % fNconvex]))
137 return kFALSE;
138 if (!fDaughters)
139 return kTRUE;
141 for (i = 0; i < nd; i++) {
143 if (poly->Contains(point))
144 return kFALSE;
145 }
146 return kTRUE;
147}
148
149////////////////////////////////////////////////////////////////////////////////
150/// Check polygon convexity.
151
153{
154 if (fNvert == 3) {
155 SetConvex();
156 return;
157 }
158 Int_t j, k;
159 Double_t point[2];
160 for (Int_t i = 0; i < fNvert; i++) {
161 j = (i + 1) % fNvert;
162 k = (i + 2) % fNvert;
163 point[0] = fX[fInd[k]];
164 point[1] = fY[fInd[k]];
165 if (!IsRightSided(point, fInd[i], fInd[j]))
166 return;
167 }
168 SetConvex();
169}
170
171////////////////////////////////////////////////////////////////////////////////
172/// Draw the polygon.
173
175{
176 if (!gGeoManager)
177 return;
179}
180
181////////////////////////////////////////////////////////////////////////////////
182/// Decompose polygon in a convex outscribed part and a list of daughter
183/// polygons that have to be subtracted to get the actual one.
184
186{
188 // check convexity
189 ConvexCheck();
190 // find outscribed convex polygon indices
192 if (IsConvex()) {
193 // printf(" -> polygon convex -> same indices\n");
194 memcpy(fIndc, fInd, fNvert * sizeof(Int_t));
195 return;
196 }
197 // printf(" -> polygon NOT convex\n");
198 // make daughters if necessary
199 if (IsConvex())
200 return;
201 // ... algorithm here
202 if (!fDaughters)
203 fDaughters = new TObjArray();
204 TGeoPolygon *poly = nullptr;
205 Int_t indconv = 0;
207 Int_t nskip;
208 while (indconv < fNconvex) {
209 indnext = (indconv + 1) % fNconvex;
211 if (nskip < 0)
212 nskip += fNvert;
213 if (nskip == 1) {
214 indconv++;
215 continue;
216 }
217 // gap -> make polygon
218 poly = new TGeoPolygon(nskip + 1);
219 poly->SetXY(fX, fY);
220 poly->SetNextIndex(fInd[fIndc[indconv]]);
221 poly->SetNextIndex(fInd[fIndc[indnext]]);
222 indback = fIndc[indnext] - 1;
223 if (indback < 0)
224 indback += fNvert;
225 while (indback != fIndc[indconv]) {
226 poly->SetNextIndex(fInd[indback]);
227 indback--;
228 if (indback < 0)
229 indback += fNvert;
230 }
231 poly->FinishPolygon();
233 indconv++;
234 }
235 for (indconv = 0; indconv < fNconvex; indconv++)
237}
238
239////////////////////////////////////////////////////////////////////////////////
240/// Fill list of vertices into provided arrays.
241
243{
244 memcpy(x, fX, fNvert * sizeof(Double_t));
245 memcpy(y, fY, fNvert * sizeof(Double_t));
246}
247
248////////////////////////////////////////////////////////////////////////////////
249/// Fill list of vertices of the convex outscribed polygon into provided arrays.
250
252{
253 for (Int_t ic = 0; ic < fNconvex; ic++) {
254 x[ic] = fX[fIndc[ic]];
255 y[ic] = fY[fIndc[ic]];
256 }
257}
258
259////////////////////////////////////////////////////////////////////////////////
260/// Check if POINT is right-sided with respect to the segment defined by IND1 and IND2.
261
263{
264 Double_t dot = (point[0] - fX[ind1]) * (fY[ind2] - fY[ind1]) - (point[1] - fY[ind1]) * (fX[ind2] - fX[ind1]);
265 if (!IsClockwise())
266 dot = -dot;
267 if (dot < -1.E-10)
268 return kFALSE;
269 return kTRUE;
270}
271
272////////////////////////////////////////////////////////////////////////////////
273/// Check if a segment [0..fNvert-1] belongs to the outscribed convex pgon.
274
276{
277 if (i2 < 0)
278 i2 = (i1 + 1) % fNvert;
279 Double_t point[2];
280 for (Int_t i = 0; i < fNvert; i++) {
281 if (i == i1 || i == i2)
282 continue;
283 point[0] = fX[fInd[i]];
284 point[1] = fY[fInd[i]];
285 if (!IsRightSided(point, fInd[i1], fInd[i2]))
286 return kFALSE;
287 }
288 return kTRUE;
289}
290
291////////////////////////////////////////////////////////////////////////////////
292/// Check for illegal crossings between non-consecutive segments
293
295{
296 if (fNvert < 4)
297 return kFALSE;
299 Double_t x1, y1, x2, y2, x3, y3, x4, y4;
300 for (Int_t i = 0; i < fNvert - 2; i++) {
301 // Check segment i
302 for (Int_t j = i + 2; j < fNvert; j++) {
303 // Versus segment j
304 if (i == 0 && j == (fNvert - 1))
305 continue;
306 x1 = fX[i];
307 y1 = fY[i];
308 x2 = fX[i + 1];
309 y2 = fY[i + 1];
310 x3 = fX[j];
311 y3 = fY[j];
312 x4 = fX[(j + 1) % fNvert];
313 y4 = fY[(j + 1) % fNvert];
314 if (TGeoShape::IsSegCrossing(x1, y1, x2, y2, x3, y3, x4, y4)) {
315 Error("IsIllegalCheck", "Illegal crossing of segment %d vs. segment %d", i, j);
317 }
318 }
319 }
320 return is_illegal;
321}
322
323////////////////////////////////////////////////////////////////////////////////
324/// Compute indices for the outscribed convex polygon.
325
327{
328 fNconvex = 0;
329 Int_t iseg = 0;
330 Int_t ivnew;
331 Bool_t conv;
332 Int_t *indconv = new Int_t[fNvert];
333 memset(indconv, 0, fNvert * sizeof(Int_t));
334 while (iseg < fNvert) {
335 if (!IsSegConvex(iseg)) {
336 if (iseg + 2 > fNvert)
337 break;
338 ivnew = (iseg + 2) % fNvert;
339 conv = kFALSE;
340 // check iseg with next vertices
341 while (ivnew) {
342 if (IsSegConvex(iseg, ivnew)) {
343 conv = kTRUE;
344 break;
345 }
346 ivnew = (ivnew + 1) % fNvert;
347 }
348 if (!conv) {
349 // Error("OutscribedConvex","NO convex line connection to vertex %d\n", iseg);
350 iseg++;
351 continue;
352 }
353 } else {
354 ivnew = (iseg + 1) % fNvert;
355 }
356 // segment belonging to convex outscribed polygon
357 if (!fNconvex)
358 indconv[fNconvex++] = iseg;
359 else if (indconv[fNconvex - 1] != iseg)
360 indconv[fNconvex++] = iseg;
361 if (iseg < fNvert - 1)
363 if (ivnew < iseg)
364 break;
365 iseg = ivnew;
366 }
367 if (!fNconvex) {
368 delete[] indconv;
369 Fatal("OutscribedConvex", "cannot build outscribed convex");
370 return;
371 }
372 fIndc = new Int_t[fNvert];
373 memcpy(fIndc, indconv, fNconvex * sizeof(Int_t)); // does not contain real indices yet
374 delete[] indconv;
375}
376
377////////////////////////////////////////////////////////////////////////////////
378/// Compute minimum distance from POINT to any segment. Returns segment index.
379
381{
382 Int_t i1, i2;
383 Double_t p1[2], p2[3];
384 Double_t lsq, ssq, dx, dy, dpx, dpy, u;
385 Double_t safe = 1E30;
386 Int_t isegmin = 0;
387 for (i1 = 0; i1 < fNvert; i1++) {
390 return 0.;
391 }
392 i2 = (i1 + 1) % fNvert;
393 p1[0] = fX[i1];
394 p1[1] = fY[i1];
395 p2[0] = fX[i2];
396 p2[1] = fY[i2];
397
398 dx = p2[0] - p1[0];
399 dy = p2[1] - p1[1];
400 dpx = point[0] - p1[0];
401 dpy = point[1] - p1[1];
402
403 lsq = dx * dx + dy * dy;
405 ssq = dpx * dpx + dpy * dpy;
406 if (ssq < safe) {
407 safe = ssq;
408 isegmin = i1;
409 }
410 continue;
411 }
412 u = (dpx * dx + dpy * dy) / lsq;
413 if (u > 1) {
414 dpx = point[0] - p2[0];
415 dpy = point[1] - p2[1];
416 } else {
417 if (u >= 0) {
418 dpx -= u * dx;
419 dpy -= u * dy;
420 }
421 }
422 ssq = dpx * dpx + dpy * dpy;
423 if (ssq < safe) {
424 safe = ssq;
425 isegmin = i1;
426 }
427 }
430 // printf("== segment %d: (%f, %f) - (%f, %f) safe=%f\n", isegment,
431 // fX[isegment],fY[isegment],fX[(isegment+1)%fNvert],fY[(isegment+1)%fNvert],safe);
432 return safe;
433}
434
435////////////////////////////////////////////////////////////////////////////////
436/// Sets the next polygone index. If index<0 sets all indices consecutive
437/// in increasing order.
438
440{
441 Int_t i;
442 if (index < 0) {
443 for (i = 0; i < fNvert; i++)
444 fInd[i] = i;
445 return;
446 }
447 if (fNconvex >= fNvert) {
448 Error("SetNextIndex", "all indices already set");
449 return;
450 }
451 fInd[fNconvex++] = index;
452 if (fNconvex == fNvert) {
453 if (!fX || !fY)
454 return;
455 Double_t area = 0.0;
456 for (i = 0; i < fNvert; i++)
457 area += fX[fInd[i]] * fY[fInd[(i + 1) % fNvert]] - fX[fInd[(i + 1) % fNvert]] * fY[fInd[i]];
458 if (area < 0)
460 else
462 }
463}
464
465////////////////////////////////////////////////////////////////////////////////
466/// Set X/Y array pointer for the polygon and daughters.
467
469{
470 Int_t i;
471 fX = x;
472 fY = y;
473 Double_t area = 0.0;
474 for (i = 0; i < fNvert; i++)
475 area += fX[fInd[i]] * fY[fInd[(i + 1) % fNvert]] - fX[fInd[(i + 1) % fNvert]] * fY[fInd[i]];
476 if (area < 0)
478 else
480
481 if (!fDaughters)
482 return;
485 for (i = 0; i < nd; i++) {
486 poly = (TGeoPolygon *)fDaughters->At(i);
487 if (poly)
488 poly->SetXY(x, y);
489 }
490}
constexpr Bool_t kFALSE
Definition RtypesCore.h:108
constexpr Bool_t kTRUE
Definition RtypesCore.h:107
const char Option_t
Option string (const char)
Definition RtypesCore.h:80
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t TPoint TPoint const char y1
R__EXTERN TGeoManager * gGeoManager
TVirtualGeoPainter * GetGeomPainter()
Make a default painter if none present. Returns pointer to it.
An arbitrary polygon defined by vertices.
Definition TGeoPolygon.h:19
Bool_t IsConvex() const
Definition TGeoPolygon.h:59
Bool_t IsSegConvex(Int_t i1, Int_t i2=-1) const
Check if a segment [0..fNvert-1] belongs to the outscribed convex pgon.
TGeoPolygon()
Dummy constructor.
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.
Int_t fNconvex
Definition TGeoPolygon.h:26
void GetVertices(Double_t *x, Double_t *y) const
Fill list of vertices into provided arrays.
void SetNextIndex(Int_t index=-1)
Sets the next polygone index.
Int_t * fInd
Definition TGeoPolygon.h:27
Double_t * fY
pointer to list of current X coordinates of vertices
Definition TGeoPolygon.h:30
void Draw(Option_t *option="") override
Draw the polygon.
Bool_t IsClockwise() const
Definition TGeoPolygon.h:58
void OutscribedConvex()
Compute indices for the outscribed convex polygon.
Double_t * fX
Definition TGeoPolygon.h:29
~TGeoPolygon() override
Destructor.
Int_t * fIndc
Definition TGeoPolygon.h:28
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.
void ConvexCheck()
Check polygon convexity.
TObjArray * fDaughters
pointer to list of current Y coordinates of vertices
Definition TGeoPolygon.h:31
void GetConvexVertices(Double_t *x, Double_t *y) const
Fill list of vertices of the convex outscribed polygon into provided arrays.
void FinishPolygon()
Decompose polygon in a convex outscribed part and a list of daughter polygons that have to be subtrac...
void SetConvex(Bool_t flag=kTRUE)
Definition TGeoPolygon.h:63
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),...
static Bool_t IsSameWithinTolerance(Double_t a, Double_t b)
Check if two numbers differ with less than a tolerance.
An array of TObjects.
Definition TObjArray.h:31
Int_t GetEntriesFast() const
Definition TObjArray.h:58
void Delete(Option_t *option="") override
Remove all objects from the array AND delete all heap based objects.
TObject * At(Int_t idx) const override
Definition TObjArray.h:164
TObject * UncheckedAt(Int_t i) const
Definition TObjArray.h:84
void Add(TObject *obj) override
Definition TObjArray.h:68
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:864
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1071
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition TObject.cxx:1099
virtual void DrawPolygon(const TGeoPolygon *poly)=0
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:673
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:124