Main Page | Namespace List | Class Hierarchy | Class List | Directories | File List | Namespace Members | Class Members | File Members

fastjet::VoronoiDiagramGenerator Class Reference

#include <Voronoi.hh>

Collaboration diagram for fastjet::VoronoiDiagramGenerator:

Collaboration graph
[legend]
List of all members.

Public Member Functions

 VoronoiDiagramGenerator ()
 ~VoronoiDiagramGenerator ()
bool generateVoronoi (std::vector< Point > *_parent_sites, double minX, double maxX, double minY, double maxY, double minDist=0)
void resetIterator ()
bool getNext (GraphEdge **e)

Public Attributes

std::vector< Point > * parent_sites
int n_parent_sites

Private Member Functions

void cleanup ()
void cleanupEdges ()
char * getfree (Freelist *fl)
HalfedgePQfind ()
int PQempty ()
HalfedgeHEcreate ()
Halfedge ** ELleft ()
Halfedge *** ELright ()
Halfedge **** ELleftbnd ()
HalfedgeHEcreate (Edge *e, int pm)
Point PQ_min ()
HalfedgePQextractmin ()
void freeinit (Freelist *fl, int size)
void makefree (Freenode *curr, Freelist *fl)
void geominit ()
void plotinit ()
bool voronoi (int triangulate)
void ref (Site *v)
void deref (Site *v)
void endpoint (Edge *e, int lr, Site *s)
void ELdelete (Halfedge *he)
HalfedgeELleftbnd (Point *p)
HalfedgeELright (Halfedge *he)
void makevertex (Site *v)
void out_triple (Site *s1, Site *s2, Site *s3)
void PQinsert (Halfedge *he, Site *v, double offset)
void PQdelete (Halfedge *he)
bool ELinitialize ()
void ELinsert (Halfedge *lb, Halfedge *newHe)
HalfedgeELgethash (int b)
HalfedgeELleft (Halfedge *he)
Siteleftreg (Halfedge *he)
void out_site (Site *s)
bool PQinitialize ()
int PQbucket (Halfedge *he)
void clip_line (Edge *e)
char * myalloc (unsigned n)
int right_of (Halfedge *el, Point *p)
Siterightreg (Halfedge *he)
Edgebisect (Site *s1, Site *s2)
double dist (Site *s, Site *t)
Siteintersect (Halfedge *el1, Halfedge *el2, Point *p=0)
void out_bisector (Edge *e)
void out_ep (Edge *e)
void out_vertex (Site *v)
Sitenextone ()
void pushGraphEdge (double x1, double y1, double x2, double y2, Site *s1, Site *s2)
void openpl ()
void circle (double x, double y, double radius)
void range (double minX, double minY, double maxX, double maxY)

Private Attributes

Halfedge ** ELhash
Freelist hfl
HalfedgeELleftend
HalfedgeELrightend
int ELhashsize
int triangulate
int sorted
int plot
int debug
double xmin
double xmax
double ymin
double ymax
double deltax
double deltay
Sitesites
int nsites
int siteidx
int sqrt_nsites
int nvertices
Freelist sfl
Sitebottomsite
int nedges
Freelist efl
int PQhashsize
HalfedgePQhash
int PQcount
int PQmin
int ntry
int totalsearch
double pxmin
double pxmax
double pymin
double pymax
double cradius
int total_alloc
double borderMinX
double borderMaxX
double borderMinY
double borderMaxY
FreeNodeArrayListallMemoryList
FreeNodeArrayListcurrentMemoryBlock
GraphEdgeallEdges
GraphEdgeiteratorEdges
double minDistanceBetweenSites

Constructor & Destructor Documentation

fastjet::VoronoiDiagramGenerator::VoronoiDiagramGenerator  ) 
 

Definition at line 64 of file Voronoi.cc.

References allEdges, allMemoryList, currentMemoryBlock, ELhash, iteratorEdges, fastjet::FreeNodeArrayList::memory, minDistanceBetweenSites, fastjet::FreeNodeArrayList::next, parent_sites, PQhash, siteidx, and sites.

00064                                                 {
00065   siteidx = 0;
00066   sites = NULL;
00067   parent_sites = NULL;
00068 
00069   allMemoryList = new FreeNodeArrayList;
00070   allMemoryList->memory = NULL;
00071   allMemoryList->next = NULL;
00072   currentMemoryBlock = allMemoryList;
00073   allEdges = NULL;
00074   iteratorEdges = 0;
00075   minDistanceBetweenSites = 0;
00076 
00077   ELhash = NULL;
00078   PQhash = NULL;
00079 }

fastjet::VoronoiDiagramGenerator::~VoronoiDiagramGenerator  ) 
 

Definition at line 81 of file Voronoi.cc.

References allMemoryList, cleanup(), and cleanupEdges().

00081                                                  {
00082   cleanup();
00083   cleanupEdges();
00084 
00085   if (allMemoryList != NULL)
00086     delete allMemoryList;
00087 }


Member Function Documentation

Edge * fastjet::VoronoiDiagramGenerator::bisect Site s1,
Site s2
[private]
 

Definition at line 338 of file Voronoi.cc.

References fastjet::Edge::a, fastjet::Edge::b, fastjet::Edge::c, fastjet::Edge::edgenbr, efl, fastjet::Edge::ep, getfree(), nedges, ref(), and fastjet::Edge::reg.

Referenced by voronoi().

00338                                                         {
00339   double dx,dy,adx,ady;
00340   Edge *newedge;        
00341 
00342   newedge = (Edge*) getfree(&efl);
00343         
00344   newedge->reg[0] = s1; //store the sites that this edge is bisecting
00345   newedge->reg[1] = s2;
00346   ref(s1); 
00347   ref(s2);
00348 
00349   // to begin with, there are no endpoints on the bisector - it goes
00350   // to infinity
00351   newedge->ep[0] = (Site*) NULL;
00352   newedge->ep[1] = (Site*) NULL;
00353   
00354   // get the difference in x dist between the sites
00355   dx = s2->coord.x - s1->coord.x;
00356   dy = s2->coord.y - s1->coord.y;
00357 
00358   // make sure that the difference is positive
00359   adx = dx>0 ? dx : -dx;
00360   ady = dy>0 ? dy : -dy;
00361 
00362   // get the slope of the line
00363   newedge->c = (double)(s1->coord.x * dx + s1->coord.y * dy
00364                         + (dx*dx + dy*dy)*0.5);
00365 
00366   if (adx>ady){ 
00367     //set formula of line, with x fixed to 1
00368     newedge->a = 1.0; newedge->b = dy/dx; newedge->c /= dx;
00369   } else {      
00370     //set formula of line, with y fixed to 1
00371     newedge->b = 1.0; newedge->a = dx/dy; newedge->c /= dy;
00372   }
00373         
00374   newedge->edgenbr = nedges;
00375   nedges++;
00376 
00377   return newedge;
00378 }

void fastjet::VoronoiDiagramGenerator::circle double  x,
double  y,
double  radius
[private]
 

Definition at line 741 of file Voronoi.cc.

Referenced by out_site().

00741 {}

void fastjet::VoronoiDiagramGenerator::cleanup  )  [private]
 

Definition at line 660 of file Voronoi.cc.

References allMemoryList, currentMemoryBlock, ELhash, fastjet::FreeNodeArrayList::memory, fastjet::FreeNodeArrayList::next, PQhash, and sites.

Referenced by generateVoronoi(), and ~VoronoiDiagramGenerator().

00661 {
00662   if(sites != NULL){
00663     free(sites);
00664     sites = 0;
00665   }
00666 
00667   FreeNodeArrayList* current=NULL, *prev=NULL;
00668 
00669   current = prev = allMemoryList;
00670 
00671   while(current->next!=NULL){
00672     prev = current;
00673     current = current->next;
00674     free(prev->memory);
00675     delete prev;
00676     prev = NULL;
00677   }
00678 
00679   if(current!=NULL){
00680     if (current->memory!=NULL){
00681       free(current->memory);
00682     }
00683     delete current;
00684   }
00685 
00686   allMemoryList = new FreeNodeArrayList;
00687   allMemoryList->next = NULL;
00688   allMemoryList->memory = NULL;
00689   currentMemoryBlock = allMemoryList;
00690 
00691   if (ELhash!=NULL)
00692     free(ELhash);
00693 
00694   if (PQhash!=NULL)
00695     free(PQhash);
00696 }

void fastjet::VoronoiDiagramGenerator::cleanupEdges  )  [private]
 

Definition at line 698 of file Voronoi.cc.

References allEdges, and fastjet::GraphEdge::next.

Referenced by generateVoronoi(), and ~VoronoiDiagramGenerator().

00699 {
00700   GraphEdge* geCurrent = NULL, *gePrev = NULL;
00701   geCurrent = gePrev = allEdges;
00702 
00703   while(geCurrent!=NULL){
00704     gePrev = geCurrent;
00705     geCurrent = geCurrent->next;
00706     delete gePrev;
00707   }
00708 
00709   allEdges = 0;
00710 
00711 }

void fastjet::VoronoiDiagramGenerator::clip_line Edge e  )  [private]
 

Definition at line 797 of file Voronoi.cc.

References fastjet::Edge::a, fastjet::Edge::b, borderMaxX, borderMaxY, borderMinX, borderMinY, fastjet::Edge::c, fastjet::Edge::ep, pushGraphEdge(), pxmax, pxmin, pymax, pymin, and fastjet::Edge::reg.

Referenced by endpoint(), and voronoi().

00798 {
00799   Site *s1, *s2;
00800   double x1=0,x2=0,y1=0,y2=0; //, temp = 0;
00801 
00802   x1 = e->reg[0]->coord.x;
00803   x2 = e->reg[1]->coord.x;
00804   y1 = e->reg[0]->coord.y;
00805   y2 = e->reg[1]->coord.y;
00806 
00807   //if the distance between the two points this line was created from is less than 
00808   //the square root of 2, then ignore it
00809   //TODO improve/remove
00810   //if(sqrt(((x2 - x1) * (x2 - x1)) + ((y2 - y1) * (y2 - y1))) < minDistanceBetweenSites)
00811   //  {
00812   //    return;
00813   //  }
00814   pxmin = borderMinX;
00815   pxmax = borderMaxX;
00816   pymin = borderMinY;
00817   pymax = borderMaxY;
00818 
00819   if(e->a == 1.0 && e ->b >= 0.0)
00820     {   
00821       s1 = e->ep[1];
00822       s2 = e->ep[0];
00823     }
00824   else 
00825     {
00826       s1 = e->ep[0];
00827       s2 = e->ep[1];
00828     };
00829         
00830   if(e->a == 1.0)
00831     {
00832       y1 = pymin;
00833       if (s1!=(Site *)NULL && s1->coord.y > pymin)
00834         {
00835           y1 = s1->coord.y;
00836         }
00837       if(y1>pymax) 
00838         {
00839           //    printf("\nClipped (1) y1 = %f to %f",y1,pymax);
00840           y1 = pymax;
00841           //return;
00842         }
00843       x1 = e->c - e->b * y1;
00844       y2 = pymax;
00845       if (s2!=(Site *)NULL && s2->coord.y < pymax) 
00846         y2 = s2->coord.y;
00847 
00848       if(y2<pymin) 
00849         {
00850           //printf("\nClipped (2) y2 = %f to %f",y2,pymin);
00851           y2 = pymin;
00852           //return;
00853         }
00854       x2 = (e->c) - (e->b) * y2;
00855       if (((x1> pxmax) & (x2>pxmax)) | ((x1<pxmin)&(x2<pxmin))) 
00856         {
00857           //printf("\nClipLine jumping out(3), x1 = %f, pxmin = %f, pxmax = %f",x1,pxmin,pxmax);
00858           return;
00859         }
00860       if(x1> pxmax)
00861         {       x1 = pxmax; y1 = (e->c - x1)/e->b;};
00862       if(x1<pxmin)
00863         {       x1 = pxmin; y1 = (e->c - x1)/e->b;};
00864       if(x2>pxmax)
00865         {       x2 = pxmax; y2 = (e->c - x2)/e->b;};
00866       if(x2<pxmin)
00867         {       x2 = pxmin; y2 = (e->c - x2)/e->b;};
00868     }
00869   else
00870     {
00871       x1 = pxmin;
00872       if (s1!=(Site *)NULL && s1->coord.x > pxmin) 
00873         x1 = s1->coord.x;
00874       if(x1>pxmax) 
00875         {
00876           //printf("\nClipped (3) x1 = %f to %f",x1,pxmin);
00877           //return;
00878           x1 = pxmax;
00879         }
00880       y1 = e->c - e->a * x1;
00881       x2 = pxmax;
00882       if (s2!=(Site *)NULL && s2->coord.x < pxmax) 
00883         x2 = s2->coord.x;
00884       if(x2<pxmin) 
00885         {
00886           //printf("\nClipped (4) x2 = %f to %f",x2,pxmin);
00887           //return;
00888           x2 = pxmin;
00889         }
00890       y2 = e->c - e->a * x2;
00891       if (((y1> pymax) & (y2>pymax)) | ((y1<pymin)&(y2<pymin))) 
00892         {
00893           //printf("\nClipLine jumping out(6), y1 = %f, pymin = %f, pymax = %f",y2,pymin,pymax);
00894           return;
00895         }
00896       if(y1> pymax)
00897         {       y1 = pymax; x1 = (e->c - y1)/e->a;};
00898       if(y1<pymin)
00899         {       y1 = pymin; x1 = (e->c - y1)/e->a;};
00900       if(y2>pymax)
00901         {       y2 = pymax; x2 = (e->c - y2)/e->a;};
00902       if(y2<pymin)
00903         {       y2 = pymin; x2 = (e->c - y2)/e->a;};
00904     };
00905         
00906   //printf("\nPushing line (%f,%f,%f,%f)",x1,y1,x2,y2);
00907   //fprintf(stdout, "Line with vertices (%f,%f) and (%f,%f)\n", 
00908   //    e->reg[0]->coord.x, e->reg[1]->coord.x, e->reg[0]->coord.y, e->reg[1]->coord.y);
00909   pushGraphEdge(x1,y1,x2,y2,e->reg[0],e->reg[1]);
00910 }

void fastjet::VoronoiDiagramGenerator::deref Site v  )  [private]
 

Definition at line 512 of file Voronoi.cc.

References makefree(), and sfl.

Referenced by endpoint(), PQdelete(), and voronoi().

00513 {
00514   v->refcnt -= 1;
00515   if (v->refcnt == 0 ) 
00516     makefree((Freenode*)v, &sfl);
00517 }

double fastjet::VoronoiDiagramGenerator::dist Site s,
Site t
[private]
 

Definition at line 495 of file Voronoi.cc.

Referenced by voronoi().

00496 {
00497   double dx,dy;
00498   dx = s->coord.x - t->coord.x;
00499   dy = s->coord.y - t->coord.y;
00500   return (double)(sqrt(dx*dx + dy*dy));
00501 }

void fastjet::VoronoiDiagramGenerator::ELdelete Halfedge he  )  [private]
 

Definition at line 287 of file Voronoi.cc.

References DELETED, fastjet::Halfedge::ELedge, fastjet::Halfedge::ELleft, and fastjet::Halfedge::ELright.

Referenced by voronoi().

00287                                                   {
00288   (he->ELleft)->ELright = he->ELright;
00289   (he->ELright)->ELleft = he->ELleft;
00290   he->ELedge = (Edge*) DELETED;
00291 }

Halfedge * fastjet::VoronoiDiagramGenerator::ELgethash int  b  )  [private]
 

Definition at line 212 of file Voronoi.cc.

References DELETED, fastjet::Halfedge::ELedge, ELhash, ELhashsize, fastjet::Halfedge::ELrefcnt, hfl, and makefree().

Referenced by ELleftbnd().

00212                                                  {
00213   Halfedge *he;
00214         
00215   if ((b<0) || (b>=ELhashsize)) 
00216     return (Halfedge *) NULL;
00217 
00218   he = ELhash[b]; 
00219   if ((he==(Halfedge*) NULL) || (he->ELedge!=(Edge*) DELETED)) 
00220     return he;
00221         
00222   /* Hash table points to deleted half edge.  Patch as necessary. */
00223   ELhash[b] = (Halfedge*) NULL;
00224   if ((he->ELrefcnt -= 1) == 0) 
00225     makefree((Freenode*)he, &hfl);
00226   return (Halfedge*) NULL;
00227 }       

bool fastjet::VoronoiDiagramGenerator::ELinitialize  )  [private]
 

Definition at line 166 of file Voronoi.cc.

References ELhash, ELhashsize, fastjet::Halfedge::ELleft, ELleftend, fastjet::Halfedge::ELright, ELrightend, freeinit(), HEcreate(), hfl, myalloc(), and sqrt_nsites.

Referenced by voronoi().

00166                                           {
00167   int i;
00168   freeinit(&hfl, sizeof(Halfedge));
00169   ELhashsize = 2 * sqrt_nsites;
00170   //ELhash = (Halfedge **) myalloc ( sizeof *ELhash * ELhashsize);
00171   ELhash = (Halfedge **) myalloc ( sizeof(Halfedge*)*ELhashsize);
00172   
00173   if(ELhash == 0)
00174     return false;
00175   
00176   for(i=0; i<ELhashsize; i +=1) ELhash[i] = (Halfedge *)NULL;
00177   ELleftend = HEcreate((Edge*) NULL, 0);
00178   ELrightend = HEcreate((Edge*) NULL, 0);
00179   ELleftend->ELleft = (Halfedge*) NULL;
00180   ELleftend->ELright = ELrightend;
00181   ELrightend->ELleft = ELleftend;
00182   ELrightend->ELright = (Halfedge*) NULL;
00183   ELhash[0] = ELleftend;
00184   ELhash[ELhashsize-1] = ELrightend;
00185 
00186   return true;
00187 }

void fastjet::VoronoiDiagramGenerator::ELinsert Halfedge lb,
Halfedge newHe
[private]
 

Definition at line 202 of file Voronoi.cc.

References fastjet::Halfedge::ELleft, and fastjet::Halfedge::ELright.

Referenced by voronoi().

00203 {
00204   newHe->ELleft = lb;
00205   newHe->ELright = lb->ELright;
00206   (lb->ELright)->ELleft = newHe;
00207   lb->ELright = newHe;
00208 }

Halfedge * fastjet::VoronoiDiagramGenerator::ELleft Halfedge he  )  [private]
 

Definition at line 299 of file Voronoi.cc.

References fastjet::Halfedge::ELleft.

00299                                                      {
00300   return he->ELleft;
00301 }

Halfedge* * fastjet::VoronoiDiagramGenerator::ELleft  )  [private]
 

Referenced by voronoi().

Halfedge * fastjet::VoronoiDiagramGenerator::ELleftbnd Point p  )  [private]
 

Definition at line 229 of file Voronoi.cc.

References deltax, ELgethash(), ELhash, ELhashsize, fastjet::Halfedge::ELleft, ELleftend, fastjet::Halfedge::ELrefcnt, fastjet::Halfedge::ELright, ELrightend, ntry, right_of(), totalsearch, and xmin.

00229                                                      {
00230   int i, bucket;
00231   Halfedge *he;
00232         
00233   /* Use hash table to get close to desired halfedge */
00234   // use the hash function to find the place in the hash map that this
00235   // HalfEdge should be
00236   bucket = (int)((p->x - xmin)/deltax * ELhashsize);    
00237 
00238   // make sure that the bucket position in within the range of the
00239   //hash array
00240   if (bucket<0) bucket =0;
00241   if (bucket>=ELhashsize) bucket = ELhashsize - 1;
00242 
00243   he = ELgethash(bucket);
00244 
00245   // if the HE isn't found, search backwards and forwards in the hash
00246   // map for the first non-null entry
00247   if (he == (Halfedge*) NULL){   
00248     for(i=1;1;i++){     
00249       if ((he=ELgethash(bucket-i)) != (Halfedge*) NULL) 
00250         break;
00251       if ((he=ELgethash(bucket+i)) != (Halfedge*) NULL) 
00252         break;
00253     };
00254     totalsearch += i;
00255   };
00256   ntry += 1;
00257   /* Now search linear list of halfedges for the correct one */
00258   if ((he==ELleftend)  || (he != ELrightend && right_of(he,p))){
00259     do{
00260       he = he->ELright;
00261     } while (he!=ELrightend && right_of(he,p));
00262       // keep going right on the list until either the end is reached,
00263       // or you find the 1st edge which the point
00264       he = he->ELleft;  //isn't to the right of
00265   } else 
00266     // if the point is to the left of the HalfEdge, then search left
00267     // for the HE just to the left of the point
00268     do{
00269       he = he->ELleft;
00270     } while ((he!=ELleftend )&& (!right_of(he,p)));
00271                 
00272   /* Update hash table and reference counts */
00273   if ((bucket > 0) && (bucket <ELhashsize-1)){  
00274     if(ELhash[bucket] != (Halfedge *) NULL){
00275       ELhash[bucket]->ELrefcnt -= 1;
00276     }
00277     ELhash[bucket] = he;
00278     ELhash[bucket]->ELrefcnt += 1;
00279   };
00280 
00281   return he;
00282 }

Halfedge* * * * fastjet::VoronoiDiagramGenerator::ELleftbnd  )  [private]
 

Referenced by voronoi().

Halfedge * fastjet::VoronoiDiagramGenerator::ELright Halfedge he  )  [private]
 

Definition at line 294 of file Voronoi.cc.

References fastjet::Halfedge::ELright.

00294                                                       {
00295   return he->ELright;
00296 }

Halfedge* * * fastjet::VoronoiDiagramGenerator::ELright  )  [private]
 

Referenced by voronoi().

void fastjet::VoronoiDiagramGenerator::endpoint Edge e,
int  lr,
Site s
[private]
 

Definition at line 480 of file Voronoi.cc.

References clip_line(), deref(), efl, le, makefree(), re, and ref().

Referenced by voronoi().

00481 {
00482   e->ep[lr] = s;
00483   ref(s);
00484   if(e->ep[re-lr]== (Site *) NULL) 
00485     return;
00486 
00487   clip_line(e);
00488 
00489   deref(e->reg[le]);
00490   deref(e->reg[re]);
00491   makefree((Freenode*)e, &efl);
00492 }

void fastjet::VoronoiDiagramGenerator::freeinit Freelist fl,
int  size
[private]
 

Definition at line 621 of file Voronoi.cc.

References fastjet::Freelist::head, and fastjet::Freelist::nodesize.

Referenced by ELinitialize(), generateVoronoi(), and geominit().

00622 {
00623   fl->head = (Freenode *) NULL;
00624   fl->nodesize = size;
00625 }

bool fastjet::VoronoiDiagramGenerator::generateVoronoi std::vector< Point > *  _parent_sites,
double  minX,
double  maxX,
double  minY,
double  maxY,
double  minDist = 0
 

Definition at line 91 of file Voronoi.cc.

References borderMaxX, borderMaxY, borderMinX, borderMinY, cleanup(), cleanupEdges(), debug, freeinit(), geominit(), minDistanceBetweenSites, myalloc(), n_parent_sites, nsites, parent_sites, plot, fastjet::scomp(), sfl, siteidx, sites, sorted, triangulate, voronoi(), xmax, xmin, ymax, and ymin.

Referenced by fastjet::ClusterSequenceVoronoiArea::VoronoiAreaCalc::VoronoiAreaCalc().

00094                                                              {
00095   cleanup();
00096   cleanupEdges();
00097   int i;
00098   double x, y;
00099 
00100   minDistanceBetweenSites = minDist;
00101 
00102   parent_sites = _parent_sites;
00103 
00104   nsites = n_parent_sites = parent_sites->size();
00105   plot = 0;
00106   triangulate = 0;      
00107   debug = 1;
00108   sorted = 0; 
00109   freeinit(&sfl, sizeof (Site));
00110                 
00111   //sites = (Site *) myalloc(nsites*sizeof( *sites));
00112   sites = (Site *) myalloc(nsites*sizeof( Site));
00113   //parent_sites = (Site *) myalloc(nsites*sizeof( Site));
00114  
00115   if (sites == 0)
00116     return false;
00117 
00118   xmax = xmin = (*parent_sites)[0].x;
00119   ymax = ymin = (*parent_sites)[0].y;
00120   
00121   for(i=0;i<nsites;i++){
00122     x = (*parent_sites)[i].x;
00123     y = (*parent_sites)[i].y;
00124     sites[i].coord.x = x;
00125     sites[i].coord.y = y;
00126     sites[i].sitenbr = i;
00127     sites[i].refcnt  = 0;
00128     
00129     if (x<xmin)
00130       xmin=x;
00131     else if (x>xmax)
00132       xmax=x;
00133     
00134     if (y<ymin)
00135       ymin=y;
00136     else if (y>ymax)
00137       ymax=y;
00138   }
00139         
00140   qsort(sites, nsites, sizeof (*sites), scomp);
00141         
00142   siteidx = 0;
00143   geominit();
00144   double temp = 0;
00145   if(minX > maxX){
00146     temp = minX;
00147     minX = maxX;
00148     maxX = temp;
00149   }
00150   if(minY > maxY){
00151     temp = minY;
00152     minY = maxY;
00153     maxY = temp;
00154   }
00155   borderMinX = minX;
00156   borderMinY = minY;
00157   borderMaxX = maxX;
00158   borderMaxY = maxY;
00159         
00160   siteidx = 0;
00161   voronoi(triangulate);
00162 
00163   return true;
00164 }

void fastjet::VoronoiDiagramGenerator::geominit  )  [private]
 

Definition at line 325 of file Voronoi.cc.

References deltax, deltay, efl, freeinit(), nedges, nsites, nvertices, sqrt_nsites, xmax, xmin, ymax, and ymin.

Referenced by generateVoronoi().

00325                                       {
00326   double sn;
00327   
00328   freeinit(&efl, sizeof(Edge));
00329   nvertices = 0;
00330   nedges = 0;
00331   sn = (double)nsites+4;
00332   sqrt_nsites = (int)sqrt(sn);
00333   deltay = ymax - ymin;
00334   deltax = xmax - xmin;
00335 }

char * fastjet::VoronoiDiagramGenerator::getfree Freelist fl  )  [private]
 

Definition at line 627 of file Voronoi.cc.

References currentMemoryBlock, fastjet::Freelist::head, makefree(), fastjet::FreeNodeArrayList::memory, myalloc(), fastjet::FreeNodeArrayList::next, fastjet::Freenode::nextfree, fastjet::Freelist::nodesize, and sqrt_nsites.

Referenced by bisect(), HEcreate(), and intersect().

00628 {
00629   int i; 
00630   Freenode *t;
00631 
00632   if(fl->head == (Freenode *) NULL)
00633     {   
00634       t =  (Freenode *) myalloc(sqrt_nsites * fl->nodesize);
00635 
00636       if(t == 0)
00637         return 0;
00638                 
00639       currentMemoryBlock->next = new FreeNodeArrayList;
00640       currentMemoryBlock = currentMemoryBlock->next;
00641       currentMemoryBlock->memory = t;
00642       currentMemoryBlock->next = 0;
00643 
00644       for(i=0; i<sqrt_nsites; i+=1)     
00645         makefree((Freenode *)((char *)t+i*fl->nodesize), fl);           
00646     };
00647   t = fl->head;
00648   fl->head = (fl->head)->nextfree;
00649   return((char *)t);
00650 }

bool fastjet::VoronoiDiagramGenerator::getNext GraphEdge **  e  )  [inline]
 

Definition at line 207 of file Voronoi.hh.

References fastjet::GraphEdge::next.

Referenced by fastjet::ClusterSequenceVoronoiArea::VoronoiAreaCalc::VoronoiAreaCalc().

00207                              {
00208     if(iteratorEdges == 0)
00209       return false;
00210     
00211     *e = iteratorEdges;
00212     iteratorEdges = iteratorEdges->next;
00213     return true;
00214   }

Halfedge * fastjet::VoronoiDiagramGenerator::HEcreate Edge e,
int  pm
[private]
 

Definition at line 190 of file Voronoi.cc.

References fastjet::Halfedge::ELedge, fastjet::Halfedge::ELpm, fastjet::Halfedge::ELrefcnt, getfree(), hfl, fastjet::Halfedge::PQnext, and fastjet::Halfedge::vertex.

00190                                                          {
00191   Halfedge *answer;
00192   answer = (Halfedge*) getfree(&hfl);
00193   answer->ELedge = e;
00194   answer->ELpm = pm;
00195   answer->PQnext = (Halfedge*) NULL;
00196   answer->vertex = (Site*) NULL;
00197   answer->ELrefcnt = 0;
00198   return answer;
00199 }

Halfedge* fastjet::VoronoiDiagramGenerator::HEcreate  )  [private]
 

Referenced by ELinitialize(), and voronoi().

Site * fastjet::VoronoiDiagramGenerator::intersect Halfedge el1,
Halfedge el2,
Point p = 0
[private]
 

Definition at line 384 of file Voronoi.cc.

References fastjet::Edge::a, fastjet::Edge::b, fastjet::Edge::c, fastjet::Halfedge::ELedge, fastjet::Halfedge::ELpm, getfree(), le, re, fastjet::Edge::reg, and sfl.

Referenced by voronoi().

00384                                                                               {
00385   Edge *e1,*e2, *e;
00386   Halfedge *el;
00387   double d, xint, yint;
00388   int right_of_site;
00389   Site *v;
00390         
00391   e1 = el1->ELedge;
00392   e2 = el2->ELedge;
00393   if ((e1==(Edge*)NULL) || (e2==(Edge*)NULL)) 
00394     return ((Site*) NULL);
00395 
00396   // if the two edges bisect the same parent, return null
00397   if (e1->reg[1] == e2->reg[1]) 
00398     return (Site*) NULL;
00399         
00400   d = e1->a * e2->b - e1->b * e2->a;
00401   if (-1.0e-10<d && d<1.0e-10) 
00402     return (Site*) NULL;
00403         
00404   xint = (e1->c*e2->b - e2->c*e1->b)/d;
00405   yint = (e2->c*e1->a - e1->c*e2->a)/d;
00406 
00407   volatile double local_y1 = e1->reg[1]->coord.y;
00408   volatile double local_y2 = e2->reg[1]->coord.y;
00409         
00410   if( (local_y1 < local_y2) ||
00411       ((local_y1 == local_y2) &&
00412        (e1->reg[1]->coord.x <  e2->reg[1]->coord.x)) ){ 
00413     el = el1; 
00414     e = e1;
00415   } else {      
00416     el = el2; 
00417     e = e2;
00418   }
00419         
00420   right_of_site = xint >= e->reg[1]->coord.x;
00421   if ((right_of_site && el->ELpm == le) || (!right_of_site && el->ELpm == re)) 
00422     return (Site*) NULL;
00423         
00424   // create a new site at the point of intersection - this is a new
00425   // vector event waiting to happen
00426   v = (Site*) getfree(&sfl);
00427   v->refcnt = 0;
00428   v->coord.x = xint;
00429   v->coord.y = yint;
00430   return v;
00431 }

Site * fastjet::VoronoiDiagramGenerator::leftreg Halfedge he  )  [private]
 

Definition at line 304 of file Voronoi.cc.

References bottomsite, fastjet::Halfedge::ELedge, fastjet::Halfedge::ELpm, le, re, and fastjet::Edge::reg.

Referenced by voronoi().

00304                                                    {
00305   if (he->ELedge == (Edge*) NULL) 
00306     return(bottomsite);
00307   return (he->ELpm == le) 
00308     ? he->ELedge->reg[le] 
00309     : he->ELedge->reg[re];
00310 }

void fastjet::VoronoiDiagramGenerator::makefree Freenode curr,
Freelist fl
[private]
 

Definition at line 654 of file Voronoi.cc.

References fastjet::Freelist::head, and fastjet::Freenode::nextfree.

Referenced by deref(), ELgethash(), endpoint(), and getfree().

00655 {
00656   curr->nextfree = fl->head;
00657   fl->head = curr;
00658 }

void fastjet::VoronoiDiagramGenerator::makevertex Site v  )  [private]
 

Definition at line 504 of file Voronoi.cc.

References nvertices, and out_vertex().

Referenced by voronoi().

00505 {
00506   v->sitenbr = nvertices;
00507   nvertices += 1;
00508   out_vertex(v);
00509 }

char * fastjet::VoronoiDiagramGenerator::myalloc unsigned  n  )  [private]
 

Definition at line 729 of file Voronoi.cc.

References total_alloc.

Referenced by ELinitialize(), generateVoronoi(), getfree(), and PQinitialize().

00730 {
00731   char *t=0;    
00732   t=(char*)malloc(n);
00733   total_alloc += n;
00734   return(t);
00735 }

Site * fastjet::VoronoiDiagramGenerator::nextone  )  [private]
 

Definition at line 1053 of file Voronoi.cc.

References nsites, siteidx, and sites.

Referenced by voronoi().

01054 {
01055   Site *s;
01056   if(siteidx < nsites)
01057     {   
01058       s = &sites[siteidx];
01059       siteidx += 1;
01060       return(s);
01061     }
01062   else  
01063     return( (Site *)NULL);
01064 }

void fastjet::VoronoiDiagramGenerator::openpl  )  [private]
 

Definition at line 740 of file Voronoi.cc.

Referenced by plotinit().

00740 {}

void fastjet::VoronoiDiagramGenerator::out_bisector Edge e  )  [private]
 

Definition at line 746 of file Voronoi.cc.

00747 {
00748         
00749 
00750 }

void fastjet::VoronoiDiagramGenerator::out_ep Edge e  )  [private]
 

Definition at line 753 of file Voronoi.cc.

00754 {
00755         
00756         
00757 }

void fastjet::VoronoiDiagramGenerator::out_site Site s  )  [private]
 

Definition at line 765 of file Voronoi.cc.

References circle(), cradius, debug, plot, and triangulate.

Referenced by voronoi().

00766 {
00767   if(!triangulate & plot & !debug)
00768     circle (s->coord.x, s->coord.y, cradius);
00769         
00770 }

void fastjet::VoronoiDiagramGenerator::out_triple Site s1,
Site s2,
Site s3
[private]
 

Definition at line 773 of file Voronoi.cc.

Referenced by voronoi().

00774 {
00775         
00776 }

void fastjet::VoronoiDiagramGenerator::out_vertex Site v  )  [private]
 

Definition at line 759 of file Voronoi.cc.

Referenced by makevertex().

00760 {
00761         
00762 }

void fastjet::VoronoiDiagramGenerator::plotinit  )  [private]
 

Definition at line 780 of file Voronoi.cc.

References cradius, openpl(), pxmax, pxmin, pymax, pymin, range(), xmax, xmin, ymax, and ymin.

00781 {
00782   double dx,dy,d;
00783         
00784   dy = ymax - ymin;
00785   dx = xmax - xmin;
00786   d = (double)(( dx > dy ? dx : dy) * 1.1);
00787   pxmin = (double)(xmin - (d-dx)/2.0);
00788   pxmax = (double)(xmax + (d-dx)/2.0);
00789   pymin = (double)(ymin - (d-dy)/2.0);
00790   pymax = (double)(ymax + (d-dy)/2.0);
00791   cradius = (double)((pxmax - pxmin)/350.0);
00792   openpl();
00793   range(pxmin, pymin, pxmax, pymax);
00794 }

Point fastjet::VoronoiDiagramGenerator::PQ_min  )  [private]
 

Definition at line 581 of file Voronoi.cc.

References PQhash, PQmin, fastjet::Halfedge::PQnext, fastjet::Halfedge::vertex, fastjet::Point::x, fastjet::Point::y, and fastjet::Halfedge::ystar.

Referenced by voronoi().

00582 {
00583   Point answer;
00584         
00585   while(PQhash[PQmin].PQnext == (Halfedge *)NULL) {PQmin += 1;};
00586   answer.x = PQhash[PQmin].PQnext->vertex->coord.x;
00587   answer.y = PQhash[PQmin].PQnext->ystar;
00588   return (answer);
00589 }

int fastjet::VoronoiDiagramGenerator::PQbucket Halfedge he  )  [private]
 

Definition at line 562 of file Voronoi.cc.

References deltay, PQhashsize, PQmin, ymin, and fastjet::Halfedge::ystar.

Referenced by PQdelete(), and PQinsert().

00563 {
00564   int bucket;
00565         
00566   bucket = (int)((he->ystar - ymin)/deltay * PQhashsize);
00567   if (bucket<0) bucket = 0;
00568   if (bucket>=PQhashsize) bucket = PQhashsize-1 ;
00569   if (bucket < PQmin) PQmin = bucket;
00570   return(bucket);
00571 }

void fastjet::VoronoiDiagramGenerator::PQdelete Halfedge he  )  [private]
 

Definition at line 545 of file Voronoi.cc.

References deref(), PQbucket(), PQcount, PQhash, fastjet::Halfedge::PQnext, and fastjet::Halfedge::vertex.

Referenced by voronoi().

00546 {
00547   Halfedge *last;
00548         
00549   if(he->vertex != (Site *) NULL)
00550     {   
00551       last = &PQhash[PQbucket(he)];
00552       while (last->PQnext != he) 
00553         last = last->PQnext;
00554 
00555       last->PQnext = he->PQnext;
00556       PQcount -= 1;
00557       deref(he->vertex);
00558       he->vertex = (Site *) NULL;
00559     };
00560 }

int fastjet::VoronoiDiagramGenerator::PQempty  )  [private]
 

Definition at line 575 of file Voronoi.cc.

References PQcount.

Referenced by voronoi().

00576 {
00577   return(PQcount==0);
00578 }

Halfedge * fastjet::VoronoiDiagramGenerator::PQextractmin  )  [private]
 

Definition at line 591 of file Voronoi.cc.

References PQcount, PQhash, PQmin, and fastjet::Halfedge::PQnext.

Referenced by voronoi().

00592 {
00593   Halfedge *curr;
00594         
00595   curr = PQhash[PQmin].PQnext;
00596   PQhash[PQmin].PQnext = curr->PQnext;
00597   PQcount -= 1;
00598   return(curr);
00599 }

Halfedge* fastjet::VoronoiDiagramGenerator::PQfind  )  [private]
 

bool fastjet::VoronoiDiagramGenerator::PQinitialize  )  [private]
 

Definition at line 602 of file Voronoi.cc.

References myalloc(), PQcount, PQhash, PQhashsize, PQmin, fastjet::Halfedge::PQnext, and sqrt_nsites.

Referenced by voronoi().

00603 {
00604   int i; 
00605         
00606   PQcount = 0;
00607   PQmin = 0;
00608   PQhashsize = 4 * sqrt_nsites;
00609   //PQhash = (Halfedge *) myalloc(PQhashsize * sizeof *PQhash);
00610   PQhash = (Halfedge *) myalloc(PQhashsize * sizeof(Halfedge));
00611 
00612   if(PQhash == 0)
00613     return false;
00614 
00615   for(i=0; i<PQhashsize; i+=1) PQhash[i].PQnext = (Halfedge *)NULL;
00616 
00617   return true;
00618 }

void fastjet::VoronoiDiagramGenerator::PQinsert Halfedge he,
Site v,
double  offset
[private]
 

Definition at line 525 of file Voronoi.cc.

References PQbucket(), PQcount, PQhash, fastjet::Halfedge::PQnext, ref(), fastjet::Halfedge::vertex, and fastjet::Halfedge::ystar.

Referenced by voronoi().

00526 {
00527   Halfedge *last, *next;
00528         
00529   he->vertex = v;
00530   ref(v);
00531   he->ystar = (double)(v->coord.y + offset);
00532   last = &PQhash[PQbucket(he)];
00533   while ((next = last->PQnext) != (Halfedge *) NULL &&
00534          (he->ystar  > next->ystar  ||
00535           (he->ystar == next->ystar && v->coord.x > next->vertex->coord.x)))
00536     {   
00537       last = next;
00538     };
00539   he->PQnext = last->PQnext; 
00540   last->PQnext = he;
00541   PQcount += 1;
00542 }

void fastjet::VoronoiDiagramGenerator::pushGraphEdge double  x1,
double  y1,
double  x2,
double  y2,
Site s1,
Site s2
[private]
 

Definition at line 713 of file Voronoi.cc.

References allEdges, fastjet::GraphEdge::next, fastjet::GraphEdge::point1, fastjet::GraphEdge::point2, fastjet::GraphEdge::x1, fastjet::GraphEdge::x2, fastjet::GraphEdge::y1, and fastjet::GraphEdge::y2.

Referenced by clip_line().

00715 {
00716   GraphEdge* newEdge = new GraphEdge;
00717   newEdge->next = allEdges;
00718   allEdges = newEdge;
00719   newEdge->x1 = x1;
00720   newEdge->y1 = y1;
00721   newEdge->x2 = x2;
00722   newEdge->y2 = y2;
00723 
00724   newEdge->point1 = s1->sitenbr;
00725   newEdge->point2 = s2->sitenbr;
00726 }

void fastjet::VoronoiDiagramGenerator::range double  minX,
double  minY,
double  maxX,
double  maxY
[private]
 

Definition at line 742 of file Voronoi.cc.

Referenced by plotinit().

00742 {}

void fastjet::VoronoiDiagramGenerator::ref Site v  )  [private]
 

Definition at line 519 of file Voronoi.cc.

Referenced by bisect(), endpoint(), and PQinsert().

00520 {
00521   v->refcnt += 1;
00522 }

void fastjet::VoronoiDiagramGenerator::resetIterator  )  [inline]
 

Definition at line 203 of file Voronoi.hh.

Referenced by fastjet::ClusterSequenceVoronoiArea::VoronoiAreaCalc::VoronoiAreaCalc().

00203                              {
00204     iteratorEdges = allEdges;
00205   }

int fastjet::VoronoiDiagramGenerator::right_of Halfedge el,
Point p
[private]
 

Definition at line 436 of file Voronoi.cc.

References fastjet::Edge::a, fastjet::Edge::b, fastjet::Edge::c, fastjet::Halfedge::ELedge, fastjet::Halfedge::ELpm, le, re, fastjet::Edge::reg, fastjet::Point::x, and fastjet::Point::y.

Referenced by ELleftbnd().

00437 {
00438   Edge *e;
00439   Site *topsite;
00440   int right_of_site, above, fast;
00441   double dxp, dyp, dxs, t1, t2, t3, yl;
00442         
00443   e = el->ELedge;
00444   topsite = e->reg[1];
00445   right_of_site = p->x > topsite->coord.x;
00446   if(right_of_site && el->ELpm == le) return(1);
00447   if(!right_of_site && el->ELpm == re) return (0);
00448         
00449   if (e->a == 1.0)
00450     {   dyp = p->y - topsite->coord.y;
00451     dxp = p->x - topsite->coord.x;
00452     fast = 0;
00453     if ((!right_of_site & (e->b<0.0)) | (right_of_site & (e->b>=0.0)) )
00454       { above = dyp>= e->b*dxp; 
00455       fast = above;
00456       }
00457     else 
00458       { above = p->x + p->y*e->b > e-> c;
00459       if(e->b<0.0) above = !above;
00460       if (!above) fast = 1;
00461       };
00462     if (!fast)
00463       { dxs = topsite->coord.x - (e->reg[0])->coord.x;
00464       above = e->b * (dxp*dxp - dyp*dyp) <
00465         dxs*dyp*(1.0+2.0*dxp/dxs + e->b*e->b);
00466       if(e->b<0.0) above = !above;
00467       };
00468     }
00469   else  /*e->b==1.0 */
00470     {   yl = e->c - e->a*p->x;
00471     t1 = p->y - yl;
00472     t2 = p->x - topsite->coord.x;
00473     t3 = yl - topsite->coord.y;
00474     above = t1*t1 > t2*t2 + t3*t3;
00475     };
00476   return (el->ELpm==le ? above : !above);
00477 }

Site * fastjet::VoronoiDiagramGenerator::rightreg Halfedge he  )  [private]
 

Definition at line 312 of file Voronoi.cc.

References bottomsite, le, and re.

Referenced by voronoi().

00312                                                     {
00313   if (he->ELedge == (Edge*) NULL)
00314     // if this halfedge has no edge, return the bottom site (whatever
00315     // that is)
00316     return bottomsite;
00317   
00318   // if the ELpm field is zero, return the site 0 that this edge
00319   // bisects, otherwise return site number 1
00320   return (he->ELpm == le)
00321     ? he->ELedge->reg[re] 
00322     : he->ELedge->reg[le];
00323 }

bool fastjet::VoronoiDiagramGenerator::voronoi int  triangulate  )  [private]
 

Definition at line 918 of file Voronoi.cc.

References bisect(), bottomsite, clip_line(), deref(), dist(), ELdelete(), fastjet::Halfedge::ELedge, ELinitialize(), ELinsert(), ELleft(), ELleftbnd(), ELleftend, fastjet::Halfedge::ELpm, ELright(), ELrightend, endpoint(), HEcreate(), intersect(), le, leftreg(), makevertex(), nextone(), out_site(), out_triple(), PQ_min(), PQdelete(), PQempty(), PQextractmin(), PQinitialize(), PQinsert(), re, rightreg(), fastjet::Halfedge::vertex, fastjet::Point::x, and fastjet::Point::y.

Referenced by generateVoronoi().

00919 {
00920   Site *newsite, *bot, *top, *temp, *p;
00921   Site *v;
00922   Point newintstar;
00923   int pm;
00924   Halfedge *lbnd, *rbnd, *llbnd, *rrbnd, *bisector;
00925   Edge *e;
00926         
00927   PQinitialize();
00928   bottomsite = nextone();
00929   out_site(bottomsite);
00930   bool retval = ELinitialize();
00931 
00932   if(!retval)
00933     return false;
00934         
00935   newsite = nextone();
00936   while(1)
00937     {
00938 
00939       if(!PQempty()) 
00940         newintstar = PQ_min();
00941                 
00942       //if the lowest site has a smaller y value than the lowest vector intersection, process the site
00943       //otherwise process the vector intersection               
00944       if (newsite != (Site *)NULL){
00945         volatile double local_newsite_coord_y = newsite->coord.y;
00946         volatile double local_newintstar_y = newintstar.y;
00947         if (PQempty() ||  local_newsite_coord_y < local_newintstar_y
00948             || ( local_newsite_coord_y == local_newintstar_y && newsite->coord.x < newintstar.x))
00949           {/* new site is smallest - this is a site event*/
00950             out_site(newsite);                                          //output the site
00951             lbnd = ELleftbnd(&(newsite->coord));                                //get the first HalfEdge to the LEFT of the new site
00952             rbnd = ELright(lbnd);                                               //get the first HalfEdge to the RIGHT of the new site
00953             bot = rightreg(lbnd);                                               //if this halfedge has no edge, , bot = bottom site (whatever that is)
00954             e = bisect(bot, newsite);                                   //create a new edge that bisects 
00955             bisector = HEcreate(e, le);                                 //create a new HalfEdge, setting its ELpm field to 0                    
00956             ELinsert(lbnd, bisector);                                   //insert this new bisector edge between the left and right vectors in a linked list     
00957             
00958             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
00959               { 
00960                 PQdelete(lbnd);
00961                 PQinsert(lbnd, p, dist(p,newsite));
00962               };
00963             lbnd = bisector;                                            
00964             bisector = HEcreate(e, re);                                 //create a new HalfEdge, setting its ELpm field to 1
00965             ELinsert(lbnd, bisector);                                   //insert the new HE to the right of the original bisector earlier in the IF stmt
00966             
00967             if ((p = intersect(bisector, rbnd)) != (Site *) NULL)       //if this new bisector intersects with the
00968               { 
00969                 PQinsert(bisector, p, dist(p,newsite));                 //push the HE into the ordered linked list of vertices
00970               };
00971             newsite = nextone();        
00972           }
00973         else if (!PQempty()) /* intersection is smallest - this is a vector event */                    
00974           {     
00975             lbnd = PQextractmin();                                              //pop the HalfEdge with the lowest vector off the ordered list of vectors                               
00976             llbnd = ELleft(lbnd);                                               //get the HalfEdge to the left of the above HE
00977             rbnd = ELright(lbnd);                                               //get the HalfEdge to the right of the above HE
00978             rrbnd = ELright(rbnd);                                              //get the HalfEdge to the right of the HE to the right of the lowest HE 
00979             bot = leftreg(lbnd);                                                //get the Site to the left of the left HE which it bisects
00980             top = rightreg(rbnd);                                               //get the Site to the right of the right HE which it bisects
00981             
00982             out_triple(bot, top, rightreg(lbnd));               //output the triple of sites, stating that a circle goes through them
00983             
00984             v = lbnd->vertex;                                           //get the vertex that caused this event
00985             makevertex(v);                                                      //set the vertex number - couldn't do this earlier since we didn't know when it would be processed
00986             endpoint(lbnd->ELedge,lbnd->ELpm,v);        //set the endpoint of the left HalfEdge to be this vector
00987             endpoint(rbnd->ELedge,rbnd->ELpm,v);        //set the endpoint of the right HalfEdge to be this vector
00988             ELdelete(lbnd);                                                     //mark the lowest HE for deletion - can't delete yet because there might be pointers to it in Hash Map  
00989             PQdelete(rbnd);                                                     //remove all vertex events to do with the  right HE
00990             ELdelete(rbnd);                                                     //mark the right HE for deletion - can't delete yet because there might be pointers to it in Hash Map   
00991             pm = le;                                                            //set the pm variable to zero
00992             
00993             if (bot->coord.y > top->coord.y)            //if the site to the left of the event is higher than the Site
00994               {                                                                         //to the right of it, then swap them and set the 'pm' variable to 1
00995                 temp = bot; 
00996                 bot = top; 
00997                 top = temp; 
00998                 pm = re;
00999               }
01000             e = bisect(bot, top);                                       //create an Edge (or line) that is between the two Sites. This creates
01001             //the formula of the line, and assigns a line number to it
01002             bisector = HEcreate(e, pm);                         //create a HE from the Edge 'e', and make it point to that edge with its ELedge field
01003             ELinsert(llbnd, bisector);                          //insert the new bisector to the right of the left HE
01004             endpoint(e, re-pm, v);                                      //set one endpoint to the new edge to be the vector point 'v'.
01005             //If the site to the left of this bisector is higher than the right
01006             //Site, then this endpoint is put in position 0; otherwise in pos 1
01007             deref(v);                                                           //delete the vector 'v'
01008             
01009             //if left HE and the new bisector don't intersect, then delete the left HE, and reinsert it 
01010             if((p = intersect(llbnd, bisector)) != (Site *) NULL)
01011               { 
01012                 PQdelete(llbnd);
01013                 PQinsert(llbnd, p, dist(p,bot));
01014               };
01015             
01016             //if right HE and the new bisector don't intersect, then reinsert it 
01017             if ((p = intersect(bisector, rrbnd)) != (Site *) NULL)
01018               { 
01019                 PQinsert(bisector, p, dist(p,bot));
01020               };
01021           }
01022       } else break;
01023     };
01024 
01025         
01026 
01027 
01028   for(lbnd=ELright(ELleftend); lbnd != ELrightend; lbnd=ELright(lbnd))
01029     {   
01030       e = lbnd->ELedge;
01031 
01032       clip_line(e);
01033     };
01034 
01035   //cleanup();
01036 
01037   return true;
01038         
01039 }


Member Data Documentation

GraphEdge* fastjet::VoronoiDiagramGenerator::allEdges [private]
 

Definition at line 309 of file Voronoi.hh.

Referenced by cleanupEdges(), pushGraphEdge(), and VoronoiDiagramGenerator().

FreeNodeArrayList* fastjet::VoronoiDiagramGenerator::allMemoryList [private]
 

Definition at line 306 of file Voronoi.hh.

Referenced by cleanup(), VoronoiDiagramGenerator(), and ~VoronoiDiagramGenerator().

double fastjet::VoronoiDiagramGenerator::borderMaxX [private]
 

Definition at line 304 of file Voronoi.hh.

Referenced by clip_line(), and generateVoronoi().

double fastjet::VoronoiDiagramGenerator::borderMaxY [private]
 

Definition at line 304 of file Voronoi.hh.

Referenced by clip_line(), and generateVoronoi().

double fastjet::VoronoiDiagramGenerator::borderMinX [private]
 

Definition at line 304 of file Voronoi.hh.

Referenced by clip_line(), and generateVoronoi().

double fastjet::VoronoiDiagramGenerator::borderMinY [private]
 

Definition at line 304 of file Voronoi.hh.

Referenced by clip_line(), and generateVoronoi().

Site* fastjet::VoronoiDiagramGenerator::bottomsite [private]
 

Definition at line 291 of file Voronoi.hh.

Referenced by leftreg(), rightreg(), and voronoi().

double fastjet::VoronoiDiagramGenerator::cradius [private]
 

Definition at line 301 of file Voronoi.hh.

Referenced by out_site(), and plotinit().

FreeNodeArrayList* fastjet::VoronoiDiagramGenerator::currentMemoryBlock [private]
 

Definition at line 307 of file Voronoi.hh.

Referenced by cleanup(), getfree(), and VoronoiDiagramGenerator().

int fastjet::VoronoiDiagramGenerator::debug [private]
 

Definition at line 282 of file Voronoi.hh.

Referenced by generateVoronoi(), and out_site().

double fastjet::VoronoiDiagramGenerator::deltax [private]
 

Definition at line 283 of file Voronoi.hh.

Referenced by ELleftbnd(), and geominit().

double fastjet::VoronoiDiagramGenerator::deltay [private]
 

Definition at line 283 of file Voronoi.hh.

Referenced by geominit(), and PQbucket().

Freelist fastjet::VoronoiDiagramGenerator::efl [private]
 

Definition at line 294 of file Voronoi.hh.

Referenced by bisect(), endpoint(), and geominit().

Halfedge** fastjet::VoronoiDiagramGenerator::ELhash [private]
 

Definition at line 226 of file Voronoi.hh.

Referenced by cleanup(), ELgethash(), ELinitialize(), ELleftbnd(), and VoronoiDiagramGenerator().

int fastjet::VoronoiDiagramGenerator::ELhashsize [private]
 

Definition at line 280 of file Voronoi.hh.

Referenced by ELgethash(), ELinitialize(), and ELleftbnd().

Halfedge* fastjet::VoronoiDiagramGenerator::ELleftend [private]
 

Definition at line 279 of file Voronoi.hh.

Referenced by ELinitialize(), ELleftbnd(), and voronoi().

Halfedge * fastjet::VoronoiDiagramGenerator::ELrightend [private]
 

Definition at line 279 of file Voronoi.hh.

Referenced by ELinitialize(), ELleftbnd(), and voronoi().

Freelist fastjet::VoronoiDiagramGenerator::hfl [private]
 

Definition at line 278 of file Voronoi.hh.

Referenced by ELgethash(), ELinitialize(), and HEcreate().

GraphEdge* fastjet::VoronoiDiagramGenerator::iteratorEdges [private]
 

Definition at line 310 of file Voronoi.hh.

Referenced by VoronoiDiagramGenerator().

double fastjet::VoronoiDiagramGenerator::minDistanceBetweenSites [private]
 

Definition at line 312 of file Voronoi.hh.

Referenced by generateVoronoi(), and VoronoiDiagramGenerator().

int fastjet::VoronoiDiagramGenerator::n_parent_sites
 

Definition at line 217 of file Voronoi.hh.

Referenced by generateVoronoi().

int fastjet::VoronoiDiagramGenerator::nedges [private]
 

Definition at line 293 of file Voronoi.hh.

Referenced by bisect(), and geominit().

int fastjet::VoronoiDiagramGenerator::nsites [private]
 

Definition at line 286 of file Voronoi.hh.

Referenced by generateVoronoi(), geominit(), and nextone().

int fastjet::VoronoiDiagramGenerator::ntry [private]
 

Definition at line 300 of file Voronoi.hh.

Referenced by ELleftbnd().

int fastjet::VoronoiDiagramGenerator::nvertices [private]
 

Definition at line 289 of file Voronoi.hh.

Referenced by geominit(), and makevertex().

std::vector<Point>* fastjet::VoronoiDiagramGenerator::parent_sites
 

Definition at line 216 of file Voronoi.hh.

Referenced by generateVoronoi(), and VoronoiDiagramGenerator().

int fastjet::VoronoiDiagramGenerator::plot [private]
 

Definition at line 282 of file Voronoi.hh.

Referenced by generateVoronoi(), and out_site().

int fastjet::VoronoiDiagramGenerator::PQcount [private]
 

Definition at line 297 of file Voronoi.hh.

Referenced by PQdelete(), PQempty(), PQextractmin(), PQinitialize(), and PQinsert().

Halfedge* fastjet::VoronoiDiagramGenerator::PQhash [private]
 

Definition at line 296 of file Voronoi.hh.

Referenced by cleanup(), PQ_min(), PQdelete(), PQextractmin(), PQinitialize(), PQinsert(), and VoronoiDiagramGenerator().

int fastjet::VoronoiDiagramGenerator::PQhashsize [private]
 

Definition at line 295 of file Voronoi.hh.

Referenced by PQbucket(), and PQinitialize().

int fastjet::VoronoiDiagramGenerator::PQmin [private]
 

Definition at line 298 of file Voronoi.hh.

Referenced by PQ_min(), PQbucket(), PQextractmin(), and PQinitialize().

double fastjet::VoronoiDiagramGenerator::pxmax [private]
 

Definition at line 301 of file Voronoi.hh.

Referenced by clip_line(), and plotinit().

double fastjet::VoronoiDiagramGenerator::pxmin [private]
 

Definition at line 301 of file Voronoi.hh.

Referenced by clip_line(), and plotinit().

double fastjet::VoronoiDiagramGenerator::pymax [private]
 

Definition at line 301 of file Voronoi.hh.

Referenced by clip_line(), and plotinit().

double fastjet::VoronoiDiagramGenerator::pymin [private]
 

Definition at line 301 of file Voronoi.hh.

Referenced by clip_line(), and plotinit().

Freelist fastjet::VoronoiDiagramGenerator::sfl [private]
 

Definition at line 290 of file Voronoi.hh.

Referenced by deref(), generateVoronoi(), and intersect().

int fastjet::VoronoiDiagramGenerator::siteidx [private]
 

Definition at line 287 of file Voronoi.hh.

Referenced by generateVoronoi(), nextone(), and VoronoiDiagramGenerator().

Site* fastjet::VoronoiDiagramGenerator::sites [private]
 

Definition at line 285 of file Voronoi.hh.

Referenced by cleanup(), generateVoronoi(), nextone(), and VoronoiDiagramGenerator().

int fastjet::VoronoiDiagramGenerator::sorted [private]
 

Definition at line 282 of file Voronoi.hh.

Referenced by generateVoronoi().

int fastjet::VoronoiDiagramGenerator::sqrt_nsites [private]
 

Definition at line 288 of file Voronoi.hh.

Referenced by ELinitialize(), geominit(), getfree(), and PQinitialize().

int fastjet::VoronoiDiagramGenerator::total_alloc [private]
 

Definition at line 302 of file Voronoi.hh.

Referenced by myalloc().

int fastjet::VoronoiDiagramGenerator::totalsearch [private]
 

Definition at line 300 of file Voronoi.hh.

Referenced by ELleftbnd().

int fastjet::VoronoiDiagramGenerator::triangulate [private]
 

Definition at line 282 of file Voronoi.hh.

Referenced by generateVoronoi(), and out_site().

double fastjet::VoronoiDiagramGenerator::xmax [private]
 

Definition at line 283 of file Voronoi.hh.

Referenced by generateVoronoi(), geominit(), and plotinit().

double fastjet::VoronoiDiagramGenerator::xmin [private]
 

Definition at line 283 of file Voronoi.hh.

Referenced by ELleftbnd(), generateVoronoi(), geominit(), and plotinit().

double fastjet::VoronoiDiagramGenerator::ymax [private]
 

Definition at line 283 of file Voronoi.hh.

Referenced by generateVoronoi(), geominit(), and plotinit().

double fastjet::VoronoiDiagramGenerator::ymin [private]
 

Definition at line 283 of file Voronoi.hh.

Referenced by generateVoronoi(), geominit(), plotinit(), and PQbucket().


The documentation for this class was generated from the following files:
Generated on Fri Aug 15 13:45:47 2008 for fastjet by  doxygen 1.4.2