FastJet 3.4.1
ClusterSequence.hh
1#ifndef __FASTJET_CLUSTERSEQUENCE_HH__
2#define __FASTJET_CLUSTERSEQUENCE_HH__
3
4//FJSTARTHEADER
5// $Id$
6//
7// Copyright (c) 2005-2023, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
8//
9//----------------------------------------------------------------------
10// This file is part of FastJet.
11//
12// FastJet is free software; you can redistribute it and/or modify
13// it under the terms of the GNU General Public License as published by
14// the Free Software Foundation; either version 2 of the License, or
15// (at your option) any later version.
16//
17// The algorithms that underlie FastJet have required considerable
18// development. They are described in the original FastJet paper,
19// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
20// FastJet as part of work towards a scientific publication, please
21// quote the version you use and include a citation to the manual and
22// optionally also to hep-ph/0512210.
23//
24// FastJet is distributed in the hope that it will be useful,
25// but WITHOUT ANY WARRANTY; without even the implied warranty of
26// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
27// GNU General Public License for more details.
28//
29// You should have received a copy of the GNU General Public License
30// along with FastJet. If not, see <http://www.gnu.org/licenses/>.
31//----------------------------------------------------------------------
32//FJENDHEADER
33
34
35#include<vector>
36#include<map>
37#include "fastjet/PseudoJet.hh"
38#include<memory>
39#include<cassert>
40#include<iostream>
41#include<string>
42#include<set>
43#include<cmath> // needed to get double std::abs(double)
44#include "fastjet/Error.hh"
45#include "fastjet/JetDefinition.hh"
46#include "fastjet/SharedPtr.hh"
47#include "fastjet/LimitedWarning.hh"
48#include "fastjet/FunctionOfPseudoJet.hh"
49#include "fastjet/ClusterSequenceStructure.hh"
50#include "fastjet/internal/thread_safety_helpers.hh" // helpers to write code w&wo thread-safety
51#include "fastjet/internal/deprecated.hh"
52
53FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
54
55
56// forward declarations
57class ClusterSequenceStructure;
58class DynamicNearestNeighbours;
59
60/// @ingroup basic_classes
61/// \class ClusterSequence
62/// deals with clustering
64
65
66 public:
67
68 /// default constructor
69 ClusterSequence () : _deletes_self_when_unused(false) {}
70
71 /// create a ClusterSequence, starting from the supplied set
72 /// of PseudoJets and clustering them with jet definition specified
73 /// by jet_def (which also specifies the clustering strategy)
74 template<class L> ClusterSequence (
75 const std::vector<L> & pseudojets,
76 const JetDefinition & jet_def,
77 const bool & writeout_combinations = false);
78
79 /// copy constructor for a ClusterSequence
80 ClusterSequence (const ClusterSequence & cs) : _deletes_self_when_unused(false) {
81 transfer_from_sequence(cs);
82 }
83
84 /// explicit assignment operator for a ClusterSequence
85 ClusterSequence & operator=(const ClusterSequence & cs);
86
87 // virtual ClusterSequence destructor, in case any derived class
88 // thinks of needing a destructor at some point
89 virtual ~ClusterSequence (); //{}
90
91 // NB: in the routines that follow, for extracting lists of jets, a
92 // list structure might be more efficient, if sometimes a little
93 // more awkward to use (at least for old fortran hands).
94
95 /// return a vector of all jets (in the sense of the inclusive
96 /// algorithm) with pt >= ptmin. Time taken should be of the order
97 /// of the number of jets returned.
98 std::vector<PseudoJet> inclusive_jets (const double ptmin = 0.0) const;
99
100 /// return the number of jets (in the sense of the exclusive
101 /// algorithm) that would be obtained when running the algorithm
102 /// with the given dcut.
103 int n_exclusive_jets (const double dcut) const;
104
105 /// return a vector of all jets (in the sense of the exclusive
106 /// algorithm) that would be obtained when running the algorithm
107 /// with the given dcut.
108 std::vector<PseudoJet> exclusive_jets (const double dcut) const;
109
110 /// return a vector of all jets when the event is clustered (in the
111 /// exclusive sense) to exactly njets.
112 ///
113 /// If there are fewer than njets particles in the ClusterSequence
114 /// an error is thrown
115 std::vector<PseudoJet> exclusive_jets (const int njets) const;
116
117 /// return a vector of all jets when the event is clustered (in the
118 /// exclusive sense) to exactly njets.
119 ///
120 /// If there are fewer than njets particles in the ClusterSequence
121 /// the function just returns however many particles there were.
122 std::vector<PseudoJet> exclusive_jets_up_to (const int njets) const;
123
124 /// return the dmin corresponding to the recombination that went
125 /// from n+1 to n jets (sometimes known as d_{n n+1}). If the number
126 /// of particles in the event is <= njets, the function returns 0.
127 double exclusive_dmerge (const int njets) const;
128
129 /// return the maximum of the dmin encountered during all recombinations
130 /// up to the one that led to an n-jet final state; identical to
131 /// exclusive_dmerge, except in cases where the dmin do not increase
132 /// monotonically.
133 double exclusive_dmerge_max (const int njets) const;
134
135 /// return the ymin corresponding to the recombination that went from
136 /// n+1 to n jets (sometimes known as y_{n n+1}).
137 double exclusive_ymerge (int njets) const {return exclusive_dmerge(njets) / Q2();}
138
139 /// same as exclusive_dmerge_max, but normalised to squared total energy
140 double exclusive_ymerge_max (int njets) const {return exclusive_dmerge_max(njets)/Q2();}
141
142 /// the number of exclusive jets at the given ycut
143 int n_exclusive_jets_ycut (double ycut) const {return n_exclusive_jets(ycut*Q2());}
144
145 /// the exclusive jets obtained at the given ycut
146 std::vector<PseudoJet> exclusive_jets_ycut (double ycut) const {
147 int njets = n_exclusive_jets_ycut(ycut);
148 return exclusive_jets(njets);
149 }
150
151
152 //int n_exclusive_jets (const PseudoJet & jet, const double dcut) const;
153
154 /// return a vector of all subjets of the current jet (in the sense
155 /// of the exclusive algorithm) that would be obtained when running
156 /// the algorithm with the given dcut.
157 ///
158 /// Time taken is O(m ln m), where m is the number of subjets that
159 /// are found. If m gets to be of order of the total number of
160 /// constituents in the jet, this could be substantially slower than
161 /// just getting that list of constituents.
162 std::vector<PseudoJet> exclusive_subjets (const PseudoJet & jet,
163 const double dcut) const;
164
165 /// return the size of exclusive_subjets(...); still n ln n with same
166 /// coefficient, but marginally more efficient than manually taking
167 /// exclusive_subjets.size()
168 int n_exclusive_subjets(const PseudoJet & jet,
169 const double dcut) const;
170
171 /// return the list of subjets obtained by unclustering the supplied
172 /// jet down to nsub subjets. Throws an error if there are fewer than
173 /// nsub particles in the jet.
174 ///
175 /// This requires nsub ln nsub time
176 std::vector<PseudoJet> exclusive_subjets (const PseudoJet & jet,
177 int nsub) const;
178
179 /// return the list of subjets obtained by unclustering the supplied
180 /// jet down to nsub subjets (or all constituents if there are fewer
181 /// than nsub).
182 ///
183 /// This requires nsub ln nsub time
184 std::vector<PseudoJet> exclusive_subjets_up_to (const PseudoJet & jet,
185 int nsub) const;
186
187 /// returns the dij that was present in the merging nsub+1 -> nsub
188 /// subjets inside this jet.
189 ///
190 /// Returns 0 if there were nsub or fewer constituents in the jet.
191 double exclusive_subdmerge(const PseudoJet & jet, int nsub) const;
192
193 /// returns the maximum dij that occurred in the whole event at the
194 /// stage that the nsub+1 -> nsub merge of subjets occurred inside
195 /// this jet.
196 ///
197 /// Returns 0 if there were nsub or fewer constituents in the jet.
198 double exclusive_subdmerge_max(const PseudoJet & jet, int nsub) const;
199
200 //std::vector<PseudoJet> exclusive_jets (const PseudoJet & jet,
201 // const int njets) const;
202 //double exclusive_dmerge (const PseudoJet & jet, const int njets) const;
203
204 /// returns the sum of all energies in the event (relevant mainly for e+e-)
205 double Q() const {return _Qtot;}
206 /// return Q()^2
207 double Q2() const {return _Qtot*_Qtot;}
208
209 /// returns true iff the object is included in the jet.
210 ///
211 /// NB: this is only sensible if the object is already registered
212 /// within the cluster sequence, so you cannot use it with an input
213 /// particle to the CS (since the particle won't have the history
214 /// index set properly).
215 ///
216 /// For nice clustering structures it should run in O(ln(N)) time
217 /// but in worst cases (certain cone plugins) it can take O(n) time,
218 /// where n is the number of particles in the jet.
219 bool object_in_jet(const PseudoJet & object, const PseudoJet & jet) const;
220
221 /// if the jet has parents in the clustering, it returns true
222 /// and sets parent1 and parent2 equal to them.
223 ///
224 /// if it has no parents it returns false and sets parent1 and
225 /// parent2 to zero
226 bool has_parents(const PseudoJet & jet, PseudoJet & parent1,
227 PseudoJet & parent2) const;
228
229 /// if the jet has a child then return true and give the child jet
230 /// otherwise return false and set the child to zero
231 bool has_child(const PseudoJet & jet, PseudoJet & child) const;
232
233 /// Version of has_child that sets a pointer to the child if the child
234 /// exists;
235 bool has_child(const PseudoJet & jet, const PseudoJet * & childp) const;
236
237 /// if this jet has a child (and so a partner) return true
238 /// and give the partner, otherwise return false and set the
239 /// partner to zero
240 bool has_partner(const PseudoJet & jet, PseudoJet & partner) const;
241
242
243 /// return a vector of the particles that make up jet
244 std::vector<PseudoJet> constituents (const PseudoJet & jet) const;
245
246
247 /// output the supplied vector of jets in a format that can be read
248 /// by an appropriate root script; the format is:
249 /// jet-n jet-px jet-py jet-pz jet-E
250 /// particle-n particle-rap particle-phi particle-pt
251 /// particle-n particle-rap particle-phi particle-pt
252 /// ...
253 /// #END
254 /// ... [i.e. above repeated]
255 void print_jets_for_root(const std::vector<PseudoJet> & jets,
256 std::ostream & ostr = std::cout) const;
257
258 /// print jets for root to the file labelled filename, with an
259 /// optional comment at the beginning
260 void print_jets_for_root(const std::vector<PseudoJet> & jets,
261 const std::string & filename,
262 const std::string & comment = "") const;
263
264// Not yet. Perhaps in a future release.
265// /// print out all inclusive jets with pt > ptmin
266// virtual void print_jets (const double ptmin=0.0) const;
267
268 /// add on to subjet_vector the constituents of jet (for internal use mainly)
269 void add_constituents (const PseudoJet & jet,
270 std::vector<PseudoJet> & subjet_vector) const;
271
272 /// return the enum value of the strategy used to cluster the event
273 inline Strategy strategy_used () const {return _strategy;}
274
275 /// return the name of the strategy used to cluster the event
276 std::string strategy_string () const {return strategy_string(_strategy);}
277
278 /// return the name of the strategy associated with the enum strategy_in
279 std::string strategy_string (Strategy strategy_in) const;
280
281
282 /// return a reference to the jet definition
283 const JetDefinition & jet_def() const {return _jet_def;}
284
285 /// by calling this routine you tell the ClusterSequence to delete
286 /// itself when all the Pseudojets associated with it have gone out
287 /// of scope.
288 ///
289 /// At the time you call this, there must be at least one jet or
290 /// other object outside the CS that is associated with the CS
291 /// (e.g. the result of inclusive_jets()).
292 ///
293 /// NB: after having made this call, the user is still allowed to
294 /// delete the CS. Jets associated with it will then simply not be
295 /// able to access their substructure after that point.
296 void delete_self_when_unused();
297
298 /// return true if the object has been told to delete itself
299 /// when unused
300 bool will_delete_self_when_unused() const {return _deletes_self_when_unused;}
301
302//std::shared_ptr: #ifdef FASTJET_HAVE_THREAD_SAFETY
303//std::shared_ptr: /// signals that a jet will no longer use the current CS (internal use only)
304//std::shared_ptr: void release_pseudojet(PseudoJet &jet) const;
305//std::shared_ptr: #else
306 /// tell the ClusterSequence it's about to be self deleted (internal use only)
307 void signal_imminent_self_deletion() const;
308//std::shared_ptr: #endif
309
310 /// returns the scale associated with a jet as required for this
311 /// clustering algorithm (kt^2 for the kt-algorithm, 1 for the
312 /// Cambridge algorithm). Intended mainly for internal use and not
313 /// valid for plugin algorithms.
314 double jet_scale_for_algorithm(const PseudoJet & jet) const;
315
316 ///
317
318 //----- next follow functions designed specifically for plugins, which
319 // may only be called when plugin_activated() returns true
320
321 /// record the fact that there has been a recombination between
322 /// jets()[jet_i] and jets()[jet_k], with the specified dij, and
323 /// return the index (newjet_k) allocated to the new jet, whose
324 /// momentum is assumed to be the 4-vector sum of that of jet_i and
325 /// jet_j
326 void plugin_record_ij_recombination(int jet_i, int jet_j, double dij,
327 int & newjet_k) {
328 assert(plugin_activated());
329 _do_ij_recombination_step(jet_i, jet_j, dij, newjet_k);
330 }
331
332 /// as for the simpler variant of plugin_record_ij_recombination,
333 /// except that the new jet is attributed the momentum and
334 /// user_index of newjet
335 void plugin_record_ij_recombination(int jet_i, int jet_j, double dij,
336 const PseudoJet & newjet,
337 int & newjet_k);
338
339 /// record the fact that there has been a recombination between
340 /// jets()[jet_i] and the beam, with the specified diB; when looking
341 /// for inclusive jets, any iB recombination will returned to the user
342 /// as a jet.
343 void plugin_record_iB_recombination(int jet_i, double diB) {
344 assert(plugin_activated());
345 _do_iB_recombination_step(jet_i, diB);
346 }
347
348 /// return a non-const reference to the jets()[i], to allow plugins to
349 /// modify its contents.
350 ///
351 /// It can only be called when the plugin is activated.
352 ///
353 /// If you reset the jet (or set it equal to another one) you _must_
354 /// ensure that final jet is given the structure_shared_ptr of the
355 /// original jet, otherwise you will end up with an inconsistent
356 /// ClusterSequence.
357 ///
358 /// ONLY USE THIS IF YOU ARE SURE YOU KNOW WHAT YOU ARE DOING (contact
359 /// FJ authors if you think you need this but are unsure)
361 assert(plugin_activated());
362 return _jets[i];
363 }
364
365 /// @ingroup extra_info
366 /// \class Extras
367 /// base class to store extra information that plugins may provide
368 ///
369 /// a class intended to serve as a base in case a plugin needs to
370 /// associate extra information with a ClusterSequence (see
371 /// SISConePlugin.* for an example).
372 class Extras {
373 public:
374 virtual ~Extras() {}
375 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";}
376 };
377
378 /// the plugin can associate some extra information with the
379 /// ClusterSequence object by calling this function. The
380 /// ClusterSequence takes ownership of the pointer (and
381 /// responsibility for deleting it when the CS gets deleted).
382 inline void plugin_associate_extras(Extras * extras_in) {
383 _extras.reset(extras_in);
384 }
385
386 /// the plugin can associate some extra information with the
387 /// ClusterSequence object by calling this function
388 ///
389 /// As of FJ v3.1, this is deprecated, in line with the deprecation
390 /// of auto_ptr in C++11
391#ifndef SWIG
392#ifdef FASTJET_HAVE_AUTO_PTR_INTERFACE
393 FASTJET_DEPRECATED_MSG("Please use ClusterSequence::plugin_associate_extras(Extras * extras_in)) instead",
394 inline void plugin_associate_extras(std::auto_ptr<Extras> extras_in)){
395 _extras.reset(extras_in.release());
396 }
397#endif
398#endif //SWIG
399
400 /// returns true when the plugin is allowed to run the show.
401 inline bool plugin_activated() const {return _plugin_activated;}
402
403 /// returns a pointer to the extras object (may be null)
404 const Extras * extras() const {return _extras.get();}
405
406 /// allows a plugin to run a templated clustering (nearest-neighbour heuristic)
407 ///
408 /// This has N^2 behaviour on "good" distance, but a worst case behaviour
409 /// of N^3 (and many algs trigger the worst case behaviour)
410 ///
411 ///
412 /// For more details on how this works, see GenBriefJet below
413 template<class GBJ> void plugin_simple_N2_cluster () {
414 assert(plugin_activated());
415 _simple_N2_cluster<GBJ>();
416 }
417
418
419public:
420//DEP /// set the default (static) jet finder across all current and future
421//DEP /// ClusterSequence objects -- deprecated and obsolescent (i.e. may be
422//DEP /// suppressed in a future release).
423//DEP static void set_jet_algorithm (JetAlgorithm jet_algorithm) {_default_jet_algorithm = jet_algorithm;}
424//DEP /// same as above for backward compatibility
425//DEP static void set_jet_finder (JetAlgorithm jet_algorithm) {_default_jet_algorithm = jet_algorithm;}
426
427
428 /// \ingroup extra_info
429 /// \struct history_element
430 /// a single element in the clustering history
431 ///
432 /// (see vector _history below).
434 /// index in _history where first parent of this jet
435 /// was created (InexistentParent if this jet is an
436 /// original particle)
438
439 /// index in _history where second parent of this jet
440 /// was created (InexistentParent if this jet is an
441 /// original particle); BeamJet if this history entry
442 /// just labels the fact that the jet has recombined
443 /// with the beam)
445
446 /// index in _history where the current jet is
447 /// recombined with another jet to form its child. It
448 /// is Invalid if this jet does not further
449 /// recombine.
450 int child;
451
452 /// index in the _jets vector where we will find the
453 /// PseudoJet object corresponding to this jet
454 /// (i.e. the jet created at this entry of the
455 /// history). NB: if this element of the history
456 /// corresponds to a beam recombination, then
457 /// jetp_index=Invalid.
459
460 /// the distance corresponding to the recombination
461 /// at this stage of the clustering.
462 double dij;
463
464 /// the largest recombination distance seen
465 /// so far in the clustering history.
467 };
468
469 enum JetType {Invalid=-3, InexistentParent = -2, BeamJet = -1};
470
471 /// allow the user to access the internally stored _jets() array,
472 /// which contains both the initial particles and the various
473 /// intermediate and final stages of recombination.
474 ///
475 /// The first n_particles() entries are the original particles,
476 /// in the order in which they were supplied to the ClusterSequence
477 /// constructor. It can be useful to access them for example when
478 /// examining whether a given input object is part of a specific
479 /// jet, via the objects_in_jet(...) member function (which only takes
480 /// PseudoJets that are registered in the ClusterSequence).
481 ///
482 /// One of the other (internal uses) is related to the fact
483 /// because we don't seem to be able to access protected elements of
484 /// the class for an object that is not "this" (at least in case where
485 /// "this" is of a slightly different kind from the object, both
486 /// derived from ClusterSequence).
487 const std::vector<PseudoJet> & jets() const;
488
489 /// allow the user to access the raw internal history.
490 ///
491 /// This is present (as for jets()) in part so that protected
492 /// derived classes can access this information about other
493 /// ClusterSequences.
494 ///
495 /// A user who wishes to follow the details of the ClusterSequence
496 /// can also make use of this information (and should consult the
497 /// history_element documentation for more information), but should
498 /// be aware that these internal structures may evolve in future
499 /// FastJet versions.
500 const std::vector<history_element> & history() const;
501
502 /// returns the number of particles that were provided to the
503 /// clustering algorithm (helps the user find their way around the
504 /// history and jets objects if they weren't paying attention
505 /// beforehand).
506 unsigned int n_particles() const;
507
508 /// returns a vector of size n_particles() which indicates, for
509 /// each of the initial particles (in the order in which they were
510 /// supplied), which of the supplied jets it belongs to; if it does
511 /// not belong to any of the supplied jets, the index is set to -1;
512 std::vector<int> particle_jet_indices(const std::vector<PseudoJet> &) const;
513
514 /// routine that returns an order in which to read the history
515 /// such that clusterings that lead to identical jet compositions
516 /// but different histories (because of degeneracies in the
517 /// clustering order) will have matching constituents for each
518 /// matching entry in the unique_history_order.
519 ///
520 /// The order has the property that an entry's parents will always
521 /// appear prior to that entry itself.
522 ///
523 /// Roughly speaking the order is such that we first provide all
524 /// steps that lead to the final jet containing particle 1; then we
525 /// have the steps that lead to reconstruction of the jet containing
526 /// the next-lowest-numbered unclustered particle, etc...
527 /// [see GPS CCN28-12 for more info -- of course a full explanation
528 /// here would be better...]
529 std::vector<int> unique_history_order() const;
530
531 /// return the set of particles that have not been clustered. For
532 /// kt and cam/aachen algorithms this should always be null, but for
533 /// cone type algorithms it can be non-null;
534 std::vector<PseudoJet> unclustered_particles() const;
535
536 /// Return the list of pseudojets in the ClusterSequence that do not
537 /// have children (and are not among the inclusive jets). They may
538 /// result from a clustering step or may be one of the pseudojets
539 /// returned by unclustered_particles().
540 std::vector<PseudoJet> childless_pseudojets() const;
541
542 /// returns true if the object (jet or particle) is contained by (ie
543 /// belongs to) this cluster sequence.
544 ///
545 /// Tests performed: if thejet's interface is this cluster sequence
546 /// and its cluster history index is in a consistent range.
547 bool contains(const PseudoJet & object) const;
548
549 /// transfer the sequence contained in other_seq into our own;
550 /// any plugin "extras" contained in the from_seq will be lost
551 /// from there.
552 ///
553 /// It also sets the ClusterSequence pointers of the PseudoJets in
554 /// the history to point to this ClusterSequence
555 ///
556 /// When specified, the second argument is an action that will be
557 /// applied on every jets in the resulting ClusterSequence
558 void transfer_from_sequence(const ClusterSequence & from_seq,
559 const FunctionOfPseudoJet<PseudoJet> * action_on_jets = 0);
560
561 /// retrieve a shared pointer to the wrapper to this ClusterSequence
562 ///
563 /// this may turn useful if you want to track when this
564 /// ClusterSequence goes out of scope
566 return _structure_shared_ptr;
567 }
568
569 /// the structure type associated with a jet belonging to a ClusterSequence
571
572 /// This is the function that is automatically called during
573 /// clustering to print the FastJet banner. Only the first call to
574 /// this function will result in the printout of the banner. Users
575 /// may wish to call this function themselves, during the
576 /// initialization phase of their program, in order to ensure that
577 /// the banner appears before other output. This call will not
578 /// affect 3rd-party banners, e.g. those from plugins.
579 static void print_banner();
580
581 /// \cond internal_doc
582 // [this line must be left as is to hide the doxygen comment]
583 /// A call to this function modifies the stream used to print
584 /// banners (by default cout). If a null pointer is passed, banner
585 /// printout is suppressed. This affects all banners, including
586 /// those from plugins.
587 ///
588 /// Please note that if you distribute 3rd party code
589 /// that links with FastJet, that 3rd party code must not
590 /// use this call turn off the printing of FastJet banners
591 /// by default. This requirement reflects the spirit of
592 /// clause 2c of the GNU Public License (v2), under which
593 /// FastJet and its plugins are distributed.
594 static void set_fastjet_banner_stream(std::ostream * ostr) {_fastjet_banner_ostr = ostr;}
595 // [this line must be left as is to hide the doxygen comment]
596 /// \endcond
597
598 /// returns a pointer to the stream to be used to print banners
599 /// (cout by default). This function is used by plugins to determine
600 /// where to direct their banners. Plugins should properly handle
601 /// the case where the pointer is null.
602 static std::ostream * fastjet_banner_stream() {return _fastjet_banner_ostr;}
603
604private:
605
606 /// \cond internal_doc
607 /// contains the actual stream to use for banners
608#ifdef FASTJET_HAVE_LIMITED_THREAD_SAFETY
609 static std::atomic<std::ostream*> _fastjet_banner_ostr;
610#else
611 static std::ostream * _fastjet_banner_ostr;
612#endif // FASTJET_HAVE_LIMITED_THREAD_SAFETY
613 /// \endcond
614
615protected:
616//DEP static JetAlgorithm _default_jet_algorithm;
617 JetDefinition _jet_def;
618
619 /// transfer the vector<L> of input jets into our own vector<PseudoJet>
620 /// _jets (with some reserved space for future growth).
621 template<class L> void _transfer_input_jets(
622 const std::vector<L> & pseudojets);
623
624 /// This is what is called to do all the initialisation and
625 /// then run the clustering (may be called by various constructors).
626 /// It assumes _jets contains the momenta to be clustered.
627 void _initialise_and_run (const JetDefinition & jet_def,
628 const bool & writeout_combinations);
629
630 //// this performs the initialisation, minus the option-decanting
631 //// stage; for low multiplicity, initialising a few things in the
632 //// constructor, calling the decant_options_partial() and then this
633 //// is faster than going through _initialise_and_run.
634 void _initialise_and_run_no_decant();
635
636//DEP /// This is an alternative routine for initialising and running the
637//DEP /// clustering, provided for legacy purposes. The jet finder is that
638//DEP /// specified in the static member _default_jet_algorithm.
639//DEP void _initialise_and_run (const double R,
640//DEP const Strategy & strategy,
641//DEP const bool & writeout_combinations);
642
643 /// fills in the various member variables with "decanted" options from
644 /// the jet_definition and writeout_combinations variables
645 void _decant_options(const JetDefinition & jet_def,
646 const bool & writeout_combinations);
647
648 /// assuming that the jet definition, writeout_combinations and
649 /// _structure_shared_ptr have been set (e.g. in an initialiser list
650 /// in the constructor), it handles the remaining decanting of
651 /// options.
652 void _decant_options_partial();
653
654 /// fill out the history (and jet cross refs) related to the initial
655 /// set of jets (assumed already to have been "transferred"),
656 /// without any clustering
657 void _fill_initial_history();
658
659 /// carry out the recombination between the jets numbered jet_i and
660 /// jet_j, at distance scale dij; return the index newjet_k of the
661 /// result of the recombination of i and j.
662 void _do_ij_recombination_step(const int jet_i, const int jet_j,
663 const double dij, int & newjet_k);
664
665 /// carry out an recombination step in which _jets[jet_i] merges with
666 /// the beam,
667 void _do_iB_recombination_step(const int jet_i, const double diB);
668
669 /// every time a jet is added internally during clustering, this
670 /// should be called to set the jet's structure shared ptr to point
671 /// to the CS (and the count of internally associated objects is
672 /// also updated). This should not be called outside construction of
673 /// a CS object.
674 void _set_structure_shared_ptr(PseudoJet & j);
675
676 /// make sure that the CS's internal tally of the use count matches
677 /// that of the _structure_shared_ptr
678 void _update_structure_use_count();
679
680 /// returns a suggestion for the best strategy to use on event
681 /// multiplicity, algorithm, R, etc.
682 Strategy _best_strategy() const;
683
684 /// \if internal_doc
685 /// \class _Parabola
686 /// returns c*(a*R**2 + b*R + 1);
687 /// Written as a class in case we want to give names to different
688 /// parabolas
689 /// \endif
690 class _Parabola {
691 public:
692 _Parabola(double a, double b, double c) : _a(a), _b(b), _c(c) {}
693 inline double operator()(const double R) const {return _c*(_a*R*R + _b*R + 1);}
694 private:
695 double _a, _b, _c;
696 };
697
698 /// \if internal_doc
699 /// \class _Line
700 /// operator()(R) returns a*R+b;
701 /// \endif
702 class _Line {
703 public:
704 _Line(double a, double b) : _a(a), _b(b) {}
705 inline double operator()(const double R) const {return _a*R + _b;}
706 private:
707 double _a, _b;
708 };
709
710 /// This contains the physical PseudoJets; for each PseudoJet one
711 /// can find the corresponding position in the _history by looking
712 /// at _jets[i].cluster_hist_index().
713 std::vector<PseudoJet> _jets;
714
715
716 /// this vector will contain the branching history; for each stage,
717 /// _history[i].jetp_index indicates where to look in the _jets
718 /// vector to get the physical PseudoJet.
719 std::vector<history_element> _history;
720
721 /// set subhist to be a set pointers to history entries corresponding to the
722 /// subjets of this jet; one stops going working down through the
723 /// subjets either when
724 /// - there is no further to go
725 /// - one has found maxjet entries
726 /// - max_dij_so_far <= dcut
727 /// By setting maxjet=0 one can use just dcut; by setting dcut<0
728 /// one can use jet maxjet
729 void get_subhist_set(std::set<const history_element*> & subhist,
730 const PseudoJet & jet, double dcut, int maxjet) const;
731
732 bool _writeout_combinations;
733 int _initial_n;
734 double _Rparam, _R2, _invR2;
735 double _Qtot;
736 Strategy _strategy;
737 JetAlgorithm _jet_algorithm;
738
739 SharedPtr<PseudoJetStructureBase> _structure_shared_ptr; //< will actually be of type ClusterSequenceStructure
740 int _structure_use_count_after_construction; //< info of use when CS handles its own memory
741 /// if true then the CS will delete itself when the last external
742 /// object referring to it disappears. It is mutable so as to ensure
743 /// that signal_imminent_self_deletion() [const] can make relevant
744 /// changes.
745#ifdef FASTJET_HAVE_THREAD_SAFETY
746 mutable std::atomic<bool> _deletes_self_when_unused;
747#else
749#endif
750
751 private:
752
753 bool _plugin_activated;
754 SharedPtr<Extras> _extras; // things the plugin might want to add
755
756 void _really_dumb_cluster ();
757 void _delaunay_cluster ();
758 //void _simple_N2_cluster ();
759 template<class BJ> void _simple_N2_cluster ();
760 void _tiled_N2_cluster ();
761 void _faster_tiled_N2_cluster ();
762
763 //
764 void _minheap_faster_tiled_N2_cluster();
765
766 // things needed specifically for Cambridge with Chan's 2D closest
767 // pairs method
768 void _CP2DChan_cluster();
769 void _CP2DChan_cluster_2pi2R ();
770 void _CP2DChan_cluster_2piMultD ();
771 void _CP2DChan_limited_cluster(double D);
772 void _do_Cambridge_inclusive_jets();
773
774 // NSqrtN method for C/A
775 void _fast_NsqrtN_cluster();
776
777 void _add_step_to_history( //const int step_number,
778 const int parent1,
779 const int parent2, const int jetp_index,
780 const double dij);
781
782 /// internal routine associated with the construction of the unique
783 /// history order (following children in the tree)
784 void _extract_tree_children(int pos, std::valarray<bool> &,
785 const std::valarray<int> &, std::vector<int> &) const;
786
787 /// internal routine associated with the construction of the unique
788 /// history order (following parents in the tree)
789 void _extract_tree_parents (int pos, std::valarray<bool> &,
790 const std::valarray<int> &, std::vector<int> &) const;
791
792
793 // these will be useful shorthands in the Voronoi-based code
794 typedef std::pair<int,int> TwoVertices;
795 typedef std::pair<double,TwoVertices> DijEntry;
796 typedef std::multimap<double,TwoVertices> DistMap;
797
798 /// currently used only in the Voronoi based code
799 void _add_ktdistance_to_map(const int ii,
800 DistMap & DijMap,
801 const DynamicNearestNeighbours * DNN);
802
803
804 /// will be set by default to be true for the first run
805 static thread_safety_helpers::FirstTimeTrue _first_time;
806
807 /// manage warnings related to exclusive jets access
808 static LimitedWarning _exclusive_warnings;
809
810 /// the limited warning member for notification of user that
811 /// their requested strategy has been overridden (usually because
812 /// they have R>2pi and not all strategies work then)
813 static LimitedWarning _changed_strategy_warning;
814
815 //----------------------------------------------------------------------
816 /// the fundamental structure which contains the minimal info about
817 /// a jet, as needed for our plain N^2 algorithm -- the idea is to
818 /// put all info that will be accessed N^2 times into an array of
819 /// BriefJets...
820 struct BriefJet {
821 double eta, phi, kt2, NN_dist;
822 BriefJet * NN;
823 int _jets_index;
824 };
825
826 /// structure analogous to BriefJet, but with the extra information
827 /// needed for dealing with tiles
828 class TiledJet {
829 public:
830 double eta, phi, kt2, NN_dist;
831 TiledJet * NN, *previous, * next;
832 int _jets_index, tile_index, diJ_posn;
833 // routines that are useful in the minheap version of tiled
834 // clustering ("misuse" the otherwise unused diJ_posn, so as
835 // to indicate whether jets need to have their minheap entries
836 // updated).
837 inline void label_minheap_update_needed() {diJ_posn = 1;}
838 inline void label_minheap_update_done() {diJ_posn = 0;}
839 inline bool minheap_update_needed() const {return diJ_posn==1;}
840 };
841
842 //-- some of the functions that follow are templates and will work
843 //as well for briefjet and tiled jets
844
845 /// set the kinematic and labelling info for jeta so that it corresponds
846 /// to _jets[_jets_index]
847 template <class J> void _bj_set_jetinfo( J * const jet,
848 const int _jets_index) const;
849
850 /// "remove" this jet, which implies updating links of neighbours and
851 /// perhaps modifying the tile structure
852 void _bj_remove_from_tiles( TiledJet * const jet) const;
853
854 /// return the distance between two BriefJet objects
855 template <class J> double _bj_dist(const J * const jeta,
856 const J * const jetb) const;
857
858 // return the diJ (multiplied by _R2) for this jet assuming its NN
859 // info is correct
860 template <class J> double _bj_diJ(const J * const jeta) const;
861
862 /// for testing purposes only: if in the range head--tail-1 there is a
863 /// a jet which corresponds to hist_index in the history, then
864 /// return a pointer to that jet; otherwise return tail.
865 template <class J> inline J * _bj_of_hindex(
866 const int hist_index,
867 J * const head, J * const tail)
868 const {
869 J * res;
870 for(res = head; res<tail; res++) {
871 if (_jets[res->_jets_index].cluster_hist_index() == hist_index) {break;}
872 }
873 return res;
874 }
875
876
877 //-- remaining functions are different in various cases, so we
878 // will use templates but are not sure if they're useful...
879
880 /// updates (only towards smaller distances) the NN for jeta without checking
881 /// whether in the process jeta itself might be a new NN of one of
882 /// the jets being scanned -- span the range head to tail-1 with
883 /// assumption that jeta is not contained in that range
884 template <class J> void _bj_set_NN_nocross(J * const jeta,
885 J * const head, const J * const tail) const;
886
887 /// reset the NN for jeta and DO check whether in the process jeta
888 /// itself might be a new NN of one of the jets being scanned --
889 /// span the range head to tail-1 with assumption that jeta is not
890 /// contained in that range
891 template <class J> void _bj_set_NN_crosscheck(J * const jeta,
892 J * const head, const J * const tail) const;
893
894
895
896 /// number of neighbours that a tile will have (rectangular geometry
897 /// gives 9 neighbours).
898 static const int n_tile_neighbours = 9;
899 //----------------------------------------------------------------------
900 /// The fundamental structures to be used for the tiled N^2 algorithm
901 /// (see CCN27-44 for some discussion of pattern of tiling)
902 struct Tile {
903 /// pointers to neighbouring tiles, including self
904 Tile * begin_tiles[n_tile_neighbours];
905 /// neighbouring tiles, excluding self
906 Tile ** surrounding_tiles;
907 /// half of neighbouring tiles, no self
908 Tile ** RH_tiles;
909 /// just beyond end of tiles
910 Tile ** end_tiles;
911 /// start of list of BriefJets contained in this tile
912 TiledJet * head;
913 /// sometimes useful to be able to tag a tile
914 bool tagged;
915 };
916 std::vector<Tile> _tiles;
917 double _tiles_eta_min, _tiles_eta_max;
918 double _tile_size_eta, _tile_size_phi;
919 int _n_tiles_phi,_tiles_ieta_min,_tiles_ieta_max;
920
921 // reasonably robust return of tile index given ieta and iphi, in particular
922 // it works even if iphi is negative
923 inline int _tile_index (int ieta, int iphi) const {
924 // note that (-1)%n = -1 so that we have to add _n_tiles_phi
925 // before performing modulo operation
926 return (ieta-_tiles_ieta_min)*_n_tiles_phi
927 + (iphi+_n_tiles_phi) % _n_tiles_phi;
928 }
929
930 // routines for tiled case, including some overloads of the plain
931 // BriefJet cases
932 int _tile_index(const double eta, const double phi) const;
933 void _tj_set_jetinfo ( TiledJet * const jet, const int _jets_index);
934 void _bj_remove_from_tiles(TiledJet * const jet);
935 void _initialise_tiles();
936 void _print_tiles(TiledJet * briefjets ) const;
937 void _add_neighbours_to_tile_union(const int tile_index,
938 std::vector<int> & tile_union, int & n_near_tiles) const;
939 void _add_untagged_neighbours_to_tile_union(const int tile_index,
940 std::vector<int> & tile_union, int & n_near_tiles);
941
942 //----------------------------------------------------------------------
943 /// fundamental structure for e+e- clustering
944 struct EEBriefJet {
945 double NN_dist; // obligatorily present
946 double kt2; // obligatorily present == E^2 in general
947 EEBriefJet * NN; // must be present too
948 int _jets_index; // must also be present!
949 //...........................................................
950 double nx, ny, nz; // our internal storage for fast distance calcs
951 };
952
953 /// identical to EEBriefJet, but the corresponding distance
954 /// calculation (_bj_dist) is overloaded to use the more accurate
955 /// cross-product approach
956 class EEAccurateBriefJet : public EEBriefJet { };
957
958 /// to help instantiation (fj 2.4.0; did not quite work on gcc 33 and os x 10.3?)
959 //void _dummy_N2_cluster_instantiation();
960
961
962 /// to avoid issues with template instantiation (OS X 10.3, gcc 3.3)
963 void _simple_N2_cluster_BriefJet();
964 /// to avoid issues with template instantiation (OS X 10.3, gcc 3.3)
965 void _simple_N2_cluster_EEBriefJet();
966 void _simple_N2_cluster_EEAccurateBriefJet();
967};
968
969
970//**********************************************************************
971//************** START OF INLINE MATERIAL ******************
972//**********************************************************************
973
974
975//----------------------------------------------------------------------
976// Transfer the initial jets into our internal structure
977template<class L> void ClusterSequence::_transfer_input_jets(
978 const std::vector<L> & pseudojets) {
979
980 // this will ensure that we can point to jets without difficulties
981 // arising.
982 _jets.reserve(pseudojets.size()*2);
983
984 // insert initial jets this way so that any type L that can be
985 // converted to a pseudojet will work fine (basically PseudoJet
986 // and any type that has [] subscript access to the momentum
987 // components, such as CLHEP HepLorentzVector).
988 for (unsigned int i = 0; i < pseudojets.size(); i++) {
989 _jets.push_back(pseudojets[i]);}
990
991}
992
993// //----------------------------------------------------------------------
994// // initialise from some generic type... Has to be made available
995// // here in order for it the template aspect of it to work...
996// template<class L> ClusterSequence::ClusterSequence (
997// const std::vector<L> & pseudojets,
998// const double R,
999// const Strategy & strategy,
1000// const bool & writeout_combinations) {
1001//
1002// // transfer the initial jets (type L) into our own array
1003// _transfer_input_jets(pseudojets);
1004//
1005// // run the clustering
1006// _initialise_and_run(R,strategy,writeout_combinations);
1007// }
1008
1009
1010//----------------------------------------------------------------------
1011/// constructor of a jet-clustering sequence from a vector of
1012/// four-momenta, with the jet definition specified by jet_def
1013template<class L> ClusterSequence::ClusterSequence (
1014 const std::vector<L> & pseudojets,
1015 const JetDefinition & jet_def_in,
1016 const bool & writeout_combinations) :
1017 _jet_def(jet_def_in), _writeout_combinations(writeout_combinations),
1018 _structure_shared_ptr(new ClusterSequenceStructure(this))
1019{
1020
1021 // transfer the initial jets (type L) into our own array
1022 _transfer_input_jets(pseudojets);
1023
1024 // transfer the remaining options
1026
1027 // run the clustering
1028 _initialise_and_run_no_decant();
1029}
1030
1031
1032inline const std::vector<PseudoJet> & ClusterSequence::jets () const {
1033 return _jets;
1034}
1035
1036inline const std::vector<ClusterSequence::history_element> & ClusterSequence::history () const {
1037 return _history;
1038}
1039
1040inline unsigned int ClusterSequence::n_particles() const {return _initial_n;}
1041
1042//----------------------------------------------------------------------
1043// implementation of JetDefinition::operator() is here to avoid nasty
1044// issues of order of implementations and includes
1045#ifndef __CINT__
1046template<class L>
1047std::vector<PseudoJet> JetDefinition::operator()(const std::vector<L> & particles) const {
1048 // create a new cluster sequence
1049 ClusterSequence * cs = new ClusterSequence(particles, *this);
1050
1051 // get the jets, and sort them according to whether the algorithm
1052 // is spherical or not
1053 std::vector<PseudoJet> jets;
1054 if (is_spherical()) {
1055 jets = sorted_by_E(cs->inclusive_jets());
1056 } else {
1057 jets = sorted_by_pt(cs->inclusive_jets());
1058 }
1059
1060 // make sure the ClusterSequence gets deleted once it's no longer
1061 // needed
1062 if (jets.size() != 0) {
1064 } else {
1065 delete cs;
1066 }
1067
1068 return jets;
1069}
1070#endif // __CINT__
1071
1072
1073//----------------------------------------------------------------------
1074template <class J> inline void ClusterSequence::_bj_set_jetinfo(
1075 J * const jetA, const int _jets_index) const {
1076 jetA->eta = _jets[_jets_index].rap();
1077 jetA->phi = _jets[_jets_index].phi_02pi();
1078 jetA->kt2 = jet_scale_for_algorithm(_jets[_jets_index]);
1079 jetA->_jets_index = _jets_index;
1080 // initialise NN info as well
1081 jetA->NN_dist = _R2;
1082 jetA->NN = NULL;
1083}
1084
1085
1086
1087
1088//----------------------------------------------------------------------
1089template <class J> inline double ClusterSequence::_bj_dist(
1090 const J * const jetA, const J * const jetB) const {
1091 //#define FASTJET_NEW_DELTA_PHI
1092#ifndef FASTJET_NEW_DELTA_PHI
1093 //GPS+MC old version of Delta phi calculation
1094 double dphi = std::abs(jetA->phi - jetB->phi);
1095 double deta = (jetA->eta - jetB->eta);
1096 if (dphi > pi) {dphi = twopi - dphi;}
1097#else
1098 //GPS+MC testing for 2015-02-faster-deltaR2
1099 double dphi = pi-std::abs(pi-std::abs(jetA->phi - jetB->phi));
1100 double deta = (jetA->eta - jetB->eta);
1101#endif
1102 return dphi*dphi + deta*deta;
1103}
1104
1105//----------------------------------------------------------------------
1106template <class J> inline double ClusterSequence::_bj_diJ(const J * const jet) const {
1107 double kt2 = jet->kt2;
1108 if (jet->NN != NULL) {if (jet->NN->kt2 < kt2) {kt2 = jet->NN->kt2;}}
1109 return jet->NN_dist * kt2;
1110}
1111
1112
1113//----------------------------------------------------------------------
1114// set the NN for jet without checking whether in the process you might
1115// have discovered a new nearest neighbour for another jet
1116template <class J> inline void ClusterSequence::_bj_set_NN_nocross(
1117 J * const jet, J * const head, const J * const tail) const {
1118 double NN_dist = _R2;
1119 J * NN = NULL;
1120 if (head < jet) {
1121 for (J * jetB = head; jetB != jet; jetB++) {
1122 double dist = _bj_dist(jet,jetB);
1123 if (dist < NN_dist) {
1124 NN_dist = dist;
1125 NN = jetB;
1126 }
1127 }
1128 }
1129 if (tail > jet) {
1130 for (J * jetB = jet+1; jetB != tail; jetB++) {
1131 double dist = _bj_dist(jet,jetB);
1132 if (dist < NN_dist) {
1133 NN_dist = dist;
1134 NN = jetB;
1135 }
1136 }
1137 }
1138 jet->NN = NN;
1139 jet->NN_dist = NN_dist;
1140}
1141
1142
1143//----------------------------------------------------------------------
1144template <class J> inline void ClusterSequence::_bj_set_NN_crosscheck(J * const jet,
1145 J * const head, const J * const tail) const {
1146 double NN_dist = _R2;
1147 J * NN = NULL;
1148 for (J * jetB = head; jetB != tail; jetB++) {
1149 double dist = _bj_dist(jet,jetB);
1150 if (dist < NN_dist) {
1151 NN_dist = dist;
1152 NN = jetB;
1153 }
1154 if (dist < jetB->NN_dist) {
1155 jetB->NN_dist = dist;
1156 jetB->NN = jet;
1157 }
1158 }
1159 jet->NN = NN;
1160 jet->NN_dist = NN_dist;
1161}
1162
1163FASTJET_END_NAMESPACE
1164
1165#endif // __FASTJET_CLUSTERSEQUENCE_HH__
Contains any information related to the clustering that should be directly accessible to PseudoJet.
base class to store extra information that plugins may provide
deals with clustering
bool _deletes_self_when_unused
if true then the CS will delete itself when the last external object referring to it disappears.
const SharedPtr< PseudoJetStructureBase > & structure_shared_ptr() const
retrieve a shared pointer to the wrapper to this ClusterSequence
bool plugin_activated() const
the plugin can associate some extra information with the ClusterSequence object by calling this funct...
std::vector< history_element > _history
this vector will contain the branching history; for each stage, _history[i].jetp_index indicates wher...
ClusterSequenceStructure StructureType
the structure type associated with a jet belonging to a ClusterSequence
void plugin_simple_N2_cluster()
allows a plugin to run a templated clustering (nearest-neighbour heuristic)
const JetDefinition & jet_def() const
return a reference to the jet definition
int n_exclusive_jets_ycut(double ycut) const
the number of exclusive jets at the given ycut
const std::vector< history_element > & history() const
allow the user to access the raw internal history.
void plugin_record_iB_recombination(int jet_i, double diB)
record the fact that there has been a recombination between jets()[jet_i] and the beam,...
void delete_self_when_unused()
by calling this routine you tell the ClusterSequence to delete itself when all the Pseudojets associa...
unsigned int n_particles() const
returns the number of particles that were provided to the clustering algorithm (helps the user find t...
std::vector< PseudoJet > inclusive_jets(const double ptmin=0.0) const
return a vector of all jets (in the sense of the inclusive algorithm) with pt >= ptmin.
std::string strategy_string() const
return the name of the strategy used to cluster the event
void _decant_options_partial()
assuming that the jet definition, writeout_combinations and _structure_shared_ptr have been set (e....
std::vector< PseudoJet > _jets
This contains the physical PseudoJets; for each PseudoJet one can find the corresponding position in ...
bool will_delete_self_when_unused() const
return true if the object has been told to delete itself when unused
void print_jets_for_root(const std::vector< PseudoJet > &jets, std::ostream &ostr=std::cout) const
output the supplied vector of jets in a format that can be read by an appropriate root script; the fo...
const std::vector< PseudoJet > & jets() const
allow the user to access the internally stored _jets() array, which contains both the initial particl...
ClusterSequence(const ClusterSequence &cs)
copy constructor for a ClusterSequence
double exclusive_ymerge(int njets) const
return the ymin corresponding to the recombination that went from n+1 to n jets (sometimes known as y...
double exclusive_ymerge_max(int njets) const
same as exclusive_dmerge_max, but normalised to squared total energy
void plugin_associate_extras(Extras *extras_in)
the plugin can associate some extra information with the ClusterSequence object by calling this funct...
double Q() const
returns the sum of all energies in the event (relevant mainly for e+e-)
ClusterSequence()
default constructor
std::vector< PseudoJet > exclusive_jets_ycut(double ycut) const
the exclusive jets obtained at the given ycut
double Q2() const
return Q()^2
PseudoJet & plugin_non_const_jet(unsigned i)
return a non-const reference to the jets()[i], to allow plugins to modify its contents.
static std::ostream * fastjet_banner_stream()
returns a pointer to the stream to be used to print banners (cout by default).
void _transfer_input_jets(const std::vector< L > &pseudojets)
transfer the vector<L> of input jets into our own vector<PseudoJet> _jets (with some reserved space f...
double jet_scale_for_algorithm(const PseudoJet &jet) const
returns the scale associated with a jet as required for this clustering algorithm (kt^2 for the kt-al...
const Extras * extras() const
returns a pointer to the extras object (may be null)
void plugin_record_ij_recombination(int jet_i, int jet_j, double dij, int &newjet_k)
record the fact that there has been a recombination between jets()[jet_i] and jets()[jet_k],...
Strategy strategy_used() const
return the enum value of the strategy used to cluster the event
class that is intended to hold a full definition of the jet clusterer
bool is_spherical() const
returns true if the jet definition involves an algorithm intended for use on a spherical geometry (e....
std::vector< PseudoJet > operator()(const std::vector< L > &particles) const
cluster the supplied particles and returns a vector of resulting jets, sorted by pt (or energy in the...
class to provide facilities for giving warnings up to some maximum number of times and to provide glo...
Class to contain pseudojets, including minimal information of use to jet-clustering routines.
Definition: PseudoJet.hh:68
An implementation of shared pointers that is broadly similar to C++11 shared_ptr (https://en....
Definition: SharedPtr.hh:341
structure analogous to BriefJet, but with the extra information needed for dealing with tiles
provides an object wich will return "true" the first time () is called and false afterwards
Strategy
the various options for the algorithmic strategy to adopt in clustering events with kt and cambridge ...
vector< PseudoJet > sorted_by_E(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing energy
Definition: PseudoJet.cc:887
JetAlgorithm
the various families of jet-clustering algorithm
vector< PseudoJet > sorted_by_pt(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing kt2
Definition: PseudoJet.cc:871
a single element in the clustering history
int parent1
index in _history where first parent of this jet was created (InexistentParent if this jet is an orig...
int jetp_index
index in the _jets vector where we will find the PseudoJet object corresponding to this jet (i....
int parent2
index in _history where second parent of this jet was created (InexistentParent if this jet is an ori...
int child
index in _history where the current jet is
double max_dij_so_far
the largest recombination distance seen so far in the clustering history.
double dij
the distance corresponding to the recombination at this stage of the clustering.