fastjet 2.4.5
Voronoi.cc
Go to the documentation of this file.
00001 //STARTHEADER
00002 // $Id: Voronoi.cc 2295 2011-06-28 10:40:21Z soyez $
00003 //
00004 // Copyright (c) 1994 by AT&T Bell Laboratories (see below)
00005 //
00006 //
00007 //----------------------------------------------------------------------
00008 // This file is included as part of FastJet but was mostly written by
00009 // S. Fortune in C, put into C++ with memory management by S
00010 // O'Sullivan, and with further interface and memory management
00011 // modifications by Gregory Soyez.
00012 //
00013 // Permission to use, copy, modify, and distribute this software for
00014 // any purpose without fee is hereby granted, provided that this
00015 // entire notice is included in all copies of any software which is or
00016 // includes a copy or modification of this software and in all copies
00017 // of the supporting documentation for such software. THIS SOFTWARE IS
00018 // BEING PROVIDED "AS IS", WITHOUT ANY EXPRESS OR IMPLIED WARRANTY.
00019 // IN PARTICULAR, NEITHER THE AUTHORS NOR AT&T MAKE ANY REPRESENTATION
00020 // OR WARRANTY OF ANY KIND CONCERNING THE MERCHANTABILITY OF THIS
00021 // SOFTWARE OR ITS FITNESS FOR ANY PARTICULAR PURPOSE.
00022 //
00023 //----------------------------------------------------------------------
00024 //ENDHEADER
00025 
00026 
00027 /*
00028  * The author of this software is Steven Fortune.  
00029  * Copyright (c) 1994 by AT&T Bell Laboratories.
00030  * Permission to use, copy, modify, and distribute this software for any
00031  * purpose without fee is hereby granted, provided that this entire notice
00032  * is included in all copies of any software which is or includes a copy
00033  * or modification of this software and in all copies of the supporting
00034  * documentation for such software.
00035  * THIS SOFTWARE IS BEING PROVIDED "AS IS", WITHOUT ANY EXPRESS OR IMPLIED
00036  * WARRANTY.  IN PARTICULAR, NEITHER THE AUTHORS NOR AT&T MAKE ANY
00037  * REPRESENTATION OR WARRANTY OF ANY KIND CONCERNING THE MERCHANTABILITY
00038  * OF THIS SOFTWARE OR ITS FITNESS FOR ANY PARTICULAR PURPOSE.
00039  */
00040 
00041 /* 
00042  * This code was originally written by Stephan Fortune in C code.  I,
00043  * Shane O'Sullivan, have since modified it, encapsulating it in a C++
00044  * class and, fixing memory leaks and adding accessors to the Voronoi
00045  * Edges.  Permission to use, copy, modify, and distribute this
00046  * software for any purpose without fee is hereby granted, provided
00047  * that this entire notice is included in all copies of any software
00048  * which is or includes a copy or modification of this software and in
00049  * all copies of the supporting documentation for such software.  THIS
00050  * SOFTWARE IS BEING PROVIDED "AS IS", WITHOUT ANY EXPRESS OR IMPLIED
00051  * WARRANTY.  IN PARTICULAR, NEITHER THE AUTHORS NOR AT&T MAKE ANY
00052  * REPRESENTATION OR WARRANTY OF ANY KIND CONCERNING THE
00053  * MERCHANTABILITY OF THIS SOFTWARE OR ITS FITNESS FOR ANY PARTICULAR
00054  * PURPOSE.
00055  */
00056 
00057 #include <stdio.h>
00058 #include "fastjet/internal/Voronoi.hh"
00059 
00060 using namespace std;
00061 
00062 FASTJET_BEGIN_NAMESPACE
00063 
00064 LimitedWarning VoronoiDiagramGenerator::_warning_degeneracy;
00065 
00066 VoronoiDiagramGenerator::VoronoiDiagramGenerator(){
00067   siteidx = 0;
00068   sites = NULL;
00069   parent_sites = NULL;
00070 
00071   allMemoryList = new FreeNodeArrayList;
00072   allMemoryList->memory = NULL;
00073   allMemoryList->next = NULL;
00074   currentMemoryBlock = allMemoryList;
00075   allEdges = NULL;
00076   iteratorEdges = 0;
00077   minDistanceBetweenSites = 0;
00078 
00079   ELhash = NULL;
00080   PQhash = NULL;
00081 }
00082 
00083 VoronoiDiagramGenerator::~VoronoiDiagramGenerator(){
00084   cleanup();
00085   cleanupEdges();
00086 
00087   if (allMemoryList != NULL)
00088     delete allMemoryList;
00089 }
00090 
00091 
00092 
00093 bool VoronoiDiagramGenerator::generateVoronoi(vector<Point> *_parent_sites,
00094                                               double minX, double maxX, 
00095                                               double minY, double maxY, 
00096                                               double minDist){
00097   cleanup();
00098   cleanupEdges();
00099   int i;
00100   double x, y;
00101 
00102   minDistanceBetweenSites = minDist;
00103 
00104   parent_sites = _parent_sites;
00105 
00106   nsites = n_parent_sites = parent_sites->size();
00107   plot = 0;
00108   triangulate = 0;      
00109   debug = 1;
00110   sorted = 0; 
00111   freeinit(&sfl, sizeof (Site));
00112                 
00113   //sites = (Site *) myalloc(nsites*sizeof( *sites));
00114   sites = (Site *) myalloc(nsites*sizeof( Site));
00115   //parent_sites = (Site *) myalloc(nsites*sizeof( Site));
00116  
00117   if (sites == 0)
00118     return false;
00119 
00120   xmax = xmin = (*parent_sites)[0].x;
00121   ymax = ymin = (*parent_sites)[0].y;
00122   
00123   for(i=0;i<nsites;i++){
00124     x = (*parent_sites)[i].x;
00125     y = (*parent_sites)[i].y;
00126     sites[i].coord.x = x;
00127     sites[i].coord.y = y;
00128     sites[i].sitenbr = i;
00129     sites[i].refcnt  = 0;
00130     
00131     if (x<xmin)
00132       xmin=x;
00133     else if (x>xmax)
00134       xmax=x;
00135     
00136     if (y<ymin)
00137       ymin=y;
00138     else if (y>ymax)
00139       ymax=y;
00140   }
00141         
00142   qsort(sites, nsites, sizeof (*sites), scomp);
00143 
00144   // Gregory Soyez
00145   // 
00146   // check if some of the particles are degenerate to avoid a crash.
00147   //
00148   // At the moment, we work under the assumption that they will be
00149   // "clustered" later on, so we just keep the 1st one and discard the
00150   // others
00151   unsigned int offset=0;
00152   for (int i=1;i<nsites;i++){
00153     if (sites[i].coord.y==sites[i-1].coord.y && sites[i].coord.x==sites[i-1].coord.x){
00154       offset++;
00155     } else if (offset>0){
00156       sites[i-offset] = sites[i];
00157     }
00158   }
00159 
00160   if (offset>0){
00161     nsites-=offset;
00162     _warning_degeneracy.warn("VoronoiDiagramGenerator: two (or more) particles are degenerate in rapidity and azimuth, Voronoi cell assigned to the first of each set of degenerate particles.");
00163   }
00164 
00165   siteidx = 0;
00166   geominit();
00167   double temp = 0;
00168   if(minX > maxX){
00169     temp = minX;
00170     minX = maxX;
00171     maxX = temp;
00172   }
00173   if(minY > maxY){
00174     temp = minY;
00175     minY = maxY;
00176     maxY = temp;
00177   }
00178   borderMinX = minX;
00179   borderMinY = minY;
00180   borderMaxX = maxX;
00181   borderMaxY = maxY;
00182         
00183   siteidx = 0;
00184   voronoi(triangulate);
00185 
00186   return true;
00187 }
00188 
00189 bool VoronoiDiagramGenerator::ELinitialize(){
00190   int i;
00191   freeinit(&hfl, sizeof(Halfedge));
00192   ELhashsize = 2 * sqrt_nsites;
00193   //ELhash = (Halfedge **) myalloc ( sizeof *ELhash * ELhashsize);
00194   ELhash = (Halfedge **) myalloc ( sizeof(Halfedge*)*ELhashsize);
00195   
00196   if(ELhash == 0)
00197     return false;
00198   
00199   for(i=0; i<ELhashsize; i +=1) ELhash[i] = (Halfedge *)NULL;
00200   ELleftend = HEcreate((Edge*) NULL, 0);
00201   ELrightend = HEcreate((Edge*) NULL, 0);
00202   ELleftend->ELleft = (Halfedge*) NULL;
00203   ELleftend->ELright = ELrightend;
00204   ELrightend->ELleft = ELleftend;
00205   ELrightend->ELright = (Halfedge*) NULL;
00206   ELhash[0] = ELleftend;
00207   ELhash[ELhashsize-1] = ELrightend;
00208 
00209   return true;
00210 }
00211 
00212 
00213 Halfedge* VoronoiDiagramGenerator::HEcreate(Edge *e,int pm){
00214   Halfedge *answer;
00215   answer = (Halfedge*) getfree(&hfl);
00216   answer->ELedge = e;
00217   answer->ELpm = pm;
00218   answer->PQnext = (Halfedge*) NULL;
00219   answer->vertex = (Site*) NULL;
00220   answer->ELrefcnt = 0;
00221   return answer;
00222 }
00223 
00224 
00225 void VoronoiDiagramGenerator::ELinsert(Halfedge *lb, Halfedge *newHe)
00226 {
00227   newHe->ELleft = lb;
00228   newHe->ELright = lb->ELright;
00229   (lb->ELright)->ELleft = newHe;
00230   lb->ELright = newHe;
00231 }
00232 
00233 
00234 /* Get entry from hash table, pruning any deleted nodes */
00235 Halfedge* VoronoiDiagramGenerator::ELgethash(int b){
00236   Halfedge *he;
00237         
00238   if ((b<0) || (b>=ELhashsize)) 
00239     return (Halfedge *) NULL;
00240 
00241   he = ELhash[b]; 
00242   if ((he==(Halfedge*) NULL) || (he->ELedge!=(Edge*) DELETED)) 
00243     return he;
00244         
00245   /* Hash table points to deleted half edge.  Patch as necessary. */
00246   ELhash[b] = (Halfedge*) NULL;
00247   if ((he->ELrefcnt -= 1) == 0) 
00248     makefree((Freenode*)he, &hfl);
00249   return (Halfedge*) NULL;
00250 }       
00251 
00252 Halfedge * VoronoiDiagramGenerator::ELleftbnd(Point *p){
00253   int i, bucket;
00254   Halfedge *he;
00255         
00256   /* Use hash table to get close to desired halfedge */
00257   // use the hash function to find the place in the hash map that this
00258   // HalfEdge should be
00259   // Gregory Soyez: the original code was 
00260   //
00261   //   bucket = (int)((p->x - xmin)/deltax * ELhashsize);       
00262   //   // make sure that the bucket position in within the range of the
00263   //   //hash array
00264   //   if (bucket<0) bucket =0;
00265   //   if (bucket>=ELhashsize) bucket = ELhashsize - 1;
00266   //
00267   // but this runs the risk of having a overflow which would 
00268   // cause bucket to be truncated at 0 instead of ELhashsize - 1
00269   // (or vice-versa)
00270   // We fix this by performing the test immediately on the double
00271   // We put in a extra bit of margin to be sure the conversion does
00272   // not play dirty tricks on us
00273 
00274   //const double &px = p->x;
00275   if (p->x < xmin){ bucket=0;}
00276   else if (p->x >= xmax){ bucket = ELhashsize - 1;}
00277   else{
00278     bucket= (int)((p->x - xmin)/deltax * ELhashsize);
00279     if (bucket>=ELhashsize) bucket = ELhashsize - 1;  // the lower cut should be robust
00280   }
00281 
00282   he = ELgethash(bucket);
00283 
00284   // if the HE isn't found, search backwards and forwards in the hash
00285   // map for the first non-null entry
00286   if (he == (Halfedge*) NULL){   
00287     for(i=1;1;i++){     
00288       if ((he=ELgethash(bucket-i)) != (Halfedge*) NULL) 
00289         break;
00290       if ((he=ELgethash(bucket+i)) != (Halfedge*) NULL) 
00291         break;
00292     };
00293     totalsearch += i;
00294   };
00295   ntry += 1;
00296   /* Now search linear list of halfedges for the correct one */
00297   if ((he==ELleftend)  || (he != ELrightend && right_of(he,p))){
00298     do{
00299       he = he->ELright;
00300     } while (he!=ELrightend && right_of(he,p));
00301       // keep going right on the list until either the end is reached,
00302       // or you find the 1st edge which the point
00303       he = he->ELleft;  //isn't to the right of
00304   } else 
00305     // if the point is to the left of the HalfEdge, then search left
00306     // for the HE just to the left of the point
00307     do{
00308       he = he->ELleft;
00309     } while ((he!=ELleftend )&& (!right_of(he,p)));
00310                 
00311   /* Update hash table and reference counts */
00312   if ((bucket > 0) && (bucket <ELhashsize-1)){  
00313     if(ELhash[bucket] != (Halfedge *) NULL){
00314       ELhash[bucket]->ELrefcnt -= 1;
00315     }
00316     ELhash[bucket] = he;
00317     ELhash[bucket]->ELrefcnt += 1;
00318   };
00319 
00320   return he;
00321 }
00322 
00323 
00324 /* This delete routine can't reclaim node, since pointers from hash
00325    table may be present.   */
00326 void VoronoiDiagramGenerator::ELdelete(Halfedge *he){
00327   (he->ELleft)->ELright = he->ELright;
00328   (he->ELright)->ELleft = he->ELleft;
00329   he->ELedge = (Edge*) DELETED;
00330 }
00331 
00332 
00333 Halfedge* VoronoiDiagramGenerator::ELright(Halfedge *he){
00334   return he->ELright;
00335 }
00336 
00337 
00338 Halfedge* VoronoiDiagramGenerator::ELleft(Halfedge *he){
00339   return he->ELleft;
00340 }
00341 
00342 
00343 Site * VoronoiDiagramGenerator::leftreg(Halfedge *he){
00344   if (he->ELedge == (Edge*) NULL) 
00345     return(bottomsite);
00346   return (he->ELpm == le) 
00347     ? he->ELedge->reg[le] 
00348     : he->ELedge->reg[re];
00349 }
00350 
00351 Site * VoronoiDiagramGenerator::rightreg(Halfedge *he){
00352   if (he->ELedge == (Edge*) NULL)
00353     // if this halfedge has no edge, return the bottom site (whatever
00354     // that is)
00355     return bottomsite;
00356   
00357   // if the ELpm field is zero, return the site 0 that this edge
00358   // bisects, otherwise return site number 1
00359   return (he->ELpm == le)
00360     ? he->ELedge->reg[re] 
00361     : he->ELedge->reg[le];
00362 }
00363 
00364 void VoronoiDiagramGenerator::geominit(){
00365   double sn;
00366   
00367   freeinit(&efl, sizeof(Edge));
00368   nvertices = 0;
00369   nedges = 0;
00370   sn = (double)nsites+4;
00371   sqrt_nsites = (int)sqrt(sn);
00372   deltay = ymax - ymin;
00373   deltax = xmax - xmin;
00374 }
00375 
00376 
00377 Edge * VoronoiDiagramGenerator::bisect(Site *s1, Site *s2){
00378   double dx,dy,adx,ady;
00379   Edge *newedge;        
00380 
00381   newedge = (Edge*) getfree(&efl);
00382         
00383   newedge->reg[0] = s1; //store the sites that this edge is bisecting
00384   newedge->reg[1] = s2;
00385   ref(s1); 
00386   ref(s2);
00387 
00388   // to begin with, there are no endpoints on the bisector - it goes
00389   // to infinity
00390   newedge->ep[0] = (Site*) NULL;
00391   newedge->ep[1] = (Site*) NULL;
00392   
00393   // get the difference in x dist between the sites
00394   dx = s2->coord.x - s1->coord.x;
00395   dy = s2->coord.y - s1->coord.y;
00396 
00397   // make sure that the difference is positive
00398   adx = dx>0 ? dx : -dx;
00399   ady = dy>0 ? dy : -dy;
00400 
00401   // get the slope of the line
00402   newedge->c = (double)(s1->coord.x * dx + s1->coord.y * dy
00403                         + (dx*dx + dy*dy)*0.5);
00404 
00405   if (adx>ady){ 
00406     //set formula of line, with x fixed to 1
00407     newedge->a = 1.0; newedge->b = dy/dx; newedge->c /= dx;
00408   } else {      
00409     //set formula of line, with y fixed to 1
00410     newedge->b = 1.0; newedge->a = dx/dy; newedge->c /= dy;
00411   }
00412         
00413   newedge->edgenbr = nedges;
00414   nedges++;
00415 
00416   return newedge;
00417 }
00418 
00419 
00420 // create a new site where the HalfEdges el1 and el2 intersect - note
00421 // that the Point in the argument list is not used, don't know why
00422 // it's there
00423 Site* VoronoiDiagramGenerator::intersect(Halfedge *el1, Halfedge *el2, Point *p){
00424   Edge *e1,*e2, *e;
00425   Halfedge *el;
00426   double d, xint, yint;
00427   int right_of_site;
00428   Site *v;
00429         
00430   e1 = el1->ELedge;
00431   e2 = el2->ELedge;
00432   if ((e1==(Edge*)NULL) || (e2==(Edge*)NULL)) 
00433     return ((Site*) NULL);
00434 
00435   // if the two edges bisect the same parent, return null
00436   if (e1->reg[1] == e2->reg[1]) 
00437     return (Site*) NULL;
00438 
00439   // Gregory Soyez:
00440   //    
00441   // if the 2 parents are too close, the intersection is going to be
00442   // computed from the "long edges" of the triangle which could causes
00443   // large rounding errors. In this case, use the bisector of the 2
00444   // parents to find the interaction point
00445   // 
00446   // The following replaces 
00447   //   d = e1->a * e2->b - e1->b * e2->a;
00448   //   if (-1.0e-10<d && d<1.0e-10) 
00449   //     return (Site*) NULL;
00450   //    
00451   //   xint = (e1->c*e2->b - e2->c*e1->b)/d;
00452   //   yint = (e2->c*e1->a - e1->c*e2->a)/d;
00453 
00454   double dx = e2->reg[1]->coord.x - e1->reg[1]->coord.x;
00455   double dy = e2->reg[1]->coord.y - e1->reg[1]->coord.y;
00456   double dxref = e1->reg[1]->coord.x - e1->reg[0]->coord.x;
00457   double dyref = e1->reg[1]->coord.y - e1->reg[0]->coord.y;
00458 
00459   if (dx*dx + dy*dy < 1e-14*(dxref*dxref+dyref*dyref)){
00460     // make sure that the difference is positive
00461     double adx = dx>0 ? dx : -dx;
00462     double ady = dy>0 ? dy : -dy;
00463     
00464     // get the slope of the line
00465     double a,b;
00466     double c = (double)(e1->reg[1]->coord.x * dx + e1->reg[1]->coord.y * dy
00467                         + (dx*dx + dy*dy)*0.5);
00468     
00469     if (adx>ady){
00470       a = 1.0; b = dy/dx; c /= dx;
00471     } else {
00472       b = 1.0; a = dx/dy; c /= dy;
00473     }
00474 
00475     d = e1->a * b - e1->b * a;
00476     if (-1.0e-10<d && d<1.0e-10) {
00477       return (Site*) NULL;
00478     }
00479         
00480     xint = (e1->c*b - c*e1->b)/d;
00481     yint = (c*e1->a - e1->c*a)/d;
00482     
00483   } else {      
00484     d = e1->a * e2->b - e1->b * e2->a;
00485     if (-1.0e-10<d && d<1.0e-10) {
00486       return (Site*) NULL;
00487     }
00488         
00489     xint = (e1->c*e2->b - e2->c*e1->b)/d;
00490     yint = (e2->c*e1->a - e1->c*e2->a)/d;
00491   }
00492   // end of Gregory Soyez's modifications
00493 
00494   volatile double local_y1 = e1->reg[1]->coord.y;
00495   volatile double local_y2 = e2->reg[1]->coord.y;
00496         
00497   if( (local_y1 < local_y2) ||
00498       ((local_y1 == local_y2) &&
00499        (e1->reg[1]->coord.x <  e2->reg[1]->coord.x)) ){ 
00500     el = el1; 
00501     e = e1;
00502   } else {      
00503     el = el2; 
00504     e = e2;
00505   }
00506         
00507   right_of_site = xint >= e->reg[1]->coord.x;
00508   if ((right_of_site && el->ELpm == le) || (!right_of_site && el->ELpm == re)) 
00509     return (Site*) NULL;
00510         
00511   // create a new site at the point of intersection - this is a new
00512   // vector event waiting to happen
00513   v = (Site*) getfree(&sfl);
00514   v->refcnt = 0;
00515   v->coord.x = xint;
00516   v->coord.y = yint;
00517   return v;
00518 }
00519 
00520 //HERE
00521 
00522 /* returns 1 if p is to right of halfedge e */
00523 int VoronoiDiagramGenerator::right_of(Halfedge *el,Point *p)
00524 {
00525   Edge *e;
00526   Site *topsite;
00527   int right_of_site, above, fast;
00528   double dxp, dyp, dxs, t1, t2, t3, yl;
00529         
00530   e = el->ELedge;
00531   topsite = e->reg[1];
00532   right_of_site = p->x > topsite->coord.x;
00533   if(right_of_site && el->ELpm == le) return(1);
00534   if(!right_of_site && el->ELpm == re) return (0);
00535         
00536   if (e->a == 1.0)
00537     {   dyp = p->y - topsite->coord.y;
00538     dxp = p->x - topsite->coord.x;
00539     fast = 0;
00540     if ((!right_of_site & (e->b<0.0)) | (right_of_site & (e->b>=0.0)) )
00541       { above = dyp>= e->b*dxp; 
00542       fast = above;
00543       }
00544     else 
00545       { above = p->x + p->y*e->b > e-> c;
00546       if(e->b<0.0) above = !above;
00547       if (!above) fast = 1;
00548       };
00549     if (!fast)
00550       { dxs = topsite->coord.x - (e->reg[0])->coord.x;
00551       above = e->b * (dxp*dxp - dyp*dyp) <
00552         dxs*dyp*(1.0+2.0*dxp/dxs + e->b*e->b);
00553       if(e->b<0.0) above = !above;
00554       };
00555     }
00556   else  /*e->b==1.0 */
00557     {   yl = e->c - e->a*p->x;
00558     t1 = p->y - yl;
00559     t2 = p->x - topsite->coord.x;
00560     t3 = yl - topsite->coord.y;
00561     above = t1*t1 > t2*t2 + t3*t3;
00562     };
00563   return (el->ELpm==le ? above : !above);
00564 }
00565 
00566 
00567 void VoronoiDiagramGenerator::endpoint(Edge *e,int lr,Site * s)
00568 {
00569   e->ep[lr] = s;
00570   ref(s);
00571   if(e->ep[re-lr]== (Site *) NULL) 
00572     return;
00573 
00574   clip_line(e);
00575 
00576   deref(e->reg[le]);
00577   deref(e->reg[re]);
00578   makefree((Freenode*)e, &efl);
00579 }
00580 
00581 
00582 double VoronoiDiagramGenerator::dist(Site *s,Site *t)
00583 {
00584   double dx,dy;
00585   dx = s->coord.x - t->coord.x;
00586   dy = s->coord.y - t->coord.y;
00587   return (double)(sqrt(dx*dx + dy*dy));
00588 }
00589 
00590 
00591 void VoronoiDiagramGenerator::makevertex(Site *v)
00592 {
00593   v->sitenbr = nvertices;
00594   nvertices += 1;
00595   out_vertex(v);
00596 }
00597 
00598 
00599 void VoronoiDiagramGenerator::deref(Site *v)
00600 {
00601   v->refcnt -= 1;
00602   if (v->refcnt == 0 ) 
00603     makefree((Freenode*)v, &sfl);
00604 }
00605 
00606 void VoronoiDiagramGenerator::ref(Site *v)
00607 {
00608   v->refcnt += 1;
00609 }
00610 
00611 //push the HalfEdge into the ordered linked list of vertices
00612 void VoronoiDiagramGenerator::PQinsert(Halfedge *he,Site * v, double offset)
00613 {
00614   Halfedge *last, *next;
00615         
00616   he->vertex = v;
00617   ref(v);
00618   he->ystar = (double)(v->coord.y + offset);
00619   last = &PQhash[PQbucket(he)];
00620   while ((next = last->PQnext) != (Halfedge *) NULL &&
00621          (he->ystar  > next->ystar  ||
00622           (he->ystar == next->ystar && v->coord.x > next->vertex->coord.x)))
00623     {   
00624       last = next;
00625     };
00626   he->PQnext = last->PQnext; 
00627   last->PQnext = he;
00628   PQcount += 1;
00629 }
00630 
00631 //remove the HalfEdge from the list of vertices 
00632 void VoronoiDiagramGenerator::PQdelete(Halfedge *he)
00633 {
00634   Halfedge *last;
00635         
00636   if(he->vertex != (Site *) NULL)
00637     {   
00638       last = &PQhash[PQbucket(he)];
00639       while (last->PQnext != he) 
00640         last = last->PQnext;
00641 
00642       last->PQnext = he->PQnext;
00643       PQcount -= 1;
00644       deref(he->vertex);
00645       he->vertex = (Site *) NULL;
00646     };
00647 }
00648 
00649 int VoronoiDiagramGenerator::PQbucket(Halfedge *he)
00650 {
00651   // Gregory Soyez: the original code was 
00652   //
00653   //   bucket = (int)((he->ystar - ymin)/deltay * PQhashsize);
00654   //   if (bucket<0) bucket = 0;
00655   //   if (bucket>=PQhashsize) bucket = PQhashsize-1 ;
00656   //   if (bucket < PQmin) PQmin = bucket;
00657   //   return(bucket);
00658   //
00659   // but this runs the risk of having a overflow which would 
00660   // cause bucket to be truncated at 0 instead of PQhashsize-1
00661   // (or vice-versa)
00662   // We fix this by performing the test immediately on the double
00663   // We put in a extra bit of margin to be sure the conversion does
00664   // not play dirty tricks on us
00665 
00666   int bucket;
00667         
00668   double hey = he->ystar;
00669   if (hey < ymin){ bucket = 0;}
00670   else if (hey >= ymax){ bucket = PQhashsize-1;}
00671   else {
00672     bucket = (int)((hey - ymin)/deltay * PQhashsize);
00673     if (bucket>=PQhashsize) bucket = PQhashsize-1 ;
00674   }
00675 
00676   if (bucket < PQmin) PQmin = bucket;
00677   return(bucket);
00678 }
00679 
00680 
00681 
00682 int VoronoiDiagramGenerator::PQempty()
00683 {
00684   return(PQcount==0);
00685 }
00686 
00687 
00688 Point VoronoiDiagramGenerator::PQ_min()
00689 {
00690   Point answer;
00691         
00692   while(PQhash[PQmin].PQnext == (Halfedge *)NULL) {PQmin += 1;};
00693   answer.x = PQhash[PQmin].PQnext->vertex->coord.x;
00694   answer.y = PQhash[PQmin].PQnext->ystar;
00695   return (answer);
00696 }
00697 
00698 Halfedge * VoronoiDiagramGenerator::PQextractmin()
00699 {
00700   Halfedge *curr;
00701         
00702   curr = PQhash[PQmin].PQnext;
00703   PQhash[PQmin].PQnext = curr->PQnext;
00704   PQcount -= 1;
00705   return(curr);
00706 }
00707 
00708 
00709 bool VoronoiDiagramGenerator::PQinitialize()
00710 {
00711   int i; 
00712         
00713   PQcount = 0;
00714   PQmin = 0;
00715   PQhashsize = 4 * sqrt_nsites;
00716   //PQhash = (Halfedge *) myalloc(PQhashsize * sizeof *PQhash);
00717   PQhash = (Halfedge *) myalloc(PQhashsize * sizeof(Halfedge));
00718 
00719   if(PQhash == 0)
00720     return false;
00721 
00722   for(i=0; i<PQhashsize; i+=1) PQhash[i].PQnext = (Halfedge *)NULL;
00723 
00724   return true;
00725 }
00726 
00727 
00728 void VoronoiDiagramGenerator::freeinit(Freelist *fl,int size)
00729 {
00730   fl->head = (Freenode *) NULL;
00731   fl->nodesize = size;
00732 }
00733 
00734 char * VoronoiDiagramGenerator::getfree(Freelist *fl)
00735 {
00736   int i; 
00737   Freenode *t;
00738 
00739   if(fl->head == (Freenode *) NULL)
00740     {   
00741       t =  (Freenode *) myalloc(sqrt_nsites * fl->nodesize);
00742 
00743       if(t == 0)
00744         return 0;
00745                 
00746       currentMemoryBlock->next = new FreeNodeArrayList;
00747       currentMemoryBlock = currentMemoryBlock->next;
00748       currentMemoryBlock->memory = t;
00749       currentMemoryBlock->next = 0;
00750 
00751       for(i=0; i<sqrt_nsites; i+=1)     
00752         makefree((Freenode *)((char *)t+i*fl->nodesize), fl);           
00753     };
00754   t = fl->head;
00755   fl->head = (fl->head)->nextfree;
00756   return((char *)t);
00757 }
00758 
00759 
00760 
00761 void VoronoiDiagramGenerator::makefree(Freenode *curr,Freelist *fl)
00762 {
00763   curr->nextfree = fl->head;
00764   fl->head = curr;
00765 }
00766 
00767 void VoronoiDiagramGenerator::cleanup()
00768 {
00769   if(sites != NULL){
00770     free(sites);
00771     sites = 0;
00772   }
00773 
00774   FreeNodeArrayList* current=NULL, *prev=NULL;
00775 
00776   current = prev = allMemoryList;
00777 
00778   while(current->next!=NULL){
00779     prev = current;
00780     current = current->next;
00781     free(prev->memory);
00782     delete prev;
00783     prev = NULL;
00784   }
00785 
00786   if(current!=NULL){
00787     if (current->memory!=NULL){
00788       free(current->memory);
00789     }
00790     delete current;
00791   }
00792 
00793   allMemoryList = new FreeNodeArrayList;
00794   allMemoryList->next = NULL;
00795   allMemoryList->memory = NULL;
00796   currentMemoryBlock = allMemoryList;
00797 
00798   if (ELhash!=NULL)
00799     free(ELhash);
00800 
00801   if (PQhash!=NULL)
00802     free(PQhash);
00803 }
00804 
00805 void VoronoiDiagramGenerator::cleanupEdges()
00806 {
00807   GraphEdge* geCurrent = NULL, *gePrev = NULL;
00808   geCurrent = gePrev = allEdges;
00809 
00810   while(geCurrent!=NULL){
00811     gePrev = geCurrent;
00812     geCurrent = geCurrent->next;
00813     delete gePrev;
00814   }
00815 
00816   allEdges = 0;
00817 
00818 }
00819 
00820 void VoronoiDiagramGenerator::pushGraphEdge(double x1, double y1, double x2, double y2,
00821                                             Site *s1, Site *s2)
00822 {
00823   GraphEdge* newEdge = new GraphEdge;
00824   newEdge->next = allEdges;
00825   allEdges = newEdge;
00826   newEdge->x1 = x1;
00827   newEdge->y1 = y1;
00828   newEdge->x2 = x2;
00829   newEdge->y2 = y2;
00830 
00831   newEdge->point1 = s1->sitenbr;
00832   newEdge->point2 = s2->sitenbr;
00833 }
00834 
00835 
00836 char * VoronoiDiagramGenerator::myalloc(unsigned n)
00837 {
00838   char *t=0;    
00839   t=(char*)malloc(n);
00840   total_alloc += n;
00841   return(t);
00842 }
00843 
00844 
00845 /* for those who don't have Cherry's plot */
00846 /* #include <plot.h> */
00847 void VoronoiDiagramGenerator::openpl(){}
00848 void VoronoiDiagramGenerator::circle(double x, double y, double radius){}
00849 void VoronoiDiagramGenerator::range(double minX, double minY, double maxX, double maxY){}
00850 
00851 
00852 
00853 void VoronoiDiagramGenerator::out_bisector(Edge *e)
00854 {
00855         
00856 
00857 }
00858 
00859 
00860 void VoronoiDiagramGenerator::out_ep(Edge *e)
00861 {
00862         
00863         
00864 }
00865 
00866 void VoronoiDiagramGenerator::out_vertex(Site *v)
00867 {
00868         
00869 }
00870 
00871 
00872 void VoronoiDiagramGenerator::out_site(Site *s)
00873 {
00874   if(!triangulate & plot & !debug)
00875     circle (s->coord.x, s->coord.y, cradius);
00876         
00877 }
00878 
00879 
00880 void VoronoiDiagramGenerator::out_triple(Site *s1, Site *s2,Site * s3)
00881 {
00882         
00883 }
00884 
00885 
00886 
00887 void VoronoiDiagramGenerator::plotinit()
00888 {
00889   double dx,dy,d;
00890         
00891   dy = ymax - ymin;
00892   dx = xmax - xmin;
00893   d = (double)(( dx > dy ? dx : dy) * 1.1);
00894   pxmin = (double)(xmin - (d-dx)/2.0);
00895   pxmax = (double)(xmax + (d-dx)/2.0);
00896   pymin = (double)(ymin - (d-dy)/2.0);
00897   pymax = (double)(ymax + (d-dy)/2.0);
00898   cradius = (double)((pxmax - pxmin)/350.0);
00899   openpl();
00900   range(pxmin, pymin, pxmax, pymax);
00901 }
00902 
00903 
00904 void VoronoiDiagramGenerator::clip_line(Edge *e)
00905 {
00906   Site *s1, *s2;
00907   double x1=0,x2=0,y1=0,y2=0; //, temp = 0;
00908 
00909   x1 = e->reg[0]->coord.x;
00910   x2 = e->reg[1]->coord.x;
00911   y1 = e->reg[0]->coord.y;
00912   y2 = e->reg[1]->coord.y;
00913 
00914   //if the distance between the two points this line was created from is less than 
00915   //the square root of 2, then ignore it
00916   //TODO improve/remove
00917   //if(sqrt(((x2 - x1) * (x2 - x1)) + ((y2 - y1) * (y2 - y1))) < minDistanceBetweenSites)
00918   //  {
00919   //    return;
00920   //  }
00921   pxmin = borderMinX;
00922   pxmax = borderMaxX;
00923   pymin = borderMinY;
00924   pymax = borderMaxY;
00925 
00926   if(e->a == 1.0 && e ->b >= 0.0)
00927     {   
00928       s1 = e->ep[1];
00929       s2 = e->ep[0];
00930     }
00931   else 
00932     {
00933       s1 = e->ep[0];
00934       s2 = e->ep[1];
00935     };
00936         
00937   if(e->a == 1.0)
00938     {
00939       y1 = pymin;
00940       if (s1!=(Site *)NULL && s1->coord.y > pymin)
00941         {
00942           y1 = s1->coord.y;
00943         }
00944       if(y1>pymax) 
00945         {
00946           //    printf("\nClipped (1) y1 = %f to %f",y1,pymax);
00947           y1 = pymax;
00948           //return;
00949         }
00950       x1 = e->c - e->b * y1;
00951       y2 = pymax;
00952       if (s2!=(Site *)NULL && s2->coord.y < pymax) 
00953         y2 = s2->coord.y;
00954 
00955       if(y2<pymin) 
00956         {
00957           //printf("\nClipped (2) y2 = %f to %f",y2,pymin);
00958           y2 = pymin;
00959           //return;
00960         }
00961       x2 = (e->c) - (e->b) * y2;
00962       if (((x1> pxmax) & (x2>pxmax)) | ((x1<pxmin)&(x2<pxmin))) 
00963         {
00964           //printf("\nClipLine jumping out(3), x1 = %f, pxmin = %f, pxmax = %f",x1,pxmin,pxmax);
00965           return;
00966         }
00967       if(x1> pxmax)
00968         {       x1 = pxmax; y1 = (e->c - x1)/e->b;};
00969       if(x1<pxmin)
00970         {       x1 = pxmin; y1 = (e->c - x1)/e->b;};
00971       if(x2>pxmax)
00972         {       x2 = pxmax; y2 = (e->c - x2)/e->b;};
00973       if(x2<pxmin)
00974         {       x2 = pxmin; y2 = (e->c - x2)/e->b;};
00975     }
00976   else
00977     {
00978       x1 = pxmin;
00979       if (s1!=(Site *)NULL && s1->coord.x > pxmin) 
00980         x1 = s1->coord.x;
00981       if(x1>pxmax) 
00982         {
00983           //printf("\nClipped (3) x1 = %f to %f",x1,pxmin);
00984           //return;
00985           x1 = pxmax;
00986         }
00987       y1 = e->c - e->a * x1;
00988       x2 = pxmax;
00989       if (s2!=(Site *)NULL && s2->coord.x < pxmax) 
00990         x2 = s2->coord.x;
00991       if(x2<pxmin) 
00992         {
00993           //printf("\nClipped (4) x2 = %f to %f",x2,pxmin);
00994           //return;
00995           x2 = pxmin;
00996         }
00997       y2 = e->c - e->a * x2;
00998       if (((y1> pymax) & (y2>pymax)) | ((y1<pymin)&(y2<pymin))) 
00999         {
01000           //printf("\nClipLine jumping out(6), y1 = %f, pymin = %f, pymax = %f",y2,pymin,pymax);
01001           return;
01002         }
01003       if(y1> pymax)
01004         {       y1 = pymax; x1 = (e->c - y1)/e->a;};
01005       if(y1<pymin)
01006         {       y1 = pymin; x1 = (e->c - y1)/e->a;};
01007       if(y2>pymax)
01008         {       y2 = pymax; x2 = (e->c - y2)/e->a;};
01009       if(y2<pymin)
01010         {       y2 = pymin; x2 = (e->c - y2)/e->a;};
01011     };
01012         
01013   //printf("\nPushing line (%f,%f,%f,%f)",x1,y1,x2,y2);
01014   //fprintf(stdout, "Line with vertices (%f,%f) and (%f,%f)\n", 
01015   //    e->reg[0]->coord.x, e->reg[1]->coord.x, e->reg[0]->coord.y, e->reg[1]->coord.y);
01016   pushGraphEdge(x1,y1,x2,y2,e->reg[0],e->reg[1]);
01017 }
01018 
01019 
01020 /* implicit parameters: nsites, sqrt_nsites, xmin, xmax, ymin, ymax,
01021    deltax, deltay (can all be estimates).
01022    Performance suffers if they are wrong; better to make nsites,
01023    deltax, and deltay too big than too small.  (?) */
01024 
01025 bool VoronoiDiagramGenerator::voronoi(int triangulate)
01026 {
01027   Site *newsite, *bot, *top, *temp, *p;
01028   Site *v;
01029   Point newintstar;
01030   int pm;
01031   Halfedge *lbnd, *rbnd, *llbnd, *rrbnd, *bisector;
01032   Edge *e;
01033         
01034   PQinitialize();
01035   bottomsite = nextone();
01036   out_site(bottomsite);
01037   bool retval = ELinitialize();
01038 
01039   if(!retval)
01040     return false;
01041         
01042   newsite = nextone();
01043   while(1)
01044     {
01045 
01046       if(!PQempty()) 
01047         newintstar = PQ_min();
01048                 
01049       //if the lowest site has a smaller y value than the lowest vector intersection, process the site
01050       //otherwise process the vector intersection               
01051       if (newsite != (Site *)NULL  && (PQempty() || newsite->coord.y < newintstar.y
01052                                        || (newsite->coord.y == newintstar.y && newsite->coord.x < newintstar.x)))
01053         {/* new site is smallest - this is a site event*/
01054           out_site(newsite);                                            //output the site
01055           lbnd = ELleftbnd(&(newsite->coord));                          //get the first HalfEdge to the LEFT of the new site
01056           rbnd = ELright(lbnd);                                         //get the first HalfEdge to the RIGHT of the new site
01057           bot = rightreg(lbnd);                                         //if this halfedge has no edge, , bot = bottom site (whatever that is)
01058           e = bisect(bot, newsite);                                     //create a new edge that bisects 
01059           bisector = HEcreate(e, le);                                   //create a new HalfEdge, setting its ELpm field to 0                    
01060           ELinsert(lbnd, bisector);                                     //insert this new bisector edge between the left and right vectors in a linked list     
01061             
01062           if ((p = intersect(lbnd, bisector)) != (Site *) NULL)         //if the new bisector intersects with the left edge, remove the left edge's vertex, and put in the new one
01063             {   
01064               PQdelete(lbnd);
01065               PQinsert(lbnd, p, dist(p,newsite));
01066             };
01067           lbnd = bisector;                                              
01068           bisector = HEcreate(e, re);                                   //create a new HalfEdge, setting its ELpm field to 1
01069           ELinsert(lbnd, bisector);                                     //insert the new HE to the right of the original bisector earlier in the IF stmt
01070             
01071           if ((p = intersect(bisector, rbnd)) != (Site *) NULL) //if this new bisector intersects with the
01072             {   
01073               PQinsert(bisector, p, dist(p,newsite));                   //push the HE into the ordered linked list of vertices
01074             };
01075           newsite = nextone();  
01076         }
01077       else if (!PQempty()) /* intersection is smallest - this is a vector event */                      
01078         {       
01079           lbnd = PQextractmin();                                                //pop the HalfEdge with the lowest vector off the ordered list of vectors                               
01080           llbnd = ELleft(lbnd);                                         //get the HalfEdge to the left of the above HE
01081           rbnd = ELright(lbnd);                                         //get the HalfEdge to the right of the above HE
01082           rrbnd = ELright(rbnd);                                                //get the HalfEdge to the right of the HE to the right of the lowest HE 
01083           bot = leftreg(lbnd);                                          //get the Site to the left of the left HE which it bisects
01084           top = rightreg(rbnd);                                         //get the Site to the right of the right HE which it bisects
01085             
01086           out_triple(bot, top, rightreg(lbnd));         //output the triple of sites, stating that a circle goes through them
01087             
01088           v = lbnd->vertex;                                             //get the vertex that caused this event
01089           makevertex(v);                                                        //set the vertex number - couldn't do this earlier since we didn't know when it would be processed
01090           endpoint(lbnd->ELedge,lbnd->ELpm,v);  //set the endpoint of the left HalfEdge to be this vector
01091           endpoint(rbnd->ELedge,rbnd->ELpm,v);  //set the endpoint of the right HalfEdge to be this vector
01092           ELdelete(lbnd);                                                       //mark the lowest HE for deletion - can't delete yet because there might be pointers to it in Hash Map  
01093           PQdelete(rbnd);                                                       //remove all vertex events to do with the  right HE
01094           ELdelete(rbnd);                                                       //mark the right HE for deletion - can't delete yet because there might be pointers to it in Hash Map   
01095           pm = le;                                                              //set the pm variable to zero
01096             
01097           if (bot->coord.y > top->coord.y)              //if the site to the left of the event is higher than the Site
01098             {                                                                           //to the right of it, then swap them and set the 'pm' variable to 1
01099               temp = bot; 
01100               bot = top; 
01101               top = temp; 
01102               pm = re;
01103             }
01104           e = bisect(bot, top);                                 //create an Edge (or line) that is between the two Sites. This creates
01105           //the formula of the line, and assigns a line number to it
01106           bisector = HEcreate(e, pm);                           //create a HE from the Edge 'e', and make it point to that edge with its ELedge field
01107           ELinsert(llbnd, bisector);                            //insert the new bisector to the right of the left HE
01108           endpoint(e, re-pm, v);                                        //set one endpoint to the new edge to be the vector point 'v'.
01109           //If the site to the left of this bisector is higher than the right
01110           //Site, then this endpoint is put in position 0; otherwise in pos 1
01111           deref(v);                                                             //delete the vector 'v'
01112             
01113           //if left HE and the new bisector don't intersect, then delete the left HE, and reinsert it 
01114           if((p = intersect(llbnd, bisector)) != (Site *) NULL)
01115             {   
01116               PQdelete(llbnd);
01117               PQinsert(llbnd, p, dist(p,bot));
01118             };
01119             
01120           //if right HE and the new bisector don't intersect, then reinsert it 
01121           if ((p = intersect(bisector, rrbnd)) != (Site *) NULL)
01122             {   
01123               PQinsert(bisector, p, dist(p,bot));
01124             };
01125         }
01126       else break;
01127     };
01128 
01129         
01130 
01131 
01132   for(lbnd=ELright(ELleftend); lbnd != ELrightend; lbnd=ELright(lbnd))
01133     {   
01134       e = lbnd->ELedge;
01135 
01136       clip_line(e);
01137     };
01138 
01139   //cleanup();
01140 
01141   return true;
01142         
01143 }
01144 
01145 
01146 int scomp(const void *p1,const void *p2)
01147 {
01148   Point *s1 = (Point*)p1, *s2=(Point*)p2;
01149   if(s1->y < s2->y) return(-1);
01150   if(s1->y > s2->y) return(1);
01151   if(s1->x < s2->x) return(-1);
01152   if(s1->x > s2->x) return(1);
01153   return(0);
01154 }
01155 
01156 /* return a single in-storage site */
01157 Site * VoronoiDiagramGenerator::nextone()
01158 {
01159   Site *s;
01160   if(siteidx < nsites)
01161     {   
01162       s = &sites[siteidx];
01163       siteidx += 1;
01164       return(s);
01165     }
01166   else  
01167     return( (Site *)NULL);
01168 }
01169 
01170 FASTJET_END_NAMESPACE
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines