FastJet 3.0beta1
ClusterSequence_Delaunay.cc
00001 //STARTHEADER
00002 // $Id: ClusterSequence_Delaunay.cc 2329 2011-06-30 14:38:40Z salam $
00003 //
00004 // Copyright (c) 2005-2006, Matteo Cacciari and Gavin Salam
00005 //
00006 //----------------------------------------------------------------------
00007 // This file is part of FastJet.
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 
00032 #include "fastjet/Error.hh"
00033 #include "fastjet/PseudoJet.hh"
00034 #include "fastjet/ClusterSequence.hh"
00035 #include<iostream>
00036 #include<sstream>
00037 #include<cmath>
00038 #include <cstdlib>
00039 #include<cassert>
00040 #include<memory>
00041 //
00042 #ifndef DROP_CGAL // in case we do not have the code for CGAL
00043 #include "fastjet/internal/Dnn4piCylinder.hh"
00044 #include "fastjet/internal/Dnn3piCylinder.hh"
00045 #include "fastjet/internal/Dnn2piCylinder.hh"
00046 #endif //  DROP_CGAL 
00047 
00048 FASTJET_BEGIN_NAMESPACE      // defined in fastjet/internal/base.hh
00049 
00050 using namespace std;
00051 
00052 
00053 //----------------------------------------------------------------------
00054 /// Run the clustering using a Hierarchical Delaunay triangulation and
00055 /// STL maps to achieve O(N*ln N) behaviour.
00056 ///
00057 /// There may be internally asserted assumptions about absence of
00058 /// points with coincident eta-phi coordinates.
00059 void ClusterSequence::_delaunay_cluster () {
00060 
00061   int n = _jets.size();
00062 
00063   vector<EtaPhi> points(n); // recall EtaPhi is just a typedef'd pair<double>
00064   for (int i = 0; i < n; i++) {
00065     points[i] = EtaPhi(_jets[i].rap(),_jets[i].phi_02pi());
00066     points[i].sanitize(); // make sure things are in the right range
00067   }
00068 
00069   // initialise our DNN structure with the set of points
00070   auto_ptr<DynamicNearestNeighbours> DNN;
00071 #ifndef DROP_CGAL // strategy = NlnN* are not supported if we drop CGAL...
00072   bool verbose = false;
00073   bool ignore_nearest_is_mirror = (_Rparam < twopi);
00074   if (_strategy == NlnN4pi) {
00075     DNN.reset(new Dnn4piCylinder(points,verbose));
00076   } else if (_strategy == NlnN3pi) {
00077     DNN.reset(new Dnn3piCylinder(points,ignore_nearest_is_mirror,verbose));
00078   } else if (_strategy == NlnN) {
00079     DNN.reset(new Dnn2piCylinder(points,ignore_nearest_is_mirror,verbose));
00080   } else 
00081 #else
00082   if (_strategy == NlnN4pi || _strategy == NlnN3pi || _strategy == NlnN) {
00083     ostringstream err;
00084     err << "ERROR: Requested strategy "<<strategy_string()<<" but it is not"<<endl;
00085     err << "       supported because FastJet was compiled without CGAL"<<endl;
00086     throw Error(err.str());
00087     //assert(false);
00088   }
00089 #endif // DROP_CGAL
00090   {
00091     ostringstream err;
00092     err << "ERROR: Unrecognized value for strategy: "<<_strategy<<endl;
00093     assert(false);
00094     throw Error(err.str());
00095   }
00096 
00097   // We will find nearest neighbour for each vertex, and include
00098   // distance in map (NB DistMap is a typedef given in the .h file)
00099   DistMap DijMap;
00100 
00101   // fill the map with the minimal (as far as we know) subset of Dij
00102   // distances (i.e. nearest neighbour ones).
00103   for (int ii = 0; ii < n; ii++) {
00104     _add_ktdistance_to_map(ii, DijMap, DNN.get());
00105   }
00106 
00107   // run the clustering (go up to i=n-1, but then will stop half-way down,
00108   // when we reach that point -- it will be the final beam jet and there
00109   // will be no nearest neighbours to find).
00110   for (int i=0;i<n;i++) {
00111     // find nearest vertices
00112     // NB: skip cases where the point is not there anymore!
00113     TwoVertices SmallestDijPair;
00114     int jet_i, jet_j;
00115     double SmallestDij;
00116     bool Valid2;
00117     bool recombine_with_beam;
00118     do { 
00119       SmallestDij = DijMap.begin()->first;
00120       SmallestDijPair = DijMap.begin()->second;
00121       jet_i = SmallestDijPair.first;
00122       jet_j = SmallestDijPair.second;
00123       // distance is immediately removed regardless of whether or not
00124       // it is used.
00125       // Some temporary testing code relating to problems with the gcc-3.2 compiler
00126       //cout << "got here and size is "<< DijMap.size()<< " and it is "<<SmallestDij <<"\n";
00127       //cout <<  jet_i << " "<< jet_j<<"\n";
00128       DijMap.erase(DijMap.begin());
00129       //cout << "got beyond here\n";
00130 
00131       // need to "prime" the validity of jet_j in such a way that 
00132       // if it corresponds to the beam then it is automatically valid.
00133       recombine_with_beam = (jet_j == BeamJet);
00134       if (!recombine_with_beam) {Valid2 = DNN->Valid(jet_j);} 
00135       else {Valid2 = true;}
00136 
00137     } while ( !DNN->Valid(jet_i) || !Valid2);
00138 
00139 
00140     // The following part acts just on jet momenta and on the history.
00141     // The action on the nearest-neighbour structures takes place
00142     // later (only if at least 2 jets are around).
00143     if (! recombine_with_beam) {
00144       int nn; // will be index of new jet
00145       _do_ij_recombination_step(jet_i, jet_j, SmallestDij, nn);
00146       //OBS // merge the two jets, add new jet, remove old ones
00147       //OBS _jets.push_back(_jets[jet_i] + _jets[jet_j]);
00148       //OBS 
00149       //OBS int nn = _jets.size()-1;
00150       //OBS _jets[nn].set_cluster_hist_index(n+i);
00151       //OBS 
00152       //OBS // get corresponding indices in history structure
00153       //OBS int hist_i = _jets[jet_i].cluster_hist_index();
00154       //OBS int hist_j = _jets[jet_j].cluster_hist_index();
00155       //OBS 
00156       //OBS 
00157       //OBS _add_step_to_history(n+i,min(hist_i,hist_j), max(hist_i,hist_j),
00158       //OBS                   _jets.size()-1, SmallestDij);
00159 
00160       // add new point to points vector
00161       EtaPhi newpoint(_jets[nn].rap(), _jets[nn].phi_02pi());
00162       newpoint.sanitize(); // make sure it is in correct range
00163       points.push_back(newpoint);
00164     } else {
00165       // recombine the jet with the beam
00166       _do_iB_recombination_step(jet_i, SmallestDij);
00167       //OBS _add_step_to_history(n+i,_jets[jet_i].cluster_hist_index(),BeamJet,
00168       //OBS                        Invalid, SmallestDij);
00169     }
00170 
00171     // exit the loop because we do not want to look for nearest neighbours
00172     // etc. of zero partons
00173     if (i == n-1) {break;}
00174 
00175     vector<int> updated_neighbours;
00176     if (! recombine_with_beam) {
00177       // update DNN
00178       int point3;
00179       DNN->RemoveCombinedAddCombination(jet_i, jet_j, 
00180                                        points[points.size()-1], point3,
00181                                        updated_neighbours);
00182       // C++ beginners' comment: static_cast to unsigned int is necessary
00183       // to do away with warnings about type mismatch between point3 (int) 
00184       // and points.size (unsigned int)
00185       if (static_cast<unsigned int> (point3) != points.size()-1) {
00186         throw Error("INTERNAL ERROR: point3 != points.size()-1");}
00187     } else {
00188       // update DNN
00189       DNN->RemovePoint(jet_i, updated_neighbours);
00190     }
00191 
00192     // update map
00193     vector<int>::iterator it = updated_neighbours.begin();
00194     for (; it != updated_neighbours.end(); ++it) {
00195       int ii = *it;
00196       _add_ktdistance_to_map(ii, DijMap, DNN.get());
00197     }
00198       
00199   } // end clustering loop 
00200   
00201 }
00202 
00203 
00204 //----------------------------------------------------------------------
00205 /// Add the current kt distance for particle ii to the map (DijMap)
00206 /// using information from the DNN object. Work as follows:
00207 /// 
00208 /// . if the kt is zero then it's nearest neighbour is taken to be the
00209 ///   the beam jet and the distance is zero.
00210 ///
00211 /// . if cylinder distance to nearest neighbour > _Rparam then it is
00212 ///   yiB that is smallest and this is added to map.
00213 ///
00214 /// . otherwise if the nearest neighbour jj has a larger kt then add
00215 ///   dij to the map.
00216 ///
00217 /// . otherwise do nothing
00218 ///
00219 void ClusterSequence::_add_ktdistance_to_map(
00220                           const int & ii, 
00221                           DistMap & DijMap,
00222                           const DynamicNearestNeighbours * DNN) {
00223   
00224   double yiB = jet_scale_for_algorithm(_jets[ii]);
00225   if (yiB == 0.0) {
00226     // in this case convention is that we do not worry about distances
00227     // but directly state that nearest neighbour is beam
00228     DijMap.insert(DijEntry(yiB,  TwoVertices(ii,-1)));
00229   } else {
00230     double DeltaR2 = DNN->NearestNeighbourDistance(ii) * _invR2;
00231     // Logic of following bit is: only add point to map if it has
00232     // smaller kt2 than nearest neighbour j (if it has larger kt,
00233     // then: either it is j's nearest neighbour and then we will
00234     // include dij when we come to j; or it is not j's nearest
00235     // neighbour and j will recombine with someone else).
00236     
00237     // If DeltaR2 > 1.0 then in any case it will recombine with beam rather
00238     // than with any neighbours.
00239     // (put general normalisation here at some point)
00240     if (DeltaR2 > 1.0) {
00241       DijMap.insert(DijEntry(yiB,  TwoVertices(ii,-1)));
00242     } else {
00243       double kt2i = jet_scale_for_algorithm(_jets[ii]);
00244       int jj = DNN->NearestNeighbourIndex(ii);
00245       if (kt2i <= jet_scale_for_algorithm(_jets[jj])) {
00246         double dij = DeltaR2 * kt2i;
00247         DijMap.insert(DijEntry(dij, TwoVertices(ii,jj)));
00248       }
00249     }
00250   }
00251 }
00252 
00253 
00254 FASTJET_END_NAMESPACE
00255 
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends