[973b92a] | 1 | // Nsubjettiness Package
|
---|
| 2 | // Questions/Comments? jthaler@jthaler.net
|
---|
| 3 | //
|
---|
| 4 | // Copyright (c) 2011-14
|
---|
| 5 | // Jesse Thaler, Ken Van Tilburg, Christopher K. Vermilion, and TJ Wilkason
|
---|
| 6 | //
|
---|
| 7 | // $Id: MeasureFunction.hh 742 2014-08-23 15:43:29Z jthaler $
|
---|
| 8 | //----------------------------------------------------------------------
|
---|
| 9 | // This file is part of FastJet contrib.
|
---|
| 10 | //
|
---|
| 11 | // It is free software; you can redistribute it and/or modify it under
|
---|
| 12 | // the terms of the GNU General Public License as published by the
|
---|
| 13 | // Free Software Foundation; either version 2 of the License, or (at
|
---|
| 14 | // your option) any later version.
|
---|
| 15 | //
|
---|
| 16 | // It is distributed in the hope that it will be useful, but WITHOUT
|
---|
| 17 | // ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
|
---|
| 18 | // or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public
|
---|
| 19 | // License for more details.
|
---|
| 20 | //
|
---|
| 21 | // You should have received a copy of the GNU General Public License
|
---|
| 22 | // along with this code. If not, see <http://www.gnu.org/licenses/>.
|
---|
| 23 | //----------------------------------------------------------------------
|
---|
| 24 |
|
---|
| 25 | #ifndef __FASTJET_CONTRIB_TAUCOMPONENTS_HH__
|
---|
| 26 | #define __FASTJET_CONTRIB_TAUCOMPONENTS_HH__
|
---|
| 27 |
|
---|
| 28 | #include "fastjet/PseudoJet.hh"
|
---|
| 29 | #include "fastjet/ClusterSequence.hh"
|
---|
| 30 | #include "fastjet/WrappedStructure.hh"
|
---|
| 31 |
|
---|
| 32 |
|
---|
| 33 | #include <cmath>
|
---|
| 34 | #include <vector>
|
---|
| 35 | #include <list>
|
---|
| 36 | #include <limits>
|
---|
| 37 |
|
---|
| 38 |
|
---|
| 39 | FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
|
---|
| 40 |
|
---|
| 41 | namespace contrib{
|
---|
| 42 |
|
---|
| 43 | // Classes defined in this file.
|
---|
| 44 | class TauComponents;
|
---|
| 45 | class TauPartition;
|
---|
| 46 | class NjettinessExtras;
|
---|
| 47 |
|
---|
| 48 | ///------------------------------------------------------------------------
|
---|
| 49 | /// \enum TauMode
|
---|
| 50 | /// Specified whether tau value has beam region or denominators
|
---|
| 51 | ///------------------------------------------------------------------------
|
---|
| 52 | enum TauMode {
|
---|
| 53 | UNDEFINED_SHAPE = -1, // Added so that constructor would default to some value
|
---|
| 54 | UNNORMALIZED_JET_SHAPE = 0,
|
---|
| 55 | NORMALIZED_JET_SHAPE = 1,
|
---|
| 56 | UNNORMALIZED_EVENT_SHAPE = 2,
|
---|
| 57 | NORMALIZED_EVENT_SHAPE = 3,
|
---|
| 58 | };
|
---|
| 59 |
|
---|
| 60 | ///////
|
---|
| 61 | //
|
---|
| 62 | // TauComponents
|
---|
| 63 | //
|
---|
| 64 | ///////
|
---|
| 65 |
|
---|
| 66 | ///------------------------------------------------------------------------
|
---|
| 67 | /// \class TauComponents
|
---|
| 68 | /// \brief Output wrapper for supplemental N-(sub)jettiness information
|
---|
| 69 | ///
|
---|
| 70 | /// This class creates a wrapper for the various tau/subtau values calculated in Njettiness. This class allows Njettiness access to these variables
|
---|
| 71 | /// without ever having to do the calculation itself. It takes in subtau numerators and tau denominator from MeasureFunction
|
---|
| 72 | /// and outputs tau numerator, and normalized tau and subtau.
|
---|
| 73 | ///------------------------------------------------------------------------
|
---|
| 74 | class TauComponents {
|
---|
| 75 |
|
---|
| 76 | public:
|
---|
| 77 |
|
---|
| 78 | /// empty constructor necessary to initialize tau_components in Njettiness
|
---|
| 79 | /// later set correctly in Njettiness::getTau function
|
---|
| 80 | TauComponents() {}
|
---|
| 81 |
|
---|
| 82 | /// This constructor takes input vector and double and calculates all necessary tau components
|
---|
| 83 | TauComponents(TauMode tau_mode,
|
---|
| 84 | const std::vector<double> & jet_pieces_numerator,
|
---|
| 85 | double beam_piece_numerator,
|
---|
| 86 | double denominator,
|
---|
| 87 | const std::vector<PseudoJet> & jets,
|
---|
| 88 | const std::vector<PseudoJet> & axes
|
---|
| 89 | );
|
---|
| 90 |
|
---|
| 91 | /// Test for denominator
|
---|
| 92 | bool has_denominator() const;
|
---|
| 93 | /// Test for beam region
|
---|
| 94 | bool has_beam() const;
|
---|
| 95 |
|
---|
| 96 | /// Return tau value
|
---|
| 97 | double tau() const { return _tau; }
|
---|
| 98 | /// Return jet regions
|
---|
| 99 | const std::vector<double>& jet_pieces() const { return _jet_pieces; }
|
---|
| 100 | /// Return beam region
|
---|
| 101 | double beam_piece() const { return _beam_piece; }
|
---|
| 102 |
|
---|
| 103 | /// Return jet regions (no denominator)
|
---|
| 104 | std::vector<double> jet_pieces_numerator() const { return _jet_pieces_numerator; }
|
---|
| 105 | /// Return beam regions (no denominator)
|
---|
| 106 | double beam_piece_numerator() const { return _beam_piece_numerator; }
|
---|
| 107 | /// Return numerator
|
---|
| 108 | double numerator() const { return _numerator; }
|
---|
| 109 | /// Return denominator
|
---|
| 110 | double denominator() const { return _denominator; }
|
---|
| 111 |
|
---|
| 112 | /// Four-vector of total jet (sum of clustered regions)
|
---|
| 113 | PseudoJet total_jet() const { return _total_jet;}
|
---|
| 114 | /// Four-vector of jet regions
|
---|
| 115 | const std::vector<PseudoJet>& jets() const { return _jets;}
|
---|
| 116 | /// Four-vector of axes
|
---|
| 117 | const std::vector<PseudoJet>& axes() const { return _axes;}
|
---|
| 118 |
|
---|
| 119 | class StructureType;
|
---|
| 120 |
|
---|
| 121 | protected:
|
---|
| 122 |
|
---|
| 123 | /// Defines whether there is a beam or denominator
|
---|
| 124 | TauMode _tau_mode;
|
---|
| 125 |
|
---|
| 126 | std::vector<double> _jet_pieces_numerator; ///< Constructor input (jet region numerator)
|
---|
| 127 | double _beam_piece_numerator; ///< Constructor input (beam region numerator)
|
---|
| 128 | double _denominator; ///< Constructor input (denominator)
|
---|
| 129 |
|
---|
| 130 | std::vector<double> _jet_pieces; ///< Derived value (jet regions)
|
---|
| 131 | double _beam_piece; ///< Derived value (beam region)
|
---|
| 132 | double _numerator; ///< Derived value (total numerator)
|
---|
| 133 | double _tau; ///< Derived value (final value)
|
---|
| 134 |
|
---|
| 135 | PseudoJet _total_jet; ///< Total jet four-vector
|
---|
| 136 | std::vector<PseudoJet> _jets; ///< Jet four-vectors
|
---|
| 137 | std::vector<PseudoJet> _axes; ///< AXes four-vectors
|
---|
| 138 |
|
---|
| 139 | };
|
---|
| 140 |
|
---|
| 141 | ///////
|
---|
| 142 | //
|
---|
| 143 | // TauPartition
|
---|
| 144 | //
|
---|
| 145 | ///////
|
---|
| 146 |
|
---|
| 147 | ///------------------------------------------------------------------------
|
---|
| 148 | /// \class TauPartition
|
---|
| 149 | /// \brief Output wrapper for N-(sub)jettiness partitioning information
|
---|
| 150 | ///
|
---|
| 151 | /// Class for storing partitioning information.
|
---|
| 152 | ///------------------------------------------------------------------------
|
---|
| 153 | class TauPartition {
|
---|
| 154 |
|
---|
| 155 | public:
|
---|
| 156 | /// empty constructor
|
---|
| 157 | TauPartition() {}
|
---|
| 158 |
|
---|
| 159 | /// Make partition of size to hold n_jet partitions
|
---|
| 160 | TauPartition(int n_jet) {
|
---|
| 161 | _jets_list.resize(n_jet);
|
---|
| 162 | _jets_partition.resize(n_jet);
|
---|
| 163 | }
|
---|
| 164 |
|
---|
| 165 | /// add a particle to the jet
|
---|
| 166 | void push_back_jet(int jet_num, const PseudoJet& part_to_add, int part_index) {
|
---|
| 167 | _jets_list[jet_num].push_back(part_index);
|
---|
| 168 | _jets_partition[jet_num].push_back(part_to_add);
|
---|
| 169 | }
|
---|
| 170 |
|
---|
| 171 | /// add a particle to the beam
|
---|
| 172 | void push_back_beam(const PseudoJet& part_to_add, int part_index) {
|
---|
| 173 | _beam_list.push_back(part_index);
|
---|
| 174 | _beam_partition.push_back(part_to_add);
|
---|
| 175 | }
|
---|
| 176 |
|
---|
| 177 | /// return jet regions
|
---|
| 178 | PseudoJet jet(int jet_num) const { return join(_jets_partition.at(jet_num)); }
|
---|
| 179 | /// return beam region
|
---|
| 180 | PseudoJet beam() const { return join(_beam_partition);}
|
---|
| 181 |
|
---|
| 182 | /// return jets
|
---|
| 183 | std::vector<PseudoJet> jets() const {
|
---|
| 184 | std::vector<PseudoJet> jets;
|
---|
| 185 | for (unsigned int i = 0; i < _jets_partition.size(); i++) {
|
---|
| 186 | jets.push_back(jet(i));
|
---|
| 187 | }
|
---|
| 188 | return jets;
|
---|
| 189 | }
|
---|
| 190 |
|
---|
| 191 | /// jets in list form
|
---|
| 192 | const std::list<int> & jet_list(int jet_num) const { return _jets_list.at(jet_num);}
|
---|
| 193 | /// beam in list form
|
---|
| 194 | const std::list<int> & beam_list() const { return _beam_list;}
|
---|
| 195 | /// all jets in list form
|
---|
| 196 | const std::vector<std::list<int> > & jets_list() const { return _jets_list;}
|
---|
| 197 |
|
---|
| 198 | private:
|
---|
| 199 |
|
---|
| 200 | std::vector<std::list<int> > _jets_list; ///< jets in list form
|
---|
| 201 | std::list<int> _beam_list; ///< beam in list form
|
---|
| 202 |
|
---|
| 203 | std::vector<std::vector<PseudoJet> > _jets_partition; ///< Partition in jet regions
|
---|
| 204 | std::vector<PseudoJet> _beam_partition; ///< Partition in beam region
|
---|
| 205 |
|
---|
| 206 | };
|
---|
| 207 |
|
---|
| 208 |
|
---|
| 209 | ///////
|
---|
| 210 | //
|
---|
| 211 | // NjettinessExtras
|
---|
| 212 | //
|
---|
| 213 | ///////
|
---|
| 214 |
|
---|
| 215 | ///------------------------------------------------------------------------
|
---|
| 216 | /// \class NjettinessExtras
|
---|
| 217 | /// \brief ClusterSequence add on for N-jettiness information
|
---|
| 218 | ///
|
---|
| 219 | /// This class contains the same information as TauComponents, but adds additional ways of linking up
|
---|
| 220 | /// the jets found in the ClusterSequence::Extras class.
|
---|
| 221 | /// This is done in order to help improve the interface for the main NjettinessPlugin class.
|
---|
| 222 | ///------------------------------------------------------------------------
|
---|
| 223 | class NjettinessExtras : public ClusterSequence::Extras, public TauComponents {
|
---|
| 224 |
|
---|
| 225 | public:
|
---|
| 226 | /// Constructor
|
---|
| 227 | NjettinessExtras(TauComponents tau_components,
|
---|
| 228 | std::vector<int> cluster_hist_indices)
|
---|
| 229 | : TauComponents(tau_components), _cluster_hist_indices(cluster_hist_indices) {}
|
---|
| 230 |
|
---|
| 231 |
|
---|
| 232 |
|
---|
| 233 | /// Ask for tau of the whole event, but by querying a jet
|
---|
| 234 | double tau(const fastjet::PseudoJet& /*jet*/) const {return _tau;}
|
---|
| 235 |
|
---|
| 236 | /// Ask for tau of an individual jet
|
---|
| 237 | double tau_piece(const fastjet::PseudoJet& jet) const {
|
---|
| 238 | if (labelOf(jet) == -1) return std::numeric_limits<double>::quiet_NaN(); // nonsense
|
---|
| 239 | return _jet_pieces[labelOf(jet)];
|
---|
| 240 | }
|
---|
| 241 |
|
---|
| 242 | /// Find axis associated with jet
|
---|
| 243 | fastjet::PseudoJet axis(const fastjet::PseudoJet& jet) const {
|
---|
| 244 | return _axes[labelOf(jet)];
|
---|
| 245 | }
|
---|
| 246 |
|
---|
| 247 | /// Check if extra information is available.
|
---|
| 248 | bool has_njettiness_extras(const fastjet::PseudoJet& jet) const {
|
---|
| 249 | return (labelOf(jet) >= 0);
|
---|
| 250 | }
|
---|
| 251 |
|
---|
| 252 | private:
|
---|
| 253 |
|
---|
| 254 | /// Store cluster history indices to link up with ClusterSequence
|
---|
| 255 | std::vector<int> _cluster_hist_indices;
|
---|
| 256 |
|
---|
| 257 | /// Figure out which jet things belonged to
|
---|
| 258 | int labelOf(const fastjet::PseudoJet& jet) const {
|
---|
| 259 | int thisJet = -1;
|
---|
| 260 | for (unsigned int i = 0; i < _jets.size(); i++) {
|
---|
| 261 | if (_cluster_hist_indices[i] == jet.cluster_hist_index()) {
|
---|
| 262 | thisJet = i;
|
---|
| 263 | break;
|
---|
| 264 | }
|
---|
| 265 | }
|
---|
| 266 | return thisJet;
|
---|
| 267 | }
|
---|
| 268 |
|
---|
| 269 | public:
|
---|
| 270 |
|
---|
| 271 | // These are old methods for gaining this information
|
---|
| 272 | // The recommended interface is given in TauComponents
|
---|
| 273 |
|
---|
| 274 | /// Tau value
|
---|
| 275 | double totalTau() const {return _tau;}
|
---|
| 276 | /// Jet regions
|
---|
| 277 | std::vector<double> subTaus() const {return _jet_pieces;}
|
---|
| 278 |
|
---|
| 279 | /// Tau value
|
---|
| 280 | double totalTau(const fastjet::PseudoJet& /*jet*/) const {
|
---|
| 281 | return _tau;
|
---|
| 282 | }
|
---|
| 283 |
|
---|
| 284 | /// Jet region
|
---|
| 285 | double subTau(const fastjet::PseudoJet& jet) const {
|
---|
| 286 | if (labelOf(jet) == -1) return std::numeric_limits<double>::quiet_NaN(); // nonsense
|
---|
| 287 | return _jet_pieces[labelOf(jet)];
|
---|
| 288 | }
|
---|
| 289 |
|
---|
| 290 | /// beam region
|
---|
| 291 | double beamTau() const {
|
---|
| 292 | return _beam_piece;
|
---|
| 293 | }
|
---|
| 294 |
|
---|
| 295 | };
|
---|
| 296 |
|
---|
| 297 |
|
---|
| 298 | /// Helper function to find out what njettiness_extras are (from jet)
|
---|
| 299 | inline const NjettinessExtras * njettiness_extras(const fastjet::PseudoJet& jet) {
|
---|
| 300 | const ClusterSequence * myCS = jet.associated_cluster_sequence();
|
---|
| 301 | if (myCS == NULL) return NULL;
|
---|
| 302 | const NjettinessExtras* extras = dynamic_cast<const NjettinessExtras*>(myCS->extras());
|
---|
| 303 | return extras;
|
---|
| 304 | }
|
---|
| 305 |
|
---|
| 306 | /// Helper function to find out what njettiness_extras are (from ClusterSequence)
|
---|
| 307 | inline const NjettinessExtras * njettiness_extras(const fastjet::ClusterSequence& myCS) {
|
---|
| 308 | const NjettinessExtras* extras = dynamic_cast<const NjettinessExtras*>(myCS.extras());
|
---|
| 309 | return extras;
|
---|
| 310 | }
|
---|
| 311 |
|
---|
| 312 | ///////
|
---|
| 313 | //
|
---|
| 314 | // TauComponents::StructureType
|
---|
| 315 | //
|
---|
| 316 | ///////
|
---|
| 317 |
|
---|
| 318 |
|
---|
| 319 | ///------------------------------------------------------------------------
|
---|
| 320 | /// \class TauComponents::StructureType
|
---|
| 321 | /// \brief Wrapped structure for jet-based N-(sub)jettiness information
|
---|
| 322 | ///
|
---|
| 323 | /// Small wrapped structure to store tau information
|
---|
| 324 | /// TODO: Can these be auto-joined?
|
---|
| 325 | ///------------------------------------------------------------------------
|
---|
| 326 | class TauComponents::StructureType : public WrappedStructure {
|
---|
| 327 |
|
---|
| 328 | public:
|
---|
| 329 | /// Constructor
|
---|
| 330 | StructureType(const PseudoJet& j) :
|
---|
| 331 | WrappedStructure(j.structure_shared_ptr())
|
---|
| 332 | {}
|
---|
| 333 |
|
---|
| 334 | /// tau associated with jet
|
---|
| 335 | double tau_piece() const { return _tau_piece; }
|
---|
| 336 |
|
---|
| 337 | /// alternative call, though might be confusing
|
---|
| 338 | double tau() const { return _tau_piece; }
|
---|
| 339 |
|
---|
| 340 | private:
|
---|
| 341 | friend class TauComponents;
|
---|
| 342 | double _tau_piece; ///< tau value associated with jet
|
---|
| 343 | };
|
---|
| 344 |
|
---|
| 345 |
|
---|
| 346 |
|
---|
| 347 |
|
---|
| 348 | } //namespace contrib
|
---|
| 349 |
|
---|
| 350 | FASTJET_END_NAMESPACE
|
---|
| 351 |
|
---|
| 352 | #endif // __FASTJET_CONTRIB_TAUCOMPONENTS_HH__
|
---|