FastJet 3.0beta1
|
00001 #ifndef __NNH_HH__ 00002 #define __NNH_HH__ 00003 00004 //STARTHEADER 00005 // $Id: NNH.hh 1937 2011-02-08 13:36:02Z soyez $ 00006 // 00007 // Copyright (c) 2009, Matteo Cacciari, Gavin Salam and Gregory Soyez 00008 // 00009 //---------------------------------------------------------------------- 00010 // This file is part of FastJet. 00011 // 00012 // FastJet is free software; you can redistribute it and/or modify 00013 // it under the terms of the GNU General Public License as published by 00014 // the Free Software Foundation; either version 2 of the License, or 00015 // (at your option) any later version. 00016 // 00017 // The algorithms that underlie FastJet have required considerable 00018 // development and are described in hep-ph/0512210. If you use 00019 // FastJet as part of work towards a scientific publication, please 00020 // include a citation to the FastJet paper. 00021 // 00022 // FastJet is distributed in the hope that it will be useful, 00023 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00024 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00025 // GNU General Public License for more details. 00026 // 00027 // You should have received a copy of the GNU General Public License 00028 // along with FastJet; if not, write to the Free Software 00029 // Foundation, Inc.: 00030 // 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00031 //---------------------------------------------------------------------- 00032 //ENDHEADER 00033 00034 #include<fastjet/ClusterSequence.hh> 00035 00036 00037 FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh 00038 00039 /// @ingroup advanced_usage 00040 /// \class _NoInfo 00041 /// dummy class, used as a default template argument 00042 class _NoInfo {}; 00043 00044 /// @ingroup advanced_usage 00045 /// \class NNHInfo 00046 /// template that will help initialise a BJ with a PseudoJet and extra information 00047 template<class I> class NNHInfo { 00048 public: 00049 NNHInfo() : _info(NULL) {} 00050 NNHInfo(I * info) : _info(info) {} 00051 template<class NNBJ> void init_jet(NNBJ * briefjet, const fastjet::PseudoJet & jet, int index) { briefjet->init(jet, index, _info);} 00052 private: 00053 I * _info; 00054 }; 00055 00056 /// @ingroup advanced_usage 00057 /// Specialisation of NNHInfo for cases where there is no extra info 00058 template<> class NNHInfo<_NoInfo> { 00059 public: 00060 NNHInfo() {} 00061 NNHInfo(_NoInfo * info) {} 00062 template<class NNBJ> void init_jet(NNBJ * briefjet, const fastjet::PseudoJet & jet, int index) { briefjet->init(jet, index);} 00063 }; 00064 00065 00066 //---------------------------------------------------------------------- 00067 /// @ingroup advanced_usage 00068 /// \class NNH 00069 /// Help solve closest pair problems with generic interparticle and 00070 /// beam distance. 00071 /// 00072 /// Class to help solve closest pair problems with generic interparticle 00073 /// distances and a beam distance, using Anderberg's Nearest Neighbour 00074 /// Heuristic. 00075 /// 00076 /// It is templated with a BJ (brief jet) class --- BJ should 00077 /// basically cache the minimal amount of information that is needed 00078 /// to efficiently calculate interparticle distances and particle-beam 00079 /// distances. 00080 /// 00081 /// This class can be used with or without an extra "Information" template, 00082 /// i.e. NNB<BJ> or NNH<BJ,I> 00083 /// 00084 /// For the NNH<BJ> version of the class to function, BJ must provide 00085 /// three member functions 00086 /// 00087 /// - void BJ::init(const PseudoJet & jet); // initialise with a PseudoJet 00088 /// - double BJ::distance(const BJ * other_bj_jet); // distance between this and other_bj_jet 00089 /// - double BJ::beam_distance() ; // distance to the beam 00090 /// 00091 /// For the NNH<BJ,I> version to function, the BJ::init(...) member 00092 /// must accept an extra argument 00093 /// 00094 /// - void BJ::init(const PseudoJet & jet, I * info); // initialise with a PseudoJet + info 00095 /// 00096 /// where info might be a pointer to a class that contains, e.g., information 00097 /// about R, or other parameters of the jet algorithm 00098 /// 00099 /// For an example of how the NNH<BJ> class is used, see the Jade (and 00100 /// EECambridge) plugins 00101 /// 00102 /// NB: the NNH algorithm is expected N^2, but has a worst case of 00103 /// N^3. Many QCD problems tend to place one closer to the N^3 end of 00104 /// the spectrum than one would like. There is scope for further 00105 /// progress (cf Eppstein, Cardinal & Eppstein), nevertheless the 00106 /// current class is already significantly faster than standard N^3 00107 /// implementations. 00108 /// 00109 /// 00110 /// Implementation note: this class derives from NNHInfo, which deals 00111 /// with storing any global information that is needed during the clustering 00112 00113 template<class BJ, class I = _NoInfo> class NNH : public NNHInfo<I> { 00114 public: 00115 00116 /// constructor with an initial set of jets (which will be assigned indices 00117 /// 0 ... jets.size()-1 00118 NNH(const std::vector<PseudoJet> & jets) {start(jets);} 00119 NNH(const std::vector<PseudoJet> & jets, I * info) : NNHInfo<I>(info) {start(jets);} 00120 00121 void start(const std::vector<PseudoJet> & jets); 00122 00123 /// return the dij_min and indices iA, iB, for the corresponding jets. 00124 /// If iB < 0 then iA recombines with the beam 00125 double dij_min(int & iA, int & iB); 00126 00127 /// remove the jet pointed to by index iA 00128 void remove_jet(int iA); 00129 00130 /// merge the jets pointed to by indices A and B and replace them with 00131 /// jet, assigning it an index jet_index. 00132 void merge_jets(int iA, int iB, const PseudoJet & jet, int jet_index); 00133 00134 /// a destructor 00135 ~NNH() { 00136 delete[] briefjets; 00137 } 00138 00139 private: 00140 class NNBJ; // forward declaration 00141 00142 /// establish the nearest neighbour for jet, and cross check constistency 00143 /// of distances for the other jets that are encountered. Assumes 00144 /// jet not contained within begin...end 00145 void set_NN_crosscheck(NNBJ * jet, NNBJ * begin, NNBJ * end); 00146 00147 /// establish the nearest neighbour for jet; don't cross check other jets' 00148 /// distances and allow jet to be contained within begin...end 00149 void set_NN_nocross (NNBJ * jet, NNBJ * begin, NNBJ * end); 00150 00151 /// contains the briefjets 00152 NNBJ * briefjets; 00153 00154 /// semaphores for the current extent of our structure 00155 NNBJ * head, * tail; 00156 00157 /// currently active number of jets 00158 int n; 00159 00160 /// where_is[i] contains a pointer to the jet with index i 00161 std::vector<NNBJ *> where_is; 00162 00163 /// a class that wraps around the BJ, supplementing it with extra information 00164 /// such as pointers to neighbours, etc. 00165 class NNBJ : public BJ { 00166 public: 00167 void init(const PseudoJet & jet, int index) { 00168 BJ::init(jet); 00169 other_init(index); 00170 } 00171 void init(const PseudoJet & jet, int index, I * info) { 00172 BJ::init(jet, info); 00173 other_init(index); 00174 } 00175 void other_init(int index) { 00176 _index = index; 00177 NN_dist = BJ::beam_distance(); 00178 NN = NULL; 00179 } 00180 int index() const {return _index;} 00181 00182 double NN_dist; 00183 NNBJ * NN; 00184 00185 private: 00186 int _index; 00187 }; 00188 00189 }; 00190 00191 00192 00193 //---------------------------------------------------------------------- 00194 template<class BJ, class I> void NNH<BJ,I>::start(const std::vector<PseudoJet> & jets) { 00195 n = jets.size(); 00196 briefjets = new NNBJ[n]; 00197 where_is.resize(2*n); 00198 00199 NNBJ * jetA = briefjets; 00200 00201 // initialise the basic jet info 00202 for (int i = 0; i< n; i++) { 00203 //jetA->init(jets[i], i); 00204 init_jet(jetA, jets[i], i); 00205 where_is[i] = jetA; 00206 jetA++; // move on to next entry of briefjets 00207 } 00208 tail = jetA; // a semaphore for the end of briefjets 00209 head = briefjets; // a nicer way of naming start 00210 00211 // now initialise the NN distances: jetA will run from 1..n-1; and 00212 // jetB from 0..jetA-1 00213 for (jetA = head + 1; jetA != tail; jetA++) { 00214 // set NN info for jetA based on jets running from head..jetA-1, 00215 // checking in the process whether jetA itself is an undiscovered 00216 // NN of one of those jets. 00217 set_NN_crosscheck(jetA, head, jetA); 00218 } 00219 //std::cout << "OOOO " << briefjets[1].NN_dist << " " << briefjets[1].NN - briefjets << std::endl; 00220 } 00221 00222 00223 //---------------------------------------------------------------------- 00224 template<class BJ, class I> double NNH<BJ,I>::dij_min(int & iA, int & iB) { 00225 // find the minimum of the diJ on this round 00226 double diJ_min = briefjets[0].NN_dist; 00227 int diJ_min_jet = 0; 00228 for (int i = 1; i < n; i++) { 00229 if (briefjets[i].NN_dist < diJ_min) { 00230 diJ_min_jet = i; 00231 diJ_min = briefjets[i].NN_dist; 00232 } 00233 } 00234 00235 // return information to user about recombination 00236 NNBJ * jetA = & briefjets[diJ_min_jet]; 00237 //std::cout << jetA - briefjets << std::endl; 00238 iA = jetA->index(); 00239 iB = jetA->NN ? jetA->NN->index() : -1; 00240 return diJ_min; 00241 } 00242 00243 00244 //---------------------------------------------------------------------- 00245 // remove jetA from the list 00246 template<class BJ, class I> void NNH<BJ,I>::remove_jet(int iA) { 00247 NNBJ * jetA = where_is[iA]; 00248 // now update our nearest neighbour info and diJ table 00249 // first reduce size of table 00250 tail--; n--; 00251 // Copy last jet contents and diJ info into position of jetA 00252 *jetA = *tail; 00253 // update the info on where the given index is stored 00254 where_is[jetA->index()] = jetA; 00255 00256 for (NNBJ * jetI = head; jetI != tail; jetI++) { 00257 // see if jetI had jetA or jetB as a NN -- if so recalculate the NN 00258 if (jetI->NN == jetA) set_NN_nocross(jetI, head, tail); 00259 00260 // if jetI's NN is the new tail then relabel it so that it becomes jetA 00261 if (jetI->NN == tail) {jetI->NN = jetA;} 00262 } 00263 } 00264 00265 00266 //---------------------------------------------------------------------- 00267 template<class BJ, class I> void NNH<BJ,I>::merge_jets(int iA, int iB, 00268 const PseudoJet & jet, int index) { 00269 00270 NNBJ * jetA = where_is[iA]; 00271 NNBJ * jetB = where_is[iB]; 00272 00273 // If necessary relabel A & B to ensure jetB < jetA, that way if 00274 // the larger of them == newtail then that ends up being jetA and 00275 // the new jet that is added as jetB is inserted in a position that 00276 // has a future! 00277 if (jetA < jetB) std::swap(jetA,jetB); 00278 00279 // initialise jetB based on the new jet 00280 //jetB->init(jet, index); 00281 init_jet(jetB, jet, index); 00282 // and record its position (making sure we have the space) 00283 if (index >= int(where_is.size())) where_is.resize(2*index); 00284 where_is[jetB->index()] = jetB; 00285 00286 // now update our nearest neighbour info 00287 // first reduce size of table 00288 tail--; n--; 00289 // Copy last jet contents into position of jetA 00290 *jetA = *tail; 00291 // update the info on where the tail's index is stored 00292 where_is[jetA->index()] = jetA; 00293 00294 for (NNBJ * jetI = head; jetI != tail; jetI++) { 00295 // see if jetI had jetA or jetB as a NN -- if so recalculate the NN 00296 if (jetI->NN == jetA || jetI->NN == jetB) { 00297 set_NN_nocross(jetI, head, tail); 00298 } 00299 00300 // check whether new jetB is closer than jetI's current NN and 00301 // if need be update things 00302 double dist = jetI->distance(jetB); 00303 if (dist < jetI->NN_dist) { 00304 if (jetI != jetB) { 00305 jetI->NN_dist = dist; 00306 jetI->NN = jetB; 00307 } 00308 } 00309 if (dist < jetB->NN_dist) { 00310 if (jetI != jetB) { 00311 jetB->NN_dist = dist; 00312 jetB->NN = jetI; 00313 } 00314 } 00315 00316 // if jetI's NN is the new tail then relabel it so that it becomes jetA 00317 if (jetI->NN == tail) {jetI->NN = jetA;} 00318 } 00319 } 00320 00321 00322 //---------------------------------------------------------------------- 00323 // this function assumes that jet is not contained within begin...end 00324 template <class BJ, class I> void NNH<BJ,I>::set_NN_crosscheck(NNBJ * jet, 00325 NNBJ * begin, NNBJ * end) { 00326 double NN_dist = jet->beam_distance(); 00327 NNBJ * NN = NULL; 00328 for (NNBJ * jetB = begin; jetB != end; jetB++) { 00329 double dist = jet->distance(jetB); 00330 if (dist < NN_dist) { 00331 NN_dist = dist; 00332 NN = jetB; 00333 } 00334 if (dist < jetB->NN_dist) { 00335 jetB->NN_dist = dist; 00336 jetB->NN = jet; 00337 } 00338 } 00339 jet->NN = NN; 00340 jet->NN_dist = NN_dist; 00341 } 00342 00343 00344 //---------------------------------------------------------------------- 00345 // set the NN for jet without checking whether in the process you might 00346 // have discovered a new nearest neighbour for another jet 00347 template <class BJ, class I> void NNH<BJ,I>::set_NN_nocross( 00348 NNBJ * jet, NNBJ * begin, NNBJ * end) { 00349 double NN_dist = jet->beam_distance(); 00350 NNBJ * NN = NULL; 00351 if (head < jet) { 00352 for (NNBJ * jetB = head; jetB != jet; jetB++) { 00353 double dist = jet->distance(jetB); 00354 if (dist < NN_dist) { 00355 NN_dist = dist; 00356 NN = jetB; 00357 } 00358 } 00359 } 00360 if (tail > jet) { 00361 for (NNBJ * jetB = jet+1; jetB != tail; jetB++) { 00362 double dist = jet->distance (jetB); 00363 if (dist < NN_dist) { 00364 NN_dist = dist; 00365 NN = jetB; 00366 } 00367 } 00368 } 00369 jet->NN = NN; 00370 jet->NN_dist = NN_dist; 00371 } 00372 00373 00374 00375 00376 FASTJET_END_NAMESPACE // defined in fastjet/internal/base.hh 00377 00378 00379 #endif // __NNH_HH__