|
fastjet 2.4.5
|
#include <Voronoi.hh>

Definition at line 195 of file Voronoi.hh.
| fastjet::VoronoiDiagramGenerator::VoronoiDiagramGenerator | ( | ) |
Definition at line 66 of file Voronoi.cc.
References fastjet::FreeNodeArrayList::memory.
{
siteidx = 0;
sites = NULL;
parent_sites = NULL;
allMemoryList = new FreeNodeArrayList;
allMemoryList->memory = NULL;
allMemoryList->next = NULL;
currentMemoryBlock = allMemoryList;
allEdges = NULL;
iteratorEdges = 0;
minDistanceBetweenSites = 0;
ELhash = NULL;
PQhash = NULL;
}
| fastjet::VoronoiDiagramGenerator::~VoronoiDiagramGenerator | ( | ) |
Definition at line 83 of file Voronoi.cc.
{
cleanup();
cleanupEdges();
if (allMemoryList != NULL)
delete allMemoryList;
}
Definition at line 377 of file Voronoi.cc.
References fastjet::Edge::a, fastjet::Edge::b, fastjet::Edge::c, fastjet::Site::coord, fastjet::Edge::edgenbr, fastjet::Edge::ep, fastjet::Edge::reg, fastjet::Point::x, and fastjet::Point::y.
{
double dx,dy,adx,ady;
Edge *newedge;
newedge = (Edge*) getfree(&efl);
newedge->reg[0] = s1; //store the sites that this edge is bisecting
newedge->reg[1] = s2;
ref(s1);
ref(s2);
// to begin with, there are no endpoints on the bisector - it goes
// to infinity
newedge->ep[0] = (Site*) NULL;
newedge->ep[1] = (Site*) NULL;
// get the difference in x dist between the sites
dx = s2->coord.x - s1->coord.x;
dy = s2->coord.y - s1->coord.y;
// make sure that the difference is positive
adx = dx>0 ? dx : -dx;
ady = dy>0 ? dy : -dy;
// get the slope of the line
newedge->c = (double)(s1->coord.x * dx + s1->coord.y * dy
+ (dx*dx + dy*dy)*0.5);
if (adx>ady){
//set formula of line, with x fixed to 1
newedge->a = 1.0; newedge->b = dy/dx; newedge->c /= dx;
} else {
//set formula of line, with y fixed to 1
newedge->b = 1.0; newedge->a = dx/dy; newedge->c /= dy;
}
newedge->edgenbr = nedges;
nedges++;
return newedge;
}
| void fastjet::VoronoiDiagramGenerator::circle | ( | double | x, |
| double | y, | ||
| double | radius | ||
| ) | [private] |
Definition at line 848 of file Voronoi.cc.
{}
| void fastjet::VoronoiDiagramGenerator::cleanup | ( | ) | [private] |
Definition at line 767 of file Voronoi.cc.
References fastjet::FreeNodeArrayList::memory, and fastjet::FreeNodeArrayList::next.
{
if(sites != NULL){
free(sites);
sites = 0;
}
FreeNodeArrayList* current=NULL, *prev=NULL;
current = prev = allMemoryList;
while(current->next!=NULL){
prev = current;
current = current->next;
free(prev->memory);
delete prev;
prev = NULL;
}
if(current!=NULL){
if (current->memory!=NULL){
free(current->memory);
}
delete current;
}
allMemoryList = new FreeNodeArrayList;
allMemoryList->next = NULL;
allMemoryList->memory = NULL;
currentMemoryBlock = allMemoryList;
if (ELhash!=NULL)
free(ELhash);
if (PQhash!=NULL)
free(PQhash);
}
| void fastjet::VoronoiDiagramGenerator::cleanupEdges | ( | ) | [private] |
Definition at line 805 of file Voronoi.cc.
References fastjet::GraphEdge::next.
| void fastjet::VoronoiDiagramGenerator::clip_line | ( | Edge * | e | ) | [private] |
Definition at line 904 of file Voronoi.cc.
References fastjet::Edge::a, fastjet::Edge::b, fastjet::Edge::c, fastjet::Site::coord, fastjet::Edge::ep, fastjet::Edge::reg, fastjet::Point::x, and fastjet::Point::y.
{
Site *s1, *s2;
double x1=0,x2=0,y1=0,y2=0; //, temp = 0;
x1 = e->reg[0]->coord.x;
x2 = e->reg[1]->coord.x;
y1 = e->reg[0]->coord.y;
y2 = e->reg[1]->coord.y;
//if the distance between the two points this line was created from is less than
//the square root of 2, then ignore it
//TODO improve/remove
//if(sqrt(((x2 - x1) * (x2 - x1)) + ((y2 - y1) * (y2 - y1))) < minDistanceBetweenSites)
// {
// return;
// }
pxmin = borderMinX;
pxmax = borderMaxX;
pymin = borderMinY;
pymax = borderMaxY;
if(e->a == 1.0 && e ->b >= 0.0)
{
s1 = e->ep[1];
s2 = e->ep[0];
}
else
{
s1 = e->ep[0];
s2 = e->ep[1];
};
if(e->a == 1.0)
{
y1 = pymin;
if (s1!=(Site *)NULL && s1->coord.y > pymin)
{
y1 = s1->coord.y;
}
if(y1>pymax)
{
// printf("\nClipped (1) y1 = %f to %f",y1,pymax);
y1 = pymax;
//return;
}
x1 = e->c - e->b * y1;
y2 = pymax;
if (s2!=(Site *)NULL && s2->coord.y < pymax)
y2 = s2->coord.y;
if(y2<pymin)
{
//printf("\nClipped (2) y2 = %f to %f",y2,pymin);
y2 = pymin;
//return;
}
x2 = (e->c) - (e->b) * y2;
if (((x1> pxmax) & (x2>pxmax)) | ((x1<pxmin)&(x2<pxmin)))
{
//printf("\nClipLine jumping out(3), x1 = %f, pxmin = %f, pxmax = %f",x1,pxmin,pxmax);
return;
}
if(x1> pxmax)
{ x1 = pxmax; y1 = (e->c - x1)/e->b;};
if(x1<pxmin)
{ x1 = pxmin; y1 = (e->c - x1)/e->b;};
if(x2>pxmax)
{ x2 = pxmax; y2 = (e->c - x2)/e->b;};
if(x2<pxmin)
{ x2 = pxmin; y2 = (e->c - x2)/e->b;};
}
else
{
x1 = pxmin;
if (s1!=(Site *)NULL && s1->coord.x > pxmin)
x1 = s1->coord.x;
if(x1>pxmax)
{
//printf("\nClipped (3) x1 = %f to %f",x1,pxmin);
//return;
x1 = pxmax;
}
y1 = e->c - e->a * x1;
x2 = pxmax;
if (s2!=(Site *)NULL && s2->coord.x < pxmax)
x2 = s2->coord.x;
if(x2<pxmin)
{
//printf("\nClipped (4) x2 = %f to %f",x2,pxmin);
//return;
x2 = pxmin;
}
y2 = e->c - e->a * x2;
if (((y1> pymax) & (y2>pymax)) | ((y1<pymin)&(y2<pymin)))
{
//printf("\nClipLine jumping out(6), y1 = %f, pymin = %f, pymax = %f",y2,pymin,pymax);
return;
}
if(y1> pymax)
{ y1 = pymax; x1 = (e->c - y1)/e->a;};
if(y1<pymin)
{ y1 = pymin; x1 = (e->c - y1)/e->a;};
if(y2>pymax)
{ y2 = pymax; x2 = (e->c - y2)/e->a;};
if(y2<pymin)
{ y2 = pymin; x2 = (e->c - y2)/e->a;};
};
//printf("\nPushing line (%f,%f,%f,%f)",x1,y1,x2,y2);
//fprintf(stdout, "Line with vertices (%f,%f) and (%f,%f)\n",
// e->reg[0]->coord.x, e->reg[1]->coord.x, e->reg[0]->coord.y, e->reg[1]->coord.y);
pushGraphEdge(x1,y1,x2,y2,e->reg[0],e->reg[1]);
}
| void fastjet::VoronoiDiagramGenerator::deref | ( | Site * | v | ) | [private] |
Definition at line 599 of file Voronoi.cc.
References fastjet::Site::refcnt.
Definition at line 582 of file Voronoi.cc.
References fastjet::Site::coord, fastjet::Point::x, and fastjet::Point::y.
{
double dx,dy;
dx = s->coord.x - t->coord.x;
dy = s->coord.y - t->coord.y;
return (double)(sqrt(dx*dx + dy*dy));
}
| void fastjet::VoronoiDiagramGenerator::ELdelete | ( | Halfedge * | he | ) | [private] |
Definition at line 326 of file Voronoi.cc.
References DELETED, fastjet::Halfedge::ELedge, fastjet::Halfedge::ELleft, and fastjet::Halfedge::ELright.
| Halfedge * fastjet::VoronoiDiagramGenerator::ELgethash | ( | int | b | ) | [private] |
Definition at line 235 of file Voronoi.cc.
References DELETED, fastjet::Halfedge::ELedge, and fastjet::Halfedge::ELrefcnt.
{
Halfedge *he;
if ((b<0) || (b>=ELhashsize))
return (Halfedge *) NULL;
he = ELhash[b];
if ((he==(Halfedge*) NULL) || (he->ELedge!=(Edge*) DELETED))
return he;
/* Hash table points to deleted half edge. Patch as necessary. */
ELhash[b] = (Halfedge*) NULL;
if ((he->ELrefcnt -= 1) == 0)
makefree((Freenode*)he, &hfl);
return (Halfedge*) NULL;
}
| bool fastjet::VoronoiDiagramGenerator::ELinitialize | ( | ) | [private] |
Definition at line 189 of file Voronoi.cc.
{
int i;
freeinit(&hfl, sizeof(Halfedge));
ELhashsize = 2 * sqrt_nsites;
//ELhash = (Halfedge **) myalloc ( sizeof *ELhash * ELhashsize);
ELhash = (Halfedge **) myalloc ( sizeof(Halfedge*)*ELhashsize);
if(ELhash == 0)
return false;
for(i=0; i<ELhashsize; i +=1) ELhash[i] = (Halfedge *)NULL;
ELleftend = HEcreate((Edge*) NULL, 0);
ELrightend = HEcreate((Edge*) NULL, 0);
ELleftend->ELleft = (Halfedge*) NULL;
ELleftend->ELright = ELrightend;
ELrightend->ELleft = ELleftend;
ELrightend->ELright = (Halfedge*) NULL;
ELhash[0] = ELleftend;
ELhash[ELhashsize-1] = ELrightend;
return true;
}
Definition at line 225 of file Voronoi.cc.
References fastjet::Halfedge::ELleft, and fastjet::Halfedge::ELright.
{
newHe->ELleft = lb;
newHe->ELright = lb->ELright;
(lb->ELright)->ELleft = newHe;
lb->ELright = newHe;
}
Definition at line 338 of file Voronoi.cc.
References fastjet::Halfedge::ELleft.
{
return he->ELleft;
}
| Halfedge* fastjet::VoronoiDiagramGenerator::ELleft | ( | ) | [private] |
Definition at line 252 of file Voronoi.cc.
References fastjet::Halfedge::ELleft, fastjet::Halfedge::ELrefcnt, fastjet::Halfedge::ELright, and fastjet::Point::x.
{
int i, bucket;
Halfedge *he;
/* Use hash table to get close to desired halfedge */
// use the hash function to find the place in the hash map that this
// HalfEdge should be
// Gregory Soyez: the original code was
//
// bucket = (int)((p->x - xmin)/deltax * ELhashsize);
// // make sure that the bucket position in within the range of the
// //hash array
// if (bucket<0) bucket =0;
// if (bucket>=ELhashsize) bucket = ELhashsize - 1;
//
// but this runs the risk of having a overflow which would
// cause bucket to be truncated at 0 instead of ELhashsize - 1
// (or vice-versa)
// We fix this by performing the test immediately on the double
// We put in a extra bit of margin to be sure the conversion does
// not play dirty tricks on us
//const double &px = p->x;
if (p->x < xmin){ bucket=0;}
else if (p->x >= xmax){ bucket = ELhashsize - 1;}
else{
bucket= (int)((p->x - xmin)/deltax * ELhashsize);
if (bucket>=ELhashsize) bucket = ELhashsize - 1; // the lower cut should be robust
}
he = ELgethash(bucket);
// if the HE isn't found, search backwards and forwards in the hash
// map for the first non-null entry
if (he == (Halfedge*) NULL){
for(i=1;1;i++){
if ((he=ELgethash(bucket-i)) != (Halfedge*) NULL)
break;
if ((he=ELgethash(bucket+i)) != (Halfedge*) NULL)
break;
};
totalsearch += i;
};
ntry += 1;
/* Now search linear list of halfedges for the correct one */
if ((he==ELleftend) || (he != ELrightend && right_of(he,p))){
do{
he = he->ELright;
} while (he!=ELrightend && right_of(he,p));
// keep going right on the list until either the end is reached,
// or you find the 1st edge which the point
he = he->ELleft; //isn't to the right of
} else
// if the point is to the left of the HalfEdge, then search left
// for the HE just to the left of the point
do{
he = he->ELleft;
} while ((he!=ELleftend )&& (!right_of(he,p)));
/* Update hash table and reference counts */
if ((bucket > 0) && (bucket <ELhashsize-1)){
if(ELhash[bucket] != (Halfedge *) NULL){
ELhash[bucket]->ELrefcnt -= 1;
}
ELhash[bucket] = he;
ELhash[bucket]->ELrefcnt += 1;
};
return he;
}
| Halfedge* fastjet::VoronoiDiagramGenerator::ELleftbnd | ( | ) | [private] |
| Halfedge* fastjet::VoronoiDiagramGenerator::ELright | ( | ) | [private] |
Definition at line 333 of file Voronoi.cc.
References fastjet::Halfedge::ELright.
{
return he->ELright;
}
Definition at line 567 of file Voronoi.cc.
References fastjet::Edge::ep, le, re, and fastjet::Edge::reg.
| void fastjet::VoronoiDiagramGenerator::freeinit | ( | Freelist * | fl, |
| int | size | ||
| ) | [private] |
Definition at line 728 of file Voronoi.cc.
References fastjet::Freelist::head, and fastjet::Freelist::nodesize.
{
fl->head = (Freenode *) NULL;
fl->nodesize = size;
}
| bool fastjet::VoronoiDiagramGenerator::generateVoronoi | ( | std::vector< Point > * | _parent_sites, |
| double | minX, | ||
| double | maxX, | ||
| double | minY, | ||
| double | maxY, | ||
| double | minDist = 0 |
||
| ) |
Definition at line 93 of file Voronoi.cc.
References fastjet::scomp(), and fastjet::d0::inline_maths::y().
Referenced by fastjet::ClusterSequenceVoronoiArea::VoronoiAreaCalc::VoronoiAreaCalc().
{
cleanup();
cleanupEdges();
int i;
double x, y;
minDistanceBetweenSites = minDist;
parent_sites = _parent_sites;
nsites = n_parent_sites = parent_sites->size();
plot = 0;
triangulate = 0;
debug = 1;
sorted = 0;
freeinit(&sfl, sizeof (Site));
//sites = (Site *) myalloc(nsites*sizeof( *sites));
sites = (Site *) myalloc(nsites*sizeof( Site));
//parent_sites = (Site *) myalloc(nsites*sizeof( Site));
if (sites == 0)
return false;
xmax = xmin = (*parent_sites)[0].x;
ymax = ymin = (*parent_sites)[0].y;
for(i=0;i<nsites;i++){
x = (*parent_sites)[i].x;
y = (*parent_sites)[i].y;
sites[i].coord.x = x;
sites[i].coord.y = y;
sites[i].sitenbr = i;
sites[i].refcnt = 0;
if (x<xmin)
xmin=x;
else if (x>xmax)
xmax=x;
if (y<ymin)
ymin=y;
else if (y>ymax)
ymax=y;
}
qsort(sites, nsites, sizeof (*sites), scomp);
// Gregory Soyez
//
// check if some of the particles are degenerate to avoid a crash.
//
// At the moment, we work under the assumption that they will be
// "clustered" later on, so we just keep the 1st one and discard the
// others
unsigned int offset=0;
for (int i=1;i<nsites;i++){
if (sites[i].coord.y==sites[i-1].coord.y && sites[i].coord.x==sites[i-1].coord.x){
offset++;
} else if (offset>0){
sites[i-offset] = sites[i];
}
}
if (offset>0){
nsites-=offset;
_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.");
}
siteidx = 0;
geominit();
double temp = 0;
if(minX > maxX){
temp = minX;
minX = maxX;
maxX = temp;
}
if(minY > maxY){
temp = minY;
minY = maxY;
maxY = temp;
}
borderMinX = minX;
borderMinY = minY;
borderMaxX = maxX;
borderMaxY = maxY;
siteidx = 0;
voronoi(triangulate);
return true;
}
| void fastjet::VoronoiDiagramGenerator::geominit | ( | ) | [private] |
| char * fastjet::VoronoiDiagramGenerator::getfree | ( | Freelist * | fl | ) | [private] |
Definition at line 734 of file Voronoi.cc.
References fastjet::Freelist::head, fastjet::FreeNodeArrayList::memory, fastjet::FreeNodeArrayList::next, and fastjet::Freelist::nodesize.
{
int i;
Freenode *t;
if(fl->head == (Freenode *) NULL)
{
t = (Freenode *) myalloc(sqrt_nsites * fl->nodesize);
if(t == 0)
return 0;
currentMemoryBlock->next = new FreeNodeArrayList;
currentMemoryBlock = currentMemoryBlock->next;
currentMemoryBlock->memory = t;
currentMemoryBlock->next = 0;
for(i=0; i<sqrt_nsites; i+=1)
makefree((Freenode *)((char *)t+i*fl->nodesize), fl);
};
t = fl->head;
fl->head = (fl->head)->nextfree;
return((char *)t);
}
| bool fastjet::VoronoiDiagramGenerator::getNext | ( | GraphEdge ** | e | ) | [inline] |
Definition at line 208 of file Voronoi.hh.
References fastjet::GraphEdge::next.
Referenced by fastjet::ClusterSequenceVoronoiArea::VoronoiAreaCalc::VoronoiAreaCalc().
{
if(iteratorEdges == 0)
return false;
*e = iteratorEdges;
iteratorEdges = iteratorEdges->next;
return true;
}
| Halfedge* fastjet::VoronoiDiagramGenerator::HEcreate | ( | ) | [private] |
Definition at line 213 of file Voronoi.cc.
References fastjet::Halfedge::ELedge, fastjet::Halfedge::ELpm, fastjet::Halfedge::ELrefcnt, fastjet::Halfedge::PQnext, and fastjet::Halfedge::vertex.
| Site * fastjet::VoronoiDiagramGenerator::intersect | ( | Halfedge * | el1, |
| Halfedge * | el2, | ||
| Point * | p = 0 |
||
| ) | [private] |
Definition at line 423 of file Voronoi.cc.
References fastjet::Edge::a, fastjet::Edge::b, fastjet::Edge::c, fastjet::Site::coord, fastjet::Halfedge::ELedge, fastjet::Halfedge::ELpm, le, re, fastjet::Site::refcnt, fastjet::Edge::reg, fastjet::Point::x, and fastjet::Point::y.
{
Edge *e1,*e2, *e;
Halfedge *el;
double d, xint, yint;
int right_of_site;
Site *v;
e1 = el1->ELedge;
e2 = el2->ELedge;
if ((e1==(Edge*)NULL) || (e2==(Edge*)NULL))
return ((Site*) NULL);
// if the two edges bisect the same parent, return null
if (e1->reg[1] == e2->reg[1])
return (Site*) NULL;
// Gregory Soyez:
//
// if the 2 parents are too close, the intersection is going to be
// computed from the "long edges" of the triangle which could causes
// large rounding errors. In this case, use the bisector of the 2
// parents to find the interaction point
//
// The following replaces
// d = e1->a * e2->b - e1->b * e2->a;
// if (-1.0e-10<d && d<1.0e-10)
// return (Site*) NULL;
//
// xint = (e1->c*e2->b - e2->c*e1->b)/d;
// yint = (e2->c*e1->a - e1->c*e2->a)/d;
double dx = e2->reg[1]->coord.x - e1->reg[1]->coord.x;
double dy = e2->reg[1]->coord.y - e1->reg[1]->coord.y;
double dxref = e1->reg[1]->coord.x - e1->reg[0]->coord.x;
double dyref = e1->reg[1]->coord.y - e1->reg[0]->coord.y;
if (dx*dx + dy*dy < 1e-14*(dxref*dxref+dyref*dyref)){
// make sure that the difference is positive
double adx = dx>0 ? dx : -dx;
double ady = dy>0 ? dy : -dy;
// get the slope of the line
double a,b;
double c = (double)(e1->reg[1]->coord.x * dx + e1->reg[1]->coord.y * dy
+ (dx*dx + dy*dy)*0.5);
if (adx>ady){
a = 1.0; b = dy/dx; c /= dx;
} else {
b = 1.0; a = dx/dy; c /= dy;
}
d = e1->a * b - e1->b * a;
if (-1.0e-10<d && d<1.0e-10) {
return (Site*) NULL;
}
xint = (e1->c*b - c*e1->b)/d;
yint = (c*e1->a - e1->c*a)/d;
} else {
d = e1->a * e2->b - e1->b * e2->a;
if (-1.0e-10<d && d<1.0e-10) {
return (Site*) NULL;
}
xint = (e1->c*e2->b - e2->c*e1->b)/d;
yint = (e2->c*e1->a - e1->c*e2->a)/d;
}
// end of Gregory Soyez's modifications
volatile double local_y1 = e1->reg[1]->coord.y;
volatile double local_y2 = e2->reg[1]->coord.y;
if( (local_y1 < local_y2) ||
((local_y1 == local_y2) &&
(e1->reg[1]->coord.x < e2->reg[1]->coord.x)) ){
el = el1;
e = e1;
} else {
el = el2;
e = e2;
}
right_of_site = xint >= e->reg[1]->coord.x;
if ((right_of_site && el->ELpm == le) || (!right_of_site && el->ELpm == re))
return (Site*) NULL;
// create a new site at the point of intersection - this is a new
// vector event waiting to happen
v = (Site*) getfree(&sfl);
v->refcnt = 0;
v->coord.x = xint;
v->coord.y = yint;
return v;
}
Definition at line 343 of file Voronoi.cc.
References fastjet::Halfedge::ELedge, fastjet::Halfedge::ELpm, le, re, and fastjet::Edge::reg.
{
if (he->ELedge == (Edge*) NULL)
return(bottomsite);
return (he->ELpm == le)
? he->ELedge->reg[le]
: he->ELedge->reg[re];
}
Definition at line 761 of file Voronoi.cc.
References fastjet::Freelist::head, and fastjet::Freenode::nextfree.
{
curr->nextfree = fl->head;
fl->head = curr;
}
| void fastjet::VoronoiDiagramGenerator::makevertex | ( | Site * | v | ) | [private] |
Definition at line 591 of file Voronoi.cc.
References fastjet::Site::sitenbr.
{
v->sitenbr = nvertices;
nvertices += 1;
out_vertex(v);
}
| char * fastjet::VoronoiDiagramGenerator::myalloc | ( | unsigned | n | ) | [private] |
Definition at line 836 of file Voronoi.cc.
{
char *t=0;
t=(char*)malloc(n);
total_alloc += n;
return(t);
}
| Site * fastjet::VoronoiDiagramGenerator::nextone | ( | ) | [private] |
| void fastjet::VoronoiDiagramGenerator::openpl | ( | ) | [private] |
Definition at line 847 of file Voronoi.cc.
{}
| void fastjet::VoronoiDiagramGenerator::out_bisector | ( | Edge * | e | ) | [private] |
Definition at line 853 of file Voronoi.cc.
{
}
| void fastjet::VoronoiDiagramGenerator::out_ep | ( | Edge * | e | ) | [private] |
Definition at line 860 of file Voronoi.cc.
{
}
| void fastjet::VoronoiDiagramGenerator::out_site | ( | Site * | s | ) | [private] |
Definition at line 872 of file Voronoi.cc.
References fastjet::Site::coord, fastjet::Point::x, and fastjet::Point::y.
{
if(!triangulate & plot & !debug)
circle (s->coord.x, s->coord.y, cradius);
}
Definition at line 880 of file Voronoi.cc.
{
}
| void fastjet::VoronoiDiagramGenerator::out_vertex | ( | Site * | v | ) | [private] |
Definition at line 866 of file Voronoi.cc.
{
}
| void fastjet::VoronoiDiagramGenerator::plotinit | ( | ) | [private] |
Definition at line 887 of file Voronoi.cc.
{
double dx,dy,d;
dy = ymax - ymin;
dx = xmax - xmin;
d = (double)(( dx > dy ? dx : dy) * 1.1);
pxmin = (double)(xmin - (d-dx)/2.0);
pxmax = (double)(xmax + (d-dx)/2.0);
pymin = (double)(ymin - (d-dy)/2.0);
pymax = (double)(ymax + (d-dy)/2.0);
cradius = (double)((pxmax - pxmin)/350.0);
openpl();
range(pxmin, pymin, pxmax, pymax);
}
| Point fastjet::VoronoiDiagramGenerator::PQ_min | ( | ) | [private] |
Definition at line 688 of file Voronoi.cc.
References fastjet::Point::x, and fastjet::Point::y.
| int fastjet::VoronoiDiagramGenerator::PQbucket | ( | Halfedge * | he | ) | [private] |
Definition at line 649 of file Voronoi.cc.
References fastjet::Halfedge::ystar.
{
// Gregory Soyez: the original code was
//
// bucket = (int)((he->ystar - ymin)/deltay * PQhashsize);
// if (bucket<0) bucket = 0;
// if (bucket>=PQhashsize) bucket = PQhashsize-1 ;
// if (bucket < PQmin) PQmin = bucket;
// return(bucket);
//
// but this runs the risk of having a overflow which would
// cause bucket to be truncated at 0 instead of PQhashsize-1
// (or vice-versa)
// We fix this by performing the test immediately on the double
// We put in a extra bit of margin to be sure the conversion does
// not play dirty tricks on us
int bucket;
double hey = he->ystar;
if (hey < ymin){ bucket = 0;}
else if (hey >= ymax){ bucket = PQhashsize-1;}
else {
bucket = (int)((hey - ymin)/deltay * PQhashsize);
if (bucket>=PQhashsize) bucket = PQhashsize-1 ;
}
if (bucket < PQmin) PQmin = bucket;
return(bucket);
}
| void fastjet::VoronoiDiagramGenerator::PQdelete | ( | Halfedge * | he | ) | [private] |
Definition at line 632 of file Voronoi.cc.
References fastjet::Halfedge::PQnext, and fastjet::Halfedge::vertex.
| int fastjet::VoronoiDiagramGenerator::PQempty | ( | ) | [private] |
Definition at line 682 of file Voronoi.cc.
{
return(PQcount==0);
}
| Halfedge * fastjet::VoronoiDiagramGenerator::PQextractmin | ( | ) | [private] |
Definition at line 698 of file Voronoi.cc.
References fastjet::Halfedge::PQnext.
| Halfedge* fastjet::VoronoiDiagramGenerator::PQfind | ( | ) | [private] |
| bool fastjet::VoronoiDiagramGenerator::PQinitialize | ( | ) | [private] |
Definition at line 709 of file Voronoi.cc.
{
int i;
PQcount = 0;
PQmin = 0;
PQhashsize = 4 * sqrt_nsites;
//PQhash = (Halfedge *) myalloc(PQhashsize * sizeof *PQhash);
PQhash = (Halfedge *) myalloc(PQhashsize * sizeof(Halfedge));
if(PQhash == 0)
return false;
for(i=0; i<PQhashsize; i+=1) PQhash[i].PQnext = (Halfedge *)NULL;
return true;
}
| void fastjet::VoronoiDiagramGenerator::PQinsert | ( | Halfedge * | he, |
| Site * | v, | ||
| double | offset | ||
| ) | [private] |
Definition at line 612 of file Voronoi.cc.
References fastjet::Site::coord, fastjet::Halfedge::PQnext, fastjet::Halfedge::vertex, fastjet::Point::x, fastjet::Point::y, and fastjet::Halfedge::ystar.
{
Halfedge *last, *next;
he->vertex = v;
ref(v);
he->ystar = (double)(v->coord.y + offset);
last = &PQhash[PQbucket(he)];
while ((next = last->PQnext) != (Halfedge *) NULL &&
(he->ystar > next->ystar ||
(he->ystar == next->ystar && v->coord.x > next->vertex->coord.x)))
{
last = next;
};
he->PQnext = last->PQnext;
last->PQnext = he;
PQcount += 1;
}
| void fastjet::VoronoiDiagramGenerator::pushGraphEdge | ( | double | x1, |
| double | y1, | ||
| double | x2, | ||
| double | y2, | ||
| Site * | s1, | ||
| Site * | s2 | ||
| ) | [private] |
Definition at line 820 of file Voronoi.cc.
References fastjet::GraphEdge::next, fastjet::GraphEdge::point1, fastjet::GraphEdge::point2, fastjet::Site::sitenbr, fastjet::GraphEdge::x1, fastjet::GraphEdge::x2, fastjet::GraphEdge::y1, and fastjet::GraphEdge::y2.
| void fastjet::VoronoiDiagramGenerator::range | ( | double | minX, |
| double | minY, | ||
| double | maxX, | ||
| double | maxY | ||
| ) | [private] |
Definition at line 849 of file Voronoi.cc.
{}
| void fastjet::VoronoiDiagramGenerator::ref | ( | Site * | v | ) | [private] |
| void fastjet::VoronoiDiagramGenerator::resetIterator | ( | ) | [inline] |
Definition at line 204 of file Voronoi.hh.
Referenced by fastjet::ClusterSequenceVoronoiArea::VoronoiAreaCalc::VoronoiAreaCalc().
{
iteratorEdges = allEdges;
}
Definition at line 523 of file Voronoi.cc.
References fastjet::Edge::a, fastjet::Edge::b, fastjet::Edge::c, fastjet::Site::coord, fastjet::Halfedge::ELedge, fastjet::Halfedge::ELpm, le, re, fastjet::Edge::reg, fastjet::Point::x, and fastjet::Point::y.
{
Edge *e;
Site *topsite;
int right_of_site, above, fast;
double dxp, dyp, dxs, t1, t2, t3, yl;
e = el->ELedge;
topsite = e->reg[1];
right_of_site = p->x > topsite->coord.x;
if(right_of_site && el->ELpm == le) return(1);
if(!right_of_site && el->ELpm == re) return (0);
if (e->a == 1.0)
{ dyp = p->y - topsite->coord.y;
dxp = p->x - topsite->coord.x;
fast = 0;
if ((!right_of_site & (e->b<0.0)) | (right_of_site & (e->b>=0.0)) )
{ above = dyp>= e->b*dxp;
fast = above;
}
else
{ above = p->x + p->y*e->b > e-> c;
if(e->b<0.0) above = !above;
if (!above) fast = 1;
};
if (!fast)
{ dxs = topsite->coord.x - (e->reg[0])->coord.x;
above = e->b * (dxp*dxp - dyp*dyp) <
dxs*dyp*(1.0+2.0*dxp/dxs + e->b*e->b);
if(e->b<0.0) above = !above;
};
}
else /*e->b==1.0 */
{ yl = e->c - e->a*p->x;
t1 = p->y - yl;
t2 = p->x - topsite->coord.x;
t3 = yl - topsite->coord.y;
above = t1*t1 > t2*t2 + t3*t3;
};
return (el->ELpm==le ? above : !above);
}
Definition at line 351 of file Voronoi.cc.
References fastjet::Halfedge::ELedge, fastjet::Halfedge::ELpm, le, re, and fastjet::Edge::reg.
{
if (he->ELedge == (Edge*) NULL)
// if this halfedge has no edge, return the bottom site (whatever
// that is)
return bottomsite;
// if the ELpm field is zero, return the site 0 that this edge
// bisects, otherwise return site number 1
return (he->ELpm == le)
? he->ELedge->reg[re]
: he->ELedge->reg[le];
}
| bool fastjet::VoronoiDiagramGenerator::voronoi | ( | int | triangulate | ) | [private] |
Definition at line 1025 of file Voronoi.cc.
References fastjet::Site::coord, fastjet::Halfedge::ELedge, fastjet::Halfedge::ELpm, le, re, fastjet::Halfedge::vertex, fastjet::Point::x, and fastjet::Point::y.
{
Site *newsite, *bot, *top, *temp, *p;
Site *v;
Point newintstar;
int pm;
Halfedge *lbnd, *rbnd, *llbnd, *rrbnd, *bisector;
Edge *e;
PQinitialize();
bottomsite = nextone();
out_site(bottomsite);
bool retval = ELinitialize();
if(!retval)
return false;
newsite = nextone();
while(1)
{
if(!PQempty())
newintstar = PQ_min();
//if the lowest site has a smaller y value than the lowest vector intersection, process the site
//otherwise process the vector intersection
if (newsite != (Site *)NULL && (PQempty() || newsite->coord.y < newintstar.y
|| (newsite->coord.y == newintstar.y && newsite->coord.x < newintstar.x)))
{/* new site is smallest - this is a site event*/
out_site(newsite); //output the site
lbnd = ELleftbnd(&(newsite->coord)); //get the first HalfEdge to the LEFT of the new site
rbnd = ELright(lbnd); //get the first HalfEdge to the RIGHT of the new site
bot = rightreg(lbnd); //if this halfedge has no edge, , bot = bottom site (whatever that is)
e = bisect(bot, newsite); //create a new edge that bisects
bisector = HEcreate(e, le); //create a new HalfEdge, setting its ELpm field to 0
ELinsert(lbnd, bisector); //insert this new bisector edge between the left and right vectors in a linked list
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
{
PQdelete(lbnd);
PQinsert(lbnd, p, dist(p,newsite));
};
lbnd = bisector;
bisector = HEcreate(e, re); //create a new HalfEdge, setting its ELpm field to 1
ELinsert(lbnd, bisector); //insert the new HE to the right of the original bisector earlier in the IF stmt
if ((p = intersect(bisector, rbnd)) != (Site *) NULL) //if this new bisector intersects with the
{
PQinsert(bisector, p, dist(p,newsite)); //push the HE into the ordered linked list of vertices
};
newsite = nextone();
}
else if (!PQempty()) /* intersection is smallest - this is a vector event */
{
lbnd = PQextractmin(); //pop the HalfEdge with the lowest vector off the ordered list of vectors
llbnd = ELleft(lbnd); //get the HalfEdge to the left of the above HE
rbnd = ELright(lbnd); //get the HalfEdge to the right of the above HE
rrbnd = ELright(rbnd); //get the HalfEdge to the right of the HE to the right of the lowest HE
bot = leftreg(lbnd); //get the Site to the left of the left HE which it bisects
top = rightreg(rbnd); //get the Site to the right of the right HE which it bisects
out_triple(bot, top, rightreg(lbnd)); //output the triple of sites, stating that a circle goes through them
v = lbnd->vertex; //get the vertex that caused this event
makevertex(v); //set the vertex number - couldn't do this earlier since we didn't know when it would be processed
endpoint(lbnd->ELedge,lbnd->ELpm,v); //set the endpoint of the left HalfEdge to be this vector
endpoint(rbnd->ELedge,rbnd->ELpm,v); //set the endpoint of the right HalfEdge to be this vector
ELdelete(lbnd); //mark the lowest HE for deletion - can't delete yet because there might be pointers to it in Hash Map
PQdelete(rbnd); //remove all vertex events to do with the right HE
ELdelete(rbnd); //mark the right HE for deletion - can't delete yet because there might be pointers to it in Hash Map
pm = le; //set the pm variable to zero
if (bot->coord.y > top->coord.y) //if the site to the left of the event is higher than the Site
{ //to the right of it, then swap them and set the 'pm' variable to 1
temp = bot;
bot = top;
top = temp;
pm = re;
}
e = bisect(bot, top); //create an Edge (or line) that is between the two Sites. This creates
//the formula of the line, and assigns a line number to it
bisector = HEcreate(e, pm); //create a HE from the Edge 'e', and make it point to that edge with its ELedge field
ELinsert(llbnd, bisector); //insert the new bisector to the right of the left HE
endpoint(e, re-pm, v); //set one endpoint to the new edge to be the vector point 'v'.
//If the site to the left of this bisector is higher than the right
//Site, then this endpoint is put in position 0; otherwise in pos 1
deref(v); //delete the vector 'v'
//if left HE and the new bisector don't intersect, then delete the left HE, and reinsert it
if((p = intersect(llbnd, bisector)) != (Site *) NULL)
{
PQdelete(llbnd);
PQinsert(llbnd, p, dist(p,bot));
};
//if right HE and the new bisector don't intersect, then reinsert it
if ((p = intersect(bisector, rrbnd)) != (Site *) NULL)
{
PQinsert(bisector, p, dist(p,bot));
};
}
else break;
};
for(lbnd=ELright(ELleftend); lbnd != ELrightend; lbnd=ELright(lbnd))
{
e = lbnd->ELedge;
clip_line(e);
};
//cleanup();
return true;
}
LimitedWarning fastjet::VoronoiDiagramGenerator::_warning_degeneracy [static, private] |
Definition at line 315 of file Voronoi.hh.
Definition at line 310 of file Voronoi.hh.
Definition at line 307 of file Voronoi.hh.
double fastjet::VoronoiDiagramGenerator::borderMaxX [private] |
Definition at line 305 of file Voronoi.hh.
double fastjet::VoronoiDiagramGenerator::borderMaxY [private] |
Definition at line 305 of file Voronoi.hh.
double fastjet::VoronoiDiagramGenerator::borderMinX [private] |
Definition at line 305 of file Voronoi.hh.
double fastjet::VoronoiDiagramGenerator::borderMinY [private] |
Definition at line 305 of file Voronoi.hh.
Site* fastjet::VoronoiDiagramGenerator::bottomsite [private] |
Definition at line 292 of file Voronoi.hh.
double fastjet::VoronoiDiagramGenerator::cradius [private] |
Definition at line 302 of file Voronoi.hh.
Definition at line 308 of file Voronoi.hh.
int fastjet::VoronoiDiagramGenerator::debug [private] |
Definition at line 283 of file Voronoi.hh.
double fastjet::VoronoiDiagramGenerator::deltax [private] |
Definition at line 284 of file Voronoi.hh.
double fastjet::VoronoiDiagramGenerator::deltay [private] |
Definition at line 284 of file Voronoi.hh.
Definition at line 295 of file Voronoi.hh.
Halfedge** fastjet::VoronoiDiagramGenerator::ELhash [private] |
Definition at line 227 of file Voronoi.hh.
int fastjet::VoronoiDiagramGenerator::ELhashsize [private] |
Definition at line 281 of file Voronoi.hh.
Definition at line 280 of file Voronoi.hh.
Halfedge * fastjet::VoronoiDiagramGenerator::ELrightend [private] |
Definition at line 280 of file Voronoi.hh.
Definition at line 279 of file Voronoi.hh.
Definition at line 311 of file Voronoi.hh.
double fastjet::VoronoiDiagramGenerator::minDistanceBetweenSites [private] |
Definition at line 313 of file Voronoi.hh.
Definition at line 218 of file Voronoi.hh.
int fastjet::VoronoiDiagramGenerator::nedges [private] |
Definition at line 294 of file Voronoi.hh.
int fastjet::VoronoiDiagramGenerator::nsites [private] |
Definition at line 287 of file Voronoi.hh.
int fastjet::VoronoiDiagramGenerator::ntry [private] |
Definition at line 301 of file Voronoi.hh.
int fastjet::VoronoiDiagramGenerator::nvertices [private] |
Definition at line 290 of file Voronoi.hh.
| std::vector<Point>* fastjet::VoronoiDiagramGenerator::parent_sites |
Definition at line 217 of file Voronoi.hh.
int fastjet::VoronoiDiagramGenerator::plot [private] |
Definition at line 283 of file Voronoi.hh.
int fastjet::VoronoiDiagramGenerator::PQcount [private] |
Definition at line 298 of file Voronoi.hh.
Halfedge* fastjet::VoronoiDiagramGenerator::PQhash [private] |
Definition at line 297 of file Voronoi.hh.
int fastjet::VoronoiDiagramGenerator::PQhashsize [private] |
Definition at line 296 of file Voronoi.hh.
int fastjet::VoronoiDiagramGenerator::PQmin [private] |
Definition at line 299 of file Voronoi.hh.
double fastjet::VoronoiDiagramGenerator::pxmax [private] |
Definition at line 302 of file Voronoi.hh.
double fastjet::VoronoiDiagramGenerator::pxmin [private] |
Definition at line 302 of file Voronoi.hh.
double fastjet::VoronoiDiagramGenerator::pymax [private] |
Definition at line 302 of file Voronoi.hh.
double fastjet::VoronoiDiagramGenerator::pymin [private] |
Definition at line 302 of file Voronoi.hh.
Definition at line 291 of file Voronoi.hh.
int fastjet::VoronoiDiagramGenerator::siteidx [private] |
Definition at line 288 of file Voronoi.hh.
Site* fastjet::VoronoiDiagramGenerator::sites [private] |
Definition at line 286 of file Voronoi.hh.
int fastjet::VoronoiDiagramGenerator::sorted [private] |
Definition at line 283 of file Voronoi.hh.
int fastjet::VoronoiDiagramGenerator::sqrt_nsites [private] |
Definition at line 289 of file Voronoi.hh.
int fastjet::VoronoiDiagramGenerator::total_alloc [private] |
Definition at line 303 of file Voronoi.hh.
int fastjet::VoronoiDiagramGenerator::totalsearch [private] |
Definition at line 301 of file Voronoi.hh.
int fastjet::VoronoiDiagramGenerator::triangulate [private] |
Definition at line 283 of file Voronoi.hh.
double fastjet::VoronoiDiagramGenerator::xmax [private] |
Definition at line 284 of file Voronoi.hh.
double fastjet::VoronoiDiagramGenerator::xmin [private] |
Definition at line 284 of file Voronoi.hh.
double fastjet::VoronoiDiagramGenerator::ymax [private] |
Definition at line 284 of file Voronoi.hh.
double fastjet::VoronoiDiagramGenerator::ymin [private] |
Definition at line 284 of file Voronoi.hh.
1.7.4