ClusterSequence_Delaunay.cc

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

Generated on Tue Dec 18 17:05:03 2007 for fastjet by  doxygen 1.5.2