[cb80e6f] | 1 | // $Id: RecursiveSoftDrop.hh 1192 2018-10-30 16:08:36Z gsoyez $
|
---|
[b7b836a] | 2 | //
|
---|
| 3 | // Copyright (c) 2014-, Gavin P. Salam, Gregory Soyez, Jesse Thaler,
|
---|
| 4 | // Kevin Zhou, Frederic Dreyer
|
---|
| 5 | //
|
---|
| 6 | //----------------------------------------------------------------------
|
---|
| 7 | // This file is part of FastJet contrib.
|
---|
| 8 | //
|
---|
| 9 | // It is free software; you can redistribute it and/or modify it under
|
---|
| 10 | // the terms of the GNU General Public License as published by the
|
---|
| 11 | // Free Software Foundation; either version 2 of the License, or (at
|
---|
| 12 | // your option) any later version.
|
---|
| 13 | //
|
---|
| 14 | // It is distributed in the hope that it will be useful, but WITHOUT
|
---|
| 15 | // ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
|
---|
| 16 | // or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public
|
---|
| 17 | // License for more details.
|
---|
| 18 | //
|
---|
| 19 | // You should have received a copy of the GNU General Public License
|
---|
| 20 | // along with this code. If not, see <http://www.gnu.org/licenses/>.
|
---|
| 21 | //----------------------------------------------------------------------
|
---|
| 22 |
|
---|
| 23 | #ifndef __RECURSIVESOFTDROP_HH__
|
---|
| 24 | #define __RECURSIVESOFTDROP_HH__
|
---|
| 25 |
|
---|
[cb80e6f] | 26 | // we'll use the native FJ class for reculstering if available
|
---|
| 27 | #if FASTJET_VERSION_NUMBER >= 30100
|
---|
| 28 | #include "fastjet/tools/Recluster.hh"
|
---|
| 29 | #else
|
---|
[b7b836a] | 30 | #include "Recluster.hh"
|
---|
[cb80e6f] | 31 | #endif
|
---|
[b7b836a] | 32 | #include "SoftDrop.hh"
|
---|
| 33 | #include "fastjet/WrappedStructure.hh"
|
---|
| 34 |
|
---|
| 35 | #include <iostream>
|
---|
| 36 | #include <queue>
|
---|
| 37 | #include <vector>
|
---|
| 38 |
|
---|
| 39 | FASTJET_BEGIN_NAMESPACE
|
---|
| 40 |
|
---|
| 41 | namespace contrib{
|
---|
| 42 |
|
---|
| 43 | //------------------------------------------------------------------------
|
---|
| 44 | /// \class RecursiveSoftDrop
|
---|
| 45 | /// An implementation of the RecursiveSoftDrop.
|
---|
| 46 | ///
|
---|
| 47 | /// Recursive Soft Drop will recursively groom a jet, removing
|
---|
| 48 | /// particles that fail the criterion
|
---|
| 49 | /// \f[
|
---|
| 50 | /// z > z_{\rm cut} (\theta/R0)^\beta
|
---|
| 51 | /// \f]
|
---|
| 52 | /// until n subjets have been found.
|
---|
| 53 | ///
|
---|
| 54 | /// Several variants are supported:
|
---|
| 55 | /// - set_fixed_depth_mode() switches to fixed depth on all branches
|
---|
| 56 | /// of the clustering tree
|
---|
| 57 | /// - set_dynamical_R0() switches to dynamical R0 implementation of
|
---|
| 58 | /// RSD
|
---|
| 59 | /// - set_hardest_branch_only() switches to following only the
|
---|
| 60 | /// hardest branch (e.g. for Iterated Soft Drop)
|
---|
| 61 | /// - set_min_deltaR_square(val) sets a minimum angle considered for
|
---|
| 62 | /// substructure (e.g. for Iterated Soft Drop)
|
---|
| 63 | ///
|
---|
| 64 | /// Notes:
|
---|
| 65 | ///
|
---|
| 66 | /// - Even though the calls to "set_tagging_mode()" or
|
---|
| 67 | /// "set_grooming_mode(false)" are allowed, they should not be used
|
---|
| 68 | /// with n=-1, and the default grooming_mode has to remain
|
---|
| 69 | /// untouched (except for beta<0 and finite n).
|
---|
| 70 | ///
|
---|
| 71 | //----------------------------------------------------------------------
|
---|
| 72 | class RecursiveSoftDrop : public SoftDrop {
|
---|
| 73 | public:
|
---|
| 74 | /// Simplified constructor. This takes the value of the "beta"
|
---|
| 75 | /// parameter and the symmetry cut (applied by default on the
|
---|
| 76 | /// scalar_z variable, as for the mMDT). It also takes an optional
|
---|
| 77 | /// subtractor.
|
---|
| 78 | ///
|
---|
| 79 | /// n is the number of times we require the SoftDrop condition to be
|
---|
| 80 | /// satisfied. n=-1 means infinity, i.e. we recurse into the jet
|
---|
| 81 | /// until individual constituents
|
---|
| 82 | ///
|
---|
| 83 | /// If the (optional) pileup subtractor can be supplied, then see
|
---|
| 84 | /// also the documentation for the set_input_jet_is_subtracted() member
|
---|
| 85 | /// function.
|
---|
| 86 | ///
|
---|
| 87 | /// \param beta the value of the beta parameter
|
---|
| 88 | /// \param symmetry_cut the value of the cut on the symmetry measure
|
---|
| 89 | /// \param n the requested number of iterations
|
---|
| 90 | /// \param R0 the angular distance normalisation [1 by default]
|
---|
| 91 | RecursiveSoftDrop(double beta,
|
---|
| 92 | double symmetry_cut,
|
---|
| 93 | int n = -1,
|
---|
| 94 | double R0 = 1,
|
---|
| 95 | const FunctionOfPseudoJet<PseudoJet> * subtractor = 0) :
|
---|
| 96 | SoftDrop(beta, symmetry_cut, R0, subtractor), _n(n) { set_defaults(); }
|
---|
| 97 |
|
---|
| 98 | /// Full constructor, which takes the following parameters:
|
---|
| 99 | ///
|
---|
| 100 | /// \param beta the value of the beta parameter
|
---|
| 101 | /// \param symmetry_cut the value of the cut on the symmetry measure
|
---|
| 102 | /// \param symmetry_measure the choice of measure to use to estimate the symmetry
|
---|
| 103 | /// \param n the requested number of iterations
|
---|
| 104 | /// \param R0 the angular distance normalisation [1 by default]
|
---|
| 105 | /// \param mu_cut the maximal allowed value of mass drop variable mu = m_heavy/m_parent
|
---|
| 106 | /// \param recursion_choice the strategy used to decide which subjet to recurse into
|
---|
| 107 | /// \param subtractor an optional pointer to a pileup subtractor (ignored if zero)
|
---|
| 108 | RecursiveSoftDrop(double beta,
|
---|
| 109 | double symmetry_cut,
|
---|
| 110 | SymmetryMeasure symmetry_measure,
|
---|
| 111 | int n = -1,
|
---|
| 112 | double R0 = 1.0,
|
---|
| 113 | double mu_cut = std::numeric_limits<double>::infinity(),
|
---|
| 114 | RecursionChoice recursion_choice = larger_pt,
|
---|
| 115 | const FunctionOfPseudoJet<PseudoJet> * subtractor = 0) :
|
---|
| 116 | SoftDrop(beta, symmetry_cut, symmetry_measure, R0, mu_cut, recursion_choice, subtractor),
|
---|
| 117 | _n(n) { set_defaults(); }
|
---|
| 118 |
|
---|
| 119 | /// default destructor
|
---|
| 120 | virtual ~RecursiveSoftDrop(){}
|
---|
| 121 |
|
---|
| 122 | //----------------------------------------------------------------------
|
---|
| 123 | // access to class info
|
---|
| 124 | int n() const { return _n; }
|
---|
| 125 |
|
---|
| 126 | //----------------------------------------------------------------------
|
---|
| 127 | // on top of the tweaks that we inherit from SoftDrop (via
|
---|
| 128 | // RecursiveSymmetryBase):
|
---|
| 129 | // - set_verbose_structure()
|
---|
| 130 | // - set_subtractor()
|
---|
| 131 | // - set_input_jet_is_subtracted()
|
---|
| 132 | // we provide several other knobs, given below
|
---|
| 133 |
|
---|
| 134 | /// initialise all the flags below to their default value
|
---|
| 135 | void set_defaults();
|
---|
| 136 |
|
---|
| 137 | /// switch to using the "same depth" variant where instead of
|
---|
| 138 | /// recursing from large to small angles and requiring n SD
|
---|
| 139 | /// conditions to be met (our default), we recurse simultaneously in
|
---|
| 140 | /// all the branches found during the previous iteration, up to a
|
---|
| 141 | /// maximum depth of n.
|
---|
| 142 | /// default: false
|
---|
| 143 | void set_fixed_depth_mode(bool value=true) { _fixed_depth = value; }
|
---|
| 144 | bool fixed_depth_mode() const { return _fixed_depth; }
|
---|
| 145 |
|
---|
| 146 | /// switch to using a dynamical R0 (used for the normalisation of
|
---|
| 147 | /// the symmetry measure) set by the last deltaR at which some
|
---|
| 148 | /// substructure was found.
|
---|
| 149 | /// default: false
|
---|
| 150 | void set_dynamical_R0(bool value=true) { _dynamical_R0 = value; }
|
---|
| 151 | bool use_dynamical_R0() const { return _dynamical_R0; }
|
---|
| 152 |
|
---|
| 153 | /// when finding some substructure, only follow the hardest branch
|
---|
| 154 | /// for the recursion
|
---|
| 155 | /// default: false (i.e. recurse in both branches)
|
---|
| 156 | void set_hardest_branch_only(bool value=true) { _hardest_branch_only = value; }
|
---|
| 157 | bool use_hardest_branch_only() const { return _hardest_branch_only; }
|
---|
| 158 |
|
---|
| 159 | /// set the minimum angle (squared) that we should consider for
|
---|
| 160 | /// substructure
|
---|
| 161 | /// default: -1.0 (i.e. no minimum)
|
---|
| 162 | void set_min_deltaR_squared(double value=-1.0) { _min_dR2 = value; }
|
---|
| 163 | double min_deltaR_squared() const { return _min_dR2; }
|
---|
| 164 |
|
---|
| 165 | /// description of the tool
|
---|
| 166 | virtual std::string description() const;
|
---|
| 167 |
|
---|
| 168 | //----------------------------------------------------------------------
|
---|
| 169 | /// action on a single jet with RecursiveSoftDrop.
|
---|
| 170 | ///
|
---|
| 171 | /// uses "result_fixed_tags" by default (i.e. recurse from R0 to
|
---|
| 172 | /// smaller angles until n SD conditions have been met), or
|
---|
| 173 | /// "result_fixed_depth" where each of the previous SD branches are
|
---|
| 174 | /// recirsed into down to a depth of n.
|
---|
| 175 | virtual PseudoJet result(const PseudoJet &jet) const;
|
---|
| 176 |
|
---|
| 177 | /// this routine applies the Soft Drop criterion recursively on the
|
---|
| 178 | /// CA tree until we find n subjets (or until it converges), and
|
---|
| 179 | /// adds them together into a groomed PseudoJet
|
---|
| 180 | PseudoJet result_fixed_tags(const PseudoJet &jet) const;
|
---|
| 181 |
|
---|
| 182 | /// this routine applies the Soft Drop criterion recursively on the
|
---|
| 183 | /// CA tree, recursing into all the branches found during the previous iteration
|
---|
| 184 | /// until n layers have been found (or until it converges)
|
---|
| 185 | PseudoJet result_fixed_depth(const PseudoJet &jet) const;
|
---|
| 186 |
|
---|
| 187 | protected:
|
---|
| 188 | /// return false if we reached desired layer of grooming _n
|
---|
| 189 | bool continue_grooming(int current_n) const {
|
---|
| 190 | return ((_n < 0) or (current_n < _n));
|
---|
| 191 | }
|
---|
| 192 |
|
---|
| 193 | private:
|
---|
| 194 | int _n; ///< the value of n
|
---|
| 195 |
|
---|
| 196 | // behaviour tweaks
|
---|
| 197 | bool _fixed_depth; ///< look in parallel into each all branches until depth n
|
---|
| 198 | bool _dynamical_R0; ///< when true, use the last deltaR with substructure as D0
|
---|
| 199 | bool _hardest_branch_only; ///< recurse only in the hardest branch
|
---|
| 200 | /// when substructure is found
|
---|
| 201 | double _min_dR2; ///< the min allowed angle to search for substructure
|
---|
| 202 | };
|
---|
| 203 |
|
---|
| 204 | // helper to get the (linear) list of prongs inside a jet resulting
|
---|
| 205 | // from RecursiveSoftDrop. This would avoid having amnually to go
|
---|
| 206 | // through the successive pairwise compositeness
|
---|
| 207 | std::vector<PseudoJet> recursive_soft_drop_prongs(const PseudoJet & rsd_jet);
|
---|
| 208 |
|
---|
| 209 | }
|
---|
| 210 |
|
---|
| 211 | FASTJET_END_NAMESPACE // defined in fastjet/internal/base.hh
|
---|
| 212 | #endif // __RECURSIVESOFTDROP_HH__
|
---|