TGeoPolygon.cxx

Go to the documentation of this file.
00001 // @(#)root/geom:$Id: TGeoPolygon.cxx 35006 2010-08-25 11:51:08Z agheata $
00002 // Author: Mihaela Gheata   5/01/04
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers.               *
00006  * All rights reserved.                                                  *
00007  *                                                                       *
00008  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00009  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00010  *************************************************************************/
00011 
00012 //____________________________________________________________________________
00013 // TGeoPolygon - Arbitrary polygon class 
00014 //____________________________________________________________________________
00015 //
00016 // A polygon is a 2D shape defined by vertices in the XY plane. It is used by
00017 // TGeoXtru class for computing Contains() and Safety(). Only the pointers to
00018 // the actual lists of XY values are used - these are not owned by the class.
00019 // 
00020 // To check if a point in XY plane is contained by a polygon, this is splitted
00021 // into an outscribed convex polygon and the remaining polygons of its subtracton
00022 // from the outscribed one. A point is INSIDE if it is 
00023 // contained by the outscribed polygon but NOT by the remaining ones. Since these
00024 // can also be arbitrary polygons at their turn, a tree structure is formed:
00025 //
00026 //  P = Pconvex - (Pconvex-P)           where (-) means 'subtraction'
00027 //  Pconvex-P = P1 + P2 + ...           where (+) means 'union'
00028 //
00029 //  *Note that P1, P2, ... do not intersect each other and they are defined
00030 //   by subsets of the list of vertices of P. They can be splitted in the same
00031 //   way as P*
00032 //
00033 // Therefore, if C(P) represents the Boolean : 'does P contains a given point?',
00034 // then:
00035 //
00036 // C(P) = C(Pconvex) .and. not(C(P1) | C(P2) | ...)
00037 //
00038 // For creating a polygon without TGeoXtru class, one has to call the constructor
00039 // TGeoPolygon(nvert) and then SetXY(Double_t *x, Double_t *y) providing the
00040 // arrays of X and Y vertex positions (defined clockwise) that have to 'live' longer 
00041 // than the polygon they will describe. This complication is due to efficiency reasons.
00042 // At the end one has to call the FinishPolygon() method.
00043 
00044 #include "TObjArray.h"
00045 #include "TGeoPolygon.h"
00046 #include "TMath.h"
00047 #include "TGeoShape.h"
00048 
00049 ClassImp(TGeoPolygon)
00050 
00051 //_____________________________________________________________________________
00052 TGeoPolygon::TGeoPolygon()
00053 {
00054 // Dummy constructor.
00055    fNvert   = 0;
00056    fNconvex = 0;
00057    fInd     = 0;
00058    fIndc    = 0;
00059    fX       = 0;
00060    fY       = 0;
00061    fDaughters = 0;
00062    SetConvex(kFALSE);
00063    TObject::SetBit(kGeoFinishPolygon, kFALSE);
00064 }
00065 
00066 //_____________________________________________________________________________
00067 TGeoPolygon::TGeoPolygon(Int_t nvert)
00068 {
00069 // Default constructor.
00070    if (nvert<3) {
00071       Fatal("Ctor", "Invalid number of vertices %i", nvert);
00072       return;
00073    }   
00074    fNvert   = nvert;
00075    fNconvex = 0;
00076    fInd     = new Int_t[nvert];
00077    fIndc    = 0;
00078    fX       = 0;
00079    fY       = 0;
00080    fDaughters = 0;
00081    SetConvex(kFALSE);
00082    TObject::SetBit(kGeoFinishPolygon, kFALSE);
00083    SetNextIndex();
00084 }
00085 
00086 //_____________________________________________________________________________
00087 TGeoPolygon::~TGeoPolygon()
00088 {
00089 // Destructor
00090    if (fInd)  delete [] fInd;
00091    if (fIndc) delete [] fIndc;
00092    if (fDaughters) {
00093       fDaughters->Delete();
00094       delete fDaughters;
00095    }   
00096 }
00097 //_____________________________________________________________________________
00098 Double_t TGeoPolygon::Area() const
00099 {
00100 // Computes area of the polygon in [length^2].
00101    Int_t ic,i,j;
00102    Double_t area = 0;
00103    // Compute area of the convex part
00104    for (ic=0; ic<fNvert; ic++) {
00105       i = fInd[ic];
00106       j = fInd[(ic+1)%fNvert];
00107       area += 0.5*(fX[i]*fY[j]-fX[j]*fY[i]);
00108    }
00109    return TMath::Abs(area);
00110 }      
00111 
00112 //_____________________________________________________________________________
00113 Bool_t TGeoPolygon::Contains(Double_t *point) const
00114 {
00115 // Check if a point given by X = point[0], Y = point[1] is inside the polygon.
00116    Int_t i;
00117    TGeoPolygon *poly;
00118    for (i=0; i<fNconvex; i++) 
00119       if (!IsRightSided(point, fIndc[i], fIndc[(i+1)%fNconvex])) return kFALSE;
00120    if (!fDaughters) return kTRUE;
00121    Int_t nd = fDaughters->GetEntriesFast();
00122    for (i=0; i<nd; i++) {
00123       poly = (TGeoPolygon*)fDaughters->UncheckedAt(i);
00124       if (poly->Contains(point)) return kFALSE;
00125    }   
00126    return kTRUE;      
00127 }
00128 
00129 //_____________________________________________________________________________
00130 void TGeoPolygon::ConvexCheck() 
00131 {
00132 // Check polygon convexity.
00133    if (fNvert==3) {
00134       SetConvex(); 
00135       return;
00136    }    
00137    Int_t j,k;
00138    Double_t point[2];
00139    for (Int_t i=0; i<fNvert; i++) {
00140       j = (i+1)%fNvert;
00141       k = (i+2)%fNvert;
00142       point[0] = fX[fInd[k]];
00143       point[1] = fY[fInd[k]];
00144       if (!IsRightSided(point, fInd[i], fInd[j])) return;
00145    }
00146    SetConvex();  
00147 }   
00148 
00149 //_____________________________________________________________________________
00150 void TGeoPolygon::FinishPolygon()
00151 {
00152 // Decompose polygon in a convex outscribed part and a list of daughter
00153 // polygons that have to be substracted to get the actual one.
00154    TObject::SetBit(kGeoFinishPolygon);
00155    // check convexity
00156    ConvexCheck();
00157    // find outscribed convex polygon indices
00158    OutscribedConvex();
00159    if (IsConvex()) {
00160 //      printf(" -> polygon convex -> same indices\n");
00161       memcpy(fIndc, fInd, fNvert*sizeof(Int_t));
00162       return;
00163    }   
00164 //   printf(" -> polygon NOT convex\n");
00165    // make daughters if necessary
00166    if (IsConvex()) return;
00167    // ... algorithm here
00168    if (!fDaughters) fDaughters = new TObjArray();
00169    TGeoPolygon *poly = 0;
00170    Int_t indconv = 0;
00171    Int_t indnext, indback;
00172    Int_t nskip;
00173    while (indconv < fNconvex) {
00174       indnext = (indconv+1)%fNconvex;
00175       nskip = fIndc[indnext]-fIndc[indconv];
00176       if (nskip<0) nskip+=fNvert;
00177       if (nskip==1) {
00178          indconv++;
00179          continue;
00180       }
00181       // gap -> make polygon
00182       poly = new TGeoPolygon(nskip+1);
00183       poly->SetXY(fX,fY);
00184       poly->SetNextIndex(fInd[fIndc[indconv]]);   
00185       poly->SetNextIndex(fInd[fIndc[indnext]]);
00186       indback = fIndc[indnext]-1;
00187       if (indback < 0) indback+=fNvert;
00188       while (indback != fIndc[indconv]) {
00189          poly->SetNextIndex(fInd[indback]); 
00190          indback--;
00191          if (indback < 0) indback+=fNvert;
00192       }
00193       poly->FinishPolygon();
00194       fDaughters->Add(poly);
00195       indconv++;
00196    }   
00197    for (indconv=0; indconv<fNconvex; indconv++) fIndc[indconv] = fInd[fIndc[indconv]];
00198 }
00199 
00200 //_____________________________________________________________________________
00201 Bool_t TGeoPolygon::IsRightSided(Double_t *point, Int_t ind1, Int_t ind2) const
00202 {
00203 // Check if POINT is right-sided with respect to the segment defined by IND1 and IND2.
00204    Double_t dot = (point[0]-fX[ind1])*(fY[ind2]-fY[ind1]) -
00205                   (point[1]-fY[ind1])*(fX[ind2]-fX[ind1]);
00206    if (!IsClockwise()) dot = -dot;
00207    if (dot<0) return kFALSE;
00208    return kTRUE;
00209 }
00210 
00211 //_____________________________________________________________________________
00212 Bool_t TGeoPolygon::IsSegConvex(Int_t i1, Int_t i2) const
00213 {
00214 // Check if a segment [0..fNvert-1] belongs to the outscribed convex pgon.
00215    if (i2<0) i2=(i1+1)%fNvert;
00216    Double_t point[2];
00217    for (Int_t i=0; i<fNvert; i++) {
00218       if (i==i1 || i==i2) continue;
00219       point[0] = fX[fInd[i]];
00220       point[1] = fY[fInd[i]];
00221       if (!IsRightSided(point, fInd[i1], fInd[i2])) return kFALSE;
00222    }
00223    return kTRUE;
00224 }      
00225 
00226 //_____________________________________________________________________________
00227 Bool_t TGeoPolygon::IsIllegalCheck() const
00228 {
00229 // Check for illegal crossings between non-consecutive segments
00230    if (fNvert<4) return kFALSE;
00231    Bool_t is_illegal = kFALSE;
00232    Double_t x1,y1,x2,y2,x3,y3,x4,y4;
00233    for (Int_t i=0; i<fNvert-2; i++) {
00234       // Check segment i
00235       for (Int_t j=i+2; j<fNvert; j++) {
00236          // Versus segment j
00237          if (i==0 && j==(fNvert-1)) continue;
00238          x1 = fX[i];
00239          y1 = fY[i];
00240          x2 = fX[i+1];
00241          y2 = fY[i+1];
00242          x3 = fX[j];
00243          y3 = fY[j];
00244          x4 = fX[(j+1)%fNvert];
00245          y4 = fY[(j+1)%fNvert];
00246          if (TGeoShape::IsSegCrossing(x1,y1,x2,y2,x3,y3,x4,y4)) {
00247             Error("IsIllegalCheck", "Illegal crossing of segment %d vs. segment %d", i,j);
00248             is_illegal = kTRUE;
00249          }
00250       }
00251    }
00252    return is_illegal;
00253 }
00254    
00255 //_____________________________________________________________________________
00256 void TGeoPolygon::OutscribedConvex()
00257 {
00258 // Compute indices for the outscribed convex polygon.
00259    fNconvex = 0;
00260    Int_t iseg = 0;
00261    Int_t ivnew;
00262    Bool_t conv;
00263    Int_t *indconv = new Int_t[fNvert];
00264    memset(indconv, 0, fNvert*sizeof(Int_t));
00265    while (iseg<fNvert) {
00266       if (!IsSegConvex(iseg)) {
00267          if (iseg+2 > fNvert) break;
00268          ivnew = (iseg+2)%fNvert;
00269          conv = kFALSE;
00270          // check iseg with next vertices
00271          while (ivnew) {
00272             if (IsSegConvex(iseg, ivnew)) {
00273                conv = kTRUE;
00274                break;
00275             } 
00276             ivnew = (ivnew+1)%fNvert;  
00277          } 
00278          if (!conv) {
00279             iseg++;
00280             continue;
00281          }   
00282       } else {
00283          ivnew = (iseg+1)%fNvert;
00284       }   
00285       // segment belonging to convex outscribed poligon
00286       if (!fNconvex) indconv[fNconvex++] = iseg;
00287       else if (indconv[fNconvex-1] != iseg) indconv[fNconvex++] = iseg;
00288       if (iseg<fNvert-1) indconv[fNconvex++] = ivnew;
00289       if (ivnew<iseg) break;
00290       iseg = ivnew;
00291    }    
00292    if (!fNconvex) {
00293       delete [] indconv;
00294       Fatal("OutscribedConvex","cannot build outscribed convex");
00295       return;
00296    }
00297    fIndc = new Int_t[fNconvex];
00298    memcpy(fIndc, indconv, fNconvex*sizeof(Int_t)); // does not contain real indices yet
00299    delete [] indconv;
00300 }
00301 
00302 //_____________________________________________________________________________
00303 Double_t TGeoPolygon::Safety(Double_t *point, Int_t &isegment) const
00304 {
00305 // Compute minimum distance from POINT to any segment. Returns segment index.
00306    Int_t i1, i2;
00307    Double_t p1[2], p2[3];
00308    Double_t lsq, ssq, dx, dy, dpx, dpy, u;
00309    Double_t safe=1E30;
00310    Int_t isegmin=0;
00311    for (i1=0; i1<fNvert; i1++) {
00312       if (TGeoShape::IsSameWithinTolerance(safe,0)) {
00313          isegment = isegmin;
00314          return 0.;
00315       }   
00316       i2 = (i1+1)%fNvert;
00317       p1[0] = fX[i1];
00318       p1[1] = fY[i1];
00319       p2[0] = fX[i2];
00320       p2[1] = fY[i2];
00321       
00322       dx = p2[0] - p1[0];
00323       dy = p2[1] - p1[1];
00324       dpx = point[0] - p1[0];
00325       dpy = point[1] - p1[1];
00326       
00327       lsq = dx*dx + dy*dy;
00328       if (TGeoShape::IsSameWithinTolerance(lsq,0)) {
00329          ssq = dpx*dpx + dpy*dpy;
00330          if (ssq < safe) {
00331             safe = ssq;
00332             isegmin = i1;
00333          }
00334          continue;
00335       } 
00336       u = (dpx*dx + dpy*dy)/lsq;
00337       if (u>1) {
00338          dpx = point[0]-p2[0];
00339          dpy = point[1]-p2[1];
00340       } else {
00341          if (u>=0) {
00342             dpx -= u*dx;
00343             dpy -= u*dy;
00344          }
00345       }
00346       ssq = dpx*dpx + dpy*dpy;      
00347       if (ssq < safe) {
00348          safe = ssq;
00349          isegmin = i1;
00350       }
00351    }
00352    isegment = isegmin;
00353    safe = TMath::Sqrt(safe);
00354 //   printf("== segment %d: (%f, %f) - (%f, %f) safe=%f\n", isegment, fX[isegment],fY[isegment],fX[(isegment+1)%fNvert],fY[(isegment+1)%fNvert],safe);
00355    return safe;
00356 }
00357 
00358 //_____________________________________________________________________________
00359 void TGeoPolygon::SetNextIndex(Int_t index)
00360 {
00361 // Sets the next polygone index. If index<0 sets all indices consecutive
00362 // in increasing order.
00363    Int_t i;
00364    if (index <0) {
00365       for (i=0; i<fNvert; i++) fInd[i] = i;
00366       return;
00367    }
00368    if (fNconvex >= fNvert) {
00369       Error("SetNextIndex", "all indices already set");
00370       return;
00371    }
00372    fInd[fNconvex++] = index;  
00373    if (fNconvex == fNvert) {
00374       if (!fX || !fY) return;
00375       Double_t area = 0.0;
00376       for (i=0; i<fNvert; i++) area += fX[fInd[i]]*fY[fInd[(i+1)%fNvert]]-fX[fInd[(i+1)%fNvert]]*fY[fInd[i]];
00377       if (area<0) TObject::SetBit(kGeoACW, kFALSE);
00378       else        TObject::SetBit(kGeoACW, kTRUE);
00379    }
00380 }
00381 
00382 //_____________________________________________________________________________
00383 void TGeoPolygon::SetXY(Double_t *x, Double_t *y)
00384 {
00385 // Set X/Y array pointer for the polygon and daughters.
00386    Int_t i;
00387    fX = x;
00388    fY = y;
00389    Double_t area = 0.0;
00390    for (i=0; i<fNvert; i++) area += fX[fInd[i]]*fY[fInd[(i+1)%fNvert]]-fX[fInd[(i+1)%fNvert]]*fY[fInd[i]];
00391    if (area<0) TObject::SetBit(kGeoACW, kFALSE);
00392    else        TObject::SetBit(kGeoACW, kTRUE);
00393    
00394    if (!fDaughters) return;
00395    TGeoPolygon *poly;
00396    Int_t nd = fDaughters->GetEntriesFast();
00397    for (i=0; i<nd; i++) {
00398       poly = (TGeoPolygon*)fDaughters->At(i);
00399       if (poly) poly->SetXY(x,y);
00400    }
00401 }

Generated on Tue Jul 5 14:12:57 2011 for ROOT_528-00b_version by  doxygen 1.5.1