Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:57:13

0001 //  Nsubjettiness Package
0002 //  Questions/Comments?  jthaler@jthaler.net
0003 //
0004 //  Copyright (c) 2011-14
0005 //  Jesse Thaler, Ken Van Tilburg, Christopher K. Vermilion, and TJ Wilkason
0006 //
0007 //  $Id: AxesDefinition.hh 1129 2018-06-06 12:08:58Z jthaler $
0008 //----------------------------------------------------------------------
0009 // This file is part of FastJet contrib.
0010 //
0011 // It is free software; you can redistribute it and/or modify it under
0012 // the terms of the GNU General Public License as published by the
0013 // Free Software Foundation; either version 2 of the License, or (at
0014 // your option) any later version.
0015 //
0016 // It is distributed in the hope that it will be useful, but WITHOUT
0017 // ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
0018 // or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public
0019 // License for more details.
0020 //
0021 // You should have received a copy of the GNU General Public License
0022 // along with this code. If not, see <http://www.gnu.org/licenses/>.
0023 //----------------------------------------------------------------------
0024 
0025 #ifndef __FASTJET_CONTRIB_AXES_DEFINITION_HH__
0026 #define __FASTJET_CONTRIB_AXES_DEFINITION_HH__
0027 
0028 
0029 #include "MeasureDefinition.hh"
0030 #include "ExtraRecombiners.hh"
0031 
0032 #include "fastjet/PseudoJet.hh"
0033 #include <fastjet/LimitedWarning.hh>
0034 
0035 #include <iomanip>
0036 #include <cmath>
0037 #include <vector>
0038 #include <list>
0039 
0040 FASTJET_BEGIN_NAMESPACE      // defined in fastjet/internal/base.hh
0041 
0042 namespace contrib {
0043    
0044 // The following AxesDefinitions are currently available (and the relevant arguments, if needed)
0045 class KT_Axes;
0046 class CA_Axes;
0047 class AntiKT_Axes;         // (R0)
0048 class WTA_KT_Axes;
0049 class WTA_CA_Axes;
0050 class GenKT_Axes;          // (p, R0 = infinity)
0051 class WTA_GenKT_Axes;      // (p, R0 = infinity)
0052 class GenET_GenKT_Axes;    // (delta, p, R0 = infinity)
0053 class Manual_Axes;
0054    
0055 class OnePass_KT_Axes;
0056 class OnePass_CA_Axes;
0057 class OnePass_AntiKT_Axes;       // (R0)
0058 class OnePass_WTA_KT_Axes;
0059 class OnePass_WTA_CA_Axes;
0060 class OnePass_GenKT_Axes;        // (p, R0 = infinity)
0061 class OnePass_WTA_GenKT_Axes;    // (p, R0 = infinity)
0062 class OnePass_GenET_GenKT_Axes;  // (delta, p, R0 = infinity)
0063 class OnePass_Manual_Axes;
0064    
0065 class MultiPass_Axes;            // (NPass) (currently only defined for KT_Axes)
0066 class MultiPass_Manual_Axes;     // (NPass)
0067 
0068 class Comb_GenKT_Axes;           // (nExtra, p, R0 = infinity)
0069 class Comb_WTA_GenKT_Axes;       // (nExtra, p, R0 = infinity)
0070 class Comb_GenET_GenKT_Axes;     // (nExtra, delta, p, R0 = infinity)
0071 
0072 ///////
0073 //
0074 // AxesDefinition
0075 //
0076 ///////
0077 
0078 ///------------------------------------------------------------------------
0079 /// \class AxesDefinition
0080 /// \brief Base class for axes definitions
0081 ///
0082 /// A generic AxesDefinition first finds a set of seed axes.
0083 /// Then, if desired, uses measure information
0084 /// (from MeasureDefinition) to refine those axes starting from those seed axes.
0085 /// The AxesDefinitions are typically based on sequential jet algorithms.
0086 ///------------------------------------------------------------------------
0087 class AxesDefinition {
0088    
0089 public:
0090    
0091    /// This function should be overloaded in all derived classes, and defines how to find the seed axes.
0092    /// If desired, the measure information (which might be NULL) can be used to test multiple axes choices, but should
0093    /// not be used for iterative refining (since that is the job of MeasureDefinition).
0094    virtual std::vector<fastjet::PseudoJet> get_starting_axes(int n_jets,
0095                                                              const std::vector<fastjet::PseudoJet>& inputs,
0096                                                              const MeasureDefinition * measure) const = 0;
0097    
0098    /// Short description of AxesDefinitions (and any parameters)
0099    virtual std::string short_description() const = 0;
0100 
0101    /// Long description of AxesDefinitions (and any parameters)
0102    virtual std::string description() const = 0;
0103    
0104    /// This has to be defined in all derived classes, and allows these to be copied around.
0105    virtual AxesDefinition* create() const = 0;
0106    
0107 public:
0108    
0109    /// Starting from seeds, refine axes using one or more passes.
0110    /// Note that in order to do >0 passes, we need information from the MeasureDefinition about how to do the appropriate minimization.
0111    std::vector<fastjet::PseudoJet> get_refined_axes(int n_jets,
0112                                                     const std::vector<fastjet::PseudoJet>& inputs,
0113                                                     const std::vector<fastjet::PseudoJet>& seedAxes,
0114                                                     const MeasureDefinition * measure = NULL) const {
0115 
0116       assert(n_jets == (int)seedAxes.size()); //added int casting to get rid of compiler warning
0117             
0118       if (_Npass == 0) {
0119          // no refining, just use seeds
0120          return seedAxes;
0121       } else if (_Npass == 1) {
0122          if (measure == NULL) throw Error("AxesDefinition:  One-pass minimization requires specifying a MeasureDefinition.");
0123          
0124          // do one pass minimum using measure definition
0125          return measure->get_one_pass_axes(n_jets, inputs, seedAxes,_nAttempts,_accuracy);
0126       } else {
0127          if (measure == NULL) throw Error("AxesDefinition:  Multi-pass minimization requires specifying a MeasureDefinition.");
0128          return get_multi_pass_axes(n_jets, inputs, seedAxes, measure);
0129       }
0130    }
0131    
0132    /// Combines get_starting_axes with get_refined_axes.
0133    /// In the Njettiness class, these two steps are done separately in order to store seed axes information.
0134    std::vector<fastjet::PseudoJet> get_axes(int n_jets,
0135                                             const std::vector<fastjet::PseudoJet>& inputs,
0136                                             const MeasureDefinition * measure = NULL) const {
0137       std::vector<fastjet::PseudoJet> seedAxes = get_starting_axes(n_jets, inputs, measure);
0138       return get_refined_axes(n_jets,inputs,seedAxes,measure);
0139    }
0140 
0141    
0142    /// Short-hand for the get_axes function.  Useful when trying to write terse code.
0143    inline std::vector<fastjet::PseudoJet> operator() (int n_jets,
0144                                                const std::vector<fastjet::PseudoJet>& inputs,
0145                                                const MeasureDefinition * measure = NULL) const {
0146       return get_axes(n_jets,inputs,measure);
0147    }
0148    
0149    /// \enum AxesRefiningEnum
0150    /// Defines the cases of zero pass and one pass for convenience
0151    enum AxesRefiningEnum {
0152       UNDEFINED_REFINE = -1, // added to create a default value
0153       NO_REFINING = 0,
0154       ONE_PASS = 1,
0155       MULTI_PASS = 100,
0156    };
0157    
0158    /// A integer that is used externally to decide how to do multi-pass minimization
0159    int nPass() const { return _Npass; }
0160 
0161    /// A flag that indicates whether results are deterministics.
0162    bool givesRandomizedResults() const {
0163       return (_Npass > 1);
0164    }
0165    
0166    /// A flag that indicates whether manual axes are being used.
0167    bool needsManualAxes() const {
0168       return _needsManualAxes; // if there is no starting axes finder
0169    }
0170    
0171    /// Allows user to change number of passes.  Also used internally to set nPass.
0172    /// Can also specify details of one/multi pass minimziation
0173    void setNPass(int nPass,
0174                  int nAttempts = 1000,
0175                  double accuracy  = 0.0001,
0176                  double noise_range = 1.0 // only needed for MultiPass minimization
0177                  )
0178    {
0179       _Npass = nPass;
0180       _nAttempts = nAttempts;
0181       _accuracy = accuracy;
0182       _noise_range = noise_range;
0183       if (nPass < 0) throw Error("AxesDefinition requires a nPass >= 0");
0184    }
0185    
0186    /// Destructor
0187    virtual ~AxesDefinition() {};
0188    
0189 protected:
0190    
0191    /// Default constructor contains no information.  Number of passes has to be set
0192    /// manually by derived classes using setNPass function.
0193    AxesDefinition() : _Npass(UNDEFINED_REFINE),
0194                      _nAttempts(0),
0195                      _accuracy(0.0),
0196                      _noise_range(0.0),
0197                      _needsManualAxes(false) {}
0198 
0199    /// Does multi-pass minimization by randomly jiggling the axes within _noise_range
0200    std::vector<fastjet::PseudoJet> get_multi_pass_axes(int n_jets,
0201                                                        const std::vector<fastjet::PseudoJet>& inputs,
0202                                                        const std::vector<fastjet::PseudoJet>& seedAxes,
0203                                                        const MeasureDefinition* measure) const;
0204    
0205    /// Function to jiggle axes within _noise_range
0206    PseudoJet jiggle(const PseudoJet& axis) const;
0207    
0208    int _Npass;            ///< Number of passes (0 = no refining, 1 = one-pass, >1 multi-pass)
0209    int _nAttempts;        ///< Number of attempts per pass
0210    double _accuracy;      ///< Accuracy goal per pass
0211    double _noise_range;   ///< Noise in rapidity/phi (for multi-pass minimization only)
0212    bool _needsManualAxes; ///< Flag to indicate special case of manual axes
0213 };
0214   
0215 ///------------------------------------------------------------------------
0216 /// \class ExclusiveJetAxes
0217 /// \brief Base class for axes defined from exclusive jet algorithm
0218 ///
0219 /// This class finds axes by clustering particles with an exclusive jet definition.
0220 /// This can be implemented with different jet algorithms.  The user can call this directly
0221 /// using their favorite fastjet::JetDefinition
0222 ///------------------------------------------------------------------------
0223 class ExclusiveJetAxes : public AxesDefinition {
0224    
0225 public:
0226    /// Constructor takes JetDefinition as an argument
0227    ExclusiveJetAxes(fastjet::JetDefinition def)
0228    : AxesDefinition(), _def(def) {
0229       setNPass(NO_REFINING);    // default to no minimization
0230    }
0231    
0232    /// Starting axes obtained by creating a cluster sequenence and running exclusive_jets.
0233    virtual std::vector<fastjet::PseudoJet> get_starting_axes(int n_jets,
0234                                                              const std::vector <fastjet::PseudoJet> & inputs,
0235                                                              const MeasureDefinition * ) const {
0236       fastjet::ClusterSequence jet_clust_seq(inputs, _def);
0237       
0238       std::vector<fastjet::PseudoJet> axes = jet_clust_seq.exclusive_jets_up_to(n_jets);
0239       
0240       if ((int)axes.size() < n_jets) {
0241          _too_few_axes_warning.warn("ExclusiveJetAxes::get_starting_axes:  Fewer than N axes found; results are unpredictable.");
0242          axes.resize(n_jets);  // resize to make sure there are enough axes to not yield an error elsewhere
0243       }
0244       
0245       return axes;
0246    }
0247    
0248    /// Short description
0249    virtual std::string short_description() const { return "ExclAxes";}
0250    /// Long description
0251    virtual std::string description() const { return "ExclAxes: " + _def.description();}
0252    
0253    /// To make it possible to copy around.
0254    virtual ExclusiveJetAxes* create() const {return new ExclusiveJetAxes(*this);}
0255 
0256 private:
0257    fastjet::JetDefinition _def; ///< Jet definition to use.
0258    static LimitedWarning _too_few_axes_warning;
0259 };
0260 
0261 ///------------------------------------------------------------------------
0262 /// \class ExclusiveCombinatorialJetAxes
0263 /// \brief Base class for axes defined from exclusive jet algorithm, checking combinatorial options
0264 ///
0265 /// This class finds axes by clustering particles with an exclusive jet definition.
0266 /// It takes an extra number of jets (specificed by the user via nExtra), and then finds the set of N that minimizes N-jettiness.
0267 /// WARNING: If one wants to be guarenteed that results improve by increasing nExtra, then one should use
0268 /// winner-take-all-style recombination schemes
0269 ///------------------------------------------------------------------------
0270 class ExclusiveCombinatorialJetAxes : public AxesDefinition {
0271    
0272 public:
0273    /// Constructor takes JetDefinition and nExtra as options (nExtra=0 acts the same as ExclusiveJetAxes)
0274    ExclusiveCombinatorialJetAxes(fastjet::JetDefinition def, int nExtra = 0)
0275    : AxesDefinition(), _def(def), _nExtra(nExtra) {
0276       if (nExtra < 0) throw Error("Need nExtra >= 0");
0277       setNPass(NO_REFINING);   // default to no minimization
0278    }
0279    
0280     /// Find n_jets + _nExtra axes, and then choose the n_jets subset with the smallest N-(sub)jettiness value.
0281     virtual std::vector<fastjet::PseudoJet> get_starting_axes(int n_jets, 
0282                                                            const std::vector<fastjet::PseudoJet> & inputs,
0283                                                            const MeasureDefinition *measure) const {
0284       int starting_number = n_jets + _nExtra;
0285       fastjet::ClusterSequence jet_clust_seq(inputs, _def);
0286       std::vector<fastjet::PseudoJet> starting_axes = jet_clust_seq.exclusive_jets_up_to(starting_number);
0287        
0288        if ((int)starting_axes.size() < n_jets) {
0289           _too_few_axes_warning.warn("ExclusiveCombinatorialJetAxes::get_starting_axes:  Fewer than N + nExtra axes found; results are unpredictable.");
0290           starting_axes.resize(n_jets);  // resize to make sure there are enough axes to not yield an error elsewhere
0291        }
0292        
0293       std::vector<fastjet::PseudoJet> final_axes;
0294 
0295       // check so that no computation time is wasted if there are no extra axes
0296       if (_nExtra == 0) final_axes = starting_axes;
0297 
0298       else {
0299  
0300         // define string of 1's based on number of desired jets
0301         std::string bitmask(n_jets, 1);
0302         // expand the array size to the total number of jets with extra 0's at the end, makes string easy to permute
0303         bitmask.resize(starting_number, 0); 
0304  
0305         double min_tau = std::numeric_limits<double>::max();
0306         std::vector<fastjet::PseudoJet> temp_axes;
0307          
0308         do {
0309 
0310           temp_axes.clear();
0311  
0312           // only take an axis if it is listed as true (1) in the string
0313           for (int i = 0; i < (int)starting_axes.size(); ++i) {
0314             if (bitmask[i]) temp_axes.push_back(starting_axes[i]);
0315           }
0316  
0317           double temp_tau = measure->result(inputs, temp_axes);
0318           if (temp_tau < min_tau) {
0319             min_tau = temp_tau;
0320             final_axes = temp_axes;
0321           }
0322  
0323           // permutes string of 1's and 0's according to next lexicographic ordering and returns true
0324           // continues to loop through all possible lexicographic orderings
0325           // returns false and breaks the loop when there are no more possible orderings
0326         } while (std::prev_permutation(bitmask.begin(), bitmask.end()));
0327       }
0328 
0329       return final_axes;
0330     }
0331 
0332    /// Short description
0333    virtual std::string short_description() const { return "ExclCombAxes";}
0334    /// Long description
0335    virtual std::string description() const { return "ExclCombAxes: " + _def.description();}
0336    /// To make it possible to copy around.
0337    virtual ExclusiveCombinatorialJetAxes* create() const {return new ExclusiveCombinatorialJetAxes(*this);}
0338 
0339 private:
0340    fastjet::JetDefinition _def;   ///< Jet definition to use
0341    int _nExtra;                   ///< Extra axes to find
0342    static LimitedWarning _too_few_axes_warning;
0343 };
0344    
0345 ///------------------------------------------------------------------------
0346 /// \class HardestJetAxes
0347 /// \brief Base class for axes defined from an inclusive jet algorithm
0348 ///
0349 /// This class finds axes by running an inclusive algorithm and then finding the n hardest jets.
0350 /// This can be implemented with different jet algorithms, and can be called by the user.
0351 ///------------------------------------------------------------------------
0352 class HardestJetAxes : public AxesDefinition {
0353 public:
0354    /// Constructor takes JetDefinition
0355    HardestJetAxes(fastjet::JetDefinition def)
0356    : AxesDefinition(), _def(def) {
0357       setNPass(NO_REFINING);    // default to no minimization
0358    }
0359    
0360    /// Finds seed axes by running a ClusterSequence, running inclusive_jets, and finding the N hardest
0361    virtual std::vector<fastjet::PseudoJet> get_starting_axes(int n_jets,
0362                                                              const std::vector <fastjet::PseudoJet> & inputs,
0363                                                              const MeasureDefinition * ) const {
0364       fastjet::ClusterSequence jet_clust_seq(inputs, _def);
0365       std::vector<fastjet::PseudoJet> axes = sorted_by_pt(jet_clust_seq.inclusive_jets());
0366       
0367       if ((int)axes.size() < n_jets) {
0368          _too_few_axes_warning.warn("HardestJetAxes::get_starting_axes:  Fewer than N axes found; results are unpredictable.");
0369       }
0370       
0371       axes.resize(n_jets);  // only keep n hardest
0372       return axes;
0373    }
0374    
0375    /// Short description
0376    virtual std::string short_description() const { return "HardAxes";}
0377    /// Long description
0378    virtual std::string description() const { return "HardAxes: " + _def.description();}
0379    /// To make it possible to copy around.
0380    virtual HardestJetAxes* create() const {return new HardestJetAxes(*this);}
0381    
0382 private:
0383    fastjet::JetDefinition _def;  ///< Jet Definition to use.
0384    
0385    static LimitedWarning _too_few_axes_warning;
0386 
0387 };
0388    
0389 ///------------------------------------------------------------------------
0390 /// \class KT_Axes
0391 /// \brief Axes from exclusive kT
0392 ///
0393 /// Axes from kT algorithm with E_scheme recombination.
0394 ///------------------------------------------------------------------------
0395 class KT_Axes : public ExclusiveJetAxes {
0396 public:
0397    /// Constructor
0398    KT_Axes()
0399    : ExclusiveJetAxes(fastjet::JetDefinition(fastjet::kt_algorithm,
0400                                              fastjet::JetDefinition::max_allowable_R, //maximum jet radius constant
0401                                              fastjet::E_scheme,
0402                                              fastjet::Best)
0403                       ) {
0404       setNPass(NO_REFINING);
0405    }
0406 
0407    /// Short description
0408    virtual std::string short_description() const {
0409       return "KT";
0410    };
0411    
0412    /// Long description
0413    virtual std::string description() const {
0414       std::stringstream stream;
0415       stream << std::fixed << std::setprecision(2)
0416       << "KT Axes";
0417       return stream.str();
0418    };
0419    
0420    /// For copying purposes
0421    virtual KT_Axes* create() const {return new KT_Axes(*this);}
0422 
0423 };
0424 
0425 ///------------------------------------------------------------------------
0426 /// \class CA_Axes
0427 /// \brief Axes from exclusive CA
0428 ///
0429 /// Axes from CA algorithm with E_scheme recombination.
0430 ///------------------------------------------------------------------------
0431 class CA_Axes : public ExclusiveJetAxes {
0432 public:
0433    /// Constructor
0434    CA_Axes()
0435    : ExclusiveJetAxes(fastjet::JetDefinition(fastjet::cambridge_algorithm,
0436                                              fastjet::JetDefinition::max_allowable_R, //maximum jet radius constant
0437                                              fastjet::E_scheme,
0438                                              fastjet::Best)
0439                       ) {
0440       setNPass(NO_REFINING);
0441    }
0442 
0443    /// Short description
0444    virtual std::string short_description() const {
0445       return "CA";
0446    };
0447    
0448    /// Long description
0449    virtual std::string description() const {
0450       std::stringstream stream;
0451       stream << std::fixed << std::setprecision(2)
0452       << "CA Axes";
0453       return stream.str();
0454    };
0455    
0456    /// For copying purposes
0457    virtual CA_Axes* create() const {return new CA_Axes(*this);}
0458    
0459 };
0460 
0461    
0462 ///------------------------------------------------------------------------
0463 /// \class AntiKT_Axes
0464 /// \brief Axes from inclusive anti-kT
0465 ///
0466 /// Axes from anti-kT algorithm and E_scheme.
0467 /// The one parameter R0 is subjet radius
0468 ///------------------------------------------------------------------------
0469 class AntiKT_Axes : public HardestJetAxes {
0470 
0471 public:
0472    /// Constructor.  Takes jet radius as argument
0473    AntiKT_Axes(double R0)
0474    : HardestJetAxes(fastjet::JetDefinition(fastjet::antikt_algorithm,
0475                                            R0,
0476                                            fastjet::E_scheme,
0477                                            fastjet::Best)
0478                     ), _R0(R0) {
0479       setNPass(NO_REFINING);
0480    }
0481 
0482    /// Short description
0483    virtual std::string short_description() const {
0484       std::stringstream stream;
0485       stream << std::fixed << std::setprecision(2)
0486       << "AKT" << _R0;
0487       return stream.str();
0488    };
0489    
0490    /// Long description
0491    virtual std::string description() const {
0492       std::stringstream stream;
0493       stream << std::fixed << std::setprecision(2)
0494       << "Anti-KT Axes (R0 = " << _R0 << ")";
0495       return stream.str();
0496    };
0497    
0498    /// For copying purposes
0499    virtual AntiKT_Axes* create() const {return new AntiKT_Axes(*this);}
0500    
0501 protected:
0502    double _R0;  ///<  AKT jet radius
0503 
0504 };
0505 
0506 ///------------------------------------------------------------------------
0507 /// \class JetDefinitionWrapper
0508 /// \brief Wrapper for jet definitions (for memory management)
0509 ///
0510 /// This class is used by all AxesDefinition with a manual recombiner to
0511 /// ensure that the delete_recombiner_when_unused function is always called
0512 ///------------------------------------------------------------------------
0513 class JetDefinitionWrapper {
0514 
0515 public: 
0516    
0517    /// Default Constructor
0518    JetDefinitionWrapper(JetAlgorithm jet_algorithm_in, double R_in, double xtra_param_in, const JetDefinition::Recombiner *recombiner) {
0519       jet_def = fastjet::JetDefinition(jet_algorithm_in, R_in, xtra_param_in);
0520       jet_def.set_recombiner(recombiner);
0521       jet_def.delete_recombiner_when_unused();        // added to prevent memory leaks
0522    }
0523 
0524    /// Additional constructor so that build-in FastJet algorithms can also be called
0525    JetDefinitionWrapper(JetAlgorithm jet_algorithm_in, double R_in, const JetDefinition::Recombiner *recombiner, fastjet::Strategy strategy_in) {
0526       jet_def = fastjet::JetDefinition(jet_algorithm_in, R_in, recombiner, strategy_in);
0527       jet_def.delete_recombiner_when_unused();
0528    }
0529 
0530    /// Return jet definition
0531    JetDefinition getJetDef() {
0532       return jet_def;
0533    }
0534 
0535 private:
0536    JetDefinition jet_def;  ///< my jet definition
0537 };
0538 
0539 ///------------------------------------------------------------------------
0540 /// \class WTA_KT_Axes
0541 /// \brief Axes from exclusive kT, winner-take-all recombination
0542 ///
0543 /// Axes from kT algorithm and winner-take-all recombination
0544 ///------------------------------------------------------------------------
0545 class WTA_KT_Axes : public ExclusiveJetAxes {
0546 public:
0547    /// Constructor
0548    WTA_KT_Axes()
0549    : ExclusiveJetAxes(JetDefinitionWrapper(fastjet::kt_algorithm,
0550                                           fastjet::JetDefinition::max_allowable_R, //maximum jet radius constant
0551                                           new WinnerTakeAllRecombiner(), // Needs to be explicitly declared (this will be deleted by JetDefinitionWrapper)
0552                                           fastjet::Best).getJetDef()
0553                       ) {
0554       setNPass(NO_REFINING);
0555     }
0556 
0557    /// Short description
0558    virtual std::string short_description() const {
0559       return "WTA KT";
0560    };
0561    
0562    /// Long description
0563    virtual std::string description() const {
0564       std::stringstream stream;
0565       stream << std::fixed << std::setprecision(2)
0566       << "Winner-Take-All KT Axes";
0567       return stream.str();
0568    };
0569    
0570    /// For copying purposes
0571    virtual WTA_KT_Axes* create() const {return new WTA_KT_Axes(*this);}
0572 
0573 };
0574    
0575 ///------------------------------------------------------------------------
0576 /// \class WTA_CA_Axes
0577 /// \brief Axes from exclusive CA, winner-take-all recombination
0578 ///
0579 /// Axes from CA algorithm and winner-take-all recombination
0580 ///------------------------------------------------------------------------
0581 class WTA_CA_Axes : public ExclusiveJetAxes {
0582 public:
0583    /// Constructor
0584    WTA_CA_Axes()
0585    : ExclusiveJetAxes(JetDefinitionWrapper(fastjet::cambridge_algorithm,
0586                                              fastjet::JetDefinition::max_allowable_R, //maximum jet radius constant
0587                                              new WinnerTakeAllRecombiner(), // Needs to be explicitly declared (this will be deleted by JetDefinitionWrapper)
0588                                              fastjet::Best).getJetDef()) {
0589     setNPass(NO_REFINING);
0590   }
0591 
0592    /// Short description
0593    virtual std::string short_description() const {
0594       return "WTA CA";
0595    };
0596    
0597    /// Long descriptions
0598    virtual std::string description() const {
0599       std::stringstream stream;
0600       stream << std::fixed << std::setprecision(2)
0601       << "Winner-Take-All CA Axes";
0602       return stream.str();
0603    };
0604    
0605    /// For copying purposes
0606    virtual WTA_CA_Axes* create() const {return new WTA_CA_Axes(*this);}
0607   
0608 };
0609 
0610 
0611 ///------------------------------------------------------------------------
0612 /// \class GenKT_Axes
0613 /// \brief Axes from exclusive generalized kT
0614 ///
0615 /// Axes from a general KT algorithm (standard E-scheme recombination)
0616 /// Requires the power of the KT algorithm to be used and the radius parameter
0617 ///------------------------------------------------------------------------
0618 class GenKT_Axes : public ExclusiveJetAxes {
0619    
0620 public:
0621    /// Constructor
0622    GenKT_Axes(double p, double R0 = fastjet::JetDefinition::max_allowable_R)
0623    : ExclusiveJetAxes(fastjet::JetDefinition(fastjet::genkt_algorithm,
0624                                              R0,
0625                                              p)), _p(p), _R0(R0) {
0626       if (p < 0) throw Error("GenKT_Axes:  Currently only p >=0 is supported.");
0627       setNPass(NO_REFINING);
0628    }
0629    
0630    /// Short description
0631    virtual std::string short_description() const {
0632       std::stringstream stream;
0633       stream << std::fixed << std::setprecision(2)
0634       << "GenKT Axes";
0635       return stream.str();
0636    };
0637    
0638    /// Long descriptions
0639    virtual std::string description() const {
0640       std::stringstream stream;
0641       stream << std::fixed << std::setprecision(2)
0642       << "General KT (p = " << _p << "), R0 = " << _R0;
0643       return stream.str();
0644    };
0645    
0646    /// For copying purposes
0647    virtual GenKT_Axes* create() const {return new GenKT_Axes(*this);}
0648    
0649 protected:
0650    double _p;   ///< genkT power
0651    double _R0;  ///< jet radius
0652 };
0653    
0654    
0655 ///------------------------------------------------------------------------
0656 /// \class WTA_GenKT_Axes
0657 /// \brief Axes from exclusive generalized kT, winner-take-all recombination
0658 ///
0659 /// Axes from a general KT algorithm with a Winner Take All Recombiner
0660 /// Requires the power of the KT algorithm to be used and the radius parameter
0661 ///------------------------------------------------------------------------
0662 class WTA_GenKT_Axes : public ExclusiveJetAxes {
0663 
0664 public:
0665    /// Constructor
0666    WTA_GenKT_Axes(double p, double R0 = fastjet::JetDefinition::max_allowable_R)
0667    : ExclusiveJetAxes(JetDefinitionWrapper(fastjet::genkt_algorithm,
0668                                             R0,
0669                                             p,
0670                                             new WinnerTakeAllRecombiner()
0671                                             ).getJetDef()), _p(p), _R0(R0) {
0672       if (p < 0) throw Error("WTA_GenKT_Axes:  Currently only p >=0 is supported.");
0673       setNPass(NO_REFINING);
0674    }
0675 
0676    /// Short description
0677    virtual std::string short_description() const {
0678       std::stringstream stream;
0679       stream << std::fixed << std::setprecision(2)
0680       << "WTA, GenKT Axes";
0681       return stream.str();
0682    };
0683    
0684    /// Long descriptions
0685    virtual std::string description() const {
0686       std::stringstream stream;
0687       stream << std::fixed << std::setprecision(2)
0688       << "Winner-Take-All General KT (p = " << _p << "), R0 = " << _R0;
0689       return stream.str();
0690    };
0691    
0692    /// For copying purposes
0693    virtual WTA_GenKT_Axes* create() const {return new WTA_GenKT_Axes(*this);}
0694    
0695 protected:
0696    double _p;   ///< genkT power
0697    double _R0;  ///< jet radius
0698 };
0699    
0700 ///------------------------------------------------------------------------
0701 /// \class GenET_GenKT_Axes
0702 /// \brief Axes from exclusive kT, generalized Et-scheme recombination
0703 ///
0704 /// Class using general KT algorithm with a more general recombination scheme
0705 /// Requires power of KT algorithm, power of recombination weights, and radius parameter
0706 ///------------------------------------------------------------------------
0707 class GenET_GenKT_Axes : public ExclusiveJetAxes {
0708 
0709 public:
0710    /// Constructor
0711    GenET_GenKT_Axes(double delta, double p, double R0 = fastjet::JetDefinition::max_allowable_R)
0712    : ExclusiveJetAxes((JetDefinitionWrapper(fastjet::genkt_algorithm, R0, p, new GeneralEtSchemeRecombiner(delta))).getJetDef() ),
0713     _delta(delta), _p(p), _R0(R0) {
0714        if (p < 0) throw Error("GenET_GenKT_Axes:  Currently only p >=0 is supported.");
0715        if (delta <= 0) throw Error("GenET_GenKT_Axes:  Currently only delta >0 is supported.");
0716        setNPass(NO_REFINING);
0717    }
0718 
0719    /// Short description
0720    virtual std::string short_description() const {
0721       std::stringstream stream;
0722       stream << std::fixed << std::setprecision(2)
0723       << "GenET, GenKT Axes";
0724       return stream.str();
0725    };
0726    
0727    /// Long description
0728    virtual std::string description() const {
0729       std::stringstream stream;
0730       stream << std::fixed << std::setprecision(2);
0731       // TODO: if _delta is huge, change to "WTA"
0732       if (_delta < std::numeric_limits<int>::max()) stream << "General Recombiner (delta = " << _delta << "), " << "General KT (p = " << _p << ") Axes, R0 = " << _R0;
0733       else stream << "Winner-Take-All General KT (p = " << _p << "), R0 = " << _R0;
0734 
0735       return stream.str();
0736    };
0737    
0738    /// For copying purposes
0739    virtual GenET_GenKT_Axes* create() const {return new GenET_GenKT_Axes(*this);}
0740    
0741 protected:
0742    double _delta; ///< Recombination pT weighting
0743    double _p;     ///< GenkT power
0744    double _R0;    ///< jet radius
0745 };
0746 
0747 ///------------------------------------------------------------------------
0748 /// \class OnePass_KT_Axes
0749 /// \brief Axes from exclusive kT, with one-pass minimization
0750 ///
0751 /// Onepass minimization from kt axes
0752 ///------------------------------------------------------------------------
0753 class OnePass_KT_Axes : public KT_Axes {
0754 public:
0755    /// Constructor
0756    OnePass_KT_Axes() : KT_Axes() {
0757       setNPass(ONE_PASS);
0758    }
0759    
0760    /// Short description
0761    virtual std::string short_description() const {
0762       return "OnePass KT";
0763    };
0764    
0765    /// Long description
0766    virtual std::string description() const {
0767       std::stringstream stream;
0768       stream << std::fixed << std::setprecision(2)
0769       << "One-Pass Minimization from KT Axes";
0770       return stream.str();
0771    };
0772    
0773    /// For copying purposes
0774    virtual OnePass_KT_Axes* create() const {return new OnePass_KT_Axes(*this);}
0775    
0776 
0777 };
0778 
0779 ///------------------------------------------------------------------------
0780 /// \class OnePass_CA_Axes
0781 /// \brief Axes from exclusive CA, with one-pass minimization
0782 ///
0783 /// Onepass minimization from CA axes
0784 ///------------------------------------------------------------------------
0785 class OnePass_CA_Axes : public CA_Axes {
0786 public:
0787    /// Constructor
0788    OnePass_CA_Axes() : CA_Axes() {
0789       setNPass(ONE_PASS);
0790    }
0791 
0792    /// Short description
0793    virtual std::string short_description() const {
0794       return "OnePass CA";
0795    };
0796    
0797    /// Long description
0798    virtual std::string description() const {
0799       std::stringstream stream;
0800       stream << std::fixed << std::setprecision(2)
0801       << "One-Pass Minimization from CA Axes";
0802       return stream.str();
0803    };
0804    
0805    /// For copying purposes
0806    virtual OnePass_CA_Axes* create() const {return new OnePass_CA_Axes(*this);}
0807 
0808 
0809 };
0810    
0811 ///------------------------------------------------------------------------
0812 /// \class OnePass_AntiKT_Axes
0813 /// \brief Axes from inclusive anti-kT, with one-pass minimization
0814 ///
0815 /// Onepass minimization from AntiKT axes, one parameter R0
0816 ///------------------------------------------------------------------------
0817 class OnePass_AntiKT_Axes : public AntiKT_Axes {
0818 
0819 public:
0820    /// Constructor
0821    OnePass_AntiKT_Axes(double R0) : AntiKT_Axes(R0) {
0822       setNPass(ONE_PASS);
0823    }
0824    
0825    /// Short Description
0826    virtual std::string short_description() const {
0827       std::stringstream stream;
0828       stream << std::fixed << std::setprecision(2)
0829       << "OnePassAKT" << _R0;
0830       return stream.str();
0831    };
0832    
0833    /// Long description
0834    virtual std::string description() const {
0835       std::stringstream stream;
0836       stream << std::fixed << std::setprecision(2)
0837       << "One-Pass Minimization from Anti-KT Axes (R0 = " << _R0 << ")";
0838       return stream.str();
0839    };
0840    
0841    /// For copying purposes
0842    virtual OnePass_AntiKT_Axes* create() const {return new OnePass_AntiKT_Axes(*this);}
0843 
0844 };
0845 
0846 ///------------------------------------------------------------------------
0847 /// \class OnePass_WTA_KT_Axes
0848 /// \brief Axes from exclusive kT, winner-take-all recombination, with one-pass minimization
0849 ///
0850 /// Onepass minimization from winner-take-all kt axes
0851 ///------------------------------------------------------------------------
0852 class OnePass_WTA_KT_Axes : public WTA_KT_Axes {
0853 public:
0854    /// Constructor
0855    OnePass_WTA_KT_Axes() : WTA_KT_Axes() {
0856       setNPass(ONE_PASS);
0857    }
0858    
0859    /// Short description
0860    virtual std::string short_description() const {
0861       return "OnePass WTA KT";
0862    };
0863    
0864    /// Long description
0865    virtual std::string description() const {
0866       std::stringstream stream;
0867       stream << std::fixed << std::setprecision(2)
0868       << "One-Pass Minimization from Winner-Take-All KT Axes";
0869       return stream.str();
0870    };
0871    
0872    /// For copying purposes
0873    virtual OnePass_WTA_KT_Axes* create() const {return new OnePass_WTA_KT_Axes(*this);}
0874    
0875 
0876 };
0877 
0878 ///------------------------------------------------------------------------
0879 /// \class OnePass_WTA_CA_Axes
0880 /// \brief Axes from exclusive CA, winner-take-all recombination, with one-pass minimization
0881 ///
0882 /// Onepass minimization from winner-take-all CA axes
0883 ///------------------------------------------------------------------------
0884 class OnePass_WTA_CA_Axes : public WTA_CA_Axes {
0885    
0886 public:
0887    /// Constructor
0888    OnePass_WTA_CA_Axes() : WTA_CA_Axes() {
0889       setNPass(ONE_PASS);
0890    }
0891 
0892    /// Short description
0893    virtual std::string short_description() const {
0894       return "OnePass WTA CA";
0895    };
0896    
0897    /// Long description
0898    virtual std::string description() const {
0899       std::stringstream stream;
0900       stream << std::fixed << std::setprecision(2)
0901       << "One-Pass Minimization from Winner-Take-All CA Axes";
0902       return stream.str();
0903    };
0904    
0905    /// For copying purposes
0906    virtual OnePass_WTA_CA_Axes* create() const {return new OnePass_WTA_CA_Axes(*this);}
0907    
0908 };
0909 
0910 ///------------------------------------------------------------------------
0911 /// \class OnePass_GenKT_Axes
0912 /// \brief Axes from exclusive generalized kT with one-pass minimization
0913 ///
0914 /// Onepass minimization, General KT Axes (standard E-scheme recombination)
0915 ///------------------------------------------------------------------------
0916 class OnePass_GenKT_Axes : public GenKT_Axes {
0917    
0918 public:
0919    /// Constructor
0920    OnePass_GenKT_Axes(double p, double R0 = fastjet::JetDefinition::max_allowable_R) : GenKT_Axes(p, R0) {
0921       setNPass(ONE_PASS);
0922    }
0923    
0924    /// Short description
0925    virtual std::string short_description() const {
0926       return "OnePass GenKT";
0927    };
0928    
0929    /// Long description
0930    virtual std::string description() const {
0931       std::stringstream stream;
0932       stream << std::fixed << std::setprecision(2)
0933       << "One-Pass Minimization from General KT (p = " << _p << "), R0 = " << _R0;
0934       return stream.str();
0935    };
0936    
0937    /// For copying purposes
0938    virtual OnePass_GenKT_Axes* create() const {return new OnePass_GenKT_Axes(*this);}
0939 };
0940    
0941 ///------------------------------------------------------------------------
0942 /// \class OnePass_WTA_GenKT_Axes
0943 /// \brief Axes from exclusive generalized kT, winner-take-all recombination, with one-pass minimization
0944 ///
0945 /// Onepass minimization from winner-take-all, General KT Axes
0946 ///------------------------------------------------------------------------
0947 class OnePass_WTA_GenKT_Axes : public WTA_GenKT_Axes {
0948    
0949 public:
0950    /// Constructor
0951    OnePass_WTA_GenKT_Axes(double p, double R0 = fastjet::JetDefinition::max_allowable_R) : WTA_GenKT_Axes(p, R0) {
0952       setNPass(ONE_PASS);
0953    }
0954 
0955    /// Short description
0956    virtual std::string short_description() const {
0957       return "OnePass WTA GenKT";
0958    };
0959    
0960    /// Long description
0961    virtual std::string description() const {
0962       std::stringstream stream;
0963       stream << std::fixed << std::setprecision(2)
0964       << "One-Pass Minimization from Winner-Take-All General KT (p = " << _p << "), R0 = " << _R0;
0965       return stream.str();
0966    };
0967    
0968    /// For copying purposes
0969    virtual OnePass_WTA_GenKT_Axes* create() const {return new OnePass_WTA_GenKT_Axes(*this);}
0970 };
0971 
0972 ///------------------------------------------------------------------------
0973 /// \class OnePass_GenET_GenKT_Axes
0974 /// \brief Axes from exclusive generalized kT, generalized Et-scheme recombination, with one-pass minimization
0975 ///
0976 /// Onepass minimization from General Recomb, General KT axes
0977 ///------------------------------------------------------------------------
0978 class OnePass_GenET_GenKT_Axes : public GenET_GenKT_Axes {
0979    
0980 public:
0981    /// Constructor
0982    OnePass_GenET_GenKT_Axes(double delta, double p, double R0 = fastjet::JetDefinition::max_allowable_R) : GenET_GenKT_Axes(delta, p, R0) {
0983       setNPass(ONE_PASS);
0984    }
0985 
0986    /// Short description
0987    virtual std::string short_description() const {
0988       return "OnePass GenET, GenKT";
0989    };
0990    
0991    /// Long description
0992    virtual std::string description() const {
0993       std::stringstream stream;
0994       stream << std::fixed << std::setprecision(2);
0995       if (_delta < std::numeric_limits<int>::max()) stream << "One-Pass Minimization from General Recombiner (delta = " 
0996         << _delta << "), " << "General KT (p = " << _p << ") Axes, R0 = " << _R0;
0997       else stream << "One-Pass Minimization from Winner-Take-All General KT (p = " << _p << "), R0 = " << _R0;
0998       return stream.str();
0999    };
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 ///------------------------------------------------------------------------
1012 class Manual_Axes : public AxesDefinition {
1013 public:
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 ///------------------------------------------------------------------------
1051 class OnePass_Manual_Axes : public Manual_Axes {
1052 public:
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 ///------------------------------------------------------------------------
1082 class MultiPass_Axes : public KT_Axes {
1083 
1084 public:
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 ///------------------------------------------------------------------------
1115 class MultiPass_Manual_Axes : public Manual_Axes {
1116 
1117 public:
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 ///------------------------------------------------------------------------
1150 class Comb_GenKT_Axes : public ExclusiveCombinatorialJetAxes {
1151 public:
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    
1176 private:
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 ///------------------------------------------------------------------------
1191 class Comb_WTA_GenKT_Axes : public ExclusiveCombinatorialJetAxes {
1192 public:
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 
1217 private:
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 ///------------------------------------------------------------------------
1230 class Comb_GenET_GenKT_Axes : public ExclusiveCombinatorialJetAxes {
1231 public:
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 
1259 private:
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 
1269 FASTJET_END_NAMESPACE
1270 
1271 #endif  // __FASTJET_CONTRIB_NJETTINESS_HH__
1272