Fork me on GitHub

source: svn/trunk/Utilities/Fastjet/src/ClusterSequenceVoronoiArea.cc@ 748

Last change on this file since 748 was 11, checked in by severine ovyn, 16 years ago

Fastjet added; CDFCones directory has been changed

File size: 10.2 KB
Line 
1//STARTHEADER
2// $Id: ClusterSequenceVoronoiArea.cc,v 1.1 2008-11-06 14:32:14 ovyn Exp $
3//
4// Copyright (c) 2006-2007 Matteo Cacciari, Gavin Salam and Gregory Soyez
5//
6//----------------------------------------------------------------------
7// This file is part of a simple command-line handling environment
8//
9// FastJet is free software; you can redistribute it and/or modify
10// it under the terms of the GNU General Public License as published by
11// the Free Software Foundation; either version 2 of the License, or
12// (at your option) any later version.
13//
14// The algorithms that underlie FastJet have required considerable
15// development and are described in hep-ph/0512210. If you use
16// FastJet as part of work towards a scientific publication, please
17// include a citation to the FastJet paper.
18//
19// FastJet is distributed in the hope that it will be useful,
20// but WITHOUT ANY WARRANTY; without even the implied warranty of
21// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
22// GNU General Public License for more details.
23//
24// You should have received a copy of the GNU General Public License
25// along with FastJet; if not, write to the Free Software
26// Foundation, Inc.:
27// 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
28//----------------------------------------------------------------------
29//ENDHEADER
30
31#include "../include/fastjet/ClusterSequenceVoronoiArea.hh"
32#include "../include/fastjet/internal/Voronoi.hh"
33#include <list>
34#include <cassert>
35#include <ostream>
36#include <fstream>
37#include <iterator>
38#include <cmath>
39#include <limits>
40
41using namespace std;
42
43FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
44
45typedef ClusterSequenceVoronoiArea::VoronoiAreaCalc VAC;
46
47/// class for carrying out a voronoi area calculation on a set of
48/// initial vectors
49class ClusterSequenceVoronoiArea::VoronoiAreaCalc {
50public:
51 /// constructor that takes a range of a vector together with the
52 /// effective radius for the intersection of discs with voronoi
53 /// cells
54 VoronoiAreaCalc(const vector<PseudoJet>::const_iterator &,
55 const vector<PseudoJet>::const_iterator &,
56 double effective_R);
57
58 /// return the area of the particle associated with the given
59 /// index
60 inline double area (int index) const {return _areas[index];};
61
62private:
63 std::vector<double> _areas; ///< areas, numbered as jets
64 double _effective_R; ///< effective radius
65 double _effective_R_squared; ///< effective radius squared
66
67 /**
68 * compute the intersection of one triangle with the circle
69 * the area is returned
70 */
71 double edge_circle_intersection(const Point &p0,
72 const GraphEdge &edge);
73
74 /// get the area of a circle of radius R centred on the point 0 with
75 /// 1 and 2 on each "side" of the arc. dij is the distance between
76 /// point i and point j and all distances are squared
77 inline double circle_area(const double d12_2, double d01_2, double d02_2){
78 return 0.5*_effective_R_squared
79 *acos((d01_2+d02_2-d12_2)/(2*sqrt(d01_2*d02_2)));
80 }
81};
82
83
84/**
85 * compute the intersection of one triangle with the circle
86 * the area is returned
87 */
88double VAC::edge_circle_intersection(const Point &p0,
89 const GraphEdge &edge){
90 Point p1(edge.x1-p0.x, edge.y1-p0.y);
91 Point p2(edge.x2-p0.x, edge.y2-p0.y);
92 Point pdiff = p2-p1;
93
94 //fprintf(stdout, "\tpt(%f,%f)\n", p0.x, p0.y);
95
96 double cross = vector_product(p1, p2);
97 double d12_2 = norm(pdiff);
98 double d01_2 = norm(p1);
99 double d02_2 = norm(p2);
100
101 // compute intersections between edge line and circle
102 double delta = d12_2*_effective_R_squared - cross*cross;
103
104 // if no intersection, area=area_circle
105 if (delta<=0){
106 return circle_area(d12_2, d01_2, d02_2);
107 }
108
109 // we'll only need delta's sqrt now
110 delta = sqrt(delta);
111
112 // b is the projection of 01 onto 12
113 double b = scalar_product(pdiff, p1);
114
115 // intersections with the circle:
116 // we compute the "coordinate along the line" of the intersection
117 // with t=0 (1) corresponding to p1 (p2)
118 // points with 0<t<1 are within the circle others are outside
119
120 // positive intersection
121 double tp = (delta-b)/d12_2;
122
123 // if tp is negative, tm also => inters = circle
124 if (tp<0)
125 return circle_area(d12_2, d01_2, d02_2);
126
127 // we need the second intersection
128 double tm = -(delta+b)/d12_2;
129
130 // if tp<1, it lies in the circle
131 if (tp<1){
132 // if tm<0, the segment has one intersection
133 // with the circle at p (t=tp)
134 // the area is a triangle from 1 to p
135 // then a circle from p to 2
136 // several tricks can be used:
137 // - the area of the triangle is tp*area triangle
138 // - the lenght for the circle are easily obtained
139 if (tm<0)
140 return tp*0.5*fabs(cross)
141 +circle_area((1-tp)*(1-tp)*d12_2, _effective_R_squared, d02_2);
142
143 // now, 0 < tm < tp < 1
144 // the segment intersects twice the circle
145 // area = 2 cirles at ends + a triangle in the middle
146 // again, simplifications are staightforward
147 return (tp-tm)*0.5*fabs(cross)
148 + circle_area(tm*tm*d12_2, d01_2, _effective_R_squared)
149 + circle_area((1-tp)*(1-tp)*d12_2, _effective_R_squared, d02_2);
150 }
151
152 // now, we have tp>1
153
154 // if in addition tm>1, intersectino is a circle
155 if (tm>1)
156 return circle_area(d12_2, d01_2, d02_2);
157
158 // if tm<0, the triangle is inside the circle
159 if (tm<0)
160 return 0.5*fabs(cross);
161
162 // otherwise, only the "tm point" is on the segment
163 // area = circle from 1 to m and triangle from m to 2
164
165 return (1-tm)*0.5*fabs(cross)
166 +circle_area(tm*tm*d12_2, d01_2, _effective_R_squared);
167}
168
169
170// the constructor...
171//----------------------------------------------------------------------
172VAC::VoronoiAreaCalc(const vector<PseudoJet>::const_iterator &jet_begin,
173 const vector<PseudoJet>::const_iterator &jet_end,
174 double effective_R) {
175
176 assert(effective_R < 0.5*pi);
177
178 vector<Point> voronoi_particles;
179 vector<int> voronoi_indices;
180
181 _effective_R = effective_R;
182 _effective_R_squared = effective_R*effective_R;
183
184 double minrap = numeric_limits<double>::max();
185 double maxrap = -minrap;
186
187 unsigned int n_tot = 0, n_added = 0;
188
189 // loop over jets and create the triangulation, as well as cross-referencing
190 // info
191 for (vector<PseudoJet>::const_iterator jet_it = jet_begin;
192 jet_it != jet_end; jet_it++) {
193 _areas.push_back(0.0);
194 if ((jet_it->perp2()) != 0.0 || (jet_it->E() != jet_it->pz())){
195 // generate the corresponding point
196 double rap = jet_it->rap(), phi = jet_it->phi();
197 voronoi_particles.push_back(Point(rap, phi));
198 voronoi_indices.push_back(n_tot);
199 n_added++;
200
201 // insert a copy of the point if it falls within 2*_R_effective
202 // of the 0,2pi borders (because we are interested in any
203 // voronoi edge within _R_effective of the other border)
204 if (phi < 2*_effective_R) {
205 voronoi_particles.push_back(Point(rap,phi+twopi));
206 voronoi_indices.push_back(-1);
207 n_added++;
208 } else if (twopi-phi < 2*_effective_R) {
209 voronoi_particles.push_back(Point(rap,phi-twopi));
210 voronoi_indices.push_back(-1);
211 n_added++;
212 }
213
214 // track the rapidity range
215 maxrap = max(maxrap,rap);
216 minrap = min(minrap,rap);
217 }
218 n_tot++;
219 }
220
221 assert(n_added > 0);
222
223 // add extreme cases (corner particles):
224 double max_extend = 2*max(maxrap-minrap+4*_effective_R, twopi+8*_effective_R);
225 voronoi_particles.push_back(Point(0.5*(minrap+maxrap)-max_extend, pi));
226 voronoi_particles.push_back(Point(0.5*(minrap+maxrap)+max_extend, pi));
227 voronoi_particles.push_back(Point(0.5*(minrap+maxrap), pi-max_extend));
228 voronoi_particles.push_back(Point(0.5*(minrap+maxrap), pi+max_extend));
229
230 // Build the VD
231 VoronoiDiagramGenerator vdg;
232 vdg.generateVoronoi(&voronoi_particles,
233 0.5*(minrap+maxrap)-max_extend, 0.5*(minrap+maxrap)+max_extend,
234 pi-max_extend, pi+max_extend);
235
236 vdg.resetIterator();
237 GraphEdge *e=NULL;
238 unsigned int v_index;
239 int p_index;
240 vector<PseudoJet>::const_iterator jet;
241
242 while(vdg.getNext(&e)){
243 v_index = e->point1;
244 if (v_index<n_added){ // this removes the corner particles
245 p_index = voronoi_indices[v_index];
246 if (p_index!=-1){ // this removes the copies
247 jet = jet_begin+voronoi_indices[v_index];
248 _areas[p_index]+=
249 edge_circle_intersection(voronoi_particles[v_index], *e);
250 }
251 }
252 v_index = e->point2;
253 if (v_index<n_added){ // this removes the corner particles
254 p_index = voronoi_indices[v_index];
255 if (p_index!=-1){ // this removes the copies
256 jet = jet_begin+voronoi_indices[v_index];
257 _areas[p_index]+=
258 edge_circle_intersection(voronoi_particles[v_index], *e);
259 }
260 }
261 }
262
263
264}
265
266
267//----------------------------------------------------------------------
268///
269void ClusterSequenceVoronoiArea::_initializeVA () {
270 // run the VAC on our original particles
271 _pa_calc = new VAC(_jets.begin(),
272 _jets.begin()+n_particles(),
273 _effective_Rfact*_jet_def.R());
274
275 // transfer the areas to our local structure
276 // -- first the initial ones
277 _voronoi_area.reserve(2*n_particles());
278 for (unsigned int i=0; i<n_particles(); i++) {
279 _voronoi_area.push_back(_pa_calc->area(i));
280 // make a stab at a 4-vector area
281 if (_jets[i].perp2() > 0) {
282 _voronoi_area_4vector.push_back((_pa_calc->area(i)/_jets[i].perp())
283 * _jets[i]);
284 } else {
285 // not sure what to do here -- just put zero (it won't be meaningful
286 // anyway)
287 _voronoi_area_4vector.push_back(PseudoJet(0.0,0.0,0.0,0.0));
288 }
289 }
290
291 // -- then the combined areas that arise from the clustering
292 for (unsigned int i = n_particles(); i < _history.size(); i++) {
293 double area;
294 PseudoJet area_4vect;
295 if (_history[i].parent2 >= 0) {
296 area = _voronoi_area[_history[i].parent1] +
297 _voronoi_area[_history[i].parent2];
298 area_4vect = _voronoi_area_4vector[_history[i].parent1] +
299 _voronoi_area_4vector[_history[i].parent2];
300 } else {
301 area = _voronoi_area[_history[i].parent1];
302 area_4vect = _voronoi_area_4vector[_history[i].parent1];
303 }
304 _voronoi_area.push_back(area);
305 _voronoi_area_4vector.push_back(area_4vect);
306 }
307
308}
309
310//----------------------------------------------------------------------
311ClusterSequenceVoronoiArea::~ClusterSequenceVoronoiArea() {
312 delete _pa_calc;
313}
314
315FASTJET_END_NAMESPACE
Note: See TracBrowser for help on using the repository browser.