00001 //STARTHEADER 00002 // $Id: ClusterSequence.hh 1985 2011-03-09 22:26:40Z soyez $ 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 #ifndef __FASTJET_CLUSTERSEQUENCE_HH__ 00033 #define __FASTJET_CLUSTERSEQUENCE_HH__ 00034 00035 #include<vector> 00036 #include<map> 00037 #include "fastjet/internal/DynamicNearestNeighbours.hh" 00038 #include "fastjet/PseudoJet.hh" 00039 #include<memory> 00040 #include<cassert> 00041 #include<iostream> 00042 #include<string> 00043 #include<set> 00044 #include<cmath> // needed to get double std::abs(double) 00045 #include "fastjet/Error.hh" 00046 #include "fastjet/JetDefinition.hh" 00047 #include "fastjet/SharedPtr.hh" 00048 #include "fastjet/internal/LimitedWarning.hh" 00049 00050 00051 FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh 00052 00053 //---------------------------------------------------------------------- 00054 // here's where we put the main page for fastjet (as explained in the 00055 // Doxygen faq) 00056 // We put in inside te fastjet namespace to have the links without 00057 // having to specify (fastjet::) 00058 //...................................................................... 00059 /** \mainpage FastJet code documentation 00060 * 00061 * These pages provide automatically generated documentation for the 00062 * FastJet package. 00063 * 00064 * \section useful_classes The most useful classes 00065 * 00066 * Many of the facilities of FastJet can be accessed through the three 00067 * following classes: 00068 * 00069 * - PseudoJet: the basic class for holding the 4-momentum of a 00070 * particle or a jet. 00071 * 00072 * - JetDefinition: the combination of a #JetAlgorithm and its 00073 * associated parameters. 00074 * 00075 * - ClusterSequence: constructed with a vector of input (PseudoJet) 00076 * particles and a JetDefinition, it computes and stores the 00077 * information on how the input particles are clustered into jets. 00078 * 00079 * \section advanced_classes Selected more advanced classes 00080 * 00081 * - ClusterSequenceArea: with the help of an AreaDefinition, provides 00082 * jets that also contain information about their area. 00083 * 00084 * \section further_info Further information 00085 * 00086 * - Selected classes ordered by topics can be found under the <a 00087 * href="modules.html">modules</a> tab. 00088 * 00089 * - The complete list of classes is available under the <a 00090 * href="annotated.html">classes</a> tab. 00091 * 00092 * - For non-class material (<a href="namespacefastjet.html#enum-members">enums</a>, 00093 * <a href="namespacefastjet.html#typedef-members">typedefs</a>, 00094 * <a href="namespacefastjet.html#func-members">functions</a>), see the 00095 * #fastjet documentation 00096 * 00097 * - For further information and normal documentation, see the main <a 00098 * href="http://www.lpthe.jussieu.fr/~salam/fastjet">FastJet</a> page. 00099 * 00100 * \section examples Examples 00101 * See our \subpage Examples page 00102 */ 00103 //---------------------------------------------------------------------- 00104 00105 // forward declaration 00106 class ClusterSequenceStructure; 00107 00108 /// @ingroup basic_classes 00109 /// \class ClusterSequence 00110 /// deals with clustering 00111 class ClusterSequence { 00112 00113 00114 public: 00115 00116 /// default constructor 00117 ClusterSequence () {} 00118 00119 /// create a clustersequence starting from the supplied set 00120 /// of pseudojets and clustering them with the long-invariant 00121 /// kt algorithm (E-scheme recombination) with the supplied 00122 /// value for R. 00123 /// 00124 /// If strategy=DumbN3 a very stupid N^3 algorithm is used for the 00125 /// clustering; otherwise strategy = NlnN* uses cylinders algorithms 00126 /// with some number of pi coverage. If writeout_combinations=true a 00127 /// summary of the recombination sequence is written out 00128 template<class L> ClusterSequence (const std::vector<L> & pseudojets, 00129 const double & R = 1.0, 00130 const Strategy & strategy = Best, 00131 const bool & writeout_combinations = false); 00132 00133 00134 /// create a clustersequence starting from the supplied set 00135 /// of pseudojets and clustering them with jet definition specified 00136 /// by jet_def (which also specifies the clustering strategy) 00137 template<class L> ClusterSequence ( 00138 const std::vector<L> & pseudojets, 00139 const JetDefinition & jet_def, 00140 const bool & writeout_combinations = false); 00141 00142 // virtual ClusterSequence destructor, in case any derived class 00143 // thinks of needing a destructor at some point 00144 virtual ~ClusterSequence (); //{} 00145 00146 // NB: in the routines that follow, for extracting lists of jets, a 00147 // list structure might be more efficient, if sometimes a little 00148 // more awkward to use (at least for old fortran hands). 00149 00150 /// return a vector of all jets (in the sense of the inclusive 00151 /// algorithm) with pt >= ptmin. Time taken should be of the order 00152 /// of the number of jets returned. 00153 std::vector<PseudoJet> inclusive_jets (const double & ptmin = 0.0) const; 00154 00155 /// return the number of jets (in the sense of the exclusive 00156 /// algorithm) that would be obtained when running the algorithm 00157 /// with the given dcut. 00158 int n_exclusive_jets (const double & dcut) const; 00159 00160 /// return a vector of all jets (in the sense of the exclusive 00161 /// algorithm) that would be obtained when running the algorithm 00162 /// with the given dcut. 00163 std::vector<PseudoJet> exclusive_jets (const double & dcut) const; 00164 00165 /// return a vector of all jets when the event is clustered (in the 00166 /// exclusive sense) to exactly njets. 00167 /// 00168 /// If there are fewer than njets particles in the ClusterSequence 00169 /// then FastJet crashes since it is not able to return njets 00170 /// exclusive jets. 00171 std::vector<PseudoJet> exclusive_jets (const int & njets) const; 00172 00173 /// return the dmin corresponding to the recombination that went from 00174 /// n+1 to n jets (sometimes known as d_{n n+1}). 00175 double exclusive_dmerge (const int & njets) const; 00176 00177 /// return the maximum of the dmin encountered during all recombinations 00178 /// up to the one that led to an n-jet final state; identical to 00179 /// exclusive_dmerge, except in cases where the dmin do not increase 00180 /// monotonically. 00181 double exclusive_dmerge_max (const int & njets) const; 00182 00183 /// return the ymin corresponding to the recombination that went from 00184 /// n+1 to n jets (sometimes known as y_{n n+1}). 00185 double exclusive_ymerge (int njets) const {return exclusive_dmerge(njets) / Q2();} 00186 00187 /// same as exclusive_dmerge_max, but normalised to squared total energy 00188 double exclusive_ymerge_max (int njets) const {return exclusive_dmerge_max(njets)/Q2();} 00189 00190 /// the number of exclusive jets at the given ycut 00191 int n_exclusive_jets_ycut (double ycut) const {return n_exclusive_jets(ycut*Q2());} 00192 00193 /// the exclusive jets obtained at the given ycut 00194 std::vector<PseudoJet> exclusive_jets_ycut (double ycut) const { 00195 int njets = n_exclusive_jets_ycut(ycut); 00196 return exclusive_jets(njets); 00197 } 00198 00199 00200 //int n_exclusive_jets (const PseudoJet & jet, const double & dcut) const; 00201 00202 /// return a vector of all subjets of the current jet (in the sense 00203 /// of the exclusive algorithm) that would be obtained when running 00204 /// the algorithm with the given dcut. 00205 /// 00206 /// Time taken is O(m ln m), where m is the number of subjets that 00207 /// are found. If m gets to be of order of the total number of 00208 /// constituents in the jet, this could be substantially slower than 00209 /// just getting that list of constituents. 00210 std::vector<PseudoJet> exclusive_subjets (const PseudoJet & jet, 00211 const double & dcut) const; 00212 00213 /// return the size of exclusive_subjets(...); still n ln n with same 00214 /// coefficient, but marginally more efficient than manually taking 00215 /// exclusive_subjets.size() 00216 int n_exclusive_subjets(const PseudoJet & jet, 00217 const double & dcut) const; 00218 00219 /// return the list of subjets obtained by unclustering the supplied 00220 /// jet down to n subjets (or all constituents if there are fewer 00221 /// than n). 00222 /// 00223 /// This requires n ln n time 00224 /// 00225 /// If the jet contains fewer than nsub particles (in which case it 00226 /// is not possible to return nsub subjets) then the vector of 00227 /// subjets that is returned is simply the list of particles (of 00228 /// size < nsub). Note that this behaviour differs from that of 00229 /// exclusive_jets(). 00230 std::vector<PseudoJet> exclusive_subjets (const PseudoJet & jet, 00231 int nsub) const; 00232 00233 /// return the dij that was present in the merging nsub+1 -> nsub 00234 /// subjets inside this jet. 00235 /// 00236 /// Returns 0 if there were nsub or fewer constituents in the jet. 00237 double exclusive_subdmerge(const PseudoJet & jet, int nsub) const; 00238 00239 /// return the maximum dij that occurred in the whole event at the 00240 /// stage that the nsub+1 -> nsub merge of subjets occurred inside 00241 /// this jet. 00242 /// 00243 /// Returns 0 if there were nsub or fewer constituents in the jet. 00244 double exclusive_subdmerge_max(const PseudoJet & jet, int nsub) const; 00245 00246 //std::vector<PseudoJet> exclusive_jets (const PseudoJet & jet, 00247 // const int & njets) const; 00248 //double exclusive_dmerge (const PseudoJet & jet, const int & njets) const; 00249 00250 /// returns the sum of all energies in the event (relevant mainly for e+e-) 00251 double Q() const {return _Qtot;} 00252 /// return Q()^2 00253 double Q2() const {return _Qtot*_Qtot;} 00254 00255 /// returns true iff the object is included in the jet. 00256 /// 00257 /// NB: this is only sensible if the object is already registered 00258 /// within the cluster sequence, so you cannot use it with an input 00259 /// particle to the CS (since the particle won't have the history 00260 /// index set properly). 00261 /// 00262 /// For nice clustering structures it should run in O(ln(N)) time 00263 /// but in worst cases (certain cone plugins) it can take O(n) time, 00264 /// where n is the number of particles in the jet. 00265 bool object_in_jet(const PseudoJet & object, const PseudoJet & jet) const; 00266 00267 /// if the jet has parents in the clustering, it returns true 00268 /// and sets parent1 and parent2 equal to them. 00269 /// 00270 /// if it has no parents it returns false and sets parent1 and 00271 /// parent2 to zero 00272 bool has_parents(const PseudoJet & jet, PseudoJet & parent1, 00273 PseudoJet & parent2) const; 00274 00275 /// if the jet has a child then return true and give the child jet 00276 /// otherwise return false and set the child to zero 00277 bool has_child(const PseudoJet & jet, PseudoJet & child) const; 00278 00279 /// Version of has_child that sets a pointer to the child if the child 00280 /// exists; 00281 bool has_child(const PseudoJet & jet, const PseudoJet * & childp) const; 00282 00283 /// if this jet has a child (and so a partner) return true 00284 /// and give the partner, otherwise return false and set the 00285 /// partner to zero 00286 bool has_partner(const PseudoJet & jet, PseudoJet & partner) const; 00287 00288 00289 /// return a vector of the particles that make up jet 00290 std::vector<PseudoJet> constituents (const PseudoJet & jet) const; 00291 00292 00293 /// output the supplied vector of jets in a format that can be read 00294 /// by an appropriate root script; the format is: 00295 /// jet-n jet-px jet-py jet-pz jet-E 00296 /// particle-n particle-rap particle-phi particle-pt 00297 /// particle-n particle-rap particle-phi particle-pt 00298 /// ... 00299 /// #END 00300 /// ... [i.e. above repeated] 00301 void print_jets_for_root(const std::vector<PseudoJet> & jets, 00302 std::ostream & ostr = std::cout) const; 00303 00304 /// print jets for root to the file labelled filename, with an 00305 /// optional comment at the beginning 00306 void print_jets_for_root(const std::vector<PseudoJet> & jets, 00307 const std::string & filename, 00308 const std::string & comment = "") const; 00309 00310 // Not yet. Perhaps in a future release. 00311 // /// print out all inclusive jets with pt > ptmin 00312 // virtual void print_jets (const double & ptmin=0.0) const; 00313 00314 /// add on to subjet_vector the constituents of jet (for internal use mainly) 00315 void add_constituents (const PseudoJet & jet, 00316 std::vector<PseudoJet> & subjet_vector) const; 00317 00318 /// return the enum value of the strategy used to cluster the event 00319 inline Strategy strategy_used () const {return _strategy;} 00320 00321 /// return the name of the strategy used to cluster the event 00322 std::string strategy_string () const {return strategy_string(_strategy);} 00323 00324 /// return the name of the strategy associated with the enum strategy_in 00325 std::string strategy_string (Strategy strategy_in) const; 00326 00327 00328 /// return a reference to the jet definition 00329 const JetDefinition & jet_def() const {return _jet_def;} 00330 00331 /// returns the scale associated with a jet as required for this 00332 /// clustering algorithm (kt^2 for the kt-algorithm, 1 for the 00333 /// Cambridge algorithm). [May become virtual at some point] 00334 double jet_scale_for_algorithm(const PseudoJet & jet) const; 00335 00336 //----- next follow functions designed specifically for plugins, which 00337 // may only be called when plugin_activated() returns true 00338 00339 /// record the fact that there has been a recombination between 00340 /// jets()[jet_i] and jets()[jet_k], with the specified dij, and 00341 /// return the index (newjet_k) allocated to the new jet, whose 00342 /// momentum is assumed to be the 4-vector sum of that of jet_i and 00343 /// jet_j 00344 void plugin_record_ij_recombination(int jet_i, int jet_j, double dij, 00345 int & newjet_k) { 00346 assert(plugin_activated()); 00347 _do_ij_recombination_step(jet_i, jet_j, dij, newjet_k); 00348 } 00349 00350 /// as for the simpler variant of plugin_record_ij_recombination, 00351 /// except that the new jet is attributed the momentum and 00352 /// user_index of newjet 00353 void plugin_record_ij_recombination(int jet_i, int jet_j, double dij, 00354 const PseudoJet & newjet, 00355 int & newjet_k); 00356 00357 /// record the fact that there has been a recombination between 00358 /// jets()[jet_i] and the beam, with the specified diB; when looking 00359 /// for inclusive jets, any iB recombination will returned to the user 00360 /// as a jet. 00361 void plugin_record_iB_recombination(int jet_i, double diB) { 00362 assert(plugin_activated()); 00363 _do_iB_recombination_step(jet_i, diB); 00364 } 00365 00366 /// @ingroup extra_info 00367 /// \class Extras 00368 /// base class to store extra information that plugins may provide 00369 /// 00370 /// a class intended to serve as a base in case a plugin needs to 00371 /// associate extra information with a ClusterSequence (see 00372 /// SISConePlugin.* for an example). 00373 class Extras { 00374 public: 00375 virtual ~Extras() {} 00376 virtual std::string description() const {return "This is a dummy extras class that contains no extra information! Derive from it if you want to use it to provide extra information from a plugin jet finder";} 00377 }; 00378 00379 /// the plugin can associated some extra information with the 00380 /// ClusterSequence object by calling this function 00381 inline void plugin_associate_extras(std::auto_ptr<Extras> extras_in) { 00382 _extras = extras_in; 00383 } 00384 00385 /// returns true when the plugin is allowed to run the show. 00386 inline bool plugin_activated() const {return _plugin_activated;} 00387 00388 /// returns a pointer to the extras object (may be null) 00389 const Extras * extras() const {return _extras.get();} 00390 00391 /// allows a plugin to run a templated clustering (nearest-neighbour heuristic) 00392 /// 00393 /// This has N^2 behaviour on "good" distance, but a worst case behaviour 00394 /// of N^3 (and many algs trigger the worst case behaviour) 00395 /// 00396 /// 00397 /// For more details on how this works, see GenBriefJet below 00398 template<class GBJ> void plugin_simple_N2_cluster () { 00399 assert(plugin_activated()); 00400 _simple_N2_cluster<GBJ>(); 00401 } 00402 00403 00404 public: 00405 /// set the default (static) jet finder across all current and future 00406 /// ClusterSequence objects -- deprecated and obsolescent (i.e. may be 00407 /// suppressed in a future release). 00408 static void set_jet_algorithm (JetAlgorithm jet_algorithm) {_default_jet_algorithm = jet_algorithm;} 00409 /// same as above for backward compatibility 00410 static void set_jet_finder (JetAlgorithm jet_algorithm) {_default_jet_algorithm = jet_algorithm;} 00411 00412 00413 /// \ingroup extra_info 00414 /// \struct history_element 00415 /// a single element in the clustering history 00416 /// 00417 /// (see vector _history below). 00418 struct history_element{ 00419 int parent1; /// index in _history where first parent of this jet 00420 /// was created (InexistentParent if this jet is an 00421 /// original particle) 00422 00423 int parent2; /// index in _history where second parent of this jet 00424 /// was created (InexistentParent if this jet is an 00425 /// original particle); BeamJet if this history entry 00426 /// just labels the fact that the jet has recombined 00427 /// with the beam) 00428 00429 int child; /// index in _history where the current jet is 00430 /// recombined with another jet to form its child. It 00431 /// is Invalid if this jet does not further 00432 /// recombine. 00433 00434 int jetp_index; /// index in the _jets vector where we will find the 00435 /// PseudoJet object corresponding to this jet 00436 /// (i.e. the jet created at this entry of the 00437 /// history). NB: if this element of the history 00438 /// corresponds to a beam recombination, then 00439 /// jetp_index=Invalid. 00440 00441 double dij; /// the distance corresponding to the recombination 00442 /// at this stage of the clustering. 00443 00444 double max_dij_so_far; /// the largest recombination distance seen 00445 /// so far in the clustering history. 00446 }; 00447 00448 enum JetType {Invalid=-3, InexistentParent = -2, BeamJet = -1}; 00449 00450 /// allow the user to access the internally stored _jets() array, 00451 /// which contains both the initial particles and the various 00452 /// intermediate and final stages of recombination. 00453 /// 00454 /// The first n_particles() entries are the original particles, 00455 /// in the order in which they were supplied to the ClusterSequence 00456 /// constructor. It can be useful to access them for example when 00457 /// examining whether a given input object is part of a specific 00458 /// jet, via the objects_in_jet(...) member function (which only takes 00459 /// PseudoJets that are registered in the ClusterSequence). 00460 /// 00461 /// One of the other (internal uses) is related to the fact 00462 /// because we don't seem to be able to access protected elements of 00463 /// the class for an object that is not "this" (at least in case where 00464 /// "this" is of a slightly different kind from the object, both 00465 /// derived from ClusterSequence). 00466 const std::vector<PseudoJet> & jets() const; 00467 00468 /// allow the user to access the raw internal history. 00469 /// 00470 /// This is present (as for jets()) in part so that protected 00471 /// derived classes can access this information about other 00472 /// ClusterSequences. 00473 /// 00474 /// A user who wishes to follow the details of the ClusterSequence 00475 /// can also make use of this information (and should consult the 00476 /// history_element documentation for more information), but should 00477 /// be aware that these internal structures may evolve in future 00478 /// FastJet versions. 00479 const std::vector<history_element> & history() const; 00480 00481 /// returns the number of particles that were provided to the 00482 /// clustering algorithm (helps the user find their way around the 00483 /// history and jets objects if they weren't paying attention 00484 /// beforehand). 00485 unsigned int n_particles() const; 00486 00487 /// returns a vector of size n_particles() which indicates, for 00488 /// each of the initial particles (in the order in which they were 00489 /// supplied), which of the supplied jets it belongs to; if it does 00490 /// not belong to any of the supplied jets, the index is set to -1; 00491 std::vector<int> particle_jet_indices(const std::vector<PseudoJet> &) const; 00492 00493 /// routine that returns an order in which to read the history 00494 /// such that clusterings that lead to identical jet compositions 00495 /// but different histories (because of degeneracies in the 00496 /// clustering order) will have matching constituents for each 00497 /// matching entry in the unique_history_order. 00498 /// 00499 /// The order has the property that an entry's parents will always 00500 /// appear prior to that entry itself. 00501 /// 00502 /// Roughly speaking the order is such that we first provide all 00503 /// steps that lead to the final jet containing particle 1; then we 00504 /// have the steps that lead to reconstruction of the jet containing 00505 /// the next-lowest-numbered unclustered particle, etc... 00506 /// [see GPS CCN28-12 for more info -- of course a full explanation 00507 /// here would be better...] 00508 std::vector<int> unique_history_order() const; 00509 00510 /// return the set of particles that have not been clustered. For 00511 /// kt and cam/aachen algorithms this should always be null, but for 00512 /// cone type algorithms it can be non-null; 00513 std::vector<PseudoJet> unclustered_particles() const; 00514 00515 /// returns true if the object (jet or particle) is contained by (ie 00516 /// belongs to) this cluster sequence. 00517 /// 00518 /// Tests performed: if thejet's interface is this cluster sequence 00519 /// and its cluster history index is in a consistent range. 00520 bool contains(const PseudoJet & object) const; 00521 00522 /// transfer the sequence contained in other_seq into our own; 00523 /// any plugin "extras" contained in the from_seq will be lost 00524 /// from there. 00525 /// 00526 /// If transfer_ownership is true, it also sets the ClusterSequence 00527 /// pointers of the PseudoJets in the history to point to this 00528 /// ClusterSequence (true by default) 00529 void transfer_from_sequence(ClusterSequence & from_seq, bool transfer_ownership=true); 00530 00531 /// retrieve a shared pointer to the wrapper to this ClusterSequence 00532 /// 00533 /// this may turn useful if you want to track when this 00534 /// ClusterSequence goes out of scope 00535 const SharedPtr<PseudoJetStructureBase> & structure_shared_ptr() const{ 00536 return _structure_shared_ptr; 00537 } 00538 00539 /// the structure type associated with a jet belonging to a ClusterSequence 00540 typedef ClusterSequenceStructure StructureType; 00541 00542 00543 protected: 00544 static JetAlgorithm _default_jet_algorithm; 00545 JetDefinition _jet_def; 00546 00547 /// transfer the vector<L> of input jets into our own vector<PseudoJet> 00548 /// _jets (with some reserved space for future growth). 00549 template<class L> void _transfer_input_jets( 00550 const std::vector<L> & pseudojets); 00551 00552 /// This is the routine that will do all the initialisation and 00553 /// then run the clustering (may be called by various constructors). 00554 /// It assumes _jets contains the momenta to be clustered. 00555 void _initialise_and_run (const JetDefinition & jet_def, 00556 const bool & writeout_combinations); 00557 00558 /// This is an alternative routine for initialising and running the 00559 /// clustering, provided for legacy purposes. The jet finder is that 00560 /// specified in the static member _default_jet_algorithm. 00561 void _initialise_and_run (const double & R, 00562 const Strategy & strategy, 00563 const bool & writeout_combinations); 00564 00565 /// fills in the various member variables with "decanted" options from 00566 /// the jet_definition and writeout_combinations variables 00567 void _decant_options(const JetDefinition & jet_def, 00568 const bool & writeout_combinations); 00569 00570 /// fill out the history (and jet cross refs) related to the initial 00571 /// set of jets (assumed already to have been "transferred"), 00572 /// without any clustering 00573 void _fill_initial_history(); 00574 00575 /// carry out the recombination between the jets numbered jet_i and 00576 /// jet_j, at distance scale dij; return the index newjet_k of the 00577 /// result of the recombination of i and j. 00578 void _do_ij_recombination_step(const int & jet_i, const int & jet_j, 00579 const double & dij, int & newjet_k); 00580 00581 /// carry out an recombination step in which _jets[jet_i] merges with 00582 /// the beam, 00583 void _do_iB_recombination_step(const int & jet_i, const double & diB); 00584 00585 00586 /// This contains the physical PseudoJets; for each PseudoJet one 00587 /// can find the corresponding position in the _history by looking 00588 /// at _jets[i].cluster_hist_index(). 00589 std::vector<PseudoJet> _jets; 00590 00591 00592 /// this vector will contain the branching history; for each stage, 00593 /// _history[i].jetp_index indicates where to look in the _jets 00594 /// vector to get the physical PseudoJet. 00595 std::vector<history_element> _history; 00596 00597 /// set subhist to be a set pointers to history entries corresponding to the 00598 /// subjets of this jet; one stops going working down through the 00599 /// subjets either when 00600 /// - there is no further to go 00601 /// - one has found maxjet entries 00602 /// - max_dij_so_far <= dcut 00603 /// By setting maxjet=0 one can use just dcut; by setting dcut<0 00604 /// one can use jet maxjet 00605 void get_subhist_set(std::set<const history_element*> & subhist, 00606 const PseudoJet & jet, double dcut, int maxjet) const; 00607 00608 bool _writeout_combinations; 00609 int _initial_n; 00610 double _Rparam, _R2, _invR2; 00611 double _Qtot; 00612 Strategy _strategy; 00613 JetAlgorithm _jet_algorithm; 00614 00615 SharedPtr<PseudoJetStructureBase> _structure_shared_ptr; //< will actually be of type ClusterSequenceStructure 00616 00617 private: 00618 00619 bool _plugin_activated; 00620 std::auto_ptr<Extras> _extras; // things the plugin might want to add 00621 00622 void _really_dumb_cluster (); 00623 void _delaunay_cluster (); 00624 //void _simple_N2_cluster (); 00625 template<class BJ> void _simple_N2_cluster (); 00626 void _tiled_N2_cluster (); 00627 void _faster_tiled_N2_cluster (); 00628 00629 // 00630 void _minheap_faster_tiled_N2_cluster(); 00631 00632 // things needed specifically for Cambridge with Chan's 2D closest 00633 // pairs method 00634 void _CP2DChan_cluster(); 00635 void _CP2DChan_cluster_2pi2R (); 00636 void _CP2DChan_cluster_2piMultD (); 00637 void _CP2DChan_limited_cluster(double D); 00638 void _do_Cambridge_inclusive_jets(); 00639 00640 // NSqrtN method for C/A 00641 void _fast_NsqrtN_cluster(); 00642 00643 void _add_step_to_history(const int & step_number, const int & parent1, 00644 const int & parent2, const int & jetp_index, 00645 const double & dij); 00646 00647 /// internal routine associated with the construction of the unique 00648 /// history order (following children in the tree) 00649 void _extract_tree_children(int pos, std::valarray<bool> &, 00650 const std::valarray<int> &, std::vector<int> &) const; 00651 00652 /// internal routine associated with the construction of the unique 00653 /// history order (following parents in the tree) 00654 void _extract_tree_parents (int pos, std::valarray<bool> &, 00655 const std::valarray<int> &, std::vector<int> &) const; 00656 00657 00658 // these will be useful shorthands in the Voronoi-based code 00659 typedef std::pair<int,int> TwoVertices; 00660 typedef std::pair<double,TwoVertices> DijEntry; 00661 typedef std::multimap<double,TwoVertices> DistMap; 00662 00663 /// currently used only in the Voronoi based code 00664 void _add_ktdistance_to_map(const int & ii, 00665 DistMap & DijMap, 00666 const DynamicNearestNeighbours * DNN); 00667 00668 /// for making sure the user knows what it is they're running... 00669 void _print_banner(); 00670 /// will be set by default to be true for the first run 00671 static bool _first_time; 00672 00673 /// record the number of warnings provided about the exclusive 00674 /// algorithm -- so that we don't print it out more than a few 00675 /// times. 00676 static int _n_exclusive_warnings; 00677 00678 /// the limited warning member for notification of user that 00679 /// their requested strategy has been overridden (usually because 00680 /// they have R>2pi and not all strategies work then) 00681 static LimitedWarning _changed_strategy_warning; 00682 00683 //---------------------------------------------------------------------- 00684 /// the fundamental structure which contains the minimal info about 00685 /// a jet, as needed for our plain N^2 algorithm -- the idea is to 00686 /// put all info that will be accessed N^2 times into an array of 00687 /// BriefJets... 00688 struct BriefJet { 00689 double eta, phi, kt2, NN_dist; 00690 BriefJet * NN; 00691 int _jets_index; 00692 }; 00693 00694 00695 /// structure analogous to BriefJet, but with the extra information 00696 /// needed for dealing with tiles 00697 class TiledJet { 00698 public: 00699 double eta, phi, kt2, NN_dist; 00700 TiledJet * NN, *previous, * next; 00701 int _jets_index, tile_index, diJ_posn; 00702 // routines that are useful in the minheap version of tiled 00703 // clustering ("misuse" the otherwise unused diJ_posn, so as 00704 // to indicate whether jets need to have their minheap entries 00705 // updated). 00706 inline void label_minheap_update_needed() {diJ_posn = 1;} 00707 inline void label_minheap_update_done() {diJ_posn = 0;} 00708 inline bool minheap_update_needed() const {return diJ_posn==1;} 00709 }; 00710 00711 //-- some of the functions that follow are templates and will work 00712 //as well for briefjet and tiled jets 00713 00714 /// set the kinematic and labelling info for jeta so that it corresponds 00715 /// to _jets[_jets_index] 00716 template <class J> void _bj_set_jetinfo( J * const jet, 00717 const int _jets_index) const; 00718 00719 /// "remove" this jet, which implies updating links of neighbours and 00720 /// perhaps modifying the tile structure 00721 void _bj_remove_from_tiles( TiledJet * const jet) const; 00722 00723 /// return the distance between two BriefJet objects 00724 template <class J> double _bj_dist(const J * const jeta, 00725 const J * const jetb) const; 00726 00727 // return the diJ (multiplied by _R2) for this jet assuming its NN 00728 // info is correct 00729 template <class J> double _bj_diJ(const J * const jeta) const; 00730 00731 /// for testing purposes only: if in the range head--tail-1 there is a 00732 /// a jet which corresponds to hist_index in the history, then 00733 /// return a pointer to that jet; otherwise return tail. 00734 template <class J> inline J * _bj_of_hindex( 00735 const int hist_index, 00736 J * const head, J * const tail) 00737 const { 00738 J * res; 00739 for(res = head; res<tail; res++) { 00740 if (_jets[res->_jets_index].cluster_hist_index() == hist_index) {break;} 00741 } 00742 return res; 00743 } 00744 00745 00746 //-- remaining functions are different in various cases, so we 00747 // will use templates but are not sure if they're useful... 00748 00749 /// updates (only towards smaller distances) the NN for jeta without checking 00750 /// whether in the process jeta itself might be a new NN of one of 00751 /// the jets being scanned -- span the range head to tail-1 with 00752 /// assumption that jeta is not contained in that range 00753 template <class J> void _bj_set_NN_nocross(J * const jeta, 00754 J * const head, const J * const tail) const; 00755 00756 /// reset the NN for jeta and DO check whether in the process jeta 00757 /// itself might be a new NN of one of the jets being scanned -- 00758 /// span the range head to tail-1 with assumption that jeta is not 00759 /// contained in that range 00760 template <class J> void _bj_set_NN_crosscheck(J * const jeta, 00761 J * const head, const J * const tail) const; 00762 00763 00764 00765 /// number of neighbours that a tile will have (rectangular geometry 00766 /// gives 9 neighbours). 00767 static const int n_tile_neighbours = 9; 00768 //---------------------------------------------------------------------- 00769 /// The fundamental structures to be used for the tiled N^2 algorithm 00770 /// (see CCN27-44 for some discussion of pattern of tiling) 00771 struct Tile { 00772 /// pointers to neighbouring tiles, including self 00773 Tile * begin_tiles[n_tile_neighbours]; 00774 /// neighbouring tiles, excluding self 00775 Tile ** surrounding_tiles; 00776 /// half of neighbouring tiles, no self 00777 Tile ** RH_tiles; 00778 /// just beyond end of tiles 00779 Tile ** end_tiles; 00780 /// start of list of BriefJets contained in this tile 00781 TiledJet * head; 00782 /// sometimes useful to be able to tag a tile 00783 bool tagged; 00784 }; 00785 std::vector<Tile> _tiles; 00786 double _tiles_eta_min, _tiles_eta_max; 00787 double _tile_size_eta, _tile_size_phi; 00788 int _n_tiles_phi,_tiles_ieta_min,_tiles_ieta_max; 00789 00790 // reasonably robust return of tile index given ieta and iphi, in particular 00791 // it works even if iphi is negative 00792 inline int _tile_index (int ieta, int iphi) const { 00793 // note that (-1)%n = -1 so that we have to add _n_tiles_phi 00794 // before performing modulo operation 00795 return (ieta-_tiles_ieta_min)*_n_tiles_phi 00796 + (iphi+_n_tiles_phi) % _n_tiles_phi; 00797 } 00798 00799 // routines for tiled case, including some overloads of the plain 00800 // BriefJet cases 00801 int _tile_index(const double & eta, const double & phi) const; 00802 void _tj_set_jetinfo ( TiledJet * const jet, const int _jets_index); 00803 void _bj_remove_from_tiles(TiledJet * const jet); 00804 void _initialise_tiles(); 00805 void _print_tiles(TiledJet * briefjets ) const; 00806 void _add_neighbours_to_tile_union(const int tile_index, 00807 std::vector<int> & tile_union, int & n_near_tiles) const; 00808 void _add_untagged_neighbours_to_tile_union(const int tile_index, 00809 std::vector<int> & tile_union, int & n_near_tiles); 00810 00811 00812 //---------------------------------------------------------------------- 00813 /// fundamental structure for e+e- clustering 00814 struct EEBriefJet { 00815 double NN_dist; // obligatorily present 00816 double kt2; // obligatorily present == E^2 in general 00817 EEBriefJet * NN; // must be present too 00818 int _jets_index; // must also be present! 00819 //........................................................... 00820 double nx, ny, nz; // our internal storage for fast distance calcs 00821 }; 00822 00823 /// to help instantiation (fj 2.4.0; did not quite work on gcc 33 and os x 10.3?) 00824 //void _dummy_N2_cluster_instantiation(); 00825 00826 00827 /// to avoid issues with template instantiation (OS X 10.3, gcc 3.3) 00828 void _simple_N2_cluster_BriefJet(); 00829 /// to avoid issues with template instantiation (OS X 10.3, gcc 3.3) 00830 void _simple_N2_cluster_EEBriefJet(); 00831 }; 00832 00833 00834 //********************************************************************** 00835 //************** START OF INLINE MATERIAL ****************** 00836 //********************************************************************** 00837 00838 00839 //---------------------------------------------------------------------- 00840 // Transfer the initial jets into our internal structure 00841 template<class L> void ClusterSequence::_transfer_input_jets( 00842 const std::vector<L> & pseudojets) { 00843 00844 // this will ensure that we can point to jets without difficulties 00845 // arising. 00846 _jets.reserve(pseudojets.size()*2); 00847 00848 // insert initial jets this way so that any type L that can be 00849 // converted to a pseudojet will work fine (basically PseudoJet 00850 // and any type that has [] subscript access to the momentum 00851 // components, such as CLHEP HepLorentzVector). 00852 for (unsigned int i = 0; i < pseudojets.size(); i++) { 00853 _jets.push_back(pseudojets[i]);} 00854 00855 } 00856 00857 //---------------------------------------------------------------------- 00858 // initialise from some generic type... Has to be made available 00859 // here in order for it the template aspect of it to work... 00860 template<class L> ClusterSequence::ClusterSequence ( 00861 const std::vector<L> & pseudojets, 00862 const double & R, 00863 const Strategy & strategy, 00864 const bool & writeout_combinations) { 00865 00866 // transfer the initial jets (type L) into our own array 00867 _transfer_input_jets(pseudojets); 00868 00869 // run the clustering 00870 _initialise_and_run(R,strategy,writeout_combinations); 00871 } 00872 00873 00874 //---------------------------------------------------------------------- 00875 /// constructor of a jet-clustering sequence from a vector of 00876 /// four-momenta, with the jet definition specified by jet_def 00877 template<class L> ClusterSequence::ClusterSequence ( 00878 const std::vector<L> & pseudojets, 00879 const JetDefinition & jet_def, 00880 const bool & writeout_combinations) { 00881 00882 // transfer the initial jets (type L) into our own array 00883 _transfer_input_jets(pseudojets); 00884 00885 // run the clustering 00886 _initialise_and_run(jet_def,writeout_combinations); 00887 } 00888 00889 00890 inline const std::vector<PseudoJet> & ClusterSequence::jets () const { 00891 return _jets; 00892 } 00893 00894 inline const std::vector<ClusterSequence::history_element> & ClusterSequence::history () const { 00895 return _history; 00896 } 00897 00898 inline unsigned int ClusterSequence::n_particles() const {return _initial_n;} 00899 00900 00901 00902 //---------------------------------------------------------------------- 00903 template <class J> inline void ClusterSequence::_bj_set_jetinfo( 00904 J * const jetA, const int _jets_index) const { 00905 jetA->eta = _jets[_jets_index].rap(); 00906 jetA->phi = _jets[_jets_index].phi_02pi(); 00907 jetA->kt2 = jet_scale_for_algorithm(_jets[_jets_index]); 00908 jetA->_jets_index = _jets_index; 00909 // initialise NN info as well 00910 jetA->NN_dist = _R2; 00911 jetA->NN = NULL; 00912 } 00913 00914 00915 00916 00917 //---------------------------------------------------------------------- 00918 template <class J> inline double ClusterSequence::_bj_dist( 00919 const J * const jetA, const J * const jetB) const { 00920 double dphi = std::abs(jetA->phi - jetB->phi); 00921 double deta = (jetA->eta - jetB->eta); 00922 if (dphi > pi) {dphi = twopi - dphi;} 00923 return dphi*dphi + deta*deta; 00924 } 00925 00926 //---------------------------------------------------------------------- 00927 template <class J> inline double ClusterSequence::_bj_diJ(const J * const jet) const { 00928 double kt2 = jet->kt2; 00929 if (jet->NN != NULL) {if (jet->NN->kt2 < kt2) {kt2 = jet->NN->kt2;}} 00930 return jet->NN_dist * kt2; 00931 } 00932 00933 00934 //---------------------------------------------------------------------- 00935 // set the NN for jet without checking whether in the process you might 00936 // have discovered a new nearest neighbour for another jet 00937 template <class J> inline void ClusterSequence::_bj_set_NN_nocross( 00938 J * const jet, J * const head, const J * const tail) const { 00939 double NN_dist = _R2; 00940 J * NN = NULL; 00941 if (head < jet) { 00942 for (J * jetB = head; jetB != jet; jetB++) { 00943 double dist = _bj_dist(jet,jetB); 00944 if (dist < NN_dist) { 00945 NN_dist = dist; 00946 NN = jetB; 00947 } 00948 } 00949 } 00950 if (tail > jet) { 00951 for (J * jetB = jet+1; jetB != tail; jetB++) { 00952 double dist = _bj_dist(jet,jetB); 00953 if (dist < NN_dist) { 00954 NN_dist = dist; 00955 NN = jetB; 00956 } 00957 } 00958 } 00959 jet->NN = NN; 00960 jet->NN_dist = NN_dist; 00961 } 00962 00963 00964 //---------------------------------------------------------------------- 00965 template <class J> inline void ClusterSequence::_bj_set_NN_crosscheck(J * const jet, 00966 J * const head, const J * const tail) const { 00967 double NN_dist = _R2; 00968 J * NN = NULL; 00969 for (J * jetB = head; jetB != tail; jetB++) { 00970 double dist = _bj_dist(jet,jetB); 00971 if (dist < NN_dist) { 00972 NN_dist = dist; 00973 NN = jetB; 00974 } 00975 if (dist < jetB->NN_dist) { 00976 jetB->NN_dist = dist; 00977 jetB->NN = jet; 00978 } 00979 } 00980 jet->NN = NN; 00981 jet->NN_dist = NN_dist; 00982 } 00983 00984 00985 00986 FASTJET_END_NAMESPACE 00987 00988 #endif // __FASTJET_CLUSTERSEQUENCE_HH__