fastjet 2.4.5
Voronoi.hh
Go to the documentation of this file.
00001 #ifndef __FASTJET__VORONOI_H__
00002 #define __FASTJET__VORONOI_H__
00003 
00004 //STARTHEADER
00005 // $Id: Voronoi.hh 2295 2011-06-28 10:40:21Z soyez $
00006 //
00007 // Copyright (c) 1994 by AT&T Bell Laboratories (see below)
00008 //
00009 //
00010 //----------------------------------------------------------------------
00011 // This file is included as part of FastJet but was mostly written by
00012 // S. Fortune in C, put into C++ with memory management by S
00013 // O'Sullivan, and with further interface and memeory management
00014 // modifications by Gregory Soyez.
00015 //
00016 // Permission to use, copy, modify, and distribute this software for
00017 // any purpose without fee is hereby granted, provided that this
00018 // entire notice is included in all copies of any software which is or
00019 // includes a copy or modification of this software and in all copies
00020 // of the supporting documentation for such software. THIS SOFTWARE IS
00021 // BEING PROVIDED "AS IS", WITHOUT ANY EXPRESS OR IMPLIED WARRANTY.
00022 // IN PARTICULAR, NEITHER THE AUTHORS NOR AT&T MAKE ANY REPRESENTATION
00023 // OR WARRANTY OF ANY KIND CONCERNING THE MERCHANTABILITY OF THIS
00024 // SOFTWARE OR ITS FITNESS FOR ANY PARTICULAR PURPOSE.
00025 //
00026 //----------------------------------------------------------------------
00027 //ENDHEADER
00028 
00029 
00030 /*
00031 * The author of this software is Steven Fortune.  
00032 * Copyright (c) 1994 by AT&T Bell Laboratories.
00033 * Permission to use, copy, modify, and distribute this software for any
00034 * purpose without fee is hereby granted, provided that this entire notice
00035 * is included in all copies of any software which is or includes a copy
00036 * or modification of this software and in all copies of the supporting
00037 * documentation for such software.
00038 * THIS SOFTWARE IS BEING PROVIDED "AS IS", WITHOUT ANY EXPRESS OR IMPLIED
00039 * WARRANTY.  IN PARTICULAR, NEITHER THE AUTHORS NOR AT&T MAKE ANY
00040 * REPRESENTATION OR WARRANTY OF ANY KIND CONCERNING THE MERCHANTABILITY
00041 * OF THIS SOFTWARE OR ITS FITNESS FOR ANY PARTICULAR PURPOSE.
00042 */
00043 
00044 /* 
00045 * This code was originally written by Stephan Fortune in C code.  I,
00046 * Shane O'Sullivan, have since modified it, encapsulating it in a C++
00047 * class and, fixing memory leaks and adding accessors to the Voronoi
00048 * Edges.  Permission to use, copy, modify, and distribute this
00049 * software for any purpose without fee is hereby granted, provided
00050 * that this entire notice is included in all copies of any software
00051 * which is or includes a copy or modification of this software and in
00052 * all copies of the supporting documentation for such software.  THIS
00053 * SOFTWARE IS BEING PROVIDED "AS IS", WITHOUT ANY EXPRESS OR IMPLIED
00054 * WARRANTY.  IN PARTICULAR, NEITHER THE AUTHORS NOR AT&T MAKE ANY
00055 * REPRESENTATION OR WARRANTY OF ANY KIND CONCERNING THE
00056 * MERCHANTABILITY OF THIS SOFTWARE OR ITS FITNESS FOR ANY PARTICULAR
00057 * PURPOSE.
00058 */
00059 
00060 #include "fastjet/internal/base.hh"
00061 #include "fastjet/internal/LimitedWarning.hh"
00062 #include <vector>
00063 #include <math.h>
00064 #include <stdlib.h>
00065 #include <string.h>
00066 
00067 #define DELETED -2
00068 #define le 0
00069 #define re 1
00070 
00071 //using namespace std;
00072 
00073 FASTJET_BEGIN_NAMESPACE      // defined in fastjet/internal/base.hh
00074 
00079 class Point{
00080 public:
00082   Point() : x(0.0), y(0.0) {}
00083 
00085   Point(double _x, double _y) : x(_x), y(_y) {}
00086 
00088   inline Point operator + (const Point &p) const{
00089     return Point(x+p.x, y+p.y);
00090   }
00091 
00093   inline Point operator - (const Point &p) const{
00094     return Point(x-p.x, y-p.y);
00095   }
00096 
00098   inline Point operator * (const double t) const{
00099     return Point(x*t, y*t);
00100   }
00101 
00103   double x,y;
00104 };
00105 
00106 
00108 inline double norm(const Point p){
00109   return p.x*p.x+p.y*p.y;
00110 }
00111 
00112 
00114 inline double vector_product(const Point &p1, const Point &p2){
00115   return p1.x*p2.y-p1.y*p2.x;
00116 }
00117 
00118 
00120 inline double scalar_product(const Point &p1, const Point &p2){
00121   return p1.x*p2.x+p1.y*p2.y;
00122 }
00123 
00124 
00129 class GraphEdge{
00130 public:
00132   double x1,y1,x2,y2;
00133 
00135   int point1, point2;
00136 
00138   GraphEdge* next;
00139 };
00140 
00141 
00146 class Site{
00147  public:
00148   Point coord;
00149   int sitenbr;
00150   int refcnt;
00151 };
00152 
00153 
00154 
00155 class Freenode{
00156 public:
00157   Freenode *nextfree;
00158 };
00159 
00160 
00161 class FreeNodeArrayList{
00162 public:
00163   Freenode* memory;
00164   FreeNodeArrayList* next;
00165 };
00166 
00167 
00168 class Freelist{
00169 public:
00170   Freenode *head;
00171   int nodesize;
00172 };
00173 
00174 class Edge{
00175 public:
00176   double a,b,c;
00177   Site *ep[2];
00178   Site *reg[2];
00179   int edgenbr;
00180 };
00181 
00182 
00183 class Halfedge{
00184 public:
00185   Halfedge *ELleft, *ELright;
00186   Edge *ELedge;
00187   int ELrefcnt;
00188   char ELpm;
00189   Site *vertex;
00190   volatile double ystar;
00191   Halfedge *PQnext;
00192 };
00193 
00194 
00195 class VoronoiDiagramGenerator{
00196 public:
00197   VoronoiDiagramGenerator();
00198   ~VoronoiDiagramGenerator();
00199 
00200   bool generateVoronoi(std::vector<Point> *_parent_sites,
00201                        double minX, double maxX, double minY, double maxY, 
00202                        double minDist=0);
00203 
00204   inline void resetIterator(){
00205     iteratorEdges = allEdges;
00206   }
00207 
00208   bool getNext(GraphEdge **e){
00209     if(iteratorEdges == 0)
00210       return false;
00211     
00212     *e = iteratorEdges;
00213     iteratorEdges = iteratorEdges->next;
00214     return true;
00215   }
00216   
00217   std::vector<Point> *parent_sites;
00218   int n_parent_sites;
00219 
00220 private:
00221   void cleanup();
00222   void cleanupEdges();
00223   char *getfree(Freelist *fl);  
00224   Halfedge *PQfind();
00225   int PQempty();
00226         
00227   Halfedge **ELhash;
00228   Halfedge *HEcreate(), *ELleft(), *ELright(), *ELleftbnd();
00229   Halfedge *HEcreate(Edge *e,int pm);
00230   
00231   Point PQ_min();
00232   Halfedge *PQextractmin();     
00233   void freeinit(Freelist *fl,int size);
00234   void makefree(Freenode *curr,Freelist *fl);
00235   void geominit();
00236   void plotinit();
00237   bool voronoi(int triangulate);
00238   void ref(Site *v);
00239   void deref(Site *v);
00240   void endpoint(Edge *e,int lr,Site * s);
00241 
00242   void ELdelete(Halfedge *he);
00243   Halfedge *ELleftbnd(Point *p);
00244   Halfedge *ELright(Halfedge *he);
00245   void makevertex(Site *v);
00246   void out_triple(Site *s1, Site *s2,Site * s3);
00247   
00248   void PQinsert(Halfedge *he,Site * v, double offset);
00249   void PQdelete(Halfedge *he);
00250   bool ELinitialize();
00251   void ELinsert(Halfedge *lb, Halfedge *newHe);
00252   Halfedge * ELgethash(int b);
00253   Halfedge *ELleft(Halfedge *he);
00254   Site *leftreg(Halfedge *he);
00255   void out_site(Site *s);
00256   bool PQinitialize();
00257   int PQbucket(Halfedge *he);
00258   void clip_line(Edge *e);
00259   char *myalloc(unsigned n);
00260   int right_of(Halfedge *el,Point *p);
00261 
00262   Site *rightreg(Halfedge *he);
00263   Edge *bisect(Site *s1, Site *s2);
00264   double dist(Site *s,Site *t);
00265   Site *intersect(Halfedge *el1, Halfedge *el2, Point *p=0);
00266 
00267   void out_bisector(Edge *e);
00268   void out_ep(Edge *e);
00269   void out_vertex(Site *v);
00270   Site *nextone();
00271 
00272   void pushGraphEdge(double x1, double y1, double x2, double y2, 
00273                      Site *s1, Site *s2);
00274 
00275   void openpl();
00276   void circle(double x, double y, double radius);
00277   void range(double minX, double minY, double maxX, double maxY);
00278 
00279   Freelist hfl;
00280   Halfedge *ELleftend, *ELrightend;
00281   int ELhashsize;
00282   
00283   int triangulate, sorted, plot, debug;
00284   double xmin, xmax, ymin, ymax, deltax, deltay;
00285   
00286   Site *sites;
00287   int nsites;
00288   int siteidx;
00289   int sqrt_nsites;
00290   int nvertices;
00291   Freelist sfl;
00292   Site *bottomsite;
00293   
00294   int nedges;
00295   Freelist efl;
00296   int PQhashsize;
00297   Halfedge *PQhash;
00298   int PQcount;
00299   int PQmin;
00300   
00301   int ntry, totalsearch;
00302   double pxmin, pxmax, pymin, pymax, cradius;
00303   int total_alloc;
00304   
00305   double borderMinX, borderMaxX, borderMinY, borderMaxY;
00306   
00307   FreeNodeArrayList* allMemoryList;
00308   FreeNodeArrayList* currentMemoryBlock;
00309   
00310   GraphEdge* allEdges;
00311   GraphEdge* iteratorEdges;
00312   
00313   double minDistanceBetweenSites;
00314 
00315   static LimitedWarning _warning_degeneracy;
00316 };
00317 
00318 int scomp(const void *p1,const void *p2);
00319 
00320 
00321 FASTJET_END_NAMESPACE
00322 
00323 #endif
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines