fastjet 2.4.3
|
00001 //STARTHEADER 00002 // $Id: ClusterSequenceVoronoiArea.cc 2047 2011-04-13 13:49:17Z salam $ 00003 // 00004 // Copyright (c) 2006-2007 Matteo Cacciari, Gavin Salam and Gregory Soyez 00005 // 00006 //---------------------------------------------------------------------- 00007 // This file is part of a simple command-line handling environment 00008 // 00009 // FastJet is free software; you can redistribute it and/or modify 00010 // it under the terms of the GNU General Public License as published by 00011 // the Free Software Foundation; either version 2 of the License, or 00012 // (at your option) any later version. 00013 // 00014 // The algorithms that underlie FastJet have required considerable 00015 // development and are described in hep-ph/0512210. If you use 00016 // FastJet as part of work towards a scientific publication, please 00017 // include a citation to the FastJet paper. 00018 // 00019 // FastJet is distributed in the hope that it will be useful, 00020 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00021 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00022 // GNU General Public License for more details. 00023 // 00024 // You should have received a copy of the GNU General Public License 00025 // along with FastJet; if not, write to the Free Software 00026 // Foundation, Inc.: 00027 // 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00028 //---------------------------------------------------------------------- 00029 //ENDHEADER 00030 00031 #include "fastjet/ClusterSequenceVoronoiArea.hh" 00032 #include "fastjet/internal/Voronoi.hh" 00033 #include <list> 00034 #include <cassert> 00035 #include <ostream> 00036 #include <fstream> 00037 #include <iterator> 00038 #include <cmath> 00039 #include <limits> 00040 00041 using namespace std; 00042 00043 FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh 00044 00045 typedef ClusterSequenceVoronoiArea::VoronoiAreaCalc VAC; 00046 00049 class ClusterSequenceVoronoiArea::VoronoiAreaCalc { 00050 public: 00054 VoronoiAreaCalc(const vector<PseudoJet>::const_iterator &, 00055 const vector<PseudoJet>::const_iterator &, 00056 double effective_R); 00057 00060 inline double area (int index) const {return _areas[index];}; 00061 00062 private: 00063 std::vector<double> _areas; 00064 double _effective_R; 00065 double _effective_R_squared; 00066 00071 double edge_circle_intersection(const Point &p0, 00072 const GraphEdge &edge); 00073 00077 inline double circle_area(const double d12_2, double d01_2, double d02_2){ 00078 return 0.5*_effective_R_squared 00079 *acos(min(1.0,(d01_2+d02_2-d12_2)/(2*sqrt(d01_2*d02_2)))); 00080 } 00081 }; 00082 00083 00088 double VAC::edge_circle_intersection(const Point &p0, 00089 const GraphEdge &edge){ 00090 Point p1(edge.x1-p0.x, edge.y1-p0.y); 00091 Point p2(edge.x2-p0.x, edge.y2-p0.y); 00092 Point pdiff = p2-p1; 00093 00094 //fprintf(stdout, "\tpt(%f,%f)\n", p0.x, p0.y); 00095 00096 double cross = vector_product(p1, p2); 00097 double d12_2 = norm(pdiff); 00098 double d01_2 = norm(p1); 00099 double d02_2 = norm(p2); 00100 00101 // compute intersections between edge line and circle 00102 double delta = d12_2*_effective_R_squared - cross*cross; 00103 00104 // if no intersection, area=area_circle 00105 if (delta<=0){ 00106 return circle_area(d12_2, d01_2, d02_2); 00107 } 00108 00109 // we'll only need delta's sqrt now 00110 delta = sqrt(delta); 00111 00112 // b is the projection of 01 onto 12 00113 double b = scalar_product(pdiff, p1); 00114 00115 // intersections with the circle: 00116 // we compute the "coordinate along the line" of the intersection 00117 // with t=0 (1) corresponding to p1 (p2) 00118 // points with 0<t<1 are within the circle others are outside 00119 00120 // positive intersection 00121 double tp = (delta-b)/d12_2; 00122 00123 // if tp is negative, tm also => inters = circle 00124 if (tp<0) 00125 return circle_area(d12_2, d01_2, d02_2); 00126 00127 // we need the second intersection 00128 double tm = -(delta+b)/d12_2; 00129 00130 // if tp<1, it lies in the circle 00131 if (tp<1){ 00132 // if tm<0, the segment has one intersection 00133 // with the circle at p (t=tp) 00134 // the area is a triangle from 1 to p 00135 // then a circle from p to 2 00136 // several tricks can be used: 00137 // - the area of the triangle is tp*area triangle 00138 // - the lenght for the circle are easily obtained 00139 if (tm<0) 00140 return tp*0.5*fabs(cross) 00141 +circle_area((1-tp)*(1-tp)*d12_2, _effective_R_squared, d02_2); 00142 00143 // now, 0 < tm < tp < 1 00144 // the segment intersects twice the circle 00145 // area = 2 cirles at ends + a triangle in the middle 00146 // again, simplifications are staightforward 00147 return (tp-tm)*0.5*fabs(cross) 00148 + circle_area(tm*tm*d12_2, d01_2, _effective_R_squared) 00149 + circle_area((1-tp)*(1-tp)*d12_2, _effective_R_squared, d02_2); 00150 } 00151 00152 // now, we have tp>1 00153 00154 // if in addition tm>1, intersectino is a circle 00155 if (tm>1) 00156 return circle_area(d12_2, d01_2, d02_2); 00157 00158 // if tm<0, the triangle is inside the circle 00159 if (tm<0) 00160 return 0.5*fabs(cross); 00161 00162 // otherwise, only the "tm point" is on the segment 00163 // area = circle from 1 to m and triangle from m to 2 00164 00165 return (1-tm)*0.5*fabs(cross) 00166 +circle_area(tm*tm*d12_2, d01_2, _effective_R_squared); 00167 } 00168 00169 00170 // the constructor... 00171 //---------------------------------------------------------------------- 00172 VAC::VoronoiAreaCalc(const vector<PseudoJet>::const_iterator &jet_begin, 00173 const vector<PseudoJet>::const_iterator &jet_end, 00174 double effective_R) { 00175 00176 assert(effective_R < 0.5*pi); 00177 00178 vector<Point> voronoi_particles; 00179 vector<int> voronoi_indices; 00180 00181 _effective_R = effective_R; 00182 _effective_R_squared = effective_R*effective_R; 00183 00184 double minrap = numeric_limits<double>::max(); 00185 double maxrap = -minrap; 00186 00187 unsigned int n_tot = 0, n_added = 0; 00188 00189 // loop over jets and create the triangulation, as well as cross-referencing 00190 // info 00191 for (vector<PseudoJet>::const_iterator jet_it = jet_begin; 00192 jet_it != jet_end; jet_it++) { 00193 _areas.push_back(0.0); 00194 if ((jet_it->perp2()) != 0.0 || (jet_it->E() != jet_it->pz())){ 00195 // generate the corresponding point 00196 double rap = jet_it->rap(), phi = jet_it->phi(); 00197 voronoi_particles.push_back(Point(rap, phi)); 00198 voronoi_indices.push_back(n_tot); 00199 n_added++; 00200 00201 // insert a copy of the point if it falls within 2*_R_effective 00202 // of the 0,2pi borders (because we are interested in any 00203 // voronoi edge within _R_effective of the other border) 00204 if (phi < 2*_effective_R) { 00205 voronoi_particles.push_back(Point(rap,phi+twopi)); 00206 voronoi_indices.push_back(-1); 00207 n_added++; 00208 } else if (twopi-phi < 2*_effective_R) { 00209 voronoi_particles.push_back(Point(rap,phi-twopi)); 00210 voronoi_indices.push_back(-1); 00211 n_added++; 00212 } 00213 00214 // track the rapidity range 00215 maxrap = max(maxrap,rap); 00216 minrap = min(minrap,rap); 00217 } 00218 n_tot++; 00219 } 00220 00221 // allow for 0-particle case in graceful way 00222 if (n_added == 0) return; 00223 // assert(n_added > 0); // old (pre 2.4) non-graceful exit 00224 00225 // add extreme cases (corner particles): 00226 double max_extend = 2*max(maxrap-minrap+4*_effective_R, twopi+8*_effective_R); 00227 voronoi_particles.push_back(Point(0.5*(minrap+maxrap)-max_extend, pi)); 00228 voronoi_particles.push_back(Point(0.5*(minrap+maxrap)+max_extend, pi)); 00229 voronoi_particles.push_back(Point(0.5*(minrap+maxrap), pi-max_extend)); 00230 voronoi_particles.push_back(Point(0.5*(minrap+maxrap), pi+max_extend)); 00231 00232 // Build the VD 00233 VoronoiDiagramGenerator vdg; 00234 vdg.generateVoronoi(&voronoi_particles, 00235 0.5*(minrap+maxrap)-max_extend, 0.5*(minrap+maxrap)+max_extend, 00236 pi-max_extend, pi+max_extend); 00237 00238 vdg.resetIterator(); 00239 GraphEdge *e=NULL; 00240 unsigned int v_index; 00241 int p_index; 00242 vector<PseudoJet>::const_iterator jet; 00243 00244 while(vdg.getNext(&e)){ 00245 v_index = e->point1; 00246 if (v_index<n_added){ // this removes the corner particles 00247 p_index = voronoi_indices[v_index]; 00248 if (p_index!=-1){ // this removes the copies 00249 jet = jet_begin+voronoi_indices[v_index]; 00250 _areas[p_index]+= 00251 edge_circle_intersection(voronoi_particles[v_index], *e); 00252 } 00253 } 00254 v_index = e->point2; 00255 if (v_index<n_added){ // this removes the corner particles 00256 p_index = voronoi_indices[v_index]; 00257 if (p_index!=-1){ // this removes the copies 00258 jet = jet_begin+voronoi_indices[v_index]; 00259 _areas[p_index]+= 00260 edge_circle_intersection(voronoi_particles[v_index], *e); 00261 } 00262 } 00263 } 00264 00265 00266 } 00267 00268 00269 //---------------------------------------------------------------------- 00271 void ClusterSequenceVoronoiArea::_initializeVA () { 00272 // run the VAC on our original particles 00273 _pa_calc = new VAC(_jets.begin(), 00274 _jets.begin()+n_particles(), 00275 _effective_Rfact*_jet_def.R()); 00276 00277 // transfer the areas to our local structure 00278 // -- first the initial ones 00279 _voronoi_area.reserve(2*n_particles()); 00280 for (unsigned int i=0; i<n_particles(); i++) { 00281 _voronoi_area.push_back(_pa_calc->area(i)); 00282 // make a stab at a 4-vector area 00283 if (_jets[i].perp2() > 0) { 00284 _voronoi_area_4vector.push_back((_pa_calc->area(i)/_jets[i].perp()) 00285 * _jets[i]); 00286 } else { 00287 // not sure what to do here -- just put zero (it won't be meaningful 00288 // anyway) 00289 _voronoi_area_4vector.push_back(PseudoJet(0.0,0.0,0.0,0.0)); 00290 } 00291 } 00292 00293 // -- then the combined areas that arise from the clustering 00294 for (unsigned int i = n_particles(); i < _history.size(); i++) { 00295 double area; 00296 PseudoJet area_4vect; 00297 if (_history[i].parent2 >= 0) { 00298 area = _voronoi_area[_history[i].parent1] + 00299 _voronoi_area[_history[i].parent2]; 00300 area_4vect = _voronoi_area_4vector[_history[i].parent1] + 00301 _voronoi_area_4vector[_history[i].parent2]; 00302 } else { 00303 area = _voronoi_area[_history[i].parent1]; 00304 area_4vect = _voronoi_area_4vector[_history[i].parent1]; 00305 } 00306 _voronoi_area.push_back(area); 00307 _voronoi_area_4vector.push_back(area_4vect); 00308 } 00309 00310 } 00311 00312 //---------------------------------------------------------------------- 00313 ClusterSequenceVoronoiArea::~ClusterSequenceVoronoiArea() { 00314 delete _pa_calc; 00315 } 00316 00317 FASTJET_END_NAMESPACE