Fork me on GitHub

source: git/external/fastjet/contribs/RecursiveTools/RecursiveSoftDrop.hh@ dd64cff

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