|
||||
File indexing completed on 2025-01-30 10:06:42
0001 // $Id: GenericSubtractor.hh 861 2015-09-21 10:35:22Z gsoyez $ 0002 // 0003 // Copyright (c) 2012-, Matteo Cacciari, Jihun Kim, Gavin P. Salam and Gregory Soyez 0004 // 0005 //---------------------------------------------------------------------- 0006 // This file is part of the GenericSubtractor package of FastJet Contrib. 0007 // 0008 // It is free software; you can redistribute it and/or modify it under 0009 // the terms of the GNU General Public License as published by the 0010 // Free Software Foundation; either version 2 of the License, or (at 0011 // your option) any later version. 0012 // 0013 // It is distributed in the hope that it will be useful, but WITHOUT 0014 // ANY WARRANTY; without even the implied warranty of MERCHANTABILITY 0015 // or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public 0016 // License for more details. 0017 // 0018 // You should have received a copy of the GNU General Public License 0019 // along with this code. If not, see <http://www.gnu.org/licenses/>. 0020 //---------------------------------------------------------------------- 0021 0022 #ifndef __FASTJET_CONTRIB_GENERIC_SUBTRACTOR_HH__ 0023 #define __FASTJET_CONTRIB_GENERIC_SUBTRACTOR_HH__ 0024 0025 #include "fastjet/PseudoJet.hh" 0026 #include "fastjet/FunctionOfPseudoJet.hh" 0027 #include "fastjet/LimitedWarning.hh" 0028 #include "fastjet/tools/BackgroundEstimatorBase.hh" 0029 0030 #include "ShapeWithComponents.hh" 0031 0032 FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh 0033 0034 namespace contrib{ 0035 0036 // forward declaration of assistional class(es) defined below. 0037 class GenericSubtractorInfo; 0038 0039 0040 //------------------------------------------------------------------------ 0041 /// \class GenericSubtractor 0042 /// 0043 /// A class to perform subtraction of background (eg pileup or 0044 /// underlying event) from a jet shape. 0045 /// 0046 /// This class is a tool that allows one to subtract pileup from jet shapes 0047 /// (i.e. a FunctionOfPseudoJet<double>). It implements the method 0048 /// described in: 0049 /// Gregory Soyez, Gavin P. Salam, Jihun Kim, Souvik Dutta and Matteo Cacciari, 0050 /// Phys. Rev. Lett. 110, 162001 (2013) [arXiv:1211.2811] 0051 /// 0052 /// The basic usage of this class is as follows: 0053 /// 0054 /// GenericSubtractor gensub(& some_background_estimator); 0055 /// FunctionOfPseudoJet<double> shape; 0056 /// double subtracted_shape_value = gensub(shape, jet); 0057 /// 0058 /// By default, this only subtracts the transverse momentum density 0059 /// rho (obtained from the background estimator provided to the 0060 /// class). 0061 /// 0062 /// Two options allow also for the inclusion of the "particle mass" 0063 /// contribution to the background density (the "rho_m" term): one can 0064 /// either instruct the GenericSubtractor class to compute rho_m from 0065 /// the same background estimator using 0066 /// 0067 /// gensub.set_common_bge_for_rho_and_rhom(); (1) 0068 /// 0069 /// or explicitly construct gensub with two background estimators 0070 /// 0071 /// GenericSubtractor gensub(& background_estimator_for_rho, 0072 /// & background_estimator_for_rhom); (2) 0073 /// 0074 /// Note that since FastJet 3.1, the option (1) will work with any 0075 /// background estimator that internally computes rho_m (and has that 0076 /// computation enabled). For FastJet 3.0, the first option is only 0077 /// available for JetMedianBackgroundEstimator. 0078 /// 0079 /// For the option (2), 'gensub' will obtain rho_m using 0080 /// 0081 /// background_estimator_for_rhom->rho() [NOT rho_m()!] 0082 /// 0083 /// unless one calls gensub.set_use_bge_rhom_rhom() (it requires 0084 /// FJ>=3.1) in which case it will be obtained using 0085 /// 0086 /// background_estimator_for_rhom->rhom() 0087 /// 0088 /// The right choice between these two procedures for option (2) 0089 /// depends on the details of 'background_estimator_for_rhom'. 0090 /// 0091 /// 0092 /// If the background transverse momentum density rho (and optionally the 0093 /// density of sqrt(pt^2+m^2) - pt, denoted by rhom) are known, the 0094 /// subtractor can be constructed directly as 0095 /// 0096 /// GenericSubtractor gensub(rho,rhom); 0097 /// 0098 /// Extra information can be retrieved using 0099 /// 0100 /// GenericSubtractorInfo info; 0101 /// double subtracted_shape_value = gensub(shape, jet, info); 0102 /// 0103 class GenericSubtractor{ 0104 public: 0105 /// default constructor 0106 /// leaves the object unusable 0107 GenericSubtractor() : 0108 _bge_rho(0), _bge_rhom(0), _jet_pt_fraction(0.01), 0109 _common_bge(false), _rhom_from_bge_rhom(false), 0110 _rho(_invalid_rho), _externally_supplied_rho_rhom(false) {} 0111 0112 /// Constructor that takes a pointer to a background estimator for 0113 /// rho and optionally a pointer to a background estimator for 0114 /// rho_m. If the latter is not supplied, rho_m is assumed to 0115 /// always be zero (this behaviour can be changed by calling 0116 /// set_common_bge_for_rho_and_rhom). 0117 /// See also the discussion above (in particular for the use of 0118 /// bge_rhom) 0119 GenericSubtractor(BackgroundEstimatorBase *bge_rho, 0120 BackgroundEstimatorBase *bge_rhom=0) : 0121 _bge_rho(bge_rho), _bge_rhom(bge_rhom), _jet_pt_fraction(0.01), 0122 _common_bge(false), _rhom_from_bge_rhom(false), 0123 _rho(_invalid_rho), _externally_supplied_rho_rhom(false) {} 0124 0125 /// Constructor that takes an externally supplied value for rho and, 0126 /// optionally, for rho_m. The latter defaults to zero. 0127 GenericSubtractor(double rho, double rhom=0); 0128 0129 0130 /// destructor 0131 ~GenericSubtractor(){} 0132 0133 /// a description of what this class does 0134 std::string description() const; 0135 0136 /// returns the estimate of the supplied shape, for the given jet, 0137 /// after background subtraction. 0138 double operator()(const FunctionOfPseudoJet<double> &shape, 0139 const PseudoJet &) const; 0140 0141 /// returns the estimate of the supplied shape, for the given jet, 0142 /// after background subtraction. It also sets the "info" variable 0143 /// with information about the details of the subtraction (e.g. the 0144 /// shape's derivatives wrt to the amount of pileup, and the background 0145 /// densities used). 0146 double operator()(const FunctionOfPseudoJet<double> & shape, 0147 const PseudoJet & jet, GenericSubtractorInfo & info) const; 0148 0149 /// when only one background estimator ('bge_rho') is specified, 0150 /// calling this routine with argument "true", causes rho_m to be be 0151 /// calculated from the same background estimator as rho, instead of 0152 /// being set to zero. 0153 /// 0154 /// Since FastJet 3.1, this will use the rho_m calculation from the 0155 /// background estimator if it is available [the default for both 0156 /// GridMedianBackgroundEstimator and JetMedianBackgroundEstimator]. 0157 /// In all other cases, the use of a common background estimator for 0158 /// rho and rho_m only works if the estimator is a 0159 /// JetMedianBackgroundEstimator (or derived from it), [since it 0160 /// makes use of that class's set_jet_density_class(...) facility]. 0161 /// 0162 void set_common_bge_for_rho_and_rhom(bool value=true); 0163 0164 /// equivalent to 'set_common_bge_for_rho_and_rhom' 0165 /// 0166 /// NOTE: this is DEPRECATED and kept only for backwards 0167 /// compatibility reasons. Use 'set_common_bge_for_rho_and_rhom' 0168 /// instead. 0169 void use_common_bge_for_rho_and_rhom(bool value=true); 0170 0171 /// returns true if it uses a common background estimator for both 0172 /// rho and rho_m (see set_common_bge_for_rho_and_rhom for 0173 /// details) 0174 bool common_bge_for_rho_and_rhom() const { return _common_bge; } 0175 0176 /// when the GenericSubtractor has been created with two background 0177 /// estimators (one for rho, the second for rho_m), setting this to 0178 /// true will result in rho_m being estimated using 0179 /// bge_rhom->rho_m() instead of bge_rhom->rho() 0180 /// See also the discussion at the top of the class. 0181 void set_use_bge_rhom_rhom(bool value=true); 0182 0183 /// returns true if rho_m is being estimated using 0184 /// bge_rhom->rho_m(). (s 0185 bool use_bge_rhom_rhom() const { return _rhom_from_bge_rhom; } 0186 0187 /// sets the fraction of the jet pt that will be used in the 0188 /// stability condition when computing the derivatives 0189 void set_jet_pt_fraction_for_stability(double jet_pt_fraction){ 0190 _jet_pt_fraction=jet_pt_fraction; 0191 } 0192 0193 0194 0195 protected: 0196 // tools to help compute the derivatives 0197 //---------------------------------------------------------------------- 0198 /// do the computation of the various derivatives 0199 /// 0200 /// rhoA_pt_fraction is rho_p/(rho_p + rho_m) used to know along 0201 /// which trajectory in the (rho_p, rho_m) plane one must estimate 0202 /// the derivative. 0203 void _compute_derivatives(const FunctionOfPseudoJet<double> &shape, 0204 const PseudoJet &jet, 0205 const double original_ghost_scale, 0206 const double ghost_area, 0207 const double f0, 0208 const double rho_pt_fraction, 0209 GenericSubtractorInfo &info) const; 0210 0211 /// make a sweep in stepsize to determine which step is the most precise 0212 /// 0213 /// Note that this returns the values of the function at the points 0214 /// needed to compute the derivatives 0215 double _optimize_step(const FunctionOfPseudoJet<double> &shape, 0216 const PseudoJet &jet, 0217 const double original_ghost_scale, 0218 const double ghost_area, 0219 const double x_fraction, 0220 const double f0, 0221 double cached_functions[4], 0222 const double max_step) const; 0223 0224 /// the function that does the rescaling 0225 double _shape_with_rescaled_ghosts(const FunctionOfPseudoJet<double> &shape, 0226 const PseudoJet &jet, 0227 const double original_ghost_scale, 0228 const double new_ghost_scale, 0229 const double new_dmass_scale) const; 0230 0231 /// if the shape is a ShapeWithComponents, apply the subtraction on 0232 /// each of the components 0233 double _component_subtraction(const ShapeWithComponents * shape_ptr, 0234 const PseudoJet & jet, 0235 GenericSubtractorInfo &info) const; 0236 0237 BackgroundEstimatorBase *_bge_rho, *_bge_rhom; 0238 double _jet_pt_fraction; 0239 bool _common_bge, _rhom_from_bge_rhom; 0240 double _rho, _rhom; 0241 bool _externally_supplied_rho_rhom; 0242 /// a value of rho that is used as a default to label that the stored 0243 /// rho is not valid for subtraction. 0244 // 0245 // NB: there are two reasons for not having the value written here: 0246 // 1) that it caused problems on karnak with g++ 4.0.1 and 2) that 0247 // we anyway like -infinity as a default, and since that's a function, 0248 // that's not allowed in an include file. 0249 static const double _invalid_rho; 0250 0251 /// deprecated "use_common_bge_for_rho_and_rhom" 0252 static LimitedWarning _warning_depracated_use_common_bge; 0253 static LimitedWarning _warning_unused_rhom; 0254 }; 0255 0256 //------------------------------------------------------------------------ 0257 /// \class GenericSubtractorInfo 0258 /// 0259 /// Helper class that allows to get extra information about the shape 0260 /// subtraction 0261 class GenericSubtractorInfo{ 0262 public: 0263 /// returns the unsubtracted shape 0264 double unsubtracted() const{ return _unsubtracted;} 0265 0266 /// returns the result of applying the subtraction including only the 0267 /// first-order correction 0268 double first_order_subtracted() const{ return _first_order_subtracted;} 0269 0270 /// returns the result of applying the subtraction including the 0271 /// first and second-order corrections (this is the default returned 0272 /// by the subtractions). 0273 double second_order_subtracted() const{ return _second_order_subtracted;} 0274 0275 /// returns an estimate of the 3rd order correction. Note that the 0276 /// third derivative is quite likely to be poorly evaluated in some 0277 /// cases. For this reason, it is not used by default. 0278 double third_order_subtracted() const{ return _third_order_subtracted;} 0279 0280 /// returns the first derivative (wrt to rho+rhom). Derivatives 0281 /// are taken assuming a constant value for the ratio of rhom/rho, as 0282 /// determined from the background estimates for the jet. 0283 /// 0284 /// Note: derivatives are not evaluated for shapes with components, 0285 /// and therefore set to zero. 0286 double first_derivative() const{ return _first_derivative;} 0287 0288 /// returns the second derivative (wrt to rho+rhom). It is not evaluated 0289 /// for shapes with components, and therefore set to zero. 0290 double second_derivative() const{ return _second_derivative;} 0291 0292 /// returns an estimate of the 3rd derivative (wrt to rho+rhom); 0293 /// note that the third derivative is quite likely to be poorly 0294 /// evaluated in some cases. It is not used by default. It is not 0295 /// evaluated for shapes with components, and therefore set to zero. 0296 double third_derivative() const{ return _third_derivative;} 0297 0298 /// returns the ghost scale actually used for the computation 0299 double ghost_scale_used() const{ return _ghost_scale_used;} 0300 0301 /// return the value of rho and rhom used in the subtraction 0302 double rho() const { return _rho; } 0303 double rhom() const { return _rhom; } 0304 0305 protected: 0306 double _unsubtracted; 0307 double _first_order_subtracted; 0308 double _second_order_subtracted; 0309 double _third_order_subtracted; 0310 double _first_derivative, _second_derivative, _third_derivative; 0311 double _ghost_scale_used; 0312 double _rho, _rhom; 0313 0314 friend class GenericSubtractor; 0315 }; 0316 0317 0318 } 0319 0320 FASTJET_END_NAMESPACE 0321 0322 #endif // __FASTJET_CONTRIB_GENERIC_SUBTRACTION_HH__
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |