|
||||
File indexing completed on 2025-01-18 09:57:18
0001 #ifndef __FASTJET_NNH_HH__ 0002 #define __FASTJET_NNH_HH__ 0003 0004 //FJSTARTHEADER 0005 // $Id$ 0006 // 0007 // Copyright (c) 2005-2021, Matteo Cacciari, Gavin P. Salam and Gregory Soyez 0008 // 0009 //---------------------------------------------------------------------- 0010 // This file is part of FastJet. 0011 // 0012 // FastJet is free software; you can redistribute it and/or modify 0013 // it under the terms of the GNU General Public License as published by 0014 // the Free Software Foundation; either version 2 of the License, or 0015 // (at your option) any later version. 0016 // 0017 // The algorithms that underlie FastJet have required considerable 0018 // development. They are described in the original FastJet paper, 0019 // hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use 0020 // FastJet as part of work towards a scientific publication, please 0021 // quote the version you use and include a citation to the manual and 0022 // optionally also to hep-ph/0512210. 0023 // 0024 // FastJet is distributed in the hope that it will be useful, 0025 // but WITHOUT ANY WARRANTY; without even the implied warranty of 0026 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 0027 // GNU General Public License for more details. 0028 // 0029 // You should have received a copy of the GNU General Public License 0030 // along with FastJet. If not, see <http://www.gnu.org/licenses/>. 0031 //---------------------------------------------------------------------- 0032 //FJENDHEADER 0033 0034 #include "fastjet/NNBase.hh" 0035 0036 FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh 0037 0038 //---------------------------------------------------------------------- 0039 /// @ingroup advanced_usage 0040 /// \class NNH 0041 /// Help solve closest pair problems with generic interparticle and 0042 /// beam distance (generic case) 0043 /// 0044 /// (see NNBase.hh for an introductory description) 0045 /// 0046 /// This variant provides an implementation for any distance measure. 0047 /// It is templated with a BJ (brief jet) classand can be used with or 0048 /// without an extra "Information" template, i.e. NNH<BJ> or NNH<BJ,I> 0049 /// 0050 /// For the NNH<BJ> version of the class to function, BJ must provide 0051 /// three member functions 0052 /// 0053 /// - void BJ::init(const PseudoJet & jet); // initialise with a PseudoJet 0054 /// - double BJ::distance(const BJ * other_bj_jet); // distance between this and other_bj_jet 0055 /// - double BJ::beam_distance() ; // distance to the beam 0056 /// 0057 /// For the NNH<BJ,I> version to function, the BJ::init(...) member 0058 /// must accept an extra argument 0059 /// 0060 /// - void BJ::init(const PseudoJet & jet, I * info); // initialise with a PseudoJet + info 0061 /// 0062 /// NOTE: THE DISTANCE MUST BE SYMMETRIC I.E. SATISFY 0063 /// a.distance(b) == b.distance(a) 0064 /// 0065 /// For an example of how the NNH<BJ> class is used, see the Jade (and 0066 /// EECambridge) plugins 0067 /// 0068 /// NB: the NNH algorithm is expected N^2, but has a worst case of 0069 /// N^3. Many QCD problems tend to place one closer to the N^3 end of 0070 /// the spectrum than one would like. There is scope for further 0071 /// progress (cf Eppstein, Cardinal & Eppstein), nevertheless the 0072 /// current class is already significantly faster than standard N^3 0073 /// implementations. 0074 /// 0075 template<class BJ, class I = _NoInfo> class NNH : public NNBase<I> { 0076 public: 0077 0078 /// constructor with an initial set of jets (which will be assigned indices 0079 /// 0 ... jets.size()-1 0080 NNH(const std::vector<PseudoJet> & jets) : NNBase<I>() {start(jets);} 0081 NNH(const std::vector<PseudoJet> & jets, I * info) : NNBase<I>(info) {start(jets);} 0082 0083 // initialisation from a given list of particles 0084 void start(const std::vector<PseudoJet> & jets); 0085 0086 /// return the dij_min and indices iA, iB, for the corresponding jets. 0087 /// If iB < 0 then iA recombines with the beam 0088 double dij_min(int & iA, int & iB); 0089 0090 /// remove the jet pointed to by index iA 0091 void remove_jet(int iA); 0092 0093 /// merge the jets pointed to by indices A and B and replace them with 0094 /// jet, assigning it an index jet_index. 0095 void merge_jets(int iA, int iB, const PseudoJet & jet, int jet_index); 0096 0097 /// a destructor 0098 ~NNH() { 0099 delete[] briefjets; 0100 } 0101 0102 private: 0103 class NNBJ; // forward declaration 0104 0105 /// establish the nearest neighbour for jet, and cross check consistency 0106 /// of distances for the other jets that are encountered. Assumes 0107 /// jet not contained within begin...end 0108 void set_NN_crosscheck(NNBJ * jet, NNBJ * begin, NNBJ * end); 0109 0110 /// establish the nearest neighbour for jet; don't cross check other jets' 0111 /// distances and allow jet to be contained within begin...end 0112 void set_NN_nocross (NNBJ * jet, NNBJ * begin, NNBJ * end); 0113 0114 /// contains the briefjets 0115 NNBJ * briefjets; 0116 0117 /// semaphores for the current extent of our structure 0118 NNBJ * head, * tail; 0119 0120 /// currently active number of jets 0121 int n; 0122 0123 /// where_is[i] contains a pointer to the jet with index i 0124 std::vector<NNBJ *> where_is; 0125 0126 /// a class that wraps around the BJ, supplementing it with extra information 0127 /// such as pointers to neighbours, etc. 0128 class NNBJ : public BJ { 0129 public: 0130 void init(const PseudoJet & jet, int index_in) { 0131 BJ::init(jet); 0132 other_init(index_in); 0133 } 0134 void init(const PseudoJet & jet, int index_in, I * info) { 0135 BJ::init(jet, info); 0136 other_init(index_in); 0137 } 0138 void other_init(int index_in) { 0139 _index = index_in; 0140 NN_dist = BJ::beam_distance(); 0141 NN = NULL; 0142 } 0143 int index() const {return _index;} 0144 0145 double NN_dist; 0146 NNBJ * NN; 0147 0148 private: 0149 int _index; 0150 }; 0151 0152 }; 0153 0154 0155 0156 //---------------------------------------------------------------------- 0157 template<class BJ, class I> void NNH<BJ,I>::start(const std::vector<PseudoJet> & jets) { 0158 n = jets.size(); 0159 briefjets = new NNBJ[n]; 0160 where_is.resize(2*n); 0161 0162 NNBJ * jetA = briefjets; 0163 0164 // initialise the basic jet info 0165 for (int i = 0; i< n; i++) { 0166 //jetA->init(jets[i], i); 0167 this->init_jet(jetA, jets[i], i); 0168 where_is[i] = jetA; 0169 jetA++; // move on to next entry of briefjets 0170 } 0171 tail = jetA; // a semaphore for the end of briefjets 0172 head = briefjets; // a nicer way of naming start 0173 0174 // now initialise the NN distances: jetA will run from 1..n-1; and 0175 // jetB from 0..jetA-1 0176 for (jetA = head + 1; jetA != tail; jetA++) { 0177 // set NN info for jetA based on jets running from head..jetA-1, 0178 // checking in the process whether jetA itself is an undiscovered 0179 // NN of one of those jets. 0180 set_NN_crosscheck(jetA, head, jetA); 0181 } 0182 //std::cout << "OOOO " << briefjets[1].NN_dist << " " << briefjets[1].NN - briefjets << std::endl; 0183 } 0184 0185 0186 //---------------------------------------------------------------------- 0187 template<class BJ, class I> double NNH<BJ,I>::dij_min(int & iA, int & iB) { 0188 // find the minimum of the diJ on this round 0189 double diJ_min = briefjets[0].NN_dist; 0190 int diJ_min_jet = 0; 0191 for (int i = 1; i < n; i++) { 0192 if (briefjets[i].NN_dist < diJ_min) { 0193 diJ_min_jet = i; 0194 diJ_min = briefjets[i].NN_dist; 0195 } 0196 } 0197 0198 // return information to user about recombination 0199 NNBJ * jetA = & briefjets[diJ_min_jet]; 0200 //std::cout << jetA - briefjets << std::endl; 0201 iA = jetA->index(); 0202 iB = jetA->NN ? jetA->NN->index() : -1; 0203 return diJ_min; 0204 } 0205 0206 0207 //---------------------------------------------------------------------- 0208 // remove jetA from the list 0209 template<class BJ, class I> void NNH<BJ,I>::remove_jet(int iA) { 0210 NNBJ * jetA = where_is[iA]; 0211 // now update our nearest neighbour info and diJ table 0212 // first reduce size of table 0213 tail--; n--; 0214 // Copy last jet contents and diJ info into position of jetA 0215 *jetA = *tail; 0216 // update the info on where the given index is stored 0217 where_is[jetA->index()] = jetA; 0218 0219 for (NNBJ * jetI = head; jetI != tail; jetI++) { 0220 // see if jetI had jetA or jetB as a NN -- if so recalculate the NN 0221 if (jetI->NN == jetA) set_NN_nocross(jetI, head, tail); 0222 0223 // if jetI's NN is the new tail then relabel it so that it becomes jetA 0224 if (jetI->NN == tail) {jetI->NN = jetA;} 0225 } 0226 } 0227 0228 0229 //---------------------------------------------------------------------- 0230 template<class BJ, class I> void NNH<BJ,I>::merge_jets(int iA, int iB, 0231 const PseudoJet & jet, int index) { 0232 0233 NNBJ * jetA = where_is[iA]; 0234 NNBJ * jetB = where_is[iB]; 0235 0236 // If necessary relabel A & B to ensure jetB < jetA, that way if 0237 // the larger of them == newtail then that ends up being jetA and 0238 // the new jet that is added as jetB is inserted in a position that 0239 // has a future! 0240 if (jetA < jetB) std::swap(jetA,jetB); 0241 0242 // initialise jetB based on the new jet 0243 //jetB->init(jet, index); 0244 this->init_jet(jetB, jet, index); 0245 // and record its position (making sure we have the space) 0246 if (index >= int(where_is.size())) where_is.resize(2*index); 0247 where_is[jetB->index()] = jetB; 0248 0249 // now update our nearest neighbour info 0250 // first reduce size of table 0251 tail--; n--; 0252 // Copy last jet contents into position of jetA 0253 *jetA = *tail; 0254 // update the info on where the tail's index is stored 0255 where_is[jetA->index()] = jetA; 0256 0257 for (NNBJ * jetI = head; jetI != tail; jetI++) { 0258 // see if jetI had jetA or jetB as a NN -- if so recalculate the NN 0259 if (jetI->NN == jetA || jetI->NN == jetB) { 0260 set_NN_nocross(jetI, head, tail); 0261 } 0262 0263 // check whether new jetB is closer than jetI's current NN and 0264 // if need be update things 0265 double dist = jetI->distance(jetB); 0266 if (dist < jetI->NN_dist) { 0267 if (jetI != jetB) { 0268 jetI->NN_dist = dist; 0269 jetI->NN = jetB; 0270 } 0271 } 0272 if (dist < jetB->NN_dist) { 0273 if (jetI != jetB) { 0274 jetB->NN_dist = dist; 0275 jetB->NN = jetI; 0276 } 0277 } 0278 0279 // if jetI's NN is the new tail then relabel it so that it becomes jetA 0280 if (jetI->NN == tail) {jetI->NN = jetA;} 0281 } 0282 } 0283 0284 0285 //---------------------------------------------------------------------- 0286 // this function assumes that jet is not contained within begin...end 0287 template <class BJ, class I> void NNH<BJ,I>::set_NN_crosscheck(NNBJ * jet, 0288 NNBJ * begin, NNBJ * end) { 0289 double NN_dist = jet->beam_distance(); 0290 NNBJ * NN = NULL; 0291 for (NNBJ * jetB = begin; jetB != end; jetB++) { 0292 double dist = jet->distance(jetB); 0293 if (dist < NN_dist) { 0294 NN_dist = dist; 0295 NN = jetB; 0296 } 0297 if (dist < jetB->NN_dist) { 0298 jetB->NN_dist = dist; 0299 jetB->NN = jet; 0300 } 0301 } 0302 jet->NN = NN; 0303 jet->NN_dist = NN_dist; 0304 } 0305 0306 0307 //---------------------------------------------------------------------- 0308 // set the NN for jet without checking whether in the process you might 0309 // have discovered a new nearest neighbour for another jet 0310 template <class BJ, class I> void NNH<BJ,I>::set_NN_nocross( 0311 NNBJ * jet, NNBJ * begin, NNBJ * end) { 0312 double NN_dist = jet->beam_distance(); 0313 NNBJ * NN = NULL; 0314 // if (head < jet) { 0315 // for (NNBJ * jetB = head; jetB != jet; jetB++) { 0316 if (begin < jet) { 0317 for (NNBJ * jetB = begin; jetB != jet; jetB++) { 0318 double dist = jet->distance(jetB); 0319 if (dist < NN_dist) { 0320 NN_dist = dist; 0321 NN = jetB; 0322 } 0323 } 0324 } 0325 // if (tail > jet) { 0326 // for (NNBJ * jetB = jet+1; jetB != tail; jetB++) { 0327 if (end > jet) { 0328 for (NNBJ * jetB = jet+1; jetB != end; jetB++) { 0329 double dist = jet->distance (jetB); 0330 if (dist < NN_dist) { 0331 NN_dist = dist; 0332 NN = jetB; 0333 } 0334 } 0335 } 0336 jet->NN = NN; 0337 jet->NN_dist = NN_dist; 0338 } 0339 0340 0341 0342 0343 FASTJET_END_NAMESPACE // defined in fastjet/internal/base.hh 0344 0345 0346 #endif // __FASTJET_NNH_HH__
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |