Fork me on GitHub

source: git/external/fastjet/contribs/Nsubjettiness/AxesDefinition.hh@ 8713dee

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 8713dee was 973b92a, checked in by Pavel Demin <pavel.demin@…>, 9 years ago

update FastJet library to 3.1.3 and Nsubjettiness library to 2.2.1

  • Property mode set to 100644
File size: 47.0 KB
Line 
1// Nsubjettiness Package
2// Questions/Comments? jthaler@jthaler.net
3//
4// Copyright (c) 2011-14
5// Jesse Thaler, Ken Van Tilburg, Christopher K. Vermilion, and TJ Wilkason
6//
7// $Id: AxesDefinition.hh 833 2015-07-23 14:35:23Z jthaler $
8//----------------------------------------------------------------------
9// This file is part of FastJet contrib.
10//
11// It is free software; you can redistribute it and/or modify it under
12// the terms of the GNU General Public License as published by the
13// Free Software Foundation; either version 2 of the License, or (at
14// your option) any later version.
15//
16// It is distributed in the hope that it will be useful, but WITHOUT
17// ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
18// or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public
19// License for more details.
20//
21// You should have received a copy of the GNU General Public License
22// along with this code. If not, see <http://www.gnu.org/licenses/>.
23//----------------------------------------------------------------------
24
25#ifndef __FASTJET_CONTRIB_AXES_DEFINITION_HH__
26#define __FASTJET_CONTRIB_AXES_DEFINITION_HH__
27
28
29#include "MeasureDefinition.hh"
30#include "ExtraRecombiners.hh"
31
32#include "fastjet/PseudoJet.hh"
33#include <fastjet/LimitedWarning.hh>
34
35#include <iomanip>
36#include <cmath>
37#include <vector>
38#include <list>
39
40FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
41
42namespace contrib {
43
44// The following AxesDefinitions are currently available (and the relevant arguments, if needed)
45class KT_Axes;
46class CA_Axes;
47class AntiKT_Axes; // (R0)
48class WTA_KT_Axes;
49class WTA_CA_Axes;
50class GenKT_Axes; // (p, R0 = infinity)
51class WTA_GenKT_Axes; // (p, R0 = infinity)
52class GenET_GenKT_Axes; // (delta, p, R0 = infinity)
53class Manual_Axes;
54
55class OnePass_KT_Axes;
56class OnePass_CA_Axes;
57class OnePass_AntiKT_Axes; // (R0)
58class OnePass_WTA_KT_Axes;
59class OnePass_WTA_CA_Axes;
60class OnePass_GenKT_Axes; // (p, R0 = infinity)
61class OnePass_WTA_GenKT_Axes; // (p, R0 = infinity)
62class OnePass_GenET_GenKT_Axes; // (delta, p, R0 = infinity)
63class OnePass_Manual_Axes;
64
65class MultiPass_Axes; // (NPass) (currently only defined for KT_Axes)
66class MultiPass_Manual_Axes; // (NPass)
67
68class Comb_GenKT_Axes; // (nExtra, p, R0 = infinity)
69class Comb_WTA_GenKT_Axes; // (nExtra, p, R0 = infinity)
70class Comb_GenET_GenKT_Axes; // (nExtra, delta, p, R0 = infinity)
71
72///////
73//
74// AxesDefinition
75//
76///////
77
78///------------------------------------------------------------------------
79/// \class AxesDefinition
80/// \brief Base class for axes definitions
81///
82/// A generic AxesDefinition first finds a set of seed axes.
83/// Then, if desired, uses measure information
84/// (from MeasureDefinition) to refine those axes starting from those seed axes.
85/// The AxesDefinitions are typically based on sequential jet algorithms.
86///------------------------------------------------------------------------
87class AxesDefinition {
88
89public:
90
91 /// This function should be overloaded in all derived classes, and defines how to find the seed axes.
92 /// If desired, the measure information (which might be NULL) can be used to test multiple axes choices, but should
93 /// not be used for iterative refining (since that is the job of MeasureDefinition).
94 virtual std::vector<fastjet::PseudoJet> get_starting_axes(int n_jets,
95 const std::vector<fastjet::PseudoJet>& inputs,
96 const MeasureDefinition * measure) const = 0;
97
98 /// Short description of AxesDefinitions (and any parameters)
99 virtual std::string short_description() const = 0;
100
101 /// Long description of AxesDefinitions (and any parameters)
102 virtual std::string description() const = 0;
103
104 /// This has to be defined in all derived classes, and allows these to be copied around.
105 virtual AxesDefinition* create() const = 0;
106
107public:
108
109 /// Starting from seeds, refine axes using one or more passes.
110 /// Note that in order to do >0 passes, we need information from the MeasureDefinition about how to do the appropriate minimization.
111 std::vector<fastjet::PseudoJet> get_refined_axes(int n_jets,
112 const std::vector<fastjet::PseudoJet>& inputs,
113 const std::vector<fastjet::PseudoJet>& seedAxes,
114 const MeasureDefinition * measure = NULL) const {
115
116 assert(n_jets == (int)seedAxes.size()); //added int casting to get rid of compiler warning
117
118 if (_Npass == 0) {
119 // no refining, just use seeds
120 return seedAxes;
121 } else if (_Npass == 1) {
122 if (measure == NULL) throw Error("AxesDefinition: One-pass minimization requires specifying a MeasureDefinition.");
123
124 // do one pass minimum using measure definition
125 return measure->get_one_pass_axes(n_jets, inputs, seedAxes,_nAttempts,_accuracy);
126 } else {
127 if (measure == NULL) throw Error("AxesDefinition: Multi-pass minimization requires specifying a MeasureDefinition.");
128 return get_multi_pass_axes(n_jets, inputs, seedAxes, measure);
129 }
130 }
131
132 /// Combines get_starting_axes with get_refined_axes.
133 /// In the Njettiness class, these two steps are done separately in order to store seed axes information.
134 std::vector<fastjet::PseudoJet> get_axes(int n_jets,
135 const std::vector<fastjet::PseudoJet>& inputs,
136 const MeasureDefinition * measure = NULL) const {
137 std::vector<fastjet::PseudoJet> seedAxes = get_starting_axes(n_jets, inputs, measure);
138 return get_refined_axes(n_jets,inputs,seedAxes,measure);
139 }
140
141
142 /// Short-hand for the get_axes function. Useful when trying to write terse code.
143 inline std::vector<fastjet::PseudoJet> operator() (int n_jets,
144 const std::vector<fastjet::PseudoJet>& inputs,
145 const MeasureDefinition * measure = NULL) const {
146 return get_axes(n_jets,inputs,measure);
147 }
148
149 /// \enum AxesRefiningEnum
150 /// Defines the cases of zero pass and one pass for convenience
151 enum AxesRefiningEnum {
152 UNDEFINED_REFINE = -1, // added to create a default value
153 NO_REFINING = 0,
154 ONE_PASS = 1,
155 MULTI_PASS = 100,
156 };
157
158 /// A integer that is used externally to decide how to do multi-pass minimization
159 int nPass() const { return _Npass; }
160
161 /// A flag that indicates whether results are deterministics.
162 bool givesRandomizedResults() const {
163 return (_Npass > 1);
164 }
165
166 /// A flag that indicates whether manual axes are being used.
167 bool needsManualAxes() const {
168 return _needsManualAxes; // if there is no starting axes finder
169 }
170
171 /// Allows user to change number of passes. Also used internally to set nPass.
172 /// Can also specify details of one/multi pass minimziation
173 void setNPass(int nPass,
174 int nAttempts = 1000,
175 double accuracy = 0.0001,
176 double noise_range = 1.0 // only needed for MultiPass minimization
177 )
178 {
179 _Npass = nPass;
180 _nAttempts = nAttempts;
181 _accuracy = accuracy;
182 _noise_range = noise_range;
183 if (nPass < 0) throw Error("AxesDefinition requires a nPass >= 0");
184 }
185
186 /// Destructor
187 virtual ~AxesDefinition() {};
188
189protected:
190
191 /// Default constructor contains no information. Number of passes has to be set
192 /// manually by derived classes using setNPass function.
193 AxesDefinition() : _Npass(UNDEFINED_REFINE),
194 _nAttempts(0),
195 _accuracy(0.0),
196 _noise_range(0.0),
197 _needsManualAxes(false) {}
198
199 /// Does multi-pass minimization by randomly jiggling the axes within _noise_range
200 std::vector<fastjet::PseudoJet> get_multi_pass_axes(int n_jets,
201 const std::vector<fastjet::PseudoJet>& inputs,
202 const std::vector<fastjet::PseudoJet>& seedAxes,
203 const MeasureDefinition* measure) const;
204
205 /// Function to jiggle axes within _noise_range
206 PseudoJet jiggle(const PseudoJet& axis) const;
207
208 int _Npass; ///< Number of passes (0 = no refining, 1 = one-pass, >1 multi-pass)
209 int _nAttempts; ///< Number of attempts per pass
210 double _accuracy; ///< Accuracy goal per pass
211 double _noise_range; ///< Noise in rapidity/phi (for multi-pass minimization only)
212 bool _needsManualAxes; ///< Flag to indicate special case of manual axes
213};
214
215///------------------------------------------------------------------------
216/// \class ExclusiveJetAxes
217/// \brief Base class for axes defined from exclusive jet algorithm
218///
219/// This class finds axes by clustering particles with an exclusive jet definition.
220/// This can be implemented with different jet algorithms. The user can call this directly
221/// using their favorite fastjet::JetDefinition
222///------------------------------------------------------------------------
223class ExclusiveJetAxes : public AxesDefinition {
224
225public:
226 /// Constructor takes JetDefinition as an argument
227 ExclusiveJetAxes(fastjet::JetDefinition def)
228 : AxesDefinition(), _def(def) {
229 setNPass(NO_REFINING); // default to no minimization
230 }
231
232 /// Starting axes obtained by creating a cluster sequenence and running exclusive_jets.
233 virtual std::vector<fastjet::PseudoJet> get_starting_axes(int n_jets,
234 const std::vector <fastjet::PseudoJet> & inputs,
235 const MeasureDefinition * ) const {
236 fastjet::ClusterSequence jet_clust_seq(inputs, _def);
237
238 std::vector<fastjet::PseudoJet> axes = jet_clust_seq.exclusive_jets_up_to(n_jets);
239
240 if ((int)axes.size() < n_jets) {
241 _too_few_axes_warning.warn("ExclusiveJetAxes::get_starting_axes: Fewer than N axes found; results are unpredictable.");
242 axes.resize(n_jets); // resize to make sure there are enough axes to not yield an error elsewhere
243 }
244
245 return axes;
246 }
247
248 /// Short description
249 virtual std::string short_description() const { return "ExclAxes";}
250 /// Long description
251 virtual std::string description() const { return "ExclAxes: " + _def.description();}
252
253 /// To make it possible to copy around.
254 virtual ExclusiveJetAxes* create() const {return new ExclusiveJetAxes(*this);}
255
256private:
257 fastjet::JetDefinition _def; ///< Jet definition to use.
258 static LimitedWarning _too_few_axes_warning;
259};
260
261///------------------------------------------------------------------------
262/// \class ExclusiveCombinatorialJetAxes
263/// \brief Base class for axes defined from exclusive jet algorithm, checking combinatorial options
264///
265/// This class finds axes by clustering particles with an exclusive jet definition.
266/// It takes an extra number of jets (specificed by the user via nExtra), and then finds the set of N that minimizes N-jettiness.
267/// WARNING: If one wants to be guarenteed that results improve by increasing nExtra, then one should use
268/// winner-take-all-style recombination schemes
269///------------------------------------------------------------------------
270class ExclusiveCombinatorialJetAxes : public AxesDefinition {
271
272public:
273 /// Constructor takes JetDefinition and nExtra as options (nExtra=0 acts the same as ExclusiveJetAxes)
274 ExclusiveCombinatorialJetAxes(fastjet::JetDefinition def, int nExtra = 0)
275 : AxesDefinition(), _def(def), _nExtra(nExtra) {
276 if (nExtra < 0) throw Error("Need nExtra >= 0");
277 setNPass(NO_REFINING); // default to no minimization
278 }
279
280 /// Find n_jets + _nExtra axes, and then choose the n_jets subset with the smallest N-(sub)jettiness value.
281 virtual std::vector<fastjet::PseudoJet> get_starting_axes(int n_jets,
282 const std::vector<fastjet::PseudoJet> & inputs,
283 const MeasureDefinition *measure) const {
284 int starting_number = n_jets + _nExtra;
285 fastjet::ClusterSequence jet_clust_seq(inputs, _def);
286 std::vector<fastjet::PseudoJet> starting_axes = jet_clust_seq.exclusive_jets_up_to(starting_number);
287
288 if ((int)starting_axes.size() < n_jets) {
289 _too_few_axes_warning.warn("ExclusiveCombinatorialJetAxes::get_starting_axes: Fewer than N + nExtra axes found; results are unpredictable.");
290 starting_axes.resize(n_jets); // resize to make sure there are enough axes to not yield an error elsewhere
291 }
292
293 std::vector<fastjet::PseudoJet> final_axes;
294
295 // check so that no computation time is wasted if there are no extra axes
296 if (_nExtra == 0) final_axes = starting_axes;
297
298 else {
299
300 // define string of 1's based on number of desired jets
301 std::string bitmask(n_jets, 1);
302 // expand the array size to the total number of jets with extra 0's at the end, makes string easy to permute
303 bitmask.resize(starting_number, 0);
304
305 double min_tau = std::numeric_limits<double>::max();
306 std::vector<fastjet::PseudoJet> temp_axes;
307
308 do {
309
310 temp_axes.clear();
311
312 // only take an axis if it is listed as true (1) in the string
313 for (int i = 0; i < (int)starting_axes.size(); ++i) {
314 if (bitmask[i]) temp_axes.push_back(starting_axes[i]);
315 }
316
317 double temp_tau = measure->result(inputs, temp_axes);
318 if (temp_tau < min_tau) {
319 min_tau = temp_tau;
320 final_axes = temp_axes;
321 }
322
323 // permutes string of 1's and 0's according to next lexicographic ordering and returns true
324 // continues to loop through all possible lexicographic orderings
325 // returns false and breaks the loop when there are no more possible orderings
326 } while (std::prev_permutation(bitmask.begin(), bitmask.end()));
327 }
328
329 return final_axes;
330 }
331
332 /// Short description
333 virtual std::string short_description() const { return "ExclCombAxes";}
334 /// Long description
335 virtual std::string description() const { return "ExclCombAxes: " + _def.description();}
336 /// To make it possible to copy around.
337 virtual ExclusiveCombinatorialJetAxes* create() const {return new ExclusiveCombinatorialJetAxes(*this);}
338
339private:
340 fastjet::JetDefinition _def; ///< Jet definition to use
341 int _nExtra; ///< Extra axes to find
342 static LimitedWarning _too_few_axes_warning;
343};
344
345///------------------------------------------------------------------------
346/// \class HardestJetAxes
347/// \brief Base class for axes defined from an inclusive jet algorithm
348///
349/// This class finds axes by running an inclusive algorithm and then finding the n hardest jets.
350/// This can be implemented with different jet algorithms, and can be called by the user.
351///------------------------------------------------------------------------
352class HardestJetAxes : public AxesDefinition {
353public:
354 /// Constructor takes JetDefinition
355 HardestJetAxes(fastjet::JetDefinition def)
356 : AxesDefinition(), _def(def) {
357 setNPass(NO_REFINING); // default to no minimization
358 }
359
360 /// Finds seed axes by running a ClusterSequence, running inclusive_jets, and finding the N hardest
361 virtual std::vector<fastjet::PseudoJet> get_starting_axes(int n_jets,
362 const std::vector <fastjet::PseudoJet> & inputs,
363 const MeasureDefinition * ) const {
364 fastjet::ClusterSequence jet_clust_seq(inputs, _def);
365 std::vector<fastjet::PseudoJet> axes = sorted_by_pt(jet_clust_seq.inclusive_jets());
366
367 if ((int)axes.size() < n_jets) {
368 _too_few_axes_warning.warn("HardestJetAxes::get_starting_axes: Fewer than N axes found; results are unpredictable.");
369 }
370
371 axes.resize(n_jets); // only keep n hardest
372 return axes;
373 }
374
375 /// Short description
376 virtual std::string short_description() const { return "HardAxes";}
377 /// Long description
378 virtual std::string description() const { return "HardAxes: " + _def.description();}
379 /// To make it possible to copy around.
380 virtual HardestJetAxes* create() const {return new HardestJetAxes(*this);}
381
382private:
383 fastjet::JetDefinition _def; ///< Jet Definition to use.
384
385 static LimitedWarning _too_few_axes_warning;
386
387};
388
389///------------------------------------------------------------------------
390/// \class KT_Axes
391/// \brief Axes from exclusive kT
392///
393/// Axes from kT algorithm with E_scheme recombination.
394///------------------------------------------------------------------------
395class KT_Axes : public ExclusiveJetAxes {
396public:
397 /// Constructor
398 KT_Axes()
399 : ExclusiveJetAxes(fastjet::JetDefinition(fastjet::kt_algorithm,
400 fastjet::JetDefinition::max_allowable_R, //maximum jet radius constant
401 fastjet::E_scheme,
402 fastjet::Best)
403 ) {
404 setNPass(NO_REFINING);
405 }
406
407 /// Short description
408 virtual std::string short_description() const {
409 return "KT";
410 };
411
412 /// Long description
413 virtual std::string description() const {
414 std::stringstream stream;
415 stream << std::fixed << std::setprecision(2)
416 << "KT Axes";
417 return stream.str();
418 };
419
420 /// For copying purposes
421 virtual KT_Axes* create() const {return new KT_Axes(*this);}
422
423};
424
425///------------------------------------------------------------------------
426/// \class CA_Axes
427/// \brief Axes from exclusive CA
428///
429/// Axes from CA algorithm with E_scheme recombination.
430///------------------------------------------------------------------------
431class CA_Axes : public ExclusiveJetAxes {
432public:
433 /// Constructor
434 CA_Axes()
435 : ExclusiveJetAxes(fastjet::JetDefinition(fastjet::cambridge_algorithm,
436 fastjet::JetDefinition::max_allowable_R, //maximum jet radius constant
437 fastjet::E_scheme,
438 fastjet::Best)
439 ) {
440 setNPass(NO_REFINING);
441 }
442
443 /// Short description
444 virtual std::string short_description() const {
445 return "CA";
446 };
447
448 /// Long description
449 virtual std::string description() const {
450 std::stringstream stream;
451 stream << std::fixed << std::setprecision(2)
452 << "CA Axes";
453 return stream.str();
454 };
455
456 /// For copying purposes
457 virtual CA_Axes* create() const {return new CA_Axes(*this);}
458
459};
460
461
462///------------------------------------------------------------------------
463/// \class AntiKT_Axes
464/// \brief Axes from inclusive anti-kT
465///
466/// Axes from anti-kT algorithm and E_scheme.
467/// The one parameter R0 is subjet radius
468///------------------------------------------------------------------------
469class AntiKT_Axes : public HardestJetAxes {
470
471public:
472 /// Constructor. Takes jet radius as argument
473 AntiKT_Axes(double R0)
474 : HardestJetAxes(fastjet::JetDefinition(fastjet::antikt_algorithm,
475 R0,
476 fastjet::E_scheme,
477 fastjet::Best)
478 ), _R0(R0) {
479 setNPass(NO_REFINING);
480 }
481
482 /// Short description
483 virtual std::string short_description() const {
484 std::stringstream stream;
485 stream << std::fixed << std::setprecision(2)
486 << "AKT" << _R0;
487 return stream.str();
488 };
489
490 /// Long description
491 virtual std::string description() const {
492 std::stringstream stream;
493 stream << std::fixed << std::setprecision(2)
494 << "Anti-KT Axes (R0 = " << _R0 << ")";
495 return stream.str();
496 };
497
498 /// For copying purposes
499 virtual AntiKT_Axes* create() const {return new AntiKT_Axes(*this);}
500
501protected:
502 double _R0; ///< AKT jet radius
503
504};
505
506///------------------------------------------------------------------------
507/// \class JetDefinitionWrapper
508/// \brief Wrapper for jet definitions (for memory management)
509///
510/// This class was introduced to avoid issue of a FastJet bug when using genKT clustering
511/// Now using this for all AxesDefinition with a manual recombiner to use the delete_recombiner_when_unused function
512///------------------------------------------------------------------------
513class JetDefinitionWrapper {
514
515public:
516
517 /// Default Constructor
518 JetDefinitionWrapper(JetAlgorithm jet_algorithm_in, double R_in, double xtra_param_in, const JetDefinition::Recombiner *recombiner) {
519 jet_def = fastjet::JetDefinition(jet_algorithm_in, R_in, xtra_param_in);
520 jet_def.set_recombiner(recombiner);
521 jet_def.delete_recombiner_when_unused(); // added to prevent memory leaks
522 }
523
524 /// Additional constructor so that build-in FastJet algorithms can also be called
525 JetDefinitionWrapper(JetAlgorithm jet_algorithm_in, double R_in, const JetDefinition::Recombiner *recombiner, fastjet::Strategy strategy_in) {
526 jet_def = fastjet::JetDefinition(jet_algorithm_in, R_in, recombiner, strategy_in);
527 jet_def.delete_recombiner_when_unused();
528 }
529
530 /// Return jet definition
531 JetDefinition getJetDef() {
532 return jet_def;
533 }
534
535private:
536 JetDefinition jet_def; ///< my jet definition
537};
538
539///------------------------------------------------------------------------
540/// \class WTA_KT_Axes
541/// \brief Axes from exclusive kT, winner-take-all recombination
542///
543/// Axes from kT algorithm and winner-take-all recombination
544///------------------------------------------------------------------------
545class WTA_KT_Axes : public ExclusiveJetAxes {
546public:
547 /// Constructor
548 WTA_KT_Axes()
549 : ExclusiveJetAxes(JetDefinitionWrapper(fastjet::kt_algorithm,
550 fastjet::JetDefinition::max_allowable_R, //maximum jet radius constant
551 _recomb = new WinnerTakeAllRecombiner(), // Needs to be explicitly declared (this will be deleted by JetDefinitionWrapper)
552 fastjet::Best).getJetDef()
553 ) {
554 setNPass(NO_REFINING);
555 }
556
557 /// Short description
558 virtual std::string short_description() const {
559 return "WTA KT";
560 };
561
562 /// Long description
563 virtual std::string description() const {
564 std::stringstream stream;
565 stream << std::fixed << std::setprecision(2)
566 << "Winner-Take-All KT Axes";
567 return stream.str();
568 };
569
570 /// For copying purposes
571 virtual WTA_KT_Axes* create() const {return new WTA_KT_Axes(*this);}
572
573private:
574 const WinnerTakeAllRecombiner *_recomb; ///< Internal recombiner
575
576};
577
578///------------------------------------------------------------------------
579/// \class WTA_CA_Axes
580/// \brief Axes from exclusive CA, winner-take-all recombination
581///
582/// Axes from CA algorithm and winner-take-all recombination
583///------------------------------------------------------------------------
584class WTA_CA_Axes : public ExclusiveJetAxes {
585public:
586 /// Constructor
587 WTA_CA_Axes()
588 : ExclusiveJetAxes(JetDefinitionWrapper(fastjet::cambridge_algorithm,
589 fastjet::JetDefinition::max_allowable_R, //maximum jet radius constant
590 _recomb = new WinnerTakeAllRecombiner(), // Needs to be explicitly declared (this will be deleted by JetDefinitionWrapper)
591 fastjet::Best).getJetDef()) {
592 setNPass(NO_REFINING);
593 }
594
595 /// Short description
596 virtual std::string short_description() const {
597 return "WTA CA";
598 };
599
600 /// Long descriptions
601 virtual std::string description() const {
602 std::stringstream stream;
603 stream << std::fixed << std::setprecision(2)
604 << "Winner-Take-All CA Axes";
605 return stream.str();
606 };
607
608 /// For copying purposes
609 virtual WTA_CA_Axes* create() const {return new WTA_CA_Axes(*this);}
610
611private:
612 const WinnerTakeAllRecombiner *_recomb; ///< Internal recombiner
613
614};
615
616
617///------------------------------------------------------------------------
618/// \class GenKT_Axes
619/// \brief Axes from exclusive generalized kT
620///
621/// Axes from a general KT algorithm (standard E-scheme recombination)
622/// Requires the power of the KT algorithm to be used and the radius parameter
623///------------------------------------------------------------------------
624class GenKT_Axes : public ExclusiveJetAxes {
625
626public:
627 /// Constructor
628 GenKT_Axes(double p, double R0 = fastjet::JetDefinition::max_allowable_R)
629 : ExclusiveJetAxes(fastjet::JetDefinition(fastjet::genkt_algorithm,
630 R0,
631 p)), _p(p), _R0(R0) {
632 if (p < 0) throw Error("GenKT_Axes: Currently only p >=0 is supported.");
633 setNPass(NO_REFINING);
634 }
635
636 /// Short description
637 virtual std::string short_description() const {
638 std::stringstream stream;
639 stream << std::fixed << std::setprecision(2)
640 << "GenKT Axes";
641 return stream.str();
642 };
643
644 /// Long descriptions
645 virtual std::string description() const {
646 std::stringstream stream;
647 stream << std::fixed << std::setprecision(2)
648 << "General KT (p = " << _p << "), R0 = " << _R0;
649 return stream.str();
650 };
651
652 /// For copying purposes
653 virtual GenKT_Axes* create() const {return new GenKT_Axes(*this);}
654
655protected:
656 double _p; ///< genkT power
657 double _R0; ///< jet radius
658};
659
660
661///------------------------------------------------------------------------
662/// \class WTA_GenKT_Axes
663/// \brief Axes from exclusive generalized kT, winner-take-all recombination
664///
665/// Axes from a general KT algorithm with a Winner Take All Recombiner
666/// Requires the power of the KT algorithm to be used and the radius parameter
667///------------------------------------------------------------------------
668class WTA_GenKT_Axes : public ExclusiveJetAxes {
669
670public:
671 /// Constructor
672 WTA_GenKT_Axes(double p, double R0 = fastjet::JetDefinition::max_allowable_R)
673 : ExclusiveJetAxes(JetDefinitionWrapper(fastjet::genkt_algorithm,
674 R0,
675 p,
676 _recomb = new WinnerTakeAllRecombiner()
677 ).getJetDef()), _p(p), _R0(R0) {
678 if (p < 0) throw Error("WTA_GenKT_Axes: Currently only p >=0 is supported.");
679 setNPass(NO_REFINING);
680 }
681
682 /// Short description
683 virtual std::string short_description() const {
684 std::stringstream stream;
685 stream << std::fixed << std::setprecision(2)
686 << "WTA, GenKT Axes";
687 return stream.str();
688 };
689
690 /// Long descriptions
691 virtual std::string description() const {
692 std::stringstream stream;
693 stream << std::fixed << std::setprecision(2)
694 << "Winner-Take-All General KT (p = " << _p << "), R0 = " << _R0;
695 return stream.str();
696 };
697
698 /// For copying purposes
699 virtual WTA_GenKT_Axes* create() const {return new WTA_GenKT_Axes(*this);}
700
701protected:
702 double _p; ///< genkT power
703 double _R0; ///< jet radius
704 const WinnerTakeAllRecombiner *_recomb; ///< Internal recombiner
705};
706
707///------------------------------------------------------------------------
708/// \class GenET_GenKT_Axes
709/// \brief Axes from exclusive kT, generalized Et-scheme recombination
710///
711/// Class using general KT algorithm with a more general recombination scheme
712/// Requires power of KT algorithm, power of recombination weights, and radius parameter
713///------------------------------------------------------------------------
714class GenET_GenKT_Axes : public ExclusiveJetAxes {
715
716public:
717 /// Constructor
718 GenET_GenKT_Axes(double delta, double p, double R0 = fastjet::JetDefinition::max_allowable_R)
719 : ExclusiveJetAxes((JetDefinitionWrapper(fastjet::genkt_algorithm, R0, p, _recomb = new GeneralEtSchemeRecombiner(delta))).getJetDef() ),
720 _delta(delta), _p(p), _R0(R0) {
721 if (p < 0) throw Error("GenET_GenKT_Axes: Currently only p >=0 is supported.");
722 if (delta <= 0) throw Error("GenET_GenKT_Axes: Currently only delta >0 is supported.");
723 setNPass(NO_REFINING);
724 }
725
726 /// Short description
727 virtual std::string short_description() const {
728 std::stringstream stream;
729 stream << std::fixed << std::setprecision(2)
730 << "GenET, GenKT Axes";
731 return stream.str();
732 };
733
734 /// Long description
735 virtual std::string description() const {
736 std::stringstream stream;
737 stream << std::fixed << std::setprecision(2);
738 // TODO: if _delta is huge, change to "WTA"
739 if (_delta < std::numeric_limits<int>::max()) stream << "General Recombiner (delta = " << _delta << "), " << "General KT (p = " << _p << ") Axes, R0 = " << _R0;
740 else stream << "Winner-Take-All General KT (p = " << _p << "), R0 = " << _R0;
741
742 return stream.str();
743 };
744
745 /// For copying purposes
746 virtual GenET_GenKT_Axes* create() const {return new GenET_GenKT_Axes(*this);}
747
748protected:
749 double _delta; ///< Recombination pT weighting
750 double _p; ///< GenkT power
751 double _R0; ///< jet radius
752 const GeneralEtSchemeRecombiner *_recomb; ///< Internal recombiner
753};
754
755///------------------------------------------------------------------------
756/// \class OnePass_KT_Axes
757/// \brief Axes from exclusive kT, with one-pass minimization
758///
759/// Onepass minimization from kt axes
760///------------------------------------------------------------------------
761class OnePass_KT_Axes : public KT_Axes {
762public:
763 /// Constructor
764 OnePass_KT_Axes() : KT_Axes() {
765 setNPass(ONE_PASS);
766 }
767
768 /// Short description
769 virtual std::string short_description() const {
770 return "OnePass KT";
771 };
772
773 /// Long description
774 virtual std::string description() const {
775 std::stringstream stream;
776 stream << std::fixed << std::setprecision(2)
777 << "One-Pass Minimization from KT Axes";
778 return stream.str();
779 };
780
781 /// For copying purposes
782 virtual OnePass_KT_Axes* create() const {return new OnePass_KT_Axes(*this);}
783
784
785};
786
787///------------------------------------------------------------------------
788/// \class OnePass_CA_Axes
789/// \brief Axes from exclusive CA, with one-pass minimization
790///
791/// Onepass minimization from CA axes
792///------------------------------------------------------------------------
793class OnePass_CA_Axes : public CA_Axes {
794public:
795 /// Constructor
796 OnePass_CA_Axes() : CA_Axes() {
797 setNPass(ONE_PASS);
798 }
799
800 /// Short description
801 virtual std::string short_description() const {
802 return "OnePass CA";
803 };
804
805 /// Long description
806 virtual std::string description() const {
807 std::stringstream stream;
808 stream << std::fixed << std::setprecision(2)
809 << "One-Pass Minimization from CA Axes";
810 return stream.str();
811 };
812
813 /// For copying purposes
814 virtual OnePass_CA_Axes* create() const {return new OnePass_CA_Axes(*this);}
815
816
817};
818
819///------------------------------------------------------------------------
820/// \class OnePass_AntiKT_Axes
821/// \brief Axes from inclusive anti-kT, with one-pass minimization
822///
823/// Onepass minimization from AntiKT axes, one parameter R0
824///------------------------------------------------------------------------
825class OnePass_AntiKT_Axes : public AntiKT_Axes {
826
827public:
828 /// Constructor
829 OnePass_AntiKT_Axes(double R0) : AntiKT_Axes(R0) {
830 setNPass(ONE_PASS);
831 }
832
833 /// Short Description
834 virtual std::string short_description() const {
835 std::stringstream stream;
836 stream << std::fixed << std::setprecision(2)
837 << "OnePassAKT" << _R0;
838 return stream.str();
839 };
840
841 /// Long description
842 virtual std::string description() const {
843 std::stringstream stream;
844 stream << std::fixed << std::setprecision(2)
845 << "One-Pass Minimization from Anti-KT Axes (R0 = " << _R0 << ")";
846 return stream.str();
847 };
848
849 /// For copying purposes
850 virtual OnePass_AntiKT_Axes* create() const {return new OnePass_AntiKT_Axes(*this);}
851
852};
853
854///------------------------------------------------------------------------
855/// \class OnePass_WTA_KT_Axes
856/// \brief Axes from exclusive kT, winner-take-all recombination, with one-pass minimization
857///
858/// Onepass minimization from winner-take-all kt axes
859///------------------------------------------------------------------------
860class OnePass_WTA_KT_Axes : public WTA_KT_Axes {
861public:
862 /// Constructor
863 OnePass_WTA_KT_Axes() : WTA_KT_Axes() {
864 setNPass(ONE_PASS);
865 }
866
867 /// Short description
868 virtual std::string short_description() const {
869 return "OnePass WTA KT";
870 };
871
872 /// Long description
873 virtual std::string description() const {
874 std::stringstream stream;
875 stream << std::fixed << std::setprecision(2)
876 << "One-Pass Minimization from Winner-Take-All KT Axes";
877 return stream.str();
878 };
879
880 /// For copying purposes
881 virtual OnePass_WTA_KT_Axes* create() const {return new OnePass_WTA_KT_Axes(*this);}
882
883
884};
885
886///------------------------------------------------------------------------
887/// \class OnePass_WTA_CA_Axes
888/// \brief Axes from exclusive CA, winner-take-all recombination, with one-pass minimization
889///
890/// Onepass minimization from winner-take-all CA axes
891///------------------------------------------------------------------------
892class OnePass_WTA_CA_Axes : public WTA_CA_Axes {
893
894public:
895 /// Constructor
896 OnePass_WTA_CA_Axes() : WTA_CA_Axes() {
897 setNPass(ONE_PASS);
898 }
899
900 /// Short description
901 virtual std::string short_description() const {
902 return "OnePass WTA CA";
903 };
904
905 /// Long description
906 virtual std::string description() const {
907 std::stringstream stream;
908 stream << std::fixed << std::setprecision(2)
909 << "One-Pass Minimization from Winner-Take-All CA Axes";
910 return stream.str();
911 };
912
913 /// For copying purposes
914 virtual OnePass_WTA_CA_Axes* create() const {return new OnePass_WTA_CA_Axes(*this);}
915
916};
917
918///------------------------------------------------------------------------
919/// \class OnePass_GenKT_Axes
920/// \brief Axes from exclusive generalized kT with one-pass minimization
921///
922/// Onepass minimization, General KT Axes (standard E-scheme recombination)
923///------------------------------------------------------------------------
924class OnePass_GenKT_Axes : public GenKT_Axes {
925
926public:
927 /// Constructor
928 OnePass_GenKT_Axes(double p, double R0 = fastjet::JetDefinition::max_allowable_R) : GenKT_Axes(p, R0) {
929 setNPass(ONE_PASS);
930 }
931
932 /// Short description
933 virtual std::string short_description() const {
934 return "OnePass GenKT";
935 };
936
937 /// Long description
938 virtual std::string description() const {
939 std::stringstream stream;
940 stream << std::fixed << std::setprecision(2)
941 << "One-Pass Minimization from General KT (p = " << _p << "), R0 = " << _R0;
942 return stream.str();
943 };
944
945 /// For copying purposes
946 virtual OnePass_GenKT_Axes* create() const {return new OnePass_GenKT_Axes(*this);}
947};
948
949///------------------------------------------------------------------------
950/// \class OnePass_WTA_GenKT_Axes
951/// \brief Axes from exclusive generalized kT, winner-take-all recombination, with one-pass minimization
952///
953/// Onepass minimization from winner-take-all, General KT Axes
954///------------------------------------------------------------------------
955class OnePass_WTA_GenKT_Axes : public WTA_GenKT_Axes {
956
957public:
958 /// Constructor
959 OnePass_WTA_GenKT_Axes(double p, double R0 = fastjet::JetDefinition::max_allowable_R) : WTA_GenKT_Axes(p, R0) {
960 setNPass(ONE_PASS);
961 }
962
963 /// Short description
964 virtual std::string short_description() const {
965 return "OnePass WTA GenKT";
966 };
967
968 /// Long description
969 virtual std::string description() const {
970 std::stringstream stream;
971 stream << std::fixed << std::setprecision(2)
972 << "One-Pass Minimization from Winner-Take-All General KT (p = " << _p << "), R0 = " << _R0;
973 return stream.str();
974 };
975
976 /// For copying purposes
977 virtual OnePass_WTA_GenKT_Axes* create() const {return new OnePass_WTA_GenKT_Axes(*this);}
978};
979
980///------------------------------------------------------------------------
981/// \class OnePass_GenET_GenKT_Axes
982/// \brief Axes from exclusive generalized kT, generalized Et-scheme recombination, with one-pass minimization
983///
984/// Onepass minimization from General Recomb, General KT axes
985///------------------------------------------------------------------------
986class OnePass_GenET_GenKT_Axes : public GenET_GenKT_Axes {
987
988public:
989 /// Constructor
990 OnePass_GenET_GenKT_Axes(double delta, double p, double R0 = fastjet::JetDefinition::max_allowable_R) : GenET_GenKT_Axes(delta, p, R0) {
991 setNPass(ONE_PASS);
992 }
993
994 /// Short description
995 virtual std::string short_description() const {
996 return "OnePass GenET, GenKT";
997 };
998
999 /// Long description
1000 virtual std::string description() const {
1001 std::stringstream stream;
1002 stream << std::fixed << std::setprecision(2);
1003 if (_delta < std::numeric_limits<int>::max()) stream << "One-Pass Minimization from General Recombiner (delta = "
1004 << _delta << "), " << "General KT (p = " << _p << ") Axes, R0 = " << _R0;
1005 else stream << "One-Pass Minimization from Winner-Take-All General KT (p = " << _p << "), R0 = " << _R0;
1006 return stream.str();
1007 };
1008
1009 /// For copying purposes
1010 virtual OnePass_GenET_GenKT_Axes* create() const {return new OnePass_GenET_GenKT_Axes(*this);}
1011};
1012
1013
1014///------------------------------------------------------------------------
1015/// \class Manual_Axes
1016/// \brief Manual axes finding
1017///
1018/// Allows the user to set the axes manually
1019///------------------------------------------------------------------------
1020class Manual_Axes : public AxesDefinition {
1021public:
1022 /// Constructor. Note that _needsManualAxes is set to true.
1023 Manual_Axes() : AxesDefinition() {
1024 setNPass(NO_REFINING);
1025 _needsManualAxes = true;
1026 }
1027
1028 /// This is now a dummy function since this is manual mode
1029 virtual std::vector<fastjet::PseudoJet> get_starting_axes(int,
1030 const std::vector<fastjet::PseudoJet>&,
1031 const MeasureDefinition *) const;
1032
1033
1034 /// Short description
1035 virtual std::string short_description() const {
1036 return "Manual";
1037 };
1038
1039 /// Long description
1040 virtual std::string description() const {
1041 std::stringstream stream;
1042 stream << std::fixed << std::setprecision(2)
1043 << "Manual Axes";
1044 return stream.str();
1045 };
1046
1047 /// For copying purposes
1048 virtual Manual_Axes* create() const {return new Manual_Axes(*this);}
1049
1050
1051};
1052
1053///------------------------------------------------------------------------
1054/// \class OnePass_Manual_Axes
1055/// \brief Manual axes finding, with one-pass minimization
1056///
1057/// One pass minimization from manual starting point
1058///------------------------------------------------------------------------
1059class OnePass_Manual_Axes : public Manual_Axes {
1060public:
1061 /// Constructor. Note that _needsManualAxes is set to true.
1062 OnePass_Manual_Axes() : Manual_Axes() {
1063 setNPass(ONE_PASS);
1064 }
1065
1066 /// Short description
1067 virtual std::string short_description() const {
1068 return "OnePass Manual";
1069 };
1070
1071 /// Long description
1072 virtual std::string description() const {
1073 std::stringstream stream;
1074 stream << std::fixed << std::setprecision(2)
1075 << "One-Pass Minimization from Manual Axes";
1076 return stream.str();
1077 };
1078
1079 // For copying purposes
1080 virtual OnePass_Manual_Axes* create() const {return new OnePass_Manual_Axes(*this);}
1081
1082};
1083
1084///------------------------------------------------------------------------
1085/// \class MultiPass_Axes
1086/// \brief Manual axes finding, with multi-pass (randomized) minimization
1087///
1088/// Multi-pass minimization from kT starting point
1089///------------------------------------------------------------------------
1090class MultiPass_Axes : public KT_Axes {
1091
1092public:
1093
1094 /// Constructor
1095 MultiPass_Axes(unsigned int Npass) : KT_Axes() {
1096 setNPass(Npass);
1097 }
1098
1099 /// Short description
1100 virtual std::string short_description() const {
1101 return "MultiPass";
1102 };
1103
1104 /// Long description
1105 virtual std::string description() const {
1106 std::stringstream stream;
1107 stream << std::fixed << std::setprecision(2)
1108 << "Multi-Pass Axes (Npass = " << _Npass << ")";
1109 return stream.str();
1110 };
1111
1112 /// For copying purposs
1113 virtual MultiPass_Axes* create() const {return new MultiPass_Axes(*this);}
1114
1115};
1116
1117///------------------------------------------------------------------------
1118/// \class MultiPass_Manual_Axes
1119/// \brief Axes finding from exclusive kT, with multi-pass (randomized) minimization
1120///
1121/// multi-pass minimization from kT starting point
1122///------------------------------------------------------------------------
1123class MultiPass_Manual_Axes : public Manual_Axes {
1124
1125public:
1126 /// Constructor
1127 MultiPass_Manual_Axes(unsigned int Npass) : Manual_Axes() {
1128 setNPass(Npass);
1129 }
1130
1131 /// Short Description
1132 virtual std::string short_description() const {
1133 return "MultiPass Manual";
1134 };
1135
1136
1137 /// Long description
1138 virtual std::string description() const {
1139 std::stringstream stream;
1140 stream << std::fixed << std::setprecision(2)
1141 << "Multi-Pass Manual Axes (Npass = " << _Npass << ")";
1142 return stream.str();
1143 };
1144
1145 /// For copying purposes
1146 virtual MultiPass_Manual_Axes* create() const {return new MultiPass_Manual_Axes(*this);}
1147
1148};
1149
1150///------------------------------------------------------------------------
1151/// \class Comb_GenKT_Axes
1152/// \brief Axes from exclusive generalized kT with combinatorial testing
1153///
1154/// Axes from kT algorithm (standard E-scheme recombination)
1155/// Requires nExtra parameter and returns set of N that minimizes N-jettiness
1156/// Note that this method is not guaranteed to find a deeper minimum than GenKT_Axes
1157///------------------------------------------------------------------------
1158class Comb_GenKT_Axes : public ExclusiveCombinatorialJetAxes {
1159public:
1160 /// Constructor
1161 Comb_GenKT_Axes(int nExtra, double p, double R0 = fastjet::JetDefinition::max_allowable_R)
1162 : ExclusiveCombinatorialJetAxes(fastjet::JetDefinition(fastjet::genkt_algorithm, R0, p), nExtra),
1163 _p(p), _R0(R0) {
1164 if (p < 0) throw Error("Comb_GenKT_Axes: Currently only p >=0 is supported.");
1165 setNPass(NO_REFINING);
1166 }
1167
1168 /// Short description
1169 virtual std::string short_description() const {
1170 return "N Choose M GenKT";
1171 };
1172
1173 /// Long description
1174 virtual std::string description() const {
1175 std::stringstream stream;
1176 stream << std::fixed << std::setprecision(2)
1177 << "N Choose M Minimization (nExtra = " << _nExtra << ") from General KT (p = " << _p << "), R0 = " << _R0;
1178 return stream.str();
1179 };
1180
1181 /// For copying purposes
1182 virtual Comb_GenKT_Axes* create() const {return new Comb_GenKT_Axes(*this);}
1183
1184private:
1185 double _nExtra; ///< Number of extra axes
1186 double _p; ///< GenkT power
1187 double _R0; ///< jet radius
1188};
1189
1190
1191
1192///------------------------------------------------------------------------
1193/// \class Comb_WTA_GenKT_Axes
1194/// \brief Axes from exclusive generalized kT, winner-take-all recombination, with combinatorial testing
1195///
1196/// Axes from kT algorithm and winner-take-all recombination
1197/// Requires nExtra parameter and returns set of N that minimizes N-jettiness
1198///------------------------------------------------------------------------
1199class Comb_WTA_GenKT_Axes : public ExclusiveCombinatorialJetAxes {
1200public:
1201 /// Constructor
1202 Comb_WTA_GenKT_Axes(int nExtra, double p, double R0 = fastjet::JetDefinition::max_allowable_R)
1203 : ExclusiveCombinatorialJetAxes((JetDefinitionWrapper(fastjet::genkt_algorithm, R0, p, _recomb = new WinnerTakeAllRecombiner())).getJetDef(), nExtra),
1204 _p(p), _R0(R0) {
1205 if (p < 0) throw Error("Comb_WTA_GenKT_Axes: Currently only p >=0 is supported.");
1206 setNPass(NO_REFINING);
1207 }
1208
1209 /// Short description
1210 virtual std::string short_description() const {
1211 return "N Choose M WTA GenKT";
1212 };
1213
1214 /// Long description
1215 virtual std::string description() const {
1216 std::stringstream stream;
1217 stream << std::fixed << std::setprecision(2)
1218 << "N Choose M Minimization (nExtra = " << _nExtra << ") from Winner-Take-All General KT (p = " << _p << "), R0 = " << _R0;
1219 return stream.str();
1220 };
1221
1222 /// For copying purposes
1223 virtual Comb_WTA_GenKT_Axes* create() const {return new Comb_WTA_GenKT_Axes(*this);}
1224
1225private:
1226 double _nExtra; ///< Number of extra axes
1227 double _p; ///< GenkT power
1228 double _R0; ///< jet radius
1229 const WinnerTakeAllRecombiner *_recomb; ///< Internal recombiner
1230};
1231
1232///------------------------------------------------------------------------
1233/// \class Comb_GenET_GenKT_Axes
1234/// \brief Axes from exclusive generalized kT, generalized Et-scheme recombination, with combinatorial testing
1235///
1236/// Axes from kT algorithm and General Et scheme recombination
1237/// Requires nExtra parameter and returns set of N that minimizes N-jettiness
1238///------------------------------------------------------------------------
1239class Comb_GenET_GenKT_Axes : public ExclusiveCombinatorialJetAxes {
1240public:
1241 /// Constructor
1242 Comb_GenET_GenKT_Axes(int nExtra, double delta, double p, double R0 = fastjet::JetDefinition::max_allowable_R)
1243 : ExclusiveCombinatorialJetAxes((JetDefinitionWrapper(fastjet::genkt_algorithm, R0, p, _recomb = new GeneralEtSchemeRecombiner(delta))).getJetDef(), nExtra),
1244 _delta(delta), _p(p), _R0(R0) {
1245 if (p < 0) throw Error("Comb_GenET_GenKT_Axes: Currently only p >=0 is supported.");
1246 if (delta <= 0) throw Error("Comb_GenET_GenKT_Axes: Currently only delta >=0 is supported.");
1247 setNPass(NO_REFINING);
1248 }
1249
1250 /// Short description
1251 virtual std::string short_description() const {
1252 return "N Choose M GenET GenKT";
1253 };
1254
1255 /// Long description
1256 virtual std::string description() const {
1257 std::stringstream stream;
1258 stream << std::fixed << std::setprecision(2);
1259 if (_delta < std::numeric_limits<int>::max()) stream << "N choose M Minimization (nExtra = " << _nExtra
1260 << ") from General Recombiner (delta = " << _delta << "), " << "General KT (p = " << _p << ") Axes, R0 = " << _R0;
1261 else stream << "N choose M Minimization (nExtra = " << _nExtra << ") from Winner-Take-All General KT (p = " << _p << "), R0 = " << _R0;
1262 return stream.str();
1263 };
1264
1265 /// For copying purposes
1266 virtual Comb_GenET_GenKT_Axes* create() const {return new Comb_GenET_GenKT_Axes(*this);}
1267
1268private:
1269 double _nExtra; ///< Number of extra axes
1270 double _delta; ///< Recombination pT weighting exponent
1271 double _p; ///< GenkT power
1272 double _R0; ///< jet radius
1273 const GeneralEtSchemeRecombiner *_recomb; ///< Internal recombiner
1274};
1275
1276
1277} // namespace contrib
1278
1279FASTJET_END_NAMESPACE
1280
1281#endif // __FASTJET_CONTRIB_NJETTINESS_HH__
1282
Note: See TracBrowser for help on using the repository browser.