[1f1f858] | 1 | ------------------------------------------------------------------------
|
---|
| 2 | RecursiveTools FastJet contrib
|
---|
| 3 | ------------------------------------------------------------------------
|
---|
| 4 |
|
---|
| 5 | The RecursiveTools FastJet contrib aims to provide a common contrib
|
---|
| 6 | for a number of tools that involve recursive reclustering/declustering
|
---|
| 7 | of a jet for tagging or grooming purposes.
|
---|
| 8 |
|
---|
| 9 | Currently it contains:
|
---|
| 10 |
|
---|
| 11 | - ModifiedMassDropTagger
|
---|
| 12 | This corresponds to arXiv:1307.0007 by Mrinal Dasgupta, Alessandro
|
---|
| 13 | Fregoso, Simone Marzani and Gavin P. Salam
|
---|
| 14 |
|
---|
| 15 | - SoftDrop
|
---|
| 16 | This corresponds to arXiv:1402.2657 by Andrew J. Larkoski, Simone
|
---|
| 17 | Marzani, Gregory Soyez, Jesse Thaler
|
---|
| 18 |
|
---|
[b7b836a] | 19 | - RecursiveSoftDrop
|
---|
| 20 | - BottomUpSoftDrop
|
---|
| 21 | This corresponds to arXiv:1804.03657 by Frederic Dreyer, Lina
|
---|
| 22 | Necib, Gregory Soyez and Jesse Thaler
|
---|
| 23 |
|
---|
| 24 | - IteratedSoftDrop
|
---|
| 25 | This corresponds to arXiv:1704.06266 by Christopher Frye, Andrew J.
|
---|
| 26 | Larkoski, Jesse Thaler, Kevin Zhou
|
---|
| 27 |
|
---|
[1f1f858] | 28 | - Recluster
|
---|
| 29 | A generic tool to recluster a given jet into subjets
|
---|
[b7b836a] | 30 | Note: a Recluster class is available natively in FastJet since v3.1.
|
---|
| 31 | Users are therefore encouraged to use the FastJet version
|
---|
| 32 | rather than this one which is mostly provided for
|
---|
| 33 | compatibility of this contrib with older versions of FastJet.
|
---|
[1f1f858] | 34 |
|
---|
| 35 | The interface for these tools is described in more detail below, with
|
---|
| 36 | all of the available options documented in the header files.
|
---|
| 37 |
|
---|
| 38 | One note about nomenclature. A groomer is a procedure that takes a
|
---|
| 39 | PseudoJet and always returns another (non-zero) PseudoJet. A tagger is
|
---|
| 40 | a procedure that takes a PseudoJet, and either returns another PseudoJet
|
---|
| 41 | (i.e. tags it) or returns an empty PseudoJet (i.e. doesn't tag it).
|
---|
| 42 |
|
---|
| 43 | ------------------------------------------------------------------------
|
---|
| 44 | ModifiedMassDropTagger
|
---|
| 45 | ------------------------------------------------------------------------
|
---|
| 46 |
|
---|
| 47 | The Modified Mass Drop Tagger (mMDT) recursively declusters a jet,
|
---|
| 48 | following the largest pT subjet until a pair of subjets is found that
|
---|
| 49 | satisfy the symmetry condition on the energy sharing
|
---|
| 50 |
|
---|
| 51 | z > z_cut
|
---|
| 52 |
|
---|
| 53 | where z_cut is a predetermined value. By default, z is calculated as
|
---|
| 54 | the scalar pT fraction of the softest subjet. Note that larger values
|
---|
| 55 | of z_cut correspond to a more restrictive tagging criteria.
|
---|
| 56 |
|
---|
| 57 | By default, mMDT will first recluster the jet using the CA clustering
|
---|
| 58 | algorithm, which means that mMDT can be called on any jet, regardless
|
---|
| 59 | of the original jet finding measure.
|
---|
| 60 |
|
---|
| 61 | A default mMDT can be created via
|
---|
| 62 |
|
---|
| 63 | double z_cut = 0.10;
|
---|
| 64 | ModifiedMassDropTagger mMDT(z_cut);
|
---|
| 65 |
|
---|
| 66 | More options are available in the full constructor. To apply mMDT,
|
---|
| 67 | one simply calls it on the jet of interest.
|
---|
| 68 |
|
---|
| 69 | PseudoJet tagged_jet = mMDT(original_jet);
|
---|
| 70 |
|
---|
| 71 | Note that mMDT is a tagger, such that tagged_jet will only be non-zero
|
---|
| 72 | if the symmetry cut z > z_cut is satisfied by some branching of the
|
---|
| 73 | clustering tree.
|
---|
| 74 |
|
---|
| 75 | To gain additional information about the mMDT procedure, one can use
|
---|
| 76 |
|
---|
| 77 | tagged_jet.structure_of<ModifiedMassDropTagger>()
|
---|
| 78 |
|
---|
| 79 | which gives access to information about the delta_R between the tagged
|
---|
| 80 | subjets, their z value, etc.
|
---|
| 81 |
|
---|
| 82 | ------------------------------------------------------------------------
|
---|
| 83 | SoftDrop
|
---|
| 84 | ------------------------------------------------------------------------
|
---|
| 85 |
|
---|
| 86 | The SoftDrop procedure is very similar to mMDT, albeit with a
|
---|
[b7b836a] | 87 | generalised symmetry condition:
|
---|
[1f1f858] | 88 |
|
---|
| 89 | z > z_cut * (R / R0)^beta
|
---|
| 90 |
|
---|
| 91 | Note that larger z_cut and smaller beta correspond to more aggressive
|
---|
| 92 | grooming of the jet.
|
---|
| 93 |
|
---|
| 94 | SoftDrop is intended to be used as a groomer (instead of as a tagger),
|
---|
| 95 | such that if the symmetry condition fails throughout the whole
|
---|
| 96 | clustering tree, SoftDrop will still return a single particle in the
|
---|
| 97 | end. Apart from the tagger/groomer distinction, SoftDrop with beta=0 is
|
---|
| 98 | the same as mMDT.
|
---|
| 99 |
|
---|
| 100 | A default SoftDrop groomer can be created via:
|
---|
| 101 |
|
---|
| 102 | double beta = 2.0;
|
---|
[cb80e6f] | 103 | double z_cut = 0.10;
|
---|
[1f1f858] | 104 | double R0 = 1.0; // this is the default value
|
---|
[cb80e6f] | 105 | SoftDrop sd(beta,z_cut,R0);
|
---|
[1f1f858] | 106 |
|
---|
| 107 | and acts on a desired jet as
|
---|
| 108 |
|
---|
| 109 | PseudoJet groomed_jet = sd(original_jet);
|
---|
| 110 |
|
---|
| 111 | and additional information can be obtained via
|
---|
| 112 |
|
---|
| 113 | groomed_jet.structure_of<SoftDrop>()
|
---|
| 114 |
|
---|
| 115 | SoftDrop is typically called with beta > 0, though beta < 0 is still a
|
---|
| 116 | viable option. Because beta < 0 is infrared-collinear unsafe in
|
---|
| 117 | grooming mode, one probably wants to switch to tagging mode for negative
|
---|
| 118 | beta, via set_tagging_mode().
|
---|
| 119 |
|
---|
[b7b836a] | 120 | ------------------------------------------------------------------------
|
---|
| 121 | RecursiveSoftDrop
|
---|
| 122 | ------------------------------------------------------------------------
|
---|
| 123 |
|
---|
| 124 | The RecursiveSoftDrop procedure applies the Soft Drop procedure N times
|
---|
| 125 | in a jet in order to find up to N+1 prongs. N=0 makes no modification
|
---|
| 126 | to the jet, and N=1 is equivalent to the original SoftDrop.
|
---|
| 127 |
|
---|
| 128 | Once one has more than one prong, one has to decide which will be
|
---|
| 129 | declustered next. At each step of the declustering procedure, one
|
---|
| 130 | undoes the clustering which has the largest declustering angle
|
---|
| 131 | (amongst all the branches that are searched for substructure). [see
|
---|
| 132 | "set_fixed_depth" below for an alternative]
|
---|
| 133 |
|
---|
| 134 | Compared to SoftDrop, RecursiveSoftDrop takes an extra argument N
|
---|
| 135 | specifying the number of times the SoftDrop procedure is recursively
|
---|
| 136 | applied. Negative N means that the procedure is applied until no
|
---|
| 137 | further substructure is found (i.e. corresponds to taking N=infinity).
|
---|
| 138 |
|
---|
| 139 | double beta = 2.0;
|
---|
[cb80e6f] | 140 | double z_cut = 0.10;
|
---|
[b7b836a] | 141 | double R0 = 1.0; // this is the default value
|
---|
| 142 | int N = -1;
|
---|
[cb80e6f] | 143 | RecursiveSoftDrop rsd(beta, z_cut, N, R0);
|
---|
[b7b836a] | 144 |
|
---|
| 145 | One then acts on a jet as
|
---|
| 146 |
|
---|
| 147 | PseudoJet groomed_jet = rsd(jet)
|
---|
| 148 |
|
---|
| 149 | and get additional information via
|
---|
| 150 |
|
---|
| 151 | groomed_jet.structure_of<RecursiveSoftDrop>()
|
---|
| 152 |
|
---|
| 153 | ------------------------------------------------------------------------
|
---|
| 154 | IteratedSoftDrop
|
---|
| 155 | ------------------------------------------------------------------------
|
---|
| 156 |
|
---|
| 157 | Iterated Soft Drop (ISD) is a repeated variant of SoftDrop. After
|
---|
| 158 | performing the Soft Drop procedure once, it logs the groomed symmetry
|
---|
| 159 | factor, then recursively performs Soft Drop again on the harder
|
---|
| 160 | branch. This procedure is repeated down to an (optional) angular cut
|
---|
| 161 | theta_cut, yielding a set of symmetry factors from which observables
|
---|
| 162 | can be built.
|
---|
| 163 |
|
---|
| 164 | An IteratedSoftDrop tool can be created as follows:
|
---|
| 165 |
|
---|
| 166 | double beta = -1.0;
|
---|
| 167 | double z_cut = 0.005;
|
---|
| 168 | double theta_cut = 0.0;
|
---|
| 169 | double R0 = 0.5; // characteristic radius of jet algorithm
|
---|
| 170 | IteratedSoftDrop isd(beta, z_cut, double theta_cut, R0);
|
---|
| 171 |
|
---|
| 172 | By default, ISD applied on a jet gives a result of type
|
---|
| 173 | IteratedSoftDropInfo that can then be probed to obtain physical
|
---|
| 174 | observables
|
---|
| 175 |
|
---|
| 176 | IteratedSoftDropInfo isd_info = isd(jet);
|
---|
| 177 |
|
---|
| 178 | unsigned int multiplicity = isd_info.multiplicity();
|
---|
| 179 | double kappa = 1.0; // changes angular scale of ISD angularity
|
---|
| 180 | double isd_width = isd_info.angularity(kappa);
|
---|
| 181 | vector<pair<double,double> > zg_thetags = isd_info.all_zg_thetag();
|
---|
| 182 | vector<pair<double,double> > zg_thetags = isd_info();
|
---|
| 183 | for (unsigned int i=0; i< isd_info.size(); ++i){
|
---|
| 184 | cout << "(zg, theta_g)_" << i << " = "
|
---|
| 185 | << isd_info[i].first << " " << isd_info[i].second << endl;
|
---|
| 186 | }
|
---|
| 187 |
|
---|
| 188 | Alternatively, one can directly get the multiplicity, angularity, and
|
---|
| 189 | (zg,thetag) pairs from the IteratedSoftDrop class, at the expense of
|
---|
| 190 | re-running the declustering procedure:
|
---|
| 191 |
|
---|
| 192 | unsigned int multiplicity = isd.multiplicity(jet);
|
---|
| 193 | double isd_width = isd.angularity(jet, 1.0);
|
---|
| 194 | vector<pair<double,double> > zg_thetags = isd.all_zg_thetag(jet);
|
---|
| 195 |
|
---|
| 196 |
|
---|
| 197 | Note: the iterative declustering procedure is the same as what one
|
---|
| 198 | would obtain with RecursiveSoftDrop with an (optional) angular cut
|
---|
| 199 | and recursing only in the hardest branch [see the "Changing
|
---|
| 200 | behaviour" section below for details], except that it returns some
|
---|
| 201 | information about the jet instead of a modified jet as RSD does.
|
---|
| 202 |
|
---|
| 203 |
|
---|
| 204 | ------------------------------------------------------------------------
|
---|
| 205 | BottomUpSoftDrop
|
---|
| 206 | ------------------------------------------------------------------------
|
---|
| 207 |
|
---|
| 208 | This is a bottom-up version of the RecursiveSoftDrop procedure, in a
|
---|
| 209 | similar way as Pruning can be seen as a bottom-up version of Trimming.
|
---|
| 210 |
|
---|
| 211 | In practice, the jet is reclustered and at each step of the clustering
|
---|
| 212 | one checks the SoftDrop condition
|
---|
| 213 |
|
---|
| 214 | z > z_cut * (R / R0)^beta
|
---|
| 215 |
|
---|
| 216 | If the condition is met, the pair is recombined. If the condition is
|
---|
| 217 | not met, only the hardest of the two objects is kept for further
|
---|
| 218 | clustering and the softest is rejected.
|
---|
| 219 |
|
---|
[cb80e6f] | 220 | BottomUpSoftDrop takes the same arguments as SoftDrop, and a groomer
|
---|
| 221 | can be created with:
|
---|
| 222 |
|
---|
| 223 | double beta = 2.0;
|
---|
| 224 | double z_cut = 0.10;
|
---|
| 225 | double R0 = 1.0; // this is the default value
|
---|
| 226 | BottomUpSoftDrop busd(beta,z_cut,R0);
|
---|
| 227 |
|
---|
| 228 | One then acts on a jet as
|
---|
| 229 |
|
---|
| 230 | PseudoJet groomed_jet = busd(jet)
|
---|
| 231 |
|
---|
[1f1f858] | 232 | ------------------------------------------------------------------------
|
---|
| 233 | Recluster
|
---|
| 234 | ------------------------------------------------------------------------
|
---|
| 235 |
|
---|
[b7b836a] | 236 | *** NOTE: this is provided only for backwards compatibility ***
|
---|
| 237 | *** with FastJet <3.1. For FastJet >=3.1, the native ***
|
---|
| 238 | *** fastjet::Recluster is used instead ***
|
---|
| 239 |
|
---|
[1f1f858] | 240 | The Recluster class allows the constituents of a jet to be reclustered
|
---|
| 241 | with a different recursive clustering algorithm. This is used
|
---|
[b7b836a] | 242 | internally in the mMDT/SoftDrop/RecursiveSoftDrop/IteratedSoftDrop
|
---|
| 243 | code in order to recluster the jet using the CA algorithm. This is
|
---|
| 244 | achieved via
|
---|
[1f1f858] | 245 |
|
---|
| 246 | Recluster ca_reclusterer(cambridge_algorithm,
|
---|
| 247 | JetDefinition::max_allowable_R);
|
---|
| 248 | PseudoJet reclustered_jet = ca_reclusterer(original_jet);
|
---|
| 249 |
|
---|
| 250 | Note that reclustered_jet creates a new ClusterSequence that knows to
|
---|
| 251 | delete_self_when_unused.
|
---|
| 252 |
|
---|
| 253 | ------------------------------------------------------------------------
|
---|
| 254 | Changing behaviour
|
---|
| 255 | ------------------------------------------------------------------------
|
---|
| 256 |
|
---|
[b7b836a] | 257 | The behaviour of the all the tools provided here
|
---|
| 258 | (ModifiedMassDropTagger, SoftDrop, RecursiveSoftDrop and
|
---|
| 259 | IteratedSoftDrop) can be tweaked using the following options:
|
---|
[1f1f858] | 260 |
|
---|
[b7b836a] | 261 | SymmetryMeasure = {scalar_z, vector_z, y, theta_E, cos_theta_E}
|
---|
| 262 | [constructor argument]
|
---|
[1f1f858] | 263 | : The definition of the energy sharing between subjets, with 0
|
---|
[b7b836a] | 264 | corresponding to the most asymmetric.
|
---|
| 265 | . scalar_z = min(pt1,pt2)/(pt1+pt2) [default]
|
---|
| 266 | . vector_z = min(pt1,pt2)/pt_{1+2}
|
---|
| 267 | . y = min(pt1^2,pt2^2)/m_{12}^2 (original y from MDT)
|
---|
| 268 | . theta_E = min(E1,E2)/(E1+E2),
|
---|
| 269 | with angular measure theta_{12}^2
|
---|
| 270 | . cos_theta_E = min(E1,E2)/(E1+E2),
|
---|
| 271 | with angular measure 2[1-cos(theta_{12})]
|
---|
| 272 | The last two variants are meant for use in e+e- collisions,
|
---|
| 273 | together with the "larger_E" recursion choice (see below)
|
---|
| 274 |
|
---|
| 275 | RecursionChoice = {larger_pt, larger_mt, larger_m, larger_E}
|
---|
| 276 | [constructor argument]
|
---|
[1f1f858] | 277 | : The path to recurse through the tree after the symmetry condition
|
---|
[b7b836a] | 278 | fails. Options refer to transverse momentum (pt), transverse mass
|
---|
| 279 | (mt=sqrt(pt^2+m^2), mass (m) or energy (E). the latter is meant
|
---|
| 280 | for use in e+e- collisions
|
---|
[1f1f858] | 281 |
|
---|
| 282 | mu_cut [constructor argument]
|
---|
| 283 | : An optional mass drop condition
|
---|
| 284 |
|
---|
[b7b836a] | 285 | set_subtractor(subtractor*) [or subtractor as a constructor argument]
|
---|
[1f1f858] | 286 | : provide a subtractor. When a subtractor is supplied, the
|
---|
| 287 | kinematic constraints are applied on subtracted 4-vectors. In
|
---|
| 288 | this case, the result of the ModifiedMassDropTagger/SoftDrop is a
|
---|
| 289 | subtracted PseudoJet, and it is assumed that
|
---|
| 290 | ModifiedMassDropTagger/SoftDrop is applied to an unsubtracted jet.
|
---|
| 291 | The latter default can be changed by calling
|
---|
| 292 | set_input_jet_is_subtracted().
|
---|
| 293 |
|
---|
| 294 | set_reclustering(bool, Recluster*)
|
---|
| 295 | : An optional setting to recluster a jet with a different jet
|
---|
| 296 | recursive jet algorithm. The code is only designed to give sensible
|
---|
| 297 | results with the CA algorithm, but other reclustering algorithm
|
---|
| 298 | (especially kT) may be appropriate in certain contexts.
|
---|
| 299 | Use at your own risk.
|
---|
| 300 |
|
---|
| 301 | set_grooming_mode()/set_tagging_mode()
|
---|
| 302 | : In grooming mode, the algorithm will return a single particle if the
|
---|
| 303 | symmetry condition fails for the whole tree. In tagging mode, the
|
---|
| 304 | algorithm will return an zero PseudoJet if no symmetry conditions
|
---|
| 305 | passes. Note that ModifiedMassDropTagger defaults to tagging mode
|
---|
| 306 | and SoftDrop defaults to grooming mode.
|
---|
| 307 |
|
---|
[b7b836a] | 308 | set_verbose_structure(bool)
|
---|
| 309 | : when set to true, additional information will be stored in the jet
|
---|
| 310 | structure. This includes in particular values of symmetry,
|
---|
| 311 | delta_R, and mu of dropped branches
|
---|
| 312 |
|
---|
| 313 | For the specific case of RecursiveSoftDrop, additional tweaking is
|
---|
| 314 | possible via the following methods
|
---|
| 315 |
|
---|
| 316 | set_fixed_depth_mode(bool)
|
---|
| 317 | : when this is true, RSD will recurse (N times) into all the
|
---|
| 318 | branches found during the previous iteration [instead of recursing
|
---|
| 319 | through the largest declustering angle until N prongs have been
|
---|
| 320 | found]. This yields at most 2^N prong. For infinite N, the two
|
---|
| 321 | options are equivalent.
|
---|
| 322 |
|
---|
| 323 | set_dynamical_R0(bool)
|
---|
| 324 | : By default the angles in the SD condition are normalised to the
|
---|
| 325 | parameter R0. With "dynamical R0", RSD will dynamically adjust R0
|
---|
| 326 | to be the angle between the two prongs found during the previous
|
---|
| 327 | iteration.
|
---|
| 328 |
|
---|
| 329 | set_hardest_branch_only(bool)
|
---|
| 330 | : When substructure is found, only recurse into the hardest of the
|
---|
| 331 | two branches for further substructure search. This uses the class
|
---|
| 332 | RecursionChoice.
|
---|
| 333 |
|
---|
| 334 | set_min_deltaR_squared(double):
|
---|
| 335 | : set a minimal angle (squared) at which we stop the declustering
|
---|
| 336 | procedure. This cut is ineffective for negative values of the
|
---|
| 337 | argument.
|
---|
| 338 |
|
---|
[1f1f858] | 339 | ------------------------------------------------------------------------
|
---|
| 340 | Technical Details
|
---|
| 341 | ------------------------------------------------------------------------
|
---|
| 342 |
|
---|
| 343 | Both ModifiedMassDropTagger and SoftDrop inherit from
|
---|
| 344 | RecursiveSymmetryCutBase, which provides a common codebase for recursive
|
---|
| 345 | declustering of a jet with a symmetry cut condition. A generic
|
---|
| 346 | RecursiveSymmetryCutBase depends on the following (virtual) functions
|
---|
| 347 | (see header file for exact full specs, including constness):
|
---|
| 348 |
|
---|
| 349 | double symmetry_cut_fn(PseudoJet &, PseudoJet &)
|
---|
| 350 | : The function that defines the symmetry cut. This is what actually
|
---|
| 351 | defines different recursive declustering schemes, and all classes
|
---|
| 352 | that inherit from RecursiveSymmetryCutBase must define this
|
---|
| 353 | function.
|
---|
| 354 |
|
---|
| 355 | string symmetry_cut_description()
|
---|
| 356 | : the string description of the symmetry cut.
|
---|
| 357 |
|
---|
| 358 | ------------------------------------------------------------------------
|
---|