Fork me on GitHub

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

ImprovedOutputFile Timing dual_readout
Last change on this file since a1b19ea was b7b836a, checked in by Pavel Demin <pavel-demin@…>, 6 years ago

update FastJet library to 3.3.1 and FastJet Contrib library to 1.036

  • Property mode set to 100644
File size: 46.5 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 1130 2018-06-06 12:09:46Z 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 is used by all AxesDefinition with a manual recombiner to
511/// ensure that the delete_recombiner_when_unused function is always called
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 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
573};
574
575///------------------------------------------------------------------------
576/// \class WTA_CA_Axes
577/// \brief Axes from exclusive CA, winner-take-all recombination
578///
579/// Axes from CA algorithm and winner-take-all recombination
580///------------------------------------------------------------------------
581class WTA_CA_Axes : public ExclusiveJetAxes {
582public:
583 /// Constructor
584 WTA_CA_Axes()
585 : ExclusiveJetAxes(JetDefinitionWrapper(fastjet::cambridge_algorithm,
586 fastjet::JetDefinition::max_allowable_R, //maximum jet radius constant
587 new WinnerTakeAllRecombiner(), // Needs to be explicitly declared (this will be deleted by JetDefinitionWrapper)
588 fastjet::Best).getJetDef()) {
589 setNPass(NO_REFINING);
590 }
591
592 /// Short description
593 virtual std::string short_description() const {
594 return "WTA CA";
595 };
596
597 /// Long descriptions
598 virtual std::string description() const {
599 std::stringstream stream;
600 stream << std::fixed << std::setprecision(2)
601 << "Winner-Take-All CA Axes";
602 return stream.str();
603 };
604
605 /// For copying purposes
606 virtual WTA_CA_Axes* create() const {return new WTA_CA_Axes(*this);}
607
608};
609
610
611///------------------------------------------------------------------------
612/// \class GenKT_Axes
613/// \brief Axes from exclusive generalized kT
614///
615/// Axes from a general KT algorithm (standard E-scheme recombination)
616/// Requires the power of the KT algorithm to be used and the radius parameter
617///------------------------------------------------------------------------
618class GenKT_Axes : public ExclusiveJetAxes {
619
620public:
621 /// Constructor
622 GenKT_Axes(double p, double R0 = fastjet::JetDefinition::max_allowable_R)
623 : ExclusiveJetAxes(fastjet::JetDefinition(fastjet::genkt_algorithm,
624 R0,
625 p)), _p(p), _R0(R0) {
626 if (p < 0) throw Error("GenKT_Axes: Currently only p >=0 is supported.");
627 setNPass(NO_REFINING);
628 }
629
630 /// Short description
631 virtual std::string short_description() const {
632 std::stringstream stream;
633 stream << std::fixed << std::setprecision(2)
634 << "GenKT Axes";
635 return stream.str();
636 };
637
638 /// Long descriptions
639 virtual std::string description() const {
640 std::stringstream stream;
641 stream << std::fixed << std::setprecision(2)
642 << "General KT (p = " << _p << "), R0 = " << _R0;
643 return stream.str();
644 };
645
646 /// For copying purposes
647 virtual GenKT_Axes* create() const {return new GenKT_Axes(*this);}
648
649protected:
650 double _p; ///< genkT power
651 double _R0; ///< jet radius
652};
653
654
655///------------------------------------------------------------------------
656/// \class WTA_GenKT_Axes
657/// \brief Axes from exclusive generalized kT, winner-take-all recombination
658///
659/// Axes from a general KT algorithm with a Winner Take All Recombiner
660/// Requires the power of the KT algorithm to be used and the radius parameter
661///------------------------------------------------------------------------
662class WTA_GenKT_Axes : public ExclusiveJetAxes {
663
664public:
665 /// Constructor
666 WTA_GenKT_Axes(double p, double R0 = fastjet::JetDefinition::max_allowable_R)
667 : ExclusiveJetAxes(JetDefinitionWrapper(fastjet::genkt_algorithm,
668 R0,
669 p,
670 new WinnerTakeAllRecombiner()
671 ).getJetDef()), _p(p), _R0(R0) {
672 if (p < 0) throw Error("WTA_GenKT_Axes: Currently only p >=0 is supported.");
673 setNPass(NO_REFINING);
674 }
675
676 /// Short description
677 virtual std::string short_description() const {
678 std::stringstream stream;
679 stream << std::fixed << std::setprecision(2)
680 << "WTA, GenKT Axes";
681 return stream.str();
682 };
683
684 /// Long descriptions
685 virtual std::string description() const {
686 std::stringstream stream;
687 stream << std::fixed << std::setprecision(2)
688 << "Winner-Take-All General KT (p = " << _p << "), R0 = " << _R0;
689 return stream.str();
690 };
691
692 /// For copying purposes
693 virtual WTA_GenKT_Axes* create() const {return new WTA_GenKT_Axes(*this);}
694
695protected:
696 double _p; ///< genkT power
697 double _R0; ///< jet radius
698};
699
700///------------------------------------------------------------------------
701/// \class GenET_GenKT_Axes
702/// \brief Axes from exclusive kT, generalized Et-scheme recombination
703///
704/// Class using general KT algorithm with a more general recombination scheme
705/// Requires power of KT algorithm, power of recombination weights, and radius parameter
706///------------------------------------------------------------------------
707class GenET_GenKT_Axes : public ExclusiveJetAxes {
708
709public:
710 /// Constructor
711 GenET_GenKT_Axes(double delta, double p, double R0 = fastjet::JetDefinition::max_allowable_R)
712 : ExclusiveJetAxes((JetDefinitionWrapper(fastjet::genkt_algorithm, R0, p, new GeneralEtSchemeRecombiner(delta))).getJetDef() ),
713 _delta(delta), _p(p), _R0(R0) {
714 if (p < 0) throw Error("GenET_GenKT_Axes: Currently only p >=0 is supported.");
715 if (delta <= 0) throw Error("GenET_GenKT_Axes: Currently only delta >0 is supported.");
716 setNPass(NO_REFINING);
717 }
718
719 /// Short description
720 virtual std::string short_description() const {
721 std::stringstream stream;
722 stream << std::fixed << std::setprecision(2)
723 << "GenET, GenKT Axes";
724 return stream.str();
725 };
726
727 /// Long description
728 virtual std::string description() const {
729 std::stringstream stream;
730 stream << std::fixed << std::setprecision(2);
731 // TODO: if _delta is huge, change to "WTA"
732 if (_delta < std::numeric_limits<int>::max()) stream << "General Recombiner (delta = " << _delta << "), " << "General KT (p = " << _p << ") Axes, R0 = " << _R0;
733 else stream << "Winner-Take-All General KT (p = " << _p << "), R0 = " << _R0;
734
735 return stream.str();
736 };
737
738 /// For copying purposes
739 virtual GenET_GenKT_Axes* create() const {return new GenET_GenKT_Axes(*this);}
740
741protected:
742 double _delta; ///< Recombination pT weighting
743 double _p; ///< GenkT power
744 double _R0; ///< jet radius
745};
746
747///------------------------------------------------------------------------
748/// \class OnePass_KT_Axes
749/// \brief Axes from exclusive kT, with one-pass minimization
750///
751/// Onepass minimization from kt axes
752///------------------------------------------------------------------------
753class OnePass_KT_Axes : public KT_Axes {
754public:
755 /// Constructor
756 OnePass_KT_Axes() : KT_Axes() {
757 setNPass(ONE_PASS);
758 }
759
760 /// Short description
761 virtual std::string short_description() const {
762 return "OnePass KT";
763 };
764
765 /// Long description
766 virtual std::string description() const {
767 std::stringstream stream;
768 stream << std::fixed << std::setprecision(2)
769 << "One-Pass Minimization from KT Axes";
770 return stream.str();
771 };
772
773 /// For copying purposes
774 virtual OnePass_KT_Axes* create() const {return new OnePass_KT_Axes(*this);}
775
776
777};
778
779///------------------------------------------------------------------------
780/// \class OnePass_CA_Axes
781/// \brief Axes from exclusive CA, with one-pass minimization
782///
783/// Onepass minimization from CA axes
784///------------------------------------------------------------------------
785class OnePass_CA_Axes : public CA_Axes {
786public:
787 /// Constructor
788 OnePass_CA_Axes() : CA_Axes() {
789 setNPass(ONE_PASS);
790 }
791
792 /// Short description
793 virtual std::string short_description() const {
794 return "OnePass CA";
795 };
796
797 /// Long description
798 virtual std::string description() const {
799 std::stringstream stream;
800 stream << std::fixed << std::setprecision(2)
801 << "One-Pass Minimization from CA Axes";
802 return stream.str();
803 };
804
805 /// For copying purposes
806 virtual OnePass_CA_Axes* create() const {return new OnePass_CA_Axes(*this);}
807
808
809};
810
811///------------------------------------------------------------------------
812/// \class OnePass_AntiKT_Axes
813/// \brief Axes from inclusive anti-kT, with one-pass minimization
814///
815/// Onepass minimization from AntiKT axes, one parameter R0
816///------------------------------------------------------------------------
817class OnePass_AntiKT_Axes : public AntiKT_Axes {
818
819public:
820 /// Constructor
821 OnePass_AntiKT_Axes(double R0) : AntiKT_Axes(R0) {
822 setNPass(ONE_PASS);
823 }
824
825 /// Short Description
826 virtual std::string short_description() const {
827 std::stringstream stream;
828 stream << std::fixed << std::setprecision(2)
829 << "OnePassAKT" << _R0;
830 return stream.str();
831 };
832
833 /// Long description
834 virtual std::string description() const {
835 std::stringstream stream;
836 stream << std::fixed << std::setprecision(2)
837 << "One-Pass Minimization from Anti-KT Axes (R0 = " << _R0 << ")";
838 return stream.str();
839 };
840
841 /// For copying purposes
842 virtual OnePass_AntiKT_Axes* create() const {return new OnePass_AntiKT_Axes(*this);}
843
844};
845
846///------------------------------------------------------------------------
847/// \class OnePass_WTA_KT_Axes
848/// \brief Axes from exclusive kT, winner-take-all recombination, with one-pass minimization
849///
850/// Onepass minimization from winner-take-all kt axes
851///------------------------------------------------------------------------
852class OnePass_WTA_KT_Axes : public WTA_KT_Axes {
853public:
854 /// Constructor
855 OnePass_WTA_KT_Axes() : WTA_KT_Axes() {
856 setNPass(ONE_PASS);
857 }
858
859 /// Short description
860 virtual std::string short_description() const {
861 return "OnePass WTA KT";
862 };
863
864 /// Long description
865 virtual std::string description() const {
866 std::stringstream stream;
867 stream << std::fixed << std::setprecision(2)
868 << "One-Pass Minimization from Winner-Take-All KT Axes";
869 return stream.str();
870 };
871
872 /// For copying purposes
873 virtual OnePass_WTA_KT_Axes* create() const {return new OnePass_WTA_KT_Axes(*this);}
874
875
876};
877
878///------------------------------------------------------------------------
879/// \class OnePass_WTA_CA_Axes
880/// \brief Axes from exclusive CA, winner-take-all recombination, with one-pass minimization
881///
882/// Onepass minimization from winner-take-all CA axes
883///------------------------------------------------------------------------
884class OnePass_WTA_CA_Axes : public WTA_CA_Axes {
885
886public:
887 /// Constructor
888 OnePass_WTA_CA_Axes() : WTA_CA_Axes() {
889 setNPass(ONE_PASS);
890 }
891
892 /// Short description
893 virtual std::string short_description() const {
894 return "OnePass WTA CA";
895 };
896
897 /// Long description
898 virtual std::string description() const {
899 std::stringstream stream;
900 stream << std::fixed << std::setprecision(2)
901 << "One-Pass Minimization from Winner-Take-All CA Axes";
902 return stream.str();
903 };
904
905 /// For copying purposes
906 virtual OnePass_WTA_CA_Axes* create() const {return new OnePass_WTA_CA_Axes(*this);}
907
908};
909
910///------------------------------------------------------------------------
911/// \class OnePass_GenKT_Axes
912/// \brief Axes from exclusive generalized kT with one-pass minimization
913///
914/// Onepass minimization, General KT Axes (standard E-scheme recombination)
915///------------------------------------------------------------------------
916class OnePass_GenKT_Axes : public GenKT_Axes {
917
918public:
919 /// Constructor
920 OnePass_GenKT_Axes(double p, double R0 = fastjet::JetDefinition::max_allowable_R) : GenKT_Axes(p, R0) {
921 setNPass(ONE_PASS);
922 }
923
924 /// Short description
925 virtual std::string short_description() const {
926 return "OnePass GenKT";
927 };
928
929 /// Long description
930 virtual std::string description() const {
931 std::stringstream stream;
932 stream << std::fixed << std::setprecision(2)
933 << "One-Pass Minimization from General KT (p = " << _p << "), R0 = " << _R0;
934 return stream.str();
935 };
936
937 /// For copying purposes
938 virtual OnePass_GenKT_Axes* create() const {return new OnePass_GenKT_Axes(*this);}
939};
940
941///------------------------------------------------------------------------
942/// \class OnePass_WTA_GenKT_Axes
943/// \brief Axes from exclusive generalized kT, winner-take-all recombination, with one-pass minimization
944///
945/// Onepass minimization from winner-take-all, General KT Axes
946///------------------------------------------------------------------------
947class OnePass_WTA_GenKT_Axes : public WTA_GenKT_Axes {
948
949public:
950 /// Constructor
951 OnePass_WTA_GenKT_Axes(double p, double R0 = fastjet::JetDefinition::max_allowable_R) : WTA_GenKT_Axes(p, R0) {
952 setNPass(ONE_PASS);
953 }
954
955 /// Short description
956 virtual std::string short_description() const {
957 return "OnePass WTA GenKT";
958 };
959
960 /// Long description
961 virtual std::string description() const {
962 std::stringstream stream;
963 stream << std::fixed << std::setprecision(2)
964 << "One-Pass Minimization from Winner-Take-All General KT (p = " << _p << "), R0 = " << _R0;
965 return stream.str();
966 };
967
968 /// For copying purposes
969 virtual OnePass_WTA_GenKT_Axes* create() const {return new OnePass_WTA_GenKT_Axes(*this);}
970};
971
972///------------------------------------------------------------------------
973/// \class OnePass_GenET_GenKT_Axes
974/// \brief Axes from exclusive generalized kT, generalized Et-scheme recombination, with one-pass minimization
975///
976/// Onepass minimization from General Recomb, General KT axes
977///------------------------------------------------------------------------
978class OnePass_GenET_GenKT_Axes : public GenET_GenKT_Axes {
979
980public:
981 /// Constructor
982 OnePass_GenET_GenKT_Axes(double delta, double p, double R0 = fastjet::JetDefinition::max_allowable_R) : GenET_GenKT_Axes(delta, p, R0) {
983 setNPass(ONE_PASS);
984 }
985
986 /// Short description
987 virtual std::string short_description() const {
988 return "OnePass GenET, GenKT";
989 };
990
991 /// Long description
992 virtual std::string description() const {
993 std::stringstream stream;
994 stream << std::fixed << std::setprecision(2);
995 if (_delta < std::numeric_limits<int>::max()) stream << "One-Pass Minimization from General Recombiner (delta = "
996 << _delta << "), " << "General KT (p = " << _p << ") Axes, R0 = " << _R0;
997 else stream << "One-Pass Minimization from Winner-Take-All General KT (p = " << _p << "), R0 = " << _R0;
998 return stream.str();
999 };
1000
1001 /// For copying purposes
1002 virtual OnePass_GenET_GenKT_Axes* create() const {return new OnePass_GenET_GenKT_Axes(*this);}
1003};
1004
1005
1006///------------------------------------------------------------------------
1007/// \class Manual_Axes
1008/// \brief Manual axes finding
1009///
1010/// Allows the user to set the axes manually
1011///------------------------------------------------------------------------
1012class Manual_Axes : public AxesDefinition {
1013public:
1014 /// Constructor. Note that _needsManualAxes is set to true.
1015 Manual_Axes() : AxesDefinition() {
1016 setNPass(NO_REFINING);
1017 _needsManualAxes = true;
1018 }
1019
1020 /// This is now a dummy function since this is manual mode
1021 virtual std::vector<fastjet::PseudoJet> get_starting_axes(int,
1022 const std::vector<fastjet::PseudoJet>&,
1023 const MeasureDefinition *) const;
1024
1025
1026 /// Short description
1027 virtual std::string short_description() const {
1028 return "Manual";
1029 };
1030
1031 /// Long description
1032 virtual std::string description() const {
1033 std::stringstream stream;
1034 stream << std::fixed << std::setprecision(2)
1035 << "Manual Axes";
1036 return stream.str();
1037 };
1038
1039 /// For copying purposes
1040 virtual Manual_Axes* create() const {return new Manual_Axes(*this);}
1041
1042
1043};
1044
1045///------------------------------------------------------------------------
1046/// \class OnePass_Manual_Axes
1047/// \brief Manual axes finding, with one-pass minimization
1048///
1049/// One pass minimization from manual starting point
1050///------------------------------------------------------------------------
1051class OnePass_Manual_Axes : public Manual_Axes {
1052public:
1053 /// Constructor. Note that _needsManualAxes is set to true.
1054 OnePass_Manual_Axes() : Manual_Axes() {
1055 setNPass(ONE_PASS);
1056 }
1057
1058 /// Short description
1059 virtual std::string short_description() const {
1060 return "OnePass Manual";
1061 };
1062
1063 /// Long description
1064 virtual std::string description() const {
1065 std::stringstream stream;
1066 stream << std::fixed << std::setprecision(2)
1067 << "One-Pass Minimization from Manual Axes";
1068 return stream.str();
1069 };
1070
1071 // For copying purposes
1072 virtual OnePass_Manual_Axes* create() const {return new OnePass_Manual_Axes(*this);}
1073
1074};
1075
1076///------------------------------------------------------------------------
1077/// \class MultiPass_Axes
1078/// \brief Manual axes finding, with multi-pass (randomized) minimization
1079///
1080/// Multi-pass minimization from kT starting point
1081///------------------------------------------------------------------------
1082class MultiPass_Axes : public KT_Axes {
1083
1084public:
1085
1086 /// Constructor
1087 MultiPass_Axes(unsigned int Npass) : KT_Axes() {
1088 setNPass(Npass);
1089 }
1090
1091 /// Short description
1092 virtual std::string short_description() const {
1093 return "MultiPass";
1094 };
1095
1096 /// Long description
1097 virtual std::string description() const {
1098 std::stringstream stream;
1099 stream << std::fixed << std::setprecision(2)
1100 << "Multi-Pass Axes (Npass = " << _Npass << ")";
1101 return stream.str();
1102 };
1103
1104 /// For copying purposs
1105 virtual MultiPass_Axes* create() const {return new MultiPass_Axes(*this);}
1106
1107};
1108
1109///------------------------------------------------------------------------
1110/// \class MultiPass_Manual_Axes
1111/// \brief Axes finding from exclusive kT, with multi-pass (randomized) minimization
1112///
1113/// multi-pass minimization from kT starting point
1114///------------------------------------------------------------------------
1115class MultiPass_Manual_Axes : public Manual_Axes {
1116
1117public:
1118 /// Constructor
1119 MultiPass_Manual_Axes(unsigned int Npass) : Manual_Axes() {
1120 setNPass(Npass);
1121 }
1122
1123 /// Short Description
1124 virtual std::string short_description() const {
1125 return "MultiPass Manual";
1126 };
1127
1128
1129 /// Long description
1130 virtual std::string description() const {
1131 std::stringstream stream;
1132 stream << std::fixed << std::setprecision(2)
1133 << "Multi-Pass Manual Axes (Npass = " << _Npass << ")";
1134 return stream.str();
1135 };
1136
1137 /// For copying purposes
1138 virtual MultiPass_Manual_Axes* create() const {return new MultiPass_Manual_Axes(*this);}
1139
1140};
1141
1142///------------------------------------------------------------------------
1143/// \class Comb_GenKT_Axes
1144/// \brief Axes from exclusive generalized kT with combinatorial testing
1145///
1146/// Axes from kT algorithm (standard E-scheme recombination)
1147/// Requires nExtra parameter and returns set of N that minimizes N-jettiness
1148/// Note that this method is not guaranteed to find a deeper minimum than GenKT_Axes
1149///------------------------------------------------------------------------
1150class Comb_GenKT_Axes : public ExclusiveCombinatorialJetAxes {
1151public:
1152 /// Constructor
1153 Comb_GenKT_Axes(int nExtra, double p, double R0 = fastjet::JetDefinition::max_allowable_R)
1154 : ExclusiveCombinatorialJetAxes(fastjet::JetDefinition(fastjet::genkt_algorithm, R0, p), nExtra),
1155 _p(p), _R0(R0) {
1156 if (p < 0) throw Error("Comb_GenKT_Axes: Currently only p >=0 is supported.");
1157 setNPass(NO_REFINING);
1158 }
1159
1160 /// Short description
1161 virtual std::string short_description() const {
1162 return "N Choose M GenKT";
1163 };
1164
1165 /// Long description
1166 virtual std::string description() const {
1167 std::stringstream stream;
1168 stream << std::fixed << std::setprecision(2)
1169 << "N Choose M Minimization (nExtra = " << _nExtra << ") from General KT (p = " << _p << "), R0 = " << _R0;
1170 return stream.str();
1171 };
1172
1173 /// For copying purposes
1174 virtual Comb_GenKT_Axes* create() const {return new Comb_GenKT_Axes(*this);}
1175
1176private:
1177 double _nExtra; ///< Number of extra axes
1178 double _p; ///< GenkT power
1179 double _R0; ///< jet radius
1180};
1181
1182
1183
1184///------------------------------------------------------------------------
1185/// \class Comb_WTA_GenKT_Axes
1186/// \brief Axes from exclusive generalized kT, winner-take-all recombination, with combinatorial testing
1187///
1188/// Axes from kT algorithm and winner-take-all recombination
1189/// Requires nExtra parameter and returns set of N that minimizes N-jettiness
1190///------------------------------------------------------------------------
1191class Comb_WTA_GenKT_Axes : public ExclusiveCombinatorialJetAxes {
1192public:
1193 /// Constructor
1194 Comb_WTA_GenKT_Axes(int nExtra, double p, double R0 = fastjet::JetDefinition::max_allowable_R)
1195 : ExclusiveCombinatorialJetAxes((JetDefinitionWrapper(fastjet::genkt_algorithm, R0, p, new WinnerTakeAllRecombiner())).getJetDef(), nExtra),
1196 _p(p), _R0(R0) {
1197 if (p < 0) throw Error("Comb_WTA_GenKT_Axes: Currently only p >=0 is supported.");
1198 setNPass(NO_REFINING);
1199 }
1200
1201 /// Short description
1202 virtual std::string short_description() const {
1203 return "N Choose M WTA GenKT";
1204 };
1205
1206 /// Long description
1207 virtual std::string description() const {
1208 std::stringstream stream;
1209 stream << std::fixed << std::setprecision(2)
1210 << "N Choose M Minimization (nExtra = " << _nExtra << ") from Winner-Take-All General KT (p = " << _p << "), R0 = " << _R0;
1211 return stream.str();
1212 };
1213
1214 /// For copying purposes
1215 virtual Comb_WTA_GenKT_Axes* create() const {return new Comb_WTA_GenKT_Axes(*this);}
1216
1217private:
1218 double _nExtra; ///< Number of extra axes
1219 double _p; ///< GenkT power
1220 double _R0; ///< jet radius
1221};
1222
1223///------------------------------------------------------------------------
1224/// \class Comb_GenET_GenKT_Axes
1225/// \brief Axes from exclusive generalized kT, generalized Et-scheme recombination, with combinatorial testing
1226///
1227/// Axes from kT algorithm and General Et scheme recombination
1228/// Requires nExtra parameter and returns set of N that minimizes N-jettiness
1229///------------------------------------------------------------------------
1230class Comb_GenET_GenKT_Axes : public ExclusiveCombinatorialJetAxes {
1231public:
1232 /// Constructor
1233 Comb_GenET_GenKT_Axes(int nExtra, double delta, double p, double R0 = fastjet::JetDefinition::max_allowable_R)
1234 : ExclusiveCombinatorialJetAxes((JetDefinitionWrapper(fastjet::genkt_algorithm, R0, p, new GeneralEtSchemeRecombiner(delta))).getJetDef(), nExtra),
1235 _delta(delta), _p(p), _R0(R0) {
1236 if (p < 0) throw Error("Comb_GenET_GenKT_Axes: Currently only p >=0 is supported.");
1237 if (delta <= 0) throw Error("Comb_GenET_GenKT_Axes: Currently only delta >=0 is supported.");
1238 setNPass(NO_REFINING);
1239 }
1240
1241 /// Short description
1242 virtual std::string short_description() const {
1243 return "N Choose M GenET GenKT";
1244 };
1245
1246 /// Long description
1247 virtual std::string description() const {
1248 std::stringstream stream;
1249 stream << std::fixed << std::setprecision(2);
1250 if (_delta < std::numeric_limits<int>::max()) stream << "N choose M Minimization (nExtra = " << _nExtra
1251 << ") from General Recombiner (delta = " << _delta << "), " << "General KT (p = " << _p << ") Axes, R0 = " << _R0;
1252 else stream << "N choose M Minimization (nExtra = " << _nExtra << ") from Winner-Take-All General KT (p = " << _p << "), R0 = " << _R0;
1253 return stream.str();
1254 };
1255
1256 /// For copying purposes
1257 virtual Comb_GenET_GenKT_Axes* create() const {return new Comb_GenET_GenKT_Axes(*this);}
1258
1259private:
1260 double _nExtra; ///< Number of extra axes
1261 double _delta; ///< Recombination pT weighting exponent
1262 double _p; ///< GenkT power
1263 double _R0; ///< jet radius
1264};
1265
1266
1267} // namespace contrib
1268
1269FASTJET_END_NAMESPACE
1270
1271#endif // __FASTJET_CONTRIB_NJETTINESS_HH__
1272
Note: See TracBrowser for help on using the repository browser.