FastJet 3.5.0
Loading...
Searching...
No Matches
ClusterSequence_Delaunay.cc
1//FJSTARTHEADER
2// $Id$
3//
4// Copyright (c) 2005-2025, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
5//
6//----------------------------------------------------------------------
7// This file is part of FastJet.
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. They are described in the original FastJet paper,
16// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
17// FastJet as part of work towards a scientific publication, please
18// quote the version you use and include a citation to the manual and
19// optionally also to hep-ph/0512210.
20//
21// FastJet is distributed in the hope that it will be useful,
22// but WITHOUT ANY WARRANTY; without even the implied warranty of
23// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24// GNU General Public License for more details.
25//
26// You should have received a copy of the GNU General Public License
27// along with FastJet. If not, see <http://www.gnu.org/licenses/>.
28//----------------------------------------------------------------------
29//FJENDHEADER
30
31
32#include "fastjet/config.h"
33#include "fastjet/Error.hh"
34#include "fastjet/PseudoJet.hh"
35#include "fastjet/ClusterSequence.hh"
36#include "fastjet/internal/DynamicNearestNeighbours.hh"
37#include<iostream>
38#include<sstream>
39#include<cmath>
40#include <cstdlib>
41#include<cassert>
42#include<memory>
43//
44#ifndef DROP_CGAL // in case we do not have the code for CGAL
45#include "fastjet/internal/Dnn4piCylinder.hh"
46#include "fastjet/internal/Dnn3piCylinder.hh"
47#include "fastjet/internal/Dnn2piCylinder.hh"
48#endif // DROP_CGAL
49
50FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
51
52using namespace std;
53
54
55//----------------------------------------------------------------------
56/// Run the clustering using a Hierarchical Delaunay triangulation and
57/// STL maps to achieve O(N*ln N) behaviour.
58///
59/// There may be internally asserted assumptions about absence of
60/// points with coincident eta-phi coordinates.
61void ClusterSequence::_delaunay_cluster () {
62
63 int n = _jets.size();
64
65 vector<EtaPhi> points(n); // recall EtaPhi is just a typedef'd pair<double>
66 for (int i = 0; i < n; i++) {
67 points[i] = EtaPhi(_jets[i].rap(),_jets[i].phi_02pi());
68 points[i].sanitize(); // make sure things are in the right range
69 }
70
71 // initialise our DNN structure with the set of points
72 const bool verbose = false;
73 SharedPtr<DynamicNearestNeighbours> DNN;
74#ifndef DROP_CGAL // strategy = NlnN* are not supported if we drop CGAL...
75 bool ignore_nearest_is_mirror = (_Rparam < twopi);
76 if (_strategy == NlnN4pi) {
77 DNN.reset(new Dnn4piCylinder(points,verbose));
78 } else if (_strategy == NlnN3pi) {
79 DNN.reset(new Dnn3piCylinder(points,ignore_nearest_is_mirror,verbose));
80 } else if (_strategy == NlnN) {
81 DNN.reset(new Dnn2piCylinder(points,ignore_nearest_is_mirror,verbose));
82 } else
83#else
84 if (_strategy == NlnN4pi || _strategy == NlnN3pi || _strategy == NlnN) {
85 ostringstream err;
86 err << "ERROR: Requested strategy "<<strategy_string()<<" but it is not"<<endl;
87 err << " supported because FastJet was compiled without CGAL"<<endl;
88 throw Error(err.str());
89 //assert(false);
90 } else
91#endif // DROP_CGAL
92 {
93 //ostringstream err;
94 //err << "ERROR: Unrecognized value for strategy: "<<_strategy<<endl;
95 //throw Error(err.str());
96 //-----------------------------------------------------------------
97 // The code should never reach this point, because the checks above
98 // should always handle all _strategy values for which
99 // _delaunay_cluster() is called
100 assert(false);
101 }
102
103 // We will find nearest neighbour for each vertex, and include
104 // distance in map (NB DistMap is a typedef given in the .h file)
105 DistMap DijMap;
106
107 // fill the map with the minimal (as far as we know) subset of Dij
108 // distances (i.e. nearest neighbour ones).
109 for (int ii = 0; ii < n; ii++) {
110 _add_ktdistance_to_map(ii, DijMap, DNN.get());
111 }
112
113 // run the clustering (go up to i=n-1, but then will stop half-way down,
114 // when we reach that point -- it will be the final beam jet and there
115 // will be no nearest neighbours to find).
116 for (int i=0;i<n;i++) {
117 // find nearest vertices
118 // NB: skip cases where the point is not there anymore!
119 TwoVertices SmallestDijPair;
120 int jet_i, jet_j;
121 double SmallestDij;
122 bool Valid2;
123 bool recombine_with_beam;
124 do {
125 SmallestDij = DijMap.begin()->first;
126 SmallestDijPair = DijMap.begin()->second;
127 jet_i = SmallestDijPair.first;
128 jet_j = SmallestDijPair.second;
129 if (verbose) cout << "CS_Delaunay found recombination candidate: " << jet_i << " " << jet_j << " " << SmallestDij << endl; // GPS debugging
130 // distance is immediately removed regardless of whether or not
131 // it is used.
132 // Some temporary testing code relating to problems with the gcc-3.2 compiler
133 //cout << "got here and size is "<< DijMap.size()<< " and it is "<<SmallestDij <<"\n";
134 //cout << jet_i << " "<< jet_j<<"\n";
135 DijMap.erase(DijMap.begin());
136 //cout << "got beyond here\n";
137
138 // need to "prime" the validity of jet_j in such a way that
139 // if it corresponds to the beam then it is automatically valid.
140 recombine_with_beam = (jet_j == BeamJet);
141 if (!recombine_with_beam) {Valid2 = DNN->Valid(jet_j);}
142 else {Valid2 = true;}
143 if (verbose) cout << "CS_Delaunay validities i & j: " << DNN->Valid(jet_i) << " " << Valid2 << endl;
144 } while ( !DNN->Valid(jet_i) || !Valid2);
145
146
147 // The following part acts just on jet momenta and on the history.
148 // The action on the nearest-neighbour structures takes place
149 // later (only if at least 2 jets are around).
150 if (! recombine_with_beam) {
151 int nn; // will be index of new jet
152 if (verbose) cout << "CS_Delaunay call _do_ij_recomb: " << jet_i << " " << jet_j << " " << SmallestDij << endl; // GPS debug
153 _do_ij_recombination_step(jet_i, jet_j, SmallestDij, nn);
154 //OBS // merge the two jets, add new jet, remove old ones
155 //OBS _jets.push_back(_jets[jet_i] + _jets[jet_j]);
156 //OBS
157 //OBS int nn = _jets.size()-1;
158 //OBS _jets[nn].set_cluster_hist_index(n+i);
159 //OBS
160 //OBS // get corresponding indices in history structure
161 //OBS int hist_i = _jets[jet_i].cluster_hist_index();
162 //OBS int hist_j = _jets[jet_j].cluster_hist_index();
163 //OBS
164 //OBS
165 //OBS _add_step_to_history(n+i,min(hist_i,hist_j), max(hist_i,hist_j),
166 //OBS _jets.size()-1, SmallestDij);
167
168 // add new point to points vector
169 EtaPhi newpoint(_jets[nn].rap(), _jets[nn].phi_02pi());
170 newpoint.sanitize(); // make sure it is in correct range
171 points.push_back(newpoint);
172 } else {
173 // recombine the jet with the beam
174 if (verbose) cout << "CS_Delaunay call _do_iB_recomb: " << jet_i << " " << SmallestDij << endl; // GPS debug
175 _do_iB_recombination_step(jet_i, SmallestDij);
176 //OBS _add_step_to_history(n+i,_jets[jet_i].cluster_hist_index(),BeamJet,
177 //OBS Invalid, SmallestDij);
178 }
179
180 // exit the loop because we do not want to look for nearest neighbours
181 // etc. of zero partons
182 if (i == n-1) {break;}
183
184 vector<int> updated_neighbours;
185 if (! recombine_with_beam) {
186 // update DNN
187 int point3;
188 DNN->RemoveCombinedAddCombination(jet_i, jet_j,
189 points[points.size()-1], point3,
190 updated_neighbours);
191 // C++ beginners' comment: static_cast to unsigned int is necessary
192 // to do away with warnings about type mismatch between point3 (int)
193 // and points.size (unsigned int)
194 if (static_cast<unsigned int> (point3) != points.size()-1) {
195 throw Error("INTERNAL ERROR: point3 != points.size()-1");}
196 } else {
197 // update DNN
198 DNN->RemovePoint(jet_i, updated_neighbours);
199 }
200
201 // update map
202 vector<int>::iterator it = updated_neighbours.begin();
203 for (; it != updated_neighbours.end(); ++it) {
204 int ii = *it;
205 _add_ktdistance_to_map(ii, DijMap, DNN.get());
206 }
207
208 } // end clustering loop
209
210}
211
212
213//----------------------------------------------------------------------
214/// Add the current kt distance for particle ii to the map (DijMap)
215/// using information from the DNN object. Work as follows:
216///
217/// . if the kt is zero then it's nearest neighbour is taken to be the
218/// the beam jet and the distance is zero.
219///
220/// . if cylinder distance to nearest neighbour > _Rparam then it is
221/// yiB that is smallest and this is added to map.
222///
223/// . otherwise if the nearest neighbour jj has a larger kt then add
224/// dij to the map.
225///
226/// . otherwise do nothing
227///
228void ClusterSequence::_add_ktdistance_to_map(
229 const int ii,
230 DistMap & DijMap,
231 const DynamicNearestNeighbours * DNN) {
232
233 double yiB = jet_scale_for_algorithm(_jets[ii]);
234 if (yiB == 0.0) {
235 // in this case convention is that we do not worry about distances
236 // but directly state that nearest neighbour is beam
237 DijMap.insert(DijEntry(yiB, TwoVertices(ii,-1)));
238 } else {
239 double DeltaR2 = DNN->NearestNeighbourDistance(ii) * _invR2;
240 // Logic of following bit is: only add point to map if it has
241 // smaller kt2 than nearest neighbour j (if it has larger kt,
242 // then: either it is j's nearest neighbour and then we will
243 // include dij when we come to j; or it is not j's nearest
244 // neighbour and j will recombine with someone else).
245
246 // If DeltaR2 > 1.0 then in any case it will recombine with beam rather
247 // than with any neighbours.
248 // (put general normalisation here at some point)
249 if (DeltaR2 > 1.0) {
250 DijMap.insert(DijEntry(yiB, TwoVertices(ii,-1)));
251 } else {
252 double kt2i = jet_scale_for_algorithm(_jets[ii]);
253 int jj = DNN->NearestNeighbourIndex(ii);
254 if (kt2i <= jet_scale_for_algorithm(_jets[jj])) {
255 double dij = DeltaR2 * kt2i;
256 DijMap.insert(DijEntry(dij, TwoVertices(ii,jj)));
257 }
258 }
259 }
260}
261
262
263FASTJET_END_NAMESPACE
264