FastJet 3.5.0
Loading...
Searching...
No Matches
ClusterSequence.cc
1//FJSTARTHEADER
2// $Id$
3//
4// Copyright (c) 2005-2025, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
5//
6//----------------------------------------------------------------------
7// This file is part of FastJet.
8//
9// FastJet is free software; you can redistribute it and/or modify
10// it under the terms of the GNU General Public License as published by
11// the Free Software Foundation; either version 2 of the License, or
12// (at your option) any later version.
13//
14// The algorithms that underlie FastJet have required considerable
15// development. They are described in the original FastJet paper,
16// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
17// FastJet as part of work towards a scientific publication, please
18// quote the version you use and include a citation to the manual and
19// optionally also to hep-ph/0512210.
20//
21// FastJet is distributed in the hope that it will be useful,
22// but WITHOUT ANY WARRANTY; without even the implied warranty of
23// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24// GNU General Public License for more details.
25//
26// You should have received a copy of the GNU General Public License
27// along with FastJet. If not, see <http://www.gnu.org/licenses/>.
28//----------------------------------------------------------------------
29//FJENDHEADER
30
31#include "fastjet/config.h"
32#include "fastjet/Error.hh"
33#include "fastjet/PseudoJet.hh"
34#include "fastjet/ClusterSequence.hh"
35#include "fastjet/ClusterSequenceStructure.hh"
36#include "fastjet/version.hh" // stores the current version number
37#include "fastjet/internal/LazyTiling9Alt.hh"
38#include "fastjet/internal/LazyTiling9.hh"
39#include "fastjet/internal/LazyTiling25.hh"
40#ifndef __FJCORE__
41#include "fastjet/internal/LazyTiling9SeparateGhosts.hh"
42#endif // __FJCORE__
43#include<iostream>
44#include<sstream>
45#include<fstream>
46#include<cmath>
47#include<cstdlib>
48#include<cassert>
49#include<string>
50#include<set>
51
52FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
53
54//----------------------------------------------------------------------
55// here's where we put the main page for fastjet (as explained in the
56// Doxygen FAQ)
57// We put it inside the fastjet namespace to have the links without
58// having to specify (fastjet::)
59//......................................................................
60/** \mainpage FastJet code documentation
61 *
62 * These pages provide automatically generated documentation for the
63 * FastJet package.
64 *
65 * \section useful_classes The most useful classes
66 *
67 * Many of the facilities of FastJet can be accessed through the three
68 * following classes:
69 *
70 * - PseudoJet: the basic class for holding the 4-momentum of a
71 * particle or a jet.
72 *
73 * - JetDefinition: the combination of a #JetAlgorithm and its
74 * associated parameters. Can also be initialised with a \ref plugins "plugin".
75 *
76 * - ClusterSequence: constructed with a vector of input (PseudoJet)
77 * particles and a JetDefinition, it computes and stores the
78 * information on how the input particles are clustered into jets.
79 *
80 * \section advanced_classes Selected more advanced classes
81 *
82 * - ClusterSequenceArea: with the help of an AreaDefinition, provides
83 * jets that also contain information about their area.
84 *
85 * \section Tools Selected additional tools
86 *
87 * - JetMedianBackgroundEstimator: with the help of a Selector, a JetDefinition and
88 * an AreaDefinition, allows one to estimate the background noise density in an event; for a simpler, quicker, effective alternative, use GridMedianBackgroundEstimator
89 *
90 * - Transformer: class from which are derived various tools for
91 * manipulating jets and accessing their substructure. Examples are
92 * Subtractor, Filter, Pruner and various taggers (e.g. JHTopTagger
93 * and MassDropTagger).
94 *
95 * \section further_info Further information
96 *
97 * - Selected classes ordered by topics can be found under the <a
98 * href="modules.html">modules</a> tab.
99 *
100 * - The complete list of classes is available under the <a
101 * href="annotated.html">classes</a> tab.
102 *
103 * - For non-class material (<a href="namespacefastjet.html#enum-members">enums</a>,
104 * <a href="namespacefastjet.html#typedef-members">typedefs</a>,
105 * <a href="namespacefastjet.html#func-members">functions</a>), see the
106 * #fastjet documentation
107 *
108 * - For further information and normal documentation, see the main <a
109 * href="https://fastjet.fr/">FastJet</a> page.
110 *
111 * \section examples Examples
112 * See our \subpage Examples page
113 */
114
115// define the doxygen groups
116/// \defgroup basic_classes Fundamental FastJet classes
117/// \defgroup area_classes Area-related classes
118/// \defgroup sec_area_classes Secondary area-related classes
119/// \defgroup plugins Plugins for non-native jet definitions
120/// \defgroup selectors Selectors
121/// \defgroup tools FastJet tools
122/// \{ \defgroup tools_generic Generic tools
123/// \defgroup tools_background Background subtraction
124/// \defgroup tools_taggers Taggers
125/// \}
126/// \defgroup extra_info Access to extra information
127/// \defgroup error_handling Error handling
128/// \defgroup advanced_usage Advanced usage
129/// \if internal_doc
130/// \defgroup internal
131/// \endif
132
133//----------------------------------------------------------------------
134
135
136using namespace std;
137
138
139// The following variable can be modified from within user code
140// so as to redirect banners to an ostream other than cout.
141//
142// Please note that if you distribute 3rd party code
143// that links with FastJet, that 3rd party code is NOT
144// allowed to turn off the printing of FastJet banners
145// by default. This requirement reflects the spirit of
146// clause 2c of the GNU Public License (v2), under which
147// FastJet and its plugins are distributed.
148#ifdef FASTJET_HAVE_LIMITED_THREAD_SAFETY
149atomic<ostream *> ClusterSequence::_fastjet_banner_ostr{&cout};
150#else
151ostream * ClusterSequence::_fastjet_banner_ostr = &cout;
152#endif // FASTJET_HAVE_LIMITED_THREAD_SAFETY
153
154
155// destructor that guarantees proper bookkeeping for the CS Structure
156ClusterSequence::~ClusterSequence () {
157 // set the pointer in the wrapper to this object to NULL to say that
158 // we're going out of scope
159 if (_structure_shared_ptr){
160 ClusterSequenceStructure* csi = dynamic_cast<ClusterSequenceStructure*>(_structure_shared_ptr.get());
161 // normally the csi is purely internal so it really should not be
162 // NULL i.e assert should be OK
163 // (we assert rather than throw an error, since failure here is a
164 // sign of major internal problems)
165 assert(csi != NULL);
166 csi->set_associated_cs(NULL);
167
168//std::shared_ptr: // set_count is not available with CXX11 shared pointers. There
169//std::shared_ptr: // we'll use a different mechanism for handling self-deleting
170//std::shared_ptr: // ClusterSequences
171//std::shared_ptr: #ifndef FASTJET_HAVE_THREAD_SAFETY
172 // if the user had given the CS responsibility to delete itself,
173 // but then deletes the CS themselves, the following lines of
174 // code will ensure that the structure_shared_ptr will have
175 // a proper object count (so that jets associated with the CS will
176 // throw the correct error if the user tries to access their
177 // constituents).
178 if (_deletes_self_when_unused) {
179 _structure_shared_ptr.set_count(_structure_shared_ptr.use_count()
180 + _structure_use_count_after_construction);
181 }
182//std::shared_ptr: #endif // FASTJET_HAVE_THREAD_SAFETY
183 }
184}
185
186//std::shared_ptr: //-----------
187//std::shared_ptr: #ifdef FASTJET_HAVE_THREAD_SAFETY
188//std::shared_ptr: // signals that a jet will no longer use the current CS
189//std::shared_ptr: void ClusterSequence::release_pseudojet(PseudoJet &jet) const{
190//std::shared_ptr: // this only applies to self-deleting clusteer seqences
191//std::shared_ptr: if (!_deletes_self_when_unused) return;
192//std::shared_ptr:
193//std::shared_ptr: // we "free" the jet from the CS
194//std::shared_ptr: //jet.set_structure_shared_ptr(SharedPtr<PseudoJetStructureBase>());
195//std::shared_ptr: jet.force_reset_structure();
196//std::shared_ptr:
197//std::shared_ptr: // and then we can check if we need to delete the CS
198//std::shared_ptr: if (_structure_shared_ptr.use_count() == _structure_use_count_after_construction){
199//std::shared_ptr: //CXX11_SELF_DELETE_DBG: cout << "will self-delete CS (use_count=" << _structure_use_count_after_construction << ")" << endl;
200//std::shared_ptr: // we need to set delete_self_when unused to false before
201//std::shared_ptr: // triggering the deletion
202//std::shared_ptr: //
203//std::shared_ptr: // This also serves a 2nd purpose: if several threads delete a PJ
204//std::shared_ptr: // At the same time we end up in a situation where
205//std::shared_ptr: // "release_pseudojet" frees both their structure pointers and
206//std::shared_ptr: // both could delete the CS (giving a double-free
207//std::shared_ptr: // corruption). This is prevented by the construct below (where we
208//std::shared_ptr: // have made _deletes_self_when_unused atomic)
209//std::shared_ptr: bool expected = true;
210//std::shared_ptr: if (_deletes_self_when_unused.compare_exchange_strong(expected, false))
211//std::shared_ptr: delete this;
212//std::shared_ptr: }
213//std::shared_ptr: }
214//std::shared_ptr:
215//std::shared_ptr: #else // FASTJET_HAVE_THREAD_SAFETY
216
217void ClusterSequence::signal_imminent_self_deletion() const {
218 // normally if the destructor is called when
219 // _deletes_self_when_unused is true, it assumes that it's been
220 // called by the user (and it therefore resets the shared pointer
221 // count to the true count).
222 //
223 // for self deletion (called from the destructor of the CSstructure,
224 // the shared_ptr to which has just had its pointer -> 0) you do
225 // _not_ want to reset the pointer count (otherwise you will end up
226 // with a double delete on the shared pointer once you start
227 // deleting the internal structure of the CS).
228 //
229 // the following modification ensures that the count reset will not
230 // take place in the destructor
231 assert(_deletes_self_when_unused);
232 _deletes_self_when_unused = false;
233}
234
235//std::shared_ptr: #endif // FASTJET_HAVE_THREAD_SAFETY
236
237//DEP //----------------------------------------------------------------------
238//DEP void ClusterSequence::_initialise_and_run (
239//DEP const double R,
240//DEP const Strategy & strategy,
241//DEP const bool & writeout_combinations) {
242//DEP
243//DEP JetDefinition jet_def(_default_jet_algorithm, R, strategy);
244//DEP _initialise_and_run(jet_def, writeout_combinations);
245//DEP }
246
247
248//----------------------------------------------------------------------
249void ClusterSequence::_initialise_and_run (
250 const JetDefinition & jet_def_in,
251 const bool & writeout_combinations) {
252
253 // transfer all relevant info into internal variables
254 _decant_options(jet_def_in, writeout_combinations);
255
256 // now run
257 _initialise_and_run_no_decant();
258}
259
260//----------------------------------------------------------------------
261void ClusterSequence::_initialise_and_run_no_decant () {
262
263 // set up the history entries for the initial particles (those
264 // currently in _jets)
265 _fill_initial_history();
266
267 // don't run anything if the event is empty
268 if (n_particles() == 0) return;
269
270 // ----- deal with special cases: plugins & e+e- ------
271 if (_jet_algorithm == plugin_algorithm) {
272 // allows plugin_xyz() functions to modify cluster sequence
273 _plugin_activated = true;
274 // let the plugin do its work here
275 _jet_def.plugin()->run_clustering( (*this) );
276 _plugin_activated = false;
277 _update_structure_use_count();
278 return;
279 } else if (_jet_algorithm == ee_kt_algorithm ||
280 _jet_algorithm == ee_genkt_algorithm) {
281 if (_jet_algorithm == ee_kt_algorithm) {
282 // make sure that R is large enough so that "beam" recomb only
283 // occurs when a single particle is left
284 // Normally, this should be automatically set to 4 from JetDefinition
285 assert(_Rparam > 2.0);
286 // this is used to renormalise the dij to get a "standard" form
287 // and our convention in e+e- will be different from that
288 // in long.inv case; NB: _invR2 name should be changed -> _renorm_dij?
289 _invR2 = 1.0;
290 } else {
291 // as of 2009-01-09, choose R to be an angular distance, in
292 // radians. Since the algorithm uses 2(1-cos(theta)) as its
293 // squared angular measure, make sure that the _R2 is defined
294 // in a similar way.
295 if (_Rparam > pi) {
296 // choose a value that ensures that back-to-back particles will
297 // always recombine
298 //_R2 = 4.0000000000001;
299 _R2 = 2 * ( 3.0 + cos(_Rparam) );
300 } else {
301 _R2 = 2 * ( 1.0 - cos(_Rparam) );
302 }
303 _invR2 = 1.0/_R2;
304 }
305 // choose the appropriate strategy
306 if (_strategy == N2PlainEEAccurate) {
307 _simple_N2_cluster_EEAccurateBriefJet();
308 } else {
309 // otherwise we revert to the default strategy
310 _strategy = N2Plain;
311 _simple_N2_cluster_EEBriefJet();
312 }
313 return;
314 } else if (_jet_algorithm == undefined_jet_algorithm) {
315 throw Error("A ClusterSequence cannot be created with an uninitialised JetDefinition");
316 }
317
318
319 // automatically redefine the strategy according to N if that is
320 // what the user requested -- transition points (and especially
321 // their R-dependence) are based on empirical observations for a
322 // R=0.4, 0.7 and 1.0, running on toth (3.4GHz, Pentium IV D [dual
323 // core] with 2MB of cache).
324 //-------------
325 // 2011-11-15: lowered N2Plain -> N2Tiled switchover based on some
326 // new tests on an Intel Core 2 Duo T9400 @ 2.53 GHz
327 // with 6MB cache; tests performed with lines such as
328 // ./fastjet_timing_plugins -kt -nhardest 30 -repeat 50000 -strategy -3 -R 0.5 -nev 1 < ../../data/Pythia-PtMin1000-LHC-1000ev.dat
329 if (_strategy == Best) {
330 _strategy = _best_strategy();
331#ifdef DROP_CGAL
332 // fall back strategy for large N when CGAL is missing
333 if (_strategy == NlnN) _strategy = N2MHTLazy25;
334#endif // DROP_CGAL
335 } else if (_strategy == BestFJ30) {
336 int N = _jets.size();
337 //if (N <= 55*max(0.5,min(1.0,_Rparam))) {// old empirical scaling with R
338 //----------------------
339 // 2011-11-15: new empirical scaling with R; NB: low-R N2Tiled
340 // could be significantly improved at low N by limiting the
341 // minimum size of tiles when R is small
342 if (min(1.0,max(0.1,_Rparam)*3.3)*N <= 30) {
343 _strategy = N2Plain;
344 } else if (N > 6200/pow(_Rparam,2.0) && _jet_def.jet_algorithm() == cambridge_algorithm) {
345 _strategy = NlnNCam;
346#ifndef DROP_CGAL
347 } else if ((N > 16000/pow(_Rparam,1.15) && _jet_def.jet_algorithm() != antikt_algorithm)
348 || N > 35000/pow(_Rparam,1.15)) {
349 _strategy = NlnN;
350#endif // DROP_CGAL
351 } else if (N <= 450) {
352 _strategy = N2Tiled;
353 } else {
354 _strategy = N2MinHeapTiled;
355 }
356 }
357
358 // R >= 2pi is not supported by all clustering strategies owing to
359 // periodicity issues (a particle might cluster with itself). When
360 // R>=2pi, we therefore automatically switch to a strategy that is
361 // known to work.
362 if (_Rparam >= twopi) {
363 if ( _strategy == NlnN
364 || _strategy == NlnN3pi
365 || _strategy == NlnNCam
366 || _strategy == NlnNCam2pi2R
367 || _strategy == NlnNCam4pi) {
368#ifdef DROP_CGAL
369 _strategy = N2MinHeapTiled;
370#else
371 _strategy = NlnN4pi;
372#endif
373 }
374 if (_jet_def.strategy() != Best && _strategy != _jet_def.strategy()) {
375 ostringstream oss;
376 oss << "Cluster strategy " << strategy_string(_jet_def.strategy())
377 << " automatically changed to " << strategy_string()
378 << " because the former is not supported for R = " << _Rparam
379 << " >= 2pi";
380 _changed_strategy_warning.warn(oss.str());
381 }
382 }
383
384
385 // run the code containing the selected strategy
386 //
387 // We order the strategies starting from the ones used by the Best
388 // strategy in the order of increasing N, then the remaining ones
389 // again in the order of increasing N.
390 if (_strategy == N2Plain) {
391 // BriefJet provides standard long.invariant kt alg.
392 this->_simple_N2_cluster_BriefJet();
393 } else if (_strategy == N2Tiled) {
394 this->_faster_tiled_N2_cluster();
395 } else if (_strategy == N2MinHeapTiled) {
396 this->_minheap_faster_tiled_N2_cluster();
397 } else if (_strategy == N2MHTLazy9Alt) {
398 // attempt to use an external tiling routine -- it manipulates
399 // the CS history via the plugin mechanism
400 _plugin_activated = true;
401 LazyTiling9Alt tiling(*this);
402 tiling.run();
403 _plugin_activated = false;
404
405 } else if (_strategy == N2MHTLazy25) {
406 // attempt to use an external tiling routine -- it manipulates
407 // the CS history via the plugin mechanism
408 _plugin_activated = true;
409 LazyTiling25 tiling(*this);
410 tiling.run();
411 _plugin_activated = false;
412
413 } else if (_strategy == N2MHTLazy9) {
414 // attempt to use an external tiling routine -- it manipulates
415 // the CS history via the plugin mechanism
416 _plugin_activated = true;
417 LazyTiling9 tiling(*this);
418 tiling.run();
419 _plugin_activated = false;
420
421 } else if (_strategy == N2MHTLazy9AntiKtSeparateGhosts) {
422#ifndef __FJCORE__
423 // attempt to use an external tiling routine -- it manipulates
424 // the CS history via the plugin mechanism
425 _plugin_activated = true;
426 LazyTiling9SeparateGhosts tiling(*this);
427 tiling.run();
428 _plugin_activated = false;
429#else
430 throw Error("N2MHTLazy9AntiKtSeparateGhosts strategy not supported with FJCORE");
431#endif // __FJCORE__
432
433 } else if (_strategy == NlnN) {
434 this->_delaunay_cluster();
435 } else if (_strategy == NlnNCam) {
436 this->_CP2DChan_cluster_2piMultD();
437 } else if (_strategy == NlnN3pi || _strategy == NlnN4pi ) {
438 this->_delaunay_cluster();
439 } else if (_strategy == N3Dumb ) {
440 this->_really_dumb_cluster();
441 } else if (_strategy == N2PoorTiled) {
442 this->_tiled_N2_cluster();
443 } else if (_strategy == NlnNCam4pi) {
444 this->_CP2DChan_cluster();
445 } else if (_strategy == NlnNCam2pi2R) {
446 this->_CP2DChan_cluster_2pi2R();
447 } else {
448 ostringstream err;
449 err << "Unrecognised value for strategy: "<<_strategy;
450 throw Error(err.str());
451 }
452
453}
454
455
456// these needs to be defined outside the class definition.
457thread_safety_helpers::FirstTimeTrue ClusterSequence::_first_time;
458LimitedWarning ClusterSequence::_exclusive_warnings;
459
460
461//----------------------------------------------------------------------
462// the version string
464 return "FastJet version "+string(fastjet_version);
465}
466
467
468//----------------------------------------------------------------------
469// prints a banner on the first call
470void ClusterSequence::print_banner() {
471
472 if (!_first_time()) return;
473
474 // make sure the user has not set the banner stream to NULL
475 ostream * ostr = _fastjet_banner_ostr;
476 if (!ostr) return;
477
478 (*ostr) << "#--------------------------------------------------------------------------\n";
479 (*ostr) << "# FastJet release " << fastjet_version << endl;
480 (*ostr) << "# M. Cacciari, G.P. Salam and G. Soyez \n";
481 (*ostr) << "# A software package for jet finding and analysis at colliders \n";
482 (*ostr) << "# https://fastjet.fr \n";
483 (*ostr) << "# \n";
484 (*ostr) << "# Please cite EPJC72(2012)1896 [arXiv:1111.6097] if you use this package\n";
485 (*ostr) << "# for scientific work and optionally PLB641(2006)57 [hep-ph/0512210]. \n";
486 (*ostr) << "# \n";
487 (*ostr) << "# FastJet is provided without warranty under the GNU GPL v2 or higher. \n";
488 (*ostr) << "# It uses T. Chan's closest pair algorithm, S. Fortune's Voronoi code";
489#ifndef DROP_CGAL
490 (*ostr) << ",\n# CGAL ";
491#else
492 (*ostr) << "\n# ";
493#endif // DROP_CGAL
494 (*ostr) << "and 3rd party plugin jet algorithms. See COPYING file for details.\n";
495 (*ostr) << "#--------------------------------------------------------------------------\n";
496 // make sure we really have the output done.
497 ostr->flush();
498}
499
500//----------------------------------------------------------------------
501// transfer all relevant info into internal variables
502void ClusterSequence::_decant_options(const JetDefinition & jet_def_in,
503 const bool & writeout_combinations) {
504 // make a local copy of the jet definition (for future use)
505 _jet_def = jet_def_in;
506 _writeout_combinations = writeout_combinations;
507 // initialised the wrapper to the current CS
508 _structure_shared_ptr.reset(new ClusterSequenceStructure(this));
509
510 _decant_options_partial();
511}
512
513//----------------------------------------------------------------------
514// transfer all relevant info into internal variables
515void ClusterSequence::_decant_options_partial() {
516 // let the user know what's going on
517 print_banner();
518
519 _jet_algorithm = _jet_def.jet_algorithm();
520 _Rparam = _jet_def.R(); _R2 = _Rparam*_Rparam; _invR2 = 1.0/_R2;
521 _strategy = _jet_def.strategy();
522
523 // disallow interference from the plugin
524 _plugin_activated = false;
525
526 // initialised the wrapper to the current CS
527 //_structure_shared_ptr.reset(new ClusterSequenceStructure(this));
528 _update_structure_use_count(); // make sure it's correct already here
529}
530
531
532//----------------------------------------------------------------------
533// initialise the history in a standard way
534void ClusterSequence::_fill_initial_history () {
535
536 //if (_jets.size() == 0) {throw Error("Cannot run jet-finder on empty event");}
537
538 // reserve sufficient space for everything
539 _jets.reserve(_jets.size()*2);
540 _history.reserve(_jets.size()*2);
541
542 _Qtot = 0;
543
544 for (int i = 0; i < static_cast<int>(_jets.size()) ; i++) {
545 history_element element;
546 element.parent1 = InexistentParent;
547 element.parent2 = InexistentParent;
548 element.child = Invalid;
549 element.jetp_index = i;
550 element.dij = 0.0;
551 element.max_dij_so_far = 0.0;
552
553 _history.push_back(element);
554
555 // do any momentum preprocessing needed by the recombination scheme
556 _jet_def.recombiner()->preprocess(_jets[i]);
557
558 // get cross-referencing right from PseudoJets
559 _jets[i].set_cluster_hist_index(i);
560 _set_structure_shared_ptr(_jets[i]);
561
562 // determine the total energy in the event
563 _Qtot += _jets[i].E();
564 }
565 _initial_n = _jets.size();
566 _deletes_self_when_unused = false;
567}
568
569
570//----------------------------------------------------------------------
571string ClusterSequence::strategy_string (Strategy strategy_in) const {
572 string strategy;
573 switch(strategy_in) {
574 case NlnN:
575 strategy = "NlnN"; break;
576 case NlnN3pi:
577 strategy = "NlnN3pi"; break;
578 case NlnN4pi:
579 strategy = "NlnN4pi"; break;
580 case N2Plain:
581 strategy = "N2Plain"; break;
583 strategy = "N2PlainEEAccurate"; break;
584 case N2Tiled:
585 strategy = "N2Tiled"; break;
586 case N2MinHeapTiled:
587 strategy = "N2MinHeapTiled"; break;
588 case N2PoorTiled:
589 strategy = "N2PoorTiled"; break;
590 case N2MHTLazy9:
591 strategy = "N2MHTLazy9"; break;
592 case N2MHTLazy9Alt:
593 strategy = "N2MHTLazy9Alt"; break;
594 case N2MHTLazy25:
595 strategy = "N2MHTLazy25"; break;
597 strategy = "N2MHTLazy9AntiKtSeparateGhosts"; break;
598 case N3Dumb:
599 strategy = "N3Dumb"; break;
600 case NlnNCam4pi:
601 strategy = "NlnNCam4pi"; break;
602 case NlnNCam2pi2R:
603 strategy = "NlnNCam2pi2R"; break;
604 case NlnNCam:
605 strategy = "NlnNCam"; break; // 2piMultD
606 case plugin_strategy:
607 strategy = "plugin strategy"; break;
608 default:
609 strategy = "Unrecognized";
610 }
611 return strategy;
612}
613
614
615double ClusterSequence::jet_scale_for_algorithm(
616 const PseudoJet & jet) const {
617 if (_jet_algorithm == kt_algorithm) {return jet.kt2();}
618 else if (_jet_algorithm == cambridge_algorithm) {return 1.0;}
619 else if (_jet_algorithm == antikt_algorithm) {
620 double kt2=jet.kt2();
621 return kt2 > 1e-300 ? 1.0/kt2 : 1e300;
622 } else if (_jet_algorithm == genkt_algorithm) {
623 double kt2 = jet.kt2();
624 double p = jet_def().extra_param();
625 if (p <= 0 && kt2 < 1e-300) kt2 = 1e-300; // dodgy safety check
626 return pow(kt2, p);
627 } else if (_jet_algorithm == cambridge_for_passive_algorithm) {
628 double kt2 = jet.kt2();
629 double lim = _jet_def.extra_param();
630 if (kt2 < lim*lim && kt2 != 0.0) {
631 return 1.0/kt2;
632 } else {return 1.0;}
633 } else {throw Error("Unrecognised jet algorithm");}
634}
635
636//----------------------------------------------------------------------
637// returns a suggestion for the best strategy to use on event
638// multiplicity, algorithm, R, etc.
639//
640// Some of the work to establish the best strategy is collected in
641// issue-tracker/2014-07-auto-strategy-selection;
642// transition_fit_v2.fit indicates the results of the fits that we're
643// using here. (Automatically generated by transition_fit_v2.gp).
644//
645// The transition to NlnN is always present, and it is the the
646// caller's responsibility to drop back down to N2MHTLazy25 if NlnN
647// isn't available.
648//
649// This routine should be called only if the jet alg is one of kt,
650// antikt, cam or genkt.
651Strategy ClusterSequence::_best_strategy() const {
652 int N = _jets.size();
653 // define bounded R, always above 0.1, because we don't trust any
654 // of our parametrizations below R = 0.1
655 double bounded_R = max(_Rparam, 0.1);
656
657 // the very first test thing is a quick hard-coded test to decide
658 // if we immediately opt for N2Plain
659 if (N <= 30 || N <= 39.0/(bounded_R + 0.6)) {
660 return N2Plain;
661 }
662
663 // Define objects that describe our various boundaries. A prefix N_
664 // indicates that boundary is for N, while L_ means it's for log(N).
665 //
666 // Hopefully having them static will ensure minimal overhead
667 // in creating them; collecting them in one place should
668 // help with updates?
669 //
670 const static _Parabola N_Tiled_to_MHT_lowR (-45.4947,54.3528,44.6283);
671 const static _Parabola L_MHT_to_MHTLazy9_lowR (0.677807,-1.05006,10.6994);
672 const static _Parabola L_MHTLazy9_to_MHTLazy25_akt_lowR(0.169967,-0.512589,12.1572);
673 const static _Parabola L_MHTLazy9_to_MHTLazy25_kt_lowR (0.16237,-0.484612,12.3373);
674 const static _Parabola L_MHTLazy9_to_MHTLazy25_cam_lowR = L_MHTLazy9_to_MHTLazy25_kt_lowR;
675 const static _Parabola L_MHTLazy25_to_NlnN_akt_lowR (0.0472051,-0.22043,15.9196);
676 const static _Parabola L_MHTLazy25_to_NlnN_kt_lowR (0.118609,-0.326811,14.8287);
677 const static _Parabola L_MHTLazy25_to_NlnN_cam_lowR (0.10119,-0.295748,14.3924);
678
679 const static _Line L_Tiled_to_MHTLazy9_medR (-1.31304,7.29621);
680 const static _Parabola L_MHTLazy9_to_MHTLazy25_akt_medR = L_MHTLazy9_to_MHTLazy25_akt_lowR;
681 const static _Parabola L_MHTLazy9_to_MHTLazy25_kt_medR = L_MHTLazy9_to_MHTLazy25_kt_lowR;
682 const static _Parabola L_MHTLazy9_to_MHTLazy25_cam_medR = L_MHTLazy9_to_MHTLazy25_cam_lowR;
683 const static _Parabola L_MHTLazy25_to_NlnN_akt_medR = L_MHTLazy25_to_NlnN_akt_lowR;
684 const static _Parabola L_MHTLazy25_to_NlnN_kt_medR = L_MHTLazy25_to_NlnN_kt_lowR;
685 const static _Parabola L_MHTLazy25_to_NlnN_cam_medR = L_MHTLazy25_to_NlnN_cam_lowR;
686
687 const static double N_Plain_to_MHTLazy9_largeR = 75;
688 const static double N_MHTLazy9_to_MHTLazy25_akt_largeR = 700;
689 const static double N_MHTLazy9_to_MHTLazy25_kt_largeR = 1000;
690 const static double N_MHTLazy9_to_MHTLazy25_cam_largeR = 1000;
691 const static double N_MHTLazy25_to_NlnN_akt_largeR = 100000;
692 const static double N_MHTLazy25_to_NlnN_kt_largeR = 40000;
693 const static double N_MHTLazy25_to_NlnN_cam_largeR = 15000;
694
695 // We have timing studies only for kt, cam and antikt; for other
696 // algorithms we set the local jet_algorithm variable to the one of
697 // kt,cam,antikt that we think will be closest in behaviour to the
698 // other alg.
699 JetAlgorithm jet_algorithm;
700 if (_jet_algorithm == genkt_algorithm) {
701 // for genkt, then we set the local jet_algorithm variable (used
702 // only for strategy choice) to be either kt or antikt, depending on
703 // the p value.
704 double p = jet_def().extra_param();
705 if (p < 0.0) jet_algorithm = antikt_algorithm;
706 else jet_algorithm = kt_algorithm;
707 } else if (_jet_algorithm == cambridge_for_passive_algorithm) {
708 // we assume (but haven't tested) that using the kt-alg timing
709 // transitions should be adequate for cambridge_for_passive_algorithm
710 jet_algorithm = kt_algorithm;
711 } else {
712 jet_algorithm = _jet_algorithm;
713 }
714
715 if (bounded_R < 0.65) {
716 // low R case
717 if (N < N_Tiled_to_MHT_lowR(bounded_R)) return N2Tiled;
718 double logN = log(double(N));
719 if (logN < L_MHT_to_MHTLazy9_lowR(bounded_R)) return N2MinHeapTiled;
720 else {
721 if (jet_algorithm == antikt_algorithm){
722 if (logN < L_MHTLazy9_to_MHTLazy25_akt_lowR(bounded_R)) return N2MHTLazy9;
723 else if (logN < L_MHTLazy25_to_NlnN_akt_lowR(bounded_R)) return N2MHTLazy25;
724 else return NlnN;
725 } else if (jet_algorithm == kt_algorithm){
726 if (logN < L_MHTLazy9_to_MHTLazy25_kt_lowR(bounded_R)) return N2MHTLazy9;
727 else if (logN < L_MHTLazy25_to_NlnN_kt_lowR(bounded_R)) return N2MHTLazy25;
728 else return NlnN;
729 } else if (jet_algorithm == cambridge_algorithm) {
730 if (logN < L_MHTLazy9_to_MHTLazy25_cam_lowR(bounded_R)) return N2MHTLazy9;
731 else if (logN < L_MHTLazy25_to_NlnN_cam_lowR(bounded_R)) return N2MHTLazy25;
732 else return NlnNCam;
733 }
734 }
735 } else if (bounded_R < 0.5*pi) {
736 // medium R case
737 double logN = log(double(N));
738 if (logN < L_Tiled_to_MHTLazy9_medR(bounded_R)) return N2Tiled;
739 else {
740 if (jet_algorithm == antikt_algorithm){
741 if (logN < L_MHTLazy9_to_MHTLazy25_akt_medR(bounded_R)) return N2MHTLazy9;
742 else if (logN < L_MHTLazy25_to_NlnN_akt_medR(bounded_R)) return N2MHTLazy25;
743 else return NlnN;
744 } else if (jet_algorithm == kt_algorithm){
745 if (logN < L_MHTLazy9_to_MHTLazy25_kt_medR(bounded_R)) return N2MHTLazy9;
746 else if (logN < L_MHTLazy25_to_NlnN_kt_medR(bounded_R)) return N2MHTLazy25;
747 else return NlnN;
748 } else if (jet_algorithm == cambridge_algorithm) {
749 if (logN < L_MHTLazy9_to_MHTLazy25_cam_medR(bounded_R)) return N2MHTLazy9;
750 else if (logN < L_MHTLazy25_to_NlnN_cam_medR(bounded_R)) return N2MHTLazy25;
751 else return NlnNCam;
752 }
753 }
754 } else {
755 // large R case (R > pi/2)
756 if (N < N_Plain_to_MHTLazy9_largeR) return N2Plain;
757 else {
758 if (jet_algorithm == antikt_algorithm){
759 if (N < N_MHTLazy9_to_MHTLazy25_akt_largeR) return N2MHTLazy9;
760 else if (N < N_MHTLazy25_to_NlnN_akt_largeR) return N2MHTLazy25;
761 else return NlnN;
762 } else if (jet_algorithm == kt_algorithm){
763 if (N < N_MHTLazy9_to_MHTLazy25_kt_largeR) return N2MHTLazy9;
764 else if (N < N_MHTLazy25_to_NlnN_kt_largeR) return N2MHTLazy25;
765 else return NlnN;
766 } else if (jet_algorithm == cambridge_algorithm) {
767 if (N < N_MHTLazy9_to_MHTLazy25_cam_largeR) return N2MHTLazy9;
768 else if (N < N_MHTLazy25_to_NlnN_cam_largeR) return N2MHTLazy25;
769 else return NlnNCam;
770 }
771 }
772 }
773
774 //bool code_should_never_reach_here = false;
775 //assert(code_should_never_reach_here);
776
777 assert(0 && "Code should never reach here");
778
779 return N2MHTLazy9;
780
781}
782
783
784// //----------------------------------------------------------------------
785// /// transfer the sequence contained in other_seq into our own;
786// /// any plugin "extras" contained in the from_seq will be lost
787// /// from there.
788// void ClusterSequence::transfer_from_sequence(ClusterSequence & from_seq) {
789//
790// if (will_delete_self_when_unused())
791// throw(Error("cannot use CS::transfer_from_sequence after a call to delete_self_when_unused()"));
792//
793// // the metadata
794// _jet_def = from_seq._jet_def ;
795// _writeout_combinations = from_seq._writeout_combinations ;
796// _initial_n = from_seq._initial_n ;
797// _Rparam = from_seq._Rparam ;
798// _R2 = from_seq._R2 ;
799// _invR2 = from_seq._invR2 ;
800// _strategy = from_seq._strategy ;
801// _jet_algorithm = from_seq._jet_algorithm ;
802// _plugin_activated = from_seq._plugin_activated ;
803//
804// // the data
805// _jets = from_seq._jets;
806// _history = from_seq._history;
807// // the following transfers ownership of the extras from the from_seq
808// _extras = from_seq._extras;
809//
810// // transfer of ownership
811// if (_structure_shared_ptr()) {
812// // anything that is currently associated with the cluster sequence
813// // should be told that its cluster sequence no longer exists
814// ClusterSequenceStructure* csi = dynamic_cast<ClusterSequenceStructure*>(_structure_shared_ptr());
815// assert(csi != NULL);
816// csi->set_associated_cs(NULL);
817// }
818// // create a new _structure_shared_ptr to reflect the fact that
819// // this CS is essentially a new one
820// _structure_shared_ptr.reset(new ClusterSequenceStructure(this));
821// _update_structure_use_count();
822//
823// for (vector<PseudoJet>::iterator jit = _jets.begin(); jit != _jets.end(); jit++)
824// _set_structure_shared_ptr(*jit);
825// }
826
827
828ClusterSequence & ClusterSequence::operator=(const ClusterSequence & cs) {
829 // self assignment is trivial
830 if (&cs != this) {
831 _deletes_self_when_unused = false;
832 transfer_from_sequence(cs);
833 }
834 return *this;
835}
836
837//----------------------------------------------------------------------
838// transfer the sequence contained in other_seq into our own;
839// any plugin "extras" contained in the from_seq will be lost
840// from there.
841//
842// It also sets the ClusterSequence pointers of the PseudoJets in
843// the history to point to this ClusterSequence
844//
845// The second argument is an action that will be applied on every
846// jets in the resulting ClusterSequence
847void ClusterSequence::transfer_from_sequence(const ClusterSequence & from_seq,
848 const FunctionOfPseudoJet<PseudoJet> * action_on_jets){
849
850 if (will_delete_self_when_unused())
851 throw(Error("cannot use CS::transfer_from_sequence after a call to delete_self_when_unused()"));
852
853 // the metadata
854 _jet_def = from_seq._jet_def ;
855 _writeout_combinations = from_seq._writeout_combinations ;
856 _initial_n = from_seq._initial_n ;
857 _Rparam = from_seq._Rparam ;
858 _R2 = from_seq._R2 ;
859 _invR2 = from_seq._invR2 ;
860 _strategy = from_seq._strategy ;
861 _jet_algorithm = from_seq._jet_algorithm ;
862 _plugin_activated = from_seq._plugin_activated ;
863
864 // the data
865
866 // apply the transformation on the jets if needed
867 if (action_on_jets)
868 _jets = (*action_on_jets)(from_seq._jets);
869 else
870 _jets = from_seq._jets;
871 _history = from_seq._history;
872 // the following shares ownership of the extras with the from_seq;
873 // no transformations will be applied to the extras
874 _extras = from_seq._extras;
875
876 // clean up existing structure
877 if (_structure_shared_ptr) {
878 // If there are jets associated with an old version of the CS and
879 // a new one, keeping track of when to delete the CS becomes more
880 // complex; so we don't allow this situation to occur.
881 if (_deletes_self_when_unused) throw Error("transfer_from_sequence cannot be used for a cluster sequence that deletes self when unused");
882
883 // anything that is currently associated with the cluster sequence
884 // should be told that its cluster sequence no longer exists
885 ClusterSequenceStructure* csi = dynamic_cast<ClusterSequenceStructure*>(_structure_shared_ptr.get());
886 assert(csi != NULL);
887 csi->set_associated_cs(NULL);
888 }
889 // create a new _structure_shared_ptr to reflect the fact that
890 // this CS is essentially a new one
891 _structure_shared_ptr.reset(new ClusterSequenceStructure(this));
892 _update_structure_use_count();
893
894 for (unsigned int i=0; i<_jets.size(); i++){
895 // we reset the cluster history index in case action_on_jets
896 // messed up with it
897 _jets[i].set_cluster_hist_index(from_seq._jets[i].cluster_hist_index());
898
899 // reset the structure pointer
900 _set_structure_shared_ptr(_jets[i]);
901 }
902}
903
904
905//----------------------------------------------------------------------
906// record an ij recombination and reset the _jets[newjet_k] momentum and
907// user index to be those of newjet
908void ClusterSequence::plugin_record_ij_recombination(
909 int jet_i, int jet_j, double dij,
910 const PseudoJet & newjet, int & newjet_k) {
911
912 plugin_record_ij_recombination(jet_i, jet_j, dij, newjet_k);
913
914 // now transfer newjet into place
915 int tmp_index = _jets[newjet_k].cluster_hist_index();
916 _jets[newjet_k] = newjet;
917 _jets[newjet_k].set_cluster_hist_index(tmp_index);
918 _set_structure_shared_ptr(_jets[newjet_k]);
919}
920
921
922//----------------------------------------------------------------------
923// return all inclusive jets with pt > ptmin
924vector<PseudoJet> ClusterSequence::inclusive_jets (const double ptmin) const{
925 double dcut = ptmin*ptmin;
926 int i = _history.size() - 1; // last jet
927 vector<PseudoJet> jets_local;
928 if (_jet_algorithm == kt_algorithm) {
929 while (i >= 0) {
930 // with our specific definition of dij and diB (i.e. R appears only in
931 // dij), then dij==diB is the same as the jet.perp2() and we can exploit
932 // this in selecting the jets...
933 if (_history[i].max_dij_so_far < dcut) {break;}
934 if (_history[i].parent2 == BeamJet && _history[i].dij >= dcut) {
935 // for beam jets
936 int parent1 = _history[i].parent1;
937 jets_local.push_back(_jets[_history[parent1].jetp_index]);}
938 i--;
939 }
940 } else if (_jet_algorithm == cambridge_algorithm) {
941 while (i >= 0) {
942 // inclusive jets are all at end of clustering sequence in the
943 // Cambridge algorithm -- so if we find a non-exclusive jet, then
944 // we can exit
945 if (_history[i].parent2 != BeamJet) {break;}
946 int parent1 = _history[i].parent1;
947 const PseudoJet & jet = _jets[_history[parent1].jetp_index];
948 if (jet.perp2() >= dcut) {jets_local.push_back(jet);}
949 i--;
950 }
951 } else if (_jet_algorithm == plugin_algorithm
952 || _jet_algorithm == ee_kt_algorithm
953 || _jet_algorithm == antikt_algorithm
954 || _jet_algorithm == genkt_algorithm
955 || _jet_algorithm == ee_genkt_algorithm
956 || _jet_algorithm == cambridge_for_passive_algorithm) {
957 // for inclusive jets with a plugin algorithm, we make no
958 // assumptions about anything (relation of dij to momenta,
959 // ordering of the dij, etc.)
960 while (i >= 0) {
961 if (_history[i].parent2 == BeamJet) {
962 int parent1 = _history[i].parent1;
963 const PseudoJet & jet = _jets[_history[parent1].jetp_index];
964 if (jet.perp2() >= dcut) {jets_local.push_back(jet);}
965 }
966 i--;
967 }
968 } else {throw Error("cs::inclusive_jets(...): Unrecognized jet algorithm");}
969 return jets_local;
970}
971
972
973//----------------------------------------------------------------------
974// return the number of exclusive jets that would have been obtained
975// running the algorithm in exclusive mode with the given dcut
976int ClusterSequence::n_exclusive_jets (const double dcut) const {
977
978 // first locate the point where clustering would have stopped (i.e. the
979 // first time max_dij_so_far > dcut)
980 int i = _history.size() - 1; // last jet
981 while (i >= 0) {
982 if (_history[i].max_dij_so_far <= dcut) {break;}
983 i--;
984 }
985 int stop_point = i + 1;
986 // relation between stop_point, njets assumes one extra jet disappears
987 // at each clustering.
988 int njets = 2*_initial_n - stop_point;
989 return njets;
990}
991
992//----------------------------------------------------------------------
993// return all exclusive jets that would have been obtained running
994// the algorithm in exclusive mode with the given dcut
995vector<PseudoJet> ClusterSequence::exclusive_jets (const double dcut) const {
996 int njets = n_exclusive_jets(dcut);
997 return exclusive_jets(njets);
998}
999
1000
1001//----------------------------------------------------------------------
1002// return the jets obtained by clustering the event to n jets.
1003// Throw an error if there are fewer than n particles.
1004vector<PseudoJet> ClusterSequence::exclusive_jets (const int njets) const {
1005
1006 // make sure the user does not ask for more than jets than there
1007 // were particles in the first place.
1008 if (njets > _initial_n) {
1009 ostringstream err;
1010 err << "Requested " << njets << " exclusive jets, but there were only "
1011 << _initial_n << " particles in the event";
1012 throw Error(err.str());
1013 }
1014
1015 return exclusive_jets_up_to(njets);
1016}
1017
1018//----------------------------------------------------------------------
1019// return the jets obtained by clustering the event to n jets.
1020// If there are fewer than n particles, simply return all particles
1021vector<PseudoJet> ClusterSequence::exclusive_jets_up_to (const int njets) const {
1022
1023 // provide a warning when extracting exclusive jets for algorithms
1024 // that does not support it explicitly.
1025 // Native algorithm that support it are: kt, ee_kt, Cambridge/Aachen,
1026 // genkt and ee_genkt (both with p>=0)
1027 // For plugins, we check Plugin::exclusive_sequence_meaningful()
1028 if (( _jet_def.jet_algorithm() != kt_algorithm) &&
1029 ( _jet_def.jet_algorithm() != cambridge_algorithm) &&
1030 ( _jet_def.jet_algorithm() != ee_kt_algorithm) &&
1031 (((_jet_def.jet_algorithm() != genkt_algorithm) &&
1032 (_jet_def.jet_algorithm() != ee_genkt_algorithm)) ||
1033 (_jet_def.extra_param() <0)) &&
1034 ((_jet_def.jet_algorithm() != plugin_algorithm) ||
1035 (!_jet_def.plugin()->exclusive_sequence_meaningful()))) {
1036 _exclusive_warnings.warn("dcut and exclusive jets for jet-finders other than kt, C/A or genkt with p>=0 should be interpreted with care.");
1037 }
1038
1039
1040 // calculate the point where we have to stop the clustering.
1041 // relation between stop_point, njets assumes one extra jet disappears
1042 // at each clustering.
1043 int stop_point = 2*_initial_n - njets;
1044 // make sure it's safe when more jets are requested than there are particles
1045 if (stop_point < _initial_n) stop_point = _initial_n;
1046
1047 // some sanity checking to make sure that e+e- does not give us
1048 // surprises (should we ever implement e+e-)...
1049 if (2*_initial_n != static_cast<int>(_history.size())) {
1050 ostringstream err;
1051 err << "2*_initial_n != _history.size() -- this endangers internal assumptions!\n";
1052 throw Error(err.str());
1053 //assert(false);
1054 }
1055
1056 // now go forwards and reconstitute the jets that we have --
1057 // basically for any history element, see if the parent jets to
1058 // which it refers were created before the stopping point -- if they
1059 // were then add them to the list, otherwise they are subsequent
1060 // recombinations of the jets that we are looking for.
1061 vector<PseudoJet> jets_local;
1062 for (unsigned int i = stop_point; i < _history.size(); i++) {
1063 int parent1 = _history[i].parent1;
1064 if (parent1 < stop_point) {
1065 jets_local.push_back(_jets[_history[parent1].jetp_index]);
1066 }
1067 int parent2 = _history[i].parent2;
1068 if (parent2 < stop_point && parent2 > 0) {
1069 jets_local.push_back(_jets[_history[parent2].jetp_index]);
1070 }
1071
1072 }
1073
1074 // sanity check...
1075 if (int(jets_local.size()) != min(_initial_n, njets)) {
1076 ostringstream err;
1077 err << "ClusterSequence::exclusive_jets: size of returned vector ("
1078 <<jets_local.size()<<") does not coincide with requested number of jets ("
1079 <<njets<<")";
1080 throw Error(err.str());
1081 }
1082
1083 return jets_local;
1084}
1085
1086//----------------------------------------------------------------------
1087/// return the dmin corresponding to the recombination that went from
1088/// n+1 to n jets
1089double ClusterSequence::exclusive_dmerge (const int njets) const {
1090 assert(njets >= 0);
1091 if (njets >= _initial_n) {return 0.0;}
1092 return _history[2*_initial_n-njets-1].dij;
1093}
1094
1095
1096//----------------------------------------------------------------------
1097/// return the maximum of the dmin encountered during all recombinations
1098/// up to the one that led to an n-jet final state; identical to
1099/// exclusive_dmerge, except in cases where the dmin do not increase
1100/// monotonically.
1101double ClusterSequence::exclusive_dmerge_max (const int njets) const {
1102 assert(njets >= 0);
1103 if (njets >= _initial_n) {return 0.0;}
1104 return _history[2*_initial_n-njets-1].max_dij_so_far;
1105}
1106
1107
1108//----------------------------------------------------------------------
1109/// return a vector of all subjets of the current jet (in the sense
1110/// of the exclusive algorithm) that would be obtained when running
1111/// the algorithm with the given dcut.
1112std::vector<PseudoJet> ClusterSequence::exclusive_subjets
1113 (const PseudoJet & jet, const double dcut) const {
1114
1115 set<const history_element*> subhist;
1116
1117 // get the set of history elements that correspond to subjets at
1118 // scale dcut
1119 get_subhist_set(subhist, jet, dcut, 0);
1120
1121 // now transfer this into a sequence of jets
1122 vector<PseudoJet> subjets;
1123 subjets.reserve(subhist.size());
1124 for (set<const history_element*>::iterator elem = subhist.begin();
1125 elem != subhist.end(); elem++) {
1126 subjets.push_back(_jets[(*elem)->jetp_index]);
1127 }
1128 return subjets;
1129}
1130
1131//----------------------------------------------------------------------
1132/// return the size of exclusive_subjets(...); still n ln n with same
1133/// coefficient, but marginally more efficient than manually taking
1134/// exclusive_subjets.size()
1135int ClusterSequence::n_exclusive_subjets(const PseudoJet & jet,
1136 const double dcut) const {
1137 set<const history_element*> subhist;
1138 // get the set of history elements that correspond to subjets at
1139 // scale dcut
1140 get_subhist_set(subhist, jet, dcut, 0);
1141 return subhist.size();
1142}
1143
1144//----------------------------------------------------------------------
1145/// return the list of subjets obtained by unclustering the supplied
1146/// jet down to nsub subjets. Throws an error if there are fewer than
1147/// nsub particles in the jet.
1148std::vector<PseudoJet> ClusterSequence::exclusive_subjets
1149 (const PseudoJet & jet, int nsub) const {
1150 vector<PseudoJet> subjets = exclusive_subjets_up_to(jet, nsub);
1151 if (int(subjets.size()) < nsub) {
1152 ostringstream err;
1153 err << "Requested " << nsub << " exclusive subjets, but there were only "
1154 << subjets.size() << " particles in the jet";
1155 throw Error(err.str());
1156 }
1157 return subjets;
1158
1159}
1160
1161//----------------------------------------------------------------------
1162/// return the list of subjets obtained by unclustering the supplied
1163/// jet down to nsub subjets (or all constituents if there are fewer
1164/// than nsub).
1165std::vector<PseudoJet> ClusterSequence::exclusive_subjets_up_to
1166 (const PseudoJet & jet, int nsub) const {
1167
1168 set<const history_element*> subhist;
1169
1170 // prepare the vector into which we'll put the result
1171 vector<PseudoJet> subjets;
1172 if (nsub < 0) throw Error("Requested a negative number of subjets. This is nonsensical.");
1173 if (nsub == 0) return subjets;
1174
1175 // get the set of history elements that correspond to subjets at
1176 // scale dcut
1177 get_subhist_set(subhist, jet, -1.0, nsub);
1178
1179 // now transfer this into a sequence of jets
1180 subjets.reserve(subhist.size());
1181 for (set<const history_element*>::iterator elem = subhist.begin();
1182 elem != subhist.end(); elem++) {
1183 subjets.push_back(_jets[(*elem)->jetp_index]);
1184 }
1185 return subjets;
1186}
1187
1188
1189//----------------------------------------------------------------------
1190/// return the dij that was present in the merging nsub+1 -> nsub
1191/// subjets inside this jet.
1192///
1193/// If the jet has nsub or fewer constituents, it will return 0.
1194double ClusterSequence::exclusive_subdmerge(const PseudoJet & jet, int nsub) const {
1195 set<const history_element*> subhist;
1196
1197 // get the set of history elements that correspond to subjets at
1198 // scale dcut
1199 get_subhist_set(subhist, jet, -1.0, nsub);
1200
1201 set<const history_element*>::iterator highest = subhist.end();
1202 highest--;
1203 /// will be zero if nconst <= nsub, since highest will be an original
1204 /// particle have zero dij
1205 return (*highest)->dij;
1206}
1207
1208
1209//----------------------------------------------------------------------
1210/// return the maximum dij that occurred in the whole event at the
1211/// stage that the nsub+1 -> nsub merge of subjets occurred inside
1212/// this jet.
1213///
1214/// If the jet has nsub or fewer constituents, it will return 0.
1215double ClusterSequence::exclusive_subdmerge_max(const PseudoJet & jet, int nsub) const {
1216
1217 set<const history_element*> subhist;
1218
1219 // get the set of history elements that correspond to subjets at
1220 // scale dcut
1221 get_subhist_set(subhist, jet, -1.0, nsub);
1222
1223 set<const history_element*>::iterator highest = subhist.end();
1224 highest--;
1225 /// will be zero if nconst <= nsub, since highest will be an original
1226 /// particle have zero dij
1227 return (*highest)->max_dij_so_far;
1228}
1229
1230
1231
1232//----------------------------------------------------------------------
1233/// return a set of pointers to history entries corresponding to the
1234/// subjets of this jet; one stops going working down through the
1235/// subjets either when
1236/// - there is no further to go
1237/// - one has found maxjet entries
1238/// - max_dij_so_far <= dcut
1239void ClusterSequence::get_subhist_set(set<const history_element*> & subhist,
1240 const PseudoJet & jet,
1241 double dcut, int maxjet) const {
1242 assert(contains(jet));
1243
1244 subhist.clear();
1245 subhist.insert(&(_history[jet.cluster_hist_index()]));
1246
1247 // establish the set of jets that are relevant
1248 int njet = 1;
1249 while (true) {
1250 // first find out if we need to probe deeper into jet.
1251 // Get history element closest to end of sequence
1252 set<const history_element*>::iterator highest = subhist.end();
1253 assert (highest != subhist.begin());
1254 highest--;
1255 const history_element* elem = *highest;
1256 // make sure we haven't got too many jets
1257 if (njet == maxjet) break;
1258 // make sure it has parents
1259 if (elem->parent1 < 0) break;
1260 // make sure that we still resolve it at scale dcut
1261 if (elem->max_dij_so_far <= dcut) break;
1262
1263 // then do so: replace "highest" with its two parents
1264 subhist.erase(highest);
1265 subhist.insert(&(_history[elem->parent1]));
1266 subhist.insert(&(_history[elem->parent2]));
1267 njet++;
1268 }
1269}
1270
1271//----------------------------------------------------------------------
1272// work through the object's history until
1273bool ClusterSequence::object_in_jet(const PseudoJet & object,
1274 const PseudoJet & jet) const {
1275
1276 // make sure the object conceivably belongs to this clustering
1277 // sequence
1278 assert(contains(object) && contains(jet));
1279
1280 const PseudoJet * this_object = &object;
1281 const PseudoJet * childp;
1282 while(true) {
1283 if (this_object->cluster_hist_index() == jet.cluster_hist_index()) {
1284 return true;
1285 } else if (has_child(*this_object, childp)) {
1286 this_object = childp;
1287 } else {
1288 return false;
1289 }
1290 }
1291}
1292
1293//----------------------------------------------------------------------
1294/// if the jet has parents in the clustering, it returns true
1295/// and sets parent1 and parent2 equal to them.
1296///
1297/// if it has no parents it returns false and sets parent1 and
1298/// parent2 to zero
1299bool ClusterSequence::has_parents(const PseudoJet & jet, PseudoJet & parent1,
1300 PseudoJet & parent2) const {
1301
1302 const history_element & hist = _history[jet.cluster_hist_index()];
1303
1304 // make sure we do not run into any unexpected situations --
1305 // i.e. both parents valid, or neither
1306 assert ((hist.parent1 >= 0 && hist.parent2 >= 0) ||
1307 (hist.parent1 < 0 && hist.parent2 < 0));
1308
1309 if (hist.parent1 < 0) {
1310 parent1 = PseudoJet(0.0,0.0,0.0,0.0);
1311 parent2 = parent1;
1312 return false;
1313 } else {
1314 parent1 = _jets[_history[hist.parent1].jetp_index];
1315 parent2 = _jets[_history[hist.parent2].jetp_index];
1316 // order the parents in decreasing pt
1317 if (parent1.perp2() < parent2.perp2()) std::swap(parent1,parent2);
1318 return true;
1319 }
1320}
1321
1322//----------------------------------------------------------------------
1323/// if the jet has a child then return true and give the child jet
1324/// otherwise return false and set the child to zero
1325bool ClusterSequence::has_child(const PseudoJet & jet, PseudoJet & child) const {
1326
1327 //const history_element & hist = _history[jet.cluster_hist_index()];
1328 //
1329 //if (hist.child >= 0) {
1330 // child = _jets[_history[hist.child].jetp_index];
1331 // return true;
1332 //} else {
1333 // child = PseudoJet(0.0,0.0,0.0,0.0);
1334 // return false;
1335 //}
1336 const PseudoJet * childp;
1337 bool res = has_child(jet, childp);
1338 if (res) {
1339 child = *childp;
1340 return true;
1341 } else {
1342 child = PseudoJet(0.0,0.0,0.0,0.0);
1343 return false;
1344 }
1345}
1346
1347bool ClusterSequence::has_child(const PseudoJet & jet, const PseudoJet * & childp) const {
1348
1349 const history_element & hist = _history[jet.cluster_hist_index()];
1350
1351 // check that this jet has a child and that the child corresponds to
1352 // a true jet [RETHINK-IF-CHANGE-NUMBERING: what is the right
1353 // behaviour if the child is the same jet but made inclusive...?]
1354 if (hist.child >= 0 && _history[hist.child].jetp_index >= 0) {
1355 childp = &(_jets[_history[hist.child].jetp_index]);
1356 return true;
1357 } else {
1358 childp = NULL;
1359 return false;
1360 }
1361}
1362
1363
1364//----------------------------------------------------------------------
1365/// if this jet has a child (and so a partner) return true
1366/// and give the partner, otherwise return false and set the
1367/// partner to zero
1368bool ClusterSequence::has_partner(const PseudoJet & jet,
1369 PseudoJet & partner) const {
1370
1371 const history_element & hist = _history[jet.cluster_hist_index()];
1372
1373 // make sure we have a child and that the child does not correspond
1374 // to a clustering with the beam (or some other invalid quantity)
1375 if (hist.child >= 0 && _history[hist.child].parent2 >= 0) {
1376 const history_element & child_hist = _history[hist.child];
1377 if (child_hist.parent1 == jet.cluster_hist_index()) {
1378 // partner will be child's parent2 -- for iB clustering
1379 // parent2 will not be valid
1380 partner = _jets[_history[child_hist.parent2].jetp_index];
1381 } else {
1382 // partner will be child's parent1
1383 partner = _jets[_history[child_hist.parent1].jetp_index];
1384 }
1385 return true;
1386 } else {
1387 partner = PseudoJet(0.0,0.0,0.0,0.0);
1388 return false;
1389 }
1390}
1391
1392
1393//----------------------------------------------------------------------
1394// return a vector of the particles that make up a jet
1395vector<PseudoJet> ClusterSequence::constituents (const PseudoJet & jet) const {
1396 vector<PseudoJet> subjets;
1397 add_constituents(jet, subjets);
1398 return subjets;
1399}
1400
1401//----------------------------------------------------------------------
1402/// output the supplied vector of jets in a format that can be read
1403/// by an appropriate root script; the format is:
1404/// jet-n jet-px jet-py jet-pz jet-E
1405/// particle-n particle-rap particle-phi particle-pt
1406/// particle-n particle-rap particle-phi particle-pt
1407/// ...
1408/// #END
1409/// ... [i.e. above repeated]
1410void ClusterSequence::print_jets_for_root(const std::vector<PseudoJet> & jets_in,
1411 ostream & ostr) const {
1412 for (unsigned i = 0; i < jets_in.size(); i++) {
1413 ostr << i << " "
1414 << jets_in[i].px() << " "
1415 << jets_in[i].py() << " "
1416 << jets_in[i].pz() << " "
1417 << jets_in[i].E() << endl;
1418 vector<PseudoJet> cst = constituents(jets_in[i]);
1419 for (unsigned j = 0; j < cst.size() ; j++) {
1420 ostr << " " << j << " "
1421 << cst[j].rap() << " "
1422 << cst[j].phi() << " "
1423 << cst[j].perp() << endl;
1424 }
1425 ostr << "#END" << endl;
1426 }
1427}
1428
1429void ClusterSequence::print_jets_for_root(const std::vector<PseudoJet> & jets_in,
1430 const std::string & filename,
1431 const std::string & comment ) const {
1432 std::ofstream ostr(filename.c_str());
1433 if (comment != "") ostr << "# " << comment << endl;
1434 print_jets_for_root(jets_in, ostr);
1435}
1436
1437
1438// Not yet. Perhaps in a future release
1439// //----------------------------------------------------------------------
1440// // print out all inclusive jets with pt > ptmin
1441// void ClusterSequence::print_jets (const double ptmin) const{
1442// vector<PseudoJet> jets = sorted_by_pt(inclusive_jets(ptmin));
1443//
1444// for (size_t j = 0; j < jets.size(); j++) {
1445// printf("%5u %7.3f %7.3f %9.3f\n",
1446// j,jets[j].rap(),jets[j].phi(),jets[j].perp());
1447// }
1448// }
1449
1450//----------------------------------------------------------------------
1451/// returns a vector of size n_particles() which indicates, for
1452/// each of the initial particles (in the order in which they were
1453/// supplied), which of the supplied jets it belongs to; if it does
1454/// not belong to any of the supplied jets, the index is set to -1;
1455vector<int> ClusterSequence::particle_jet_indices(
1456 const vector<PseudoJet> & jets_in) const {
1457
1458 vector<int> indices(n_particles());
1459
1460 // first label all particles as not belonging to any jets
1461 for (unsigned ipart = 0; ipart < n_particles(); ipart++)
1462 indices[ipart] = -1;
1463
1464 // then for each of the jets relabel its consituents as belonging to
1465 // that jet
1466 for (unsigned ijet = 0; ijet < jets_in.size(); ijet++) {
1467
1468 vector<PseudoJet> jet_constituents(constituents(jets_in[ijet]));
1469
1470 for (unsigned ip = 0; ip < jet_constituents.size(); ip++) {
1471 // a safe (if slightly redundant) way of getting the particle
1472 // index (for initial particles it is actually safe to assume
1473 // ipart=iclust).
1474 unsigned iclust = jet_constituents[ip].cluster_hist_index();
1475 unsigned ipart = history()[iclust].jetp_index;
1476 indices[ipart] = ijet;
1477 }
1478 }
1479
1480 return indices;
1481}
1482
1483
1484//----------------------------------------------------------------------
1485// recursive routine that adds on constituents of jet to the subjet_vector
1486void ClusterSequence::add_constituents (
1487 const PseudoJet & jet, vector<PseudoJet> & subjet_vector) const {
1488 // find out position in cluster history
1489 int i = jet.cluster_hist_index();
1490 int parent1 = _history[i].parent1;
1491 int parent2 = _history[i].parent2;
1492
1493 if (parent1 == InexistentParent) {
1494 // It is an original particle (labelled by its parent having value
1495 // InexistentParent), therefore add it on to the subjet vector
1496 // Note: we add the initial particle and not simply 'jet' so that
1497 // calling add_constituents with a subtracted jet containing
1498 // only one particle will work.
1499 subjet_vector.push_back(_jets[i]);
1500 return;
1501 }
1502
1503 // add parent 1
1504 add_constituents(_jets[_history[parent1].jetp_index], subjet_vector);
1505
1506 // see if parent2 is a real jet; if it is then add its constituents
1507 if (parent2 != BeamJet) {
1508 add_constituents(_jets[_history[parent2].jetp_index], subjet_vector);
1509 }
1510}
1511
1512
1513
1514//----------------------------------------------------------------------
1515// initialise the history in a standard way
1516void ClusterSequence::_add_step_to_history (
1517 //NO_LONGER_USED: const int step_number,
1518 const int parent1,
1519 const int parent2, const int jetp_index,
1520 const double dij) {
1521
1522 history_element element;
1523 element.parent1 = parent1;
1524 element.parent2 = parent2;
1525 element.jetp_index = jetp_index;
1526 element.child = Invalid;
1527 element.dij = dij;
1528 element.max_dij_so_far = max(dij,_history[_history.size()-1].max_dij_so_far);
1529 _history.push_back(element);
1530
1531 int local_step = _history.size()-1;
1532 //#ifndef __NO_ASSERTS__
1533 //assert(local_step == step_number);
1534 //#endif
1535
1536 // sanity check: make sure the particles have not already been recombined
1537 //
1538 // Note that good practice would make this an assert (since this is
1539 // a serious internal issue). However, we decided to throw an
1540 // InternalError so that the end user can decide to catch it and
1541 // retry the clustering with a different strategy.
1542
1543 assert(parent1 >= 0);
1544 if (_history[parent1].child != Invalid){
1545 throw InternalError("trying to recombine an object that has previously been recombined");
1546 }
1547 _history[parent1].child = local_step;
1548 if (parent2 >= 0) {
1549 if (_history[parent2].child != Invalid){
1550 throw InternalError("trying to recombine an object that has previously been recombined");
1551 }
1552 _history[parent2].child = local_step;
1553 }
1554
1555 // get cross-referencing right from PseudoJets
1556 if (jetp_index != Invalid) {
1557 assert(jetp_index >= 0);
1558 //cout << _jets.size() <<" "<<jetp_index<<"\n";
1559 _jets[jetp_index].set_cluster_hist_index(local_step);
1560 _set_structure_shared_ptr(_jets[jetp_index]);
1561 }
1562
1563 if (_writeout_combinations) {
1564 cout << local_step << ": "
1565 << parent1 << " with " << parent2
1566 << "; y = "<< dij<<endl;
1567 }
1568
1569}
1570
1571
1572
1573
1574//======================================================================
1575// Return an order in which to read the history such that _history[order[i]]
1576// will always correspond to the same set of consituent particles if
1577// two branching histories are equivalent in terms of the particles
1578// contained in any given pseudojet.
1579vector<int> ClusterSequence::unique_history_order() const {
1580
1581 // first construct an array that will tell us the lowest constituent
1582 // of a given jet -- this will always be one of the original
1583 // particles, whose order is well defined and so will help us to
1584 // follow the tree in a unique manner.
1585 valarray<int> lowest_constituent(_history.size());
1586 int hist_n = _history.size();
1587 lowest_constituent = hist_n; // give it a large number
1588 for (int i = 0; i < hist_n; i++) {
1589 // sets things up for the initial partons
1590 lowest_constituent[i] = min(lowest_constituent[i],i);
1591 // propagates them through to the children of this parton
1592 if (_history[i].child > 0) lowest_constituent[_history[i].child]
1593 = min(lowest_constituent[_history[i].child],lowest_constituent[i]);
1594 }
1595
1596 // establish an array for what we have and have not extracted so far
1597 valarray<bool> extracted(_history.size()); extracted = false;
1598 vector<int> unique_tree;
1599 unique_tree.reserve(_history.size());
1600
1601 // now work our way through the tree
1602 for (unsigned i = 0; i < n_particles(); i++) {
1603 if (!extracted[i]) {
1604 unique_tree.push_back(i);
1605 extracted[i] = true;
1606 _extract_tree_children(i, extracted, lowest_constituent, unique_tree);
1607 }
1608 }
1609
1610 return unique_tree;
1611}
1612
1613//======================================================================
1614// helper for unique_history_order
1615void ClusterSequence::_extract_tree_children(
1616 int position,
1617 valarray<bool> & extracted,
1618 const valarray<int> & lowest_constituent,
1619 vector<int> & unique_tree) const {
1620 if (!extracted[position]) {
1621 // that means we may have unidentified parents around, so go and
1622 // collect them (extracted[position]) will then be made true)
1623 _extract_tree_parents(position,extracted,lowest_constituent,unique_tree);
1624 }
1625
1626 // now look after the children...
1627 int child = _history[position].child;
1628 if (child >= 0) _extract_tree_children(child,extracted,lowest_constituent,unique_tree);
1629}
1630
1631
1632//======================================================================
1633// return the list of unclustered particles
1634vector<PseudoJet> ClusterSequence::unclustered_particles() const {
1635 vector<PseudoJet> unclustered;
1636 for (unsigned i = 0; i < n_particles() ; i++) {
1637 if (_history[i].child == Invalid)
1638 unclustered.push_back(_jets[_history[i].jetp_index]);
1639 }
1640 return unclustered;
1641}
1642
1643//======================================================================
1644/// Return the list of pseudojets in the ClusterSequence that do not
1645/// have children (and are not among the inclusive jets). They may
1646/// result from a clustering step or may be one of the pseudojets
1647/// returned by unclustered_particles().
1648vector<PseudoJet> ClusterSequence::childless_pseudojets() const {
1649 vector<PseudoJet> unclustered;
1650 for (unsigned i = 0; i < _history.size() ; i++) {
1651 if ((_history[i].child == Invalid) && (_history[i].parent2 != BeamJet))
1652 unclustered.push_back(_jets[_history[i].jetp_index]);
1653 }
1654 return unclustered;
1655}
1656
1657
1658
1659//----------------------------------------------------------------------
1660// returns true if the cluster sequence contains this jet (i.e. jet's
1661// structure is this cluster sequence's and the cluster history index
1662// is in a consistent range)
1663bool ClusterSequence::contains(const PseudoJet & jet) const {
1664 return jet.cluster_hist_index() >= 0
1665 && jet.cluster_hist_index() < int(_history.size())
1667 && jet.associated_cluster_sequence() == this;
1668}
1669
1670
1671
1672//======================================================================
1673// helper for unique_history_order
1674void ClusterSequence::_extract_tree_parents(
1675 int position,
1676 valarray<bool> & extracted,
1677 const valarray<int> & lowest_constituent,
1678 vector<int> & unique_tree) const {
1679
1680 if (!extracted[position]) {
1681 int parent1 = _history[position].parent1;
1682 int parent2 = _history[position].parent2;
1683 // where relevant order parents so that we will first treat the
1684 // one containing the smaller "lowest_constituent"
1685 if (parent1 >= 0 && parent2 >= 0) {
1686 if (lowest_constituent[parent1] > lowest_constituent[parent2])
1687 std::swap(parent1, parent2);
1688 }
1689 // then actually run through the parents to extract the constituents...
1690 if (parent1 >= 0 && !extracted[parent1])
1691 _extract_tree_parents(parent1,extracted,lowest_constituent,unique_tree);
1692 if (parent2 >= 0 && !extracted[parent2])
1693 _extract_tree_parents(parent2,extracted,lowest_constituent,unique_tree);
1694 // finally declare this position to be accounted for and push it
1695 // onto our list.
1696 unique_tree.push_back(position);
1697 extracted[position] = true;
1698 }
1699}
1700
1701
1702//======================================================================
1703/// carries out the bookkeeping associated with the step of recombining
1704/// jet_i and jet_j (assuming a distance dij) and returns the index
1705/// of the recombined jet, newjet_k.
1706void ClusterSequence::_do_ij_recombination_step(
1707 const int jet_i, const int jet_j,
1708 const double dij,
1709 int & newjet_k) {
1710
1711 // Create the new jet by recombining the first two.
1712 //
1713 // For efficiency reasons, use a ctr that initialises only the
1714 // shared pointers, since the rest of the info will anyway be dealt
1715 // with by the recombiner.
1716 PseudoJet newjet(false);
1717 _jet_def.recombiner()->recombine(_jets[jet_i], _jets[jet_j], newjet);
1718 _jets.push_back(newjet);
1719 // original version...
1720 //_jets.push_back(_jets[jet_i] + _jets[jet_j]);
1721
1722 // get its index
1723 newjet_k = _jets.size()-1;
1724
1725 // get history index
1726 int newstep_k = _history.size();
1727 // and provide jet with the info
1728 _jets[newjet_k].set_cluster_hist_index(newstep_k);
1729
1730 // finally sort out the history
1731 int hist_i = _jets[jet_i].cluster_hist_index();
1732 int hist_j = _jets[jet_j].cluster_hist_index();
1733
1734 _add_step_to_history(min(hist_i, hist_j), max(hist_i,hist_j),
1735 newjet_k, dij);
1736
1737 // _add_step_to_history(newstep_k, min(hist_i, hist_j), max(hist_i,hist_j),
1738 // newjet_k, dij);
1739
1740
1741}
1742
1743
1744//======================================================================
1745/// carries out the bookkeeping associated with the step of recombining
1746/// jet_i with the beam
1747void ClusterSequence::_do_iB_recombination_step(
1748 const int jet_i, const double diB) {
1749 // recombine the jet with the beam
1750 _add_step_to_history(_jets[jet_i].cluster_hist_index(),BeamJet,
1751 Invalid, diB);
1752
1753 // // get history index
1754 // int newstep_k = _history.size();
1755 //
1756 // _add_step_to_history(newstep_k,_jets[jet_i].cluster_hist_index(),BeamJet,
1757 // Invalid, diB);
1758
1759}
1760
1761
1762// make sure the static member _changed_strategy_warning is defined.
1763LimitedWarning ClusterSequence::_changed_strategy_warning;
1764
1765
1766//----------------------------------------------------------------------
1767void ClusterSequence::_set_structure_shared_ptr(PseudoJet & j) {
1768 j.set_structure_shared_ptr(_structure_shared_ptr);
1769 // record the use count of the structure shared point to help
1770 // in case we want to ask the CS to handle its own memory
1771 _update_structure_use_count();
1772}
1773
1774
1775//----------------------------------------------------------------------
1776void ClusterSequence::_update_structure_use_count() {
1777 // record the use count of the structure shared point to help
1778 // in case we want to ask the CS to handle its own memory
1779 _structure_use_count_after_construction = _structure_shared_ptr.use_count();
1780}
1781
1782//----------------------------------------------------------------------
1783/// by calling this routine you tell the ClusterSequence to delete
1784/// itself when all the Pseudojets associated with it have gone out
1785/// of scope.
1786void ClusterSequence::delete_self_when_unused() {
1787 // the trick we use to handle this is to modify the use count;
1788 // that way the structure will be deleted when there are no external
1789 // objects left associated the CS and the structure's destructor will then
1790 // look after deleting the cluster sequence
1791
1792 // first make sure that there is at least one other object
1793 // associated with the CS
1794 int new_count = _structure_shared_ptr.use_count() - _structure_use_count_after_construction;
1795 if (new_count <= 0) {
1796 throw Error("delete_self_when_unused may only be called if at least one object outside the CS (e.g. a jet) is already associated with the CS");
1797 }
1798
1799//std::shared_ptr: // set_count is not available with CXX11 shared pointers. There
1800//std::shared_ptr: // we'll use a different mechanism for handling self-deleting
1801//std::shared_ptr: // ClusterSequences
1802//std::shared_ptr: #ifndef FASTJET_HAVE_THREAD_SAFETY
1803 _structure_shared_ptr.set_count(new_count);
1804//std::shared_ptr: #endif // FASTJET_HAVE_THREAD_SAFETY
1805 _deletes_self_when_unused = true;
1806}
1807
1808
1809FASTJET_END_NAMESPACE
1810
Contains any information related to the clustering that should be directly accessible to PseudoJet.
virtual void set_associated_cs(const ClusterSequence *new_cs)
set the associated csw
deals with clustering
std::vector< history_element > _history
this vector will contain the branching history; for each stage, _history[i].jetp_index indicates wher...
std::vector< PseudoJet > _jets
This contains the physical PseudoJets; for each PseudoJet one can find the corresponding position in ...
base class corresponding to errors that can be thrown by FastJet
Definition Error.hh:52
base class providing interface for a generic function of a PseudoJet
class corresponding to critical internal errors
Definition Error.hh:138
class that is intended to hold a full definition of the jet clusterer
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
double perp2() const
returns the squared transverse momentum
Definition PseudoJet.hh:156
void set_cluster_hist_index(const int index)
set the cluster_hist_index, intended to be used by clustering routines.
Definition PseudoJet.hh:836
double kt2() const
returns the squared transverse momentum
Definition PseudoJet.hh:160
bool has_valid_cluster_sequence() const
returns true if this PseudoJet has an associated and still valid(ated) ClusterSequence.
Definition PseudoJet.cc:544
int cluster_hist_index() const
return the cluster_hist_index, intended to be used by clustering routines.
Definition PseudoJet.hh:834
void set_structure_shared_ptr(const SharedPtr< PseudoJetStructureBase > &structure_in)
set the associated structure
Definition PseudoJet.cc:561
const ClusterSequence * associated_cluster_sequence() const
get a (const) pointer to the parent ClusterSequence (NULL if inexistent)
Definition PseudoJet.cc:534
Strategy
the various options for the algorithmic strategy to adopt in clustering events with kt and cambridge ...
@ N2Plain
fastest below 50
@ N3Dumb
worse even than the usual N^3 algorithms
@ N2MHTLazy9AntiKtSeparateGhosts
Like N2MHTLazy9 in a number of respects, but does not calculate ghost-ghost distances and so does not...
@ N2Tiled
fastest from about 50..500
@ NlnN3pi
legacy N ln N using 3pi coverage of cylinder.
@ N2MHTLazy9Alt
Like to N2MHTLazy9 but uses slightly different optimizations, e.g.
@ N2PlainEEAccurate
a variant of N2Plain strategy for native algorithms (eekt, ee_genkt) that uses a more accurate calcu...
@ plugin_strategy
the plugin has been used...
@ N2MHTLazy9
only looks into a neighbouring tile for a particle's nearest neighbour (NN) if that particle's in-til...
@ NlnN
best of the NlnN variants – best overall for N>10^4.
@ NlnN4pi
legacy N ln N using 4pi coverage of cylinder
@ NlnNCam
Chan's closest pair method (in a variant with 2pi+minimal extra variant), for use exclusively with th...
@ N2PoorTiled
legacy
@ NlnNCam2pi2R
Chan's closest pair method (in a variant with 2pi+2R coverage), for use exclusively with the Cambridg...
@ N2MHTLazy25
Similar to N2MHTLazy9, but uses tiles of size R/2 and a 5x5 tile grid around the particle.
@ NlnNCam4pi
Chan's closest pair method (in a variant with 4pi coverage), for use exclusively with the Cambridge a...
@ N2MinHeapTiled
faster that N2Tiled above about 500 particles; differs from it by retainig the di(closest j) distance...
JetAlgorithm
the various families of jet-clustering algorithm
@ ee_genkt_algorithm
the e+e- genkt algorithm (R > 2 and p=1 gives ee_kt)
@ genkt_algorithm
like the k_t but with distance measures dij = min(kti^{2p},ktj^{2p}) Delta R_{ij}^2 / R^2 diB = 1/kti...
@ plugin_algorithm
any plugin algorithm supplied by the user
@ ee_kt_algorithm
the e+e- kt algorithm
@ cambridge_algorithm
the longitudinally invariant variant of the cambridge algorithm (aka Aachen algoithm).
@ cambridge_for_passive_algorithm
a version of cambridge with a special distance measure for particles whose pt is < extra_param(); thi...
@ antikt_algorithm
like the k_t but with distance measures dij = min(1/kti^2,1/ktj^2) Delta R_{ij}^2 / R^2 diB = 1/kti^2
@ kt_algorithm
the longitudinally invariant kt algorithm
string fastjet_version_string()
return a string containing information about the release
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 recombined with another jet to form its child.
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.