|
||||
File indexing completed on 2025-01-18 09:57:13
0001 // $Id: ConstituentSubtractor.hh 1240 2020-02-23 13:51:05Z peter.berta $ 0002 // 0003 // ConstituentSubtractor package 0004 // Questions/comments: berta@ipnp.troja.mff.cuni.cz, Martin.Spousta@cern.ch, David.W.Miller@uchicago.edu 0005 // 0006 // Copyright (c) 2014-, Peter Berta, Martin Spousta, David W. Miller, Rupert Leitner 0007 // 0008 //---------------------------------------------------------------------- 0009 // This file is part of FastJet contrib. 0010 // 0011 // It is free software; you can redistribute it and/or modify it under 0012 // the terms of the GNU General Public License as published by the 0013 // Free Software Foundation; either version 2 of the License, or (at 0014 // your option) any later version. 0015 // 0016 // It is distributed in the hope that it will be useful, but WITHOUT 0017 // ANY WARRANTY; without even the implied warranty of MERCHANTABILITY 0018 // or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public 0019 // License for more details. 0020 // 0021 // You should have received a copy of the GNU General Public License 0022 // along with this code. If not, see <http://www.gnu.org/licenses/>. 0023 //---------------------------------------------------------------------- 0024 0025 #ifndef __FASTJET_CONTRIB_CONSTITUENTSUBTRACTOR_HH__ 0026 #define __FASTJET_CONTRIB_CONSTITUENTSUBTRACTOR_HH__ 0027 0028 0029 #include <fastjet/internal/base.hh> 0030 #include <fastjet/ClusterSequenceAreaBase.hh> 0031 #include <fastjet/tools/JetMedianBackgroundEstimator.hh> 0032 #include <fastjet/tools/GridMedianBackgroundEstimator.hh> 0033 #include <fastjet/PseudoJet.hh> 0034 #include <fastjet/Selector.hh> 0035 #include <fastjet/tools/BackgroundEstimatorBase.hh> 0036 #include "fastjet/tools/Transformer.hh" // to derive Subtractor from Transformer 0037 #include "fastjet/LimitedWarning.hh" 0038 0039 #include <algorithm> 0040 0041 FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh 0042 0043 namespace contrib{ 0044 0045 0046 //------------------------------------------------------------------------ 0047 /// \class ConstituentSubtractor 0048 /// A class to perform subtraction of background, e.g. pileup, from a set of particles at particle-level. The output is a jet or the whole event with corrected constituents. 0049 /// 0050 /// This class corrects the input particles for background contamination with the algorithm described in: 0051 /// Peter Berta, Martin Spousta, David W. Miller, Rupert Leitner [arXiv:1403.3108] 0052 /// 0053 /// For individual jet background subtraction, see example_jet_by_jet.cc 0054 /// For whole event background subtraction before jet clustering, see example_event_wide.cc 0055 /// 0056 /// The distance used for matching between particle i and ghost k is defined as: 0057 /// deltaR_{i,k}=pT_i * sin(theta_i)^polarAngleExp * sqrt((y_i-y_k)^2 + (phi_i-phi_k)^2) 0058 /// 0059 /// The class accounts for position-dependent (in rapidity-azimuth plane) background densities, rho and rho_m. The user is encouraged to use them to get the maximal performance. 0060 /// 0061 /// 0062 /// For the treatment of massive particles, choose one of the following options (currently the default option is 3 which seems to be one of the most optimal in our studies): 0063 /// 1. Keep the original mass and rapidity. Not recommended since jet mass is wrongly reconstructed. One can use it by adding before initialization: 0064 /// subtractor.set_keep_original_masses(); 0065 /// 2. Keep the original mass and pseudo-rapidity. Not recommended since jet mass is wrongly reconstructed. Use these two functions for this option: 0066 /// subtractor.set_keep_original_masses(); 0067 /// subtractor.set_fix_pseudorapidity(); 0068 /// 3. Set all masses to zero. Keep rapidity unchanged. This is recommended and is the default option, since observed the best performance for high pt large-R jets. 0069 /// Nothing needs to be specified. 0070 /// 4. Set all masses to zero. Keep pseudo-rapidity unchanged. Also recommended, almost the same performance as for option 3. Use function: 0071 /// subtractor.set_fix_pseudorapidity(); 0072 /// 5. Do correction of m_delta=sqrt(p_T^2+m^2)-p_t. Not recommended. One can use it by adding before initialization: 0073 /// subtractor.set_do_mass_subtraction(); 0074 /// One must additionally provide option for the rho_m background estimation. There are several possibilities: 0075 /// a) Same background estimator as for rho. Use this function: 0076 /// subtractor.set_common_bge_for_rho_and_rhom(); // this must be used after set_background_estimator function. 0077 /// b) Use a separate background estimator - you need to specify it as an additional parameter in function set_background_estimator 0078 /// c) Use scalar background density using function set_scalar_background_density(double rho, double rhom). 0079 /// 6. Same as 5 just keep pseudo-rapidity unchanged. Not recommended. Additionally to what is written for option 5, use: 0080 /// subtractor.set_fix_pseudorapidity(); 0081 /// 7. Keep rapidity and pseudo-rapidity fixed (scale fourmomentum). Recommended - observed better performance than the mass correction. Use: 0082 /// subtractor.set_scale_fourmomentum(); 0083 0084 0085 0086 namespace ConstituentSubtractorConstants{ 0087 const double zero_pt=1e-50; // This pt value is considered as zero. It is used for corrected particles with zero pt. 0088 const double zero_mass=1e-50; // This mass value is considered as zero. It is used for corrected particles with zero mass. 0089 // By default particles with zero mass and pt are discarded in case their corrected pt and mass are zero. The user can use the function "set_remove_zero_pt_and_mass_particles(false)" to keep such partiles. 0090 } 0091 0092 0093 class ConstituentSubtractor : public fastjet::Transformer{ 0094 public: 0095 0096 enum Distance { 0097 deltaR, /// deltaR=sqrt((y_i-y_j)^2+(phi_i-phi_j)^2)), longitudinal Lorentz invariant 0098 angle /// angle between two momenta in Euclidean space 0099 }; 0100 0101 /// 0102 /// default ctor 0103 ConstituentSubtractor(); 0104 0105 /// 0106 /// Constructor that takes a pointer to a background estimator for rho and optionally a pointer to a background estimator for rho_m. If the latter is not supplied, rho_m is assumed to always be zero (this behaviour can be changed by calling use_common_bge_for_rho_and_rhom). 0107 ConstituentSubtractor(fastjet::BackgroundEstimatorBase *bge_rho, fastjet::BackgroundEstimatorBase *bge_rhom=0, double alpha=0, double max_distance=-1, Distance distance=deltaR); 0108 0109 /// 0110 /// Constructor that takes an externally supplied value for rho and rho_m. 0111 ConstituentSubtractor(double rho, double rhom=0, double alpha=0, double max_distance=-1, Distance distance=deltaR); 0112 0113 /// 0114 /// default dtor 0115 virtual ~ConstituentSubtractor(){} 0116 0117 /// 0118 /// initialization. use it before event loop. 0119 virtual void initialize(); 0120 0121 /// 0122 /// Common description for jet-by-jet, event-wide, and iterative CS 0123 void description_common(std::ostringstream &descr) const; 0124 0125 /// 0126 /// a description of what this class does 0127 virtual std::string description() const; 0128 0129 /// 0130 /// action of the correction on a given jet. The output is PseudoJet object with subtracted constituents 0131 virtual fastjet::PseudoJet result(const fastjet::PseudoJet &jet) const; 0132 0133 /// 0134 /// do the constituent subtraction for the input particles using the provided background proxies. The output is a vector with corrected particles - particles with zero corrected pt and mass are removed by default. The user can modify this behaviour by using functions set_remove_particles_with_zero_pt_and_mass. 0135 std::vector<fastjet::PseudoJet> do_subtraction(std::vector<fastjet::PseudoJet> const &particles, std::vector<fastjet::PseudoJet> const &backgroundProxies,std::vector<fastjet::PseudoJet> *remaining_backgroundProxies=0) const; 0136 0137 0138 /// 0139 /// this version should be not used. Use the version below. 0140 virtual std::vector<fastjet::PseudoJet> subtract_event(std::vector<fastjet::PseudoJet> const &particles, double max_eta); 0141 0142 /// 0143 /// do the subtraction of the whole event - more user-friendly approach. The particles with |eta|>max_eta are discarded at the beginning, i.e. they are not used, nor returned. The ghosts are added automatically inside this function up to max_eta. 0144 virtual std::vector<fastjet::PseudoJet> subtract_event(std::vector<fastjet::PseudoJet> const &particles, std::vector<fastjet::PseudoJet> const *hard_proxies=0); 0145 0146 /// 0147 /// do the subtraction of the whole event using the tracking information for charged particles, i.e. the 4-momenta of charged particles from signal vertex, and 4-momenta of charged particles from background. The user can set the scaling of charged particles from background and signal using parameters charged_background_scale (CBS) and charged_signal_scale (CSS). These scales are useful if one assumes correlation between charged and neutral particles or in case the inputs from calorimeter are miscalibrated wrt tracks. In case CBS=CSS=0, the input charged particles are not used. In case CBS=CSS=1, the input charged particles are not scaled. Recommending to try several combinations for CBS and CSS from range [0.8, 1.5]. It is no more necessary to provide background estimator. The GridMedianBackgroundEstimator is used - probably more flexibility will be added in the future. The rescaling function for background estimator can be also provided - the rescaling function will be used for the event after subtracting charged scaled particles. The particles with |eta|>max_eta are discarded at the beginning, i.e. they are not used, nor returned. 0148 std::vector<fastjet::PseudoJet> subtract_event_using_charged_info(std::vector<fastjet::PseudoJet> const &particles, double charged_background_scale, std::vector<fastjet::PseudoJet> const &charged_background, double charged_signal_scale, std::vector<fastjet::PseudoJet> const &charged_signal, double max_eta); 0149 0150 0151 void set_rescaling(fastjet::FunctionOfPseudoJet<double> *rescaling); 0152 0153 /// 0154 /// set the grid size (not area) for the background estimation with GridMedianBackgroundEstimator when using charged info 0155 void set_grid_size_background_estimator(double const &grid_size_background_estimator); 0156 0157 0158 /// 0159 /// Set the pointer to a background estimator for rho and optionally a pointer to a background estimator for rho_m. If the latter is not supplied, rho_m is assumed to always be zero (this behaviour can be changed by calling use_common_bge_for_rho_and_rhom). 0160 void set_background_estimator(fastjet::BackgroundEstimatorBase *bge_rho, fastjet::BackgroundEstimatorBase *bge_rhom=0); 0161 0162 /// 0163 /// Set the scalar background densities rho and rho_m. 0164 void set_scalar_background_density(double rho, double rhom=0); 0165 0166 /// 0167 /// When only one background estimator, bge_rho, is specified, calling this function with argument true, causes rho_m to be calculated from the same background estimator as rho, instead of being set to zero. Currently this only works if the estimator is a JetMedianBackgroundEstimator or other estimator which has such function. 0168 void set_common_bge_for_rho_and_rhom(); 0169 0170 /// 0171 /// This function should not be used anymore (instead the function above should be used without any parameter). When used with "true", then it has the same behaviour as set_common_bge_for_rho_and_rhom() above. By default, no mass subtraction is done, so it makes no sense to use this function with "false". 0172 void set_common_bge_for_rho_and_rhom(bool value); 0173 0174 0175 /// 0176 /// By using this function, the original masses of particles are kept. By default, the masses of all particles are set to zero. If you do not want to set the masses to zero but you want to do the subtraction for the mass part, use the function set_common_bge_for_rho_and_rhom or specify background estimator for rho_m. 0177 void set_keep_original_masses(); 0178 0179 /// 0180 /// When two background estimators are used (one for rho, the second for rho_m), setting this to true will result in rho_m being estimated using bge_rhom->rho_m() instead of bge_rhom->rho(). 0181 /// void set_use_bge_rhom_rhom(bool value=true); 0182 0183 /// 0184 /// With this function, the user specifies if also the mass term sqrt(pT^2+m^2)-pT should be corrected during the subtraction procedure. This is automatically set when calling functions set_common_bge_for_rho_and_rhom, set_scalar_background_density or set_use_bge_rhom_rhom 0185 void set_do_mass_subtraction(); 0186 0187 /// 0188 /// set to true, if you want to remove particles which have zero both, pt and mass. By default, these particles are removed. 0189 void set_remove_particles_with_zero_pt_and_mass(bool value=true); 0190 0191 /// 0192 /// set to true, if you want to remove all zero pt particles - this means removing also particles with zero delta_m (massive particles with zero pt). By default, the zero pt particles with non-zero delta_m are not removed. 0193 void set_remove_all_zero_pt_particles(bool value=true); 0194 0195 /// 0196 /// function to change the alpha-parameter figuring in the distance measure deltaR. The larger the alpha, the more are preferred to be corrected the low pt particles. The default value is 0, i.e. by default the standard deltaR definition is used: deltaR=sqrt(deltay^2 + deltaphi^2) 0197 void set_alpha(double alpha); 0198 0199 /// 0200 /// function to change the parameter polarAngleExp 0201 void set_polarAngleExp(double polarAngleExp); 0202 0203 0204 /// 0205 /// function to change the parameter ghost_area 0206 void set_ghost_area(double ghost_area); 0207 0208 /// 0209 /// function to change the distance type 0210 void set_distance_type(Distance distance); 0211 0212 /// 0213 /// function to change the free parameter max_distance. The distance is defined by enum Distance. The particle-ghost pairs with distance>max_distance are not used. When max_distance<=0, the max_distance parameter is not used (no upper limit on distance). The default value is -1, i.e. by default there is no upper limit for possible distance values. 0214 void set_max_distance(double max_distance); 0215 0216 /// 0217 /// function to change the free parameter max_distance. It is the same as the set_max_distance function. It is added back for backward-comptibility (since it was replaced in ConstituentSubtractor version 1.1.3 with function set_max_distance). It should be not used. 0218 void set_max_standardDeltaR(double max_distance); 0219 0220 /// 0221 /// function to return the maximal deltaR distance used in the subtraction 0222 double get_max_distance(); 0223 0224 /// 0225 /// This function returns true if the first argument is smaller than the second argument, otherwise returns false. The comparison is done only on the first element in the two pairs. This function is used to sort in ascending order the deltaR values for each pair particle-ghost while keeping track of particles and ghosts 0226 static bool _function_used_for_sorting(std::pair<double,std::pair<int,int> > const &i,std::pair<double,std::pair<int,int> > const &j); 0227 0228 /// 0229 /// function used for sorting of Pseudojets according to eta 0230 static bool _rapidity_sorting(fastjet::PseudoJet const &i,fastjet::PseudoJet const &j); 0231 0232 /// 0233 /// Function to construct massless ghosts in the whole eta-phi space up to max_eta 0234 void construct_ghosts_uniformly(double max_eta); 0235 0236 /// 0237 /// Function to return vector of ghosts 0238 std::vector<fastjet::PseudoJet> get_ghosts(); 0239 0240 /// 0241 /// Function to set the maximal eta value for subraction 0242 void set_max_eta(double max_eta); 0243 0244 /// 0245 /// Function to return vector of areas associated to ghosts 0246 std::vector<double> get_ghosts_area(); 0247 0248 0249 /// 0250 /// Set ghost selector - this allows to perform Constituent Subtraction in certain phase space when using the event-wide correction. Only ghosts which pass the selector are used in the correction. No selection is done on particles (i.e. the user needs to perform the selection before applying ConstituentSubtractor if he/she wishes). The default behaviour is to not use any selector. The user should remember that the max_eta restriction is applied in any case, e.g. when using the function "subtract_event(particles)", then ghosts are distributed to maximal eta max_eta and particles outside this region are discarded. 0251 void set_ghost_selector(fastjet::Selector* selector); 0252 0253 /// 0254 /// Set particle selector - if specified, only particles which pass this selector are corrected. The other particles are copied to the final corrected vector without modification 0255 void set_particle_selector(fastjet::Selector* selector); 0256 0257 /// 0258 /// Use this function if you want to use the information about hard particle proxies in the subtraction procedure (typically the charged tracks from primary vertex and/or high pt calorimeter clusters). The parameters of this function are the "nearby_hard_radius" and "nearby_hard_factor". If there is at least one hard particle proxy within Distance specified by the nearby_hard_radius parameter, then the CS distance is multiplied by nearby_hard_factor. 0259 void set_use_nearby_hard(double const &nearby_hard_radius, double const &nearby_hard_factor); 0260 0261 /// 0262 /// Use this function, if you want to fix the pseudo-rapidity instead of rapidity for the corrected particles. By default, the rapidity is fixed. 0263 void set_fix_pseudorapidity(); 0264 0265 /// 0266 /// Use this function if you want to fix both, pseudo-rapidity and rapidity. By default, rapidity is fixed and the mass is set to zero. 0267 void set_scale_fourmomentum(); 0268 0269 0270 protected: 0271 std::vector<fastjet::PseudoJet> get_background_proxies_from_ghosts(std::vector<fastjet::PseudoJet> const &ghosts,std::vector<double> const &ghosts_area) const; 0272 void clear_ghosts(); 0273 0274 /// 0275 /// common function to initialize for Iterative and standard CS 0276 void _initialize_common(); 0277 0278 /// 0279 /// helper functions to find the minimal and maximal indeces for the vector of particles - useful to make the CS procedure faster 0280 unsigned int _find_index_after(double const &value, std::vector<double> const &vec) const; 0281 unsigned int _find_index_before(double const &value, std::vector<double> const &vec) const; 0282 0283 /// 0284 /// function to get the trabsformed distance 0285 double _get_transformed_distance(double const &distance) const; 0286 0287 fastjet::BackgroundEstimatorBase *_bge_rho, *_bge_rhom; 0288 bool _common_bge, _rhom_from_bge_rhom; 0289 double _rho, _rhom; 0290 bool _externally_supplied_rho_rhom, _do_mass_subtraction; 0291 Distance _distance; 0292 double _alpha; 0293 double _polarAngleExp; 0294 double _max_distance; 0295 bool _remove_particles_with_zero_pt_and_mass; 0296 bool _remove_all_zero_pt_particles; 0297 bool _use_max_distance; 0298 double _ghost_area; 0299 double _grid_size_phi; 0300 double _grid_size_rap; 0301 bool _ghosts_constructed; 0302 bool _ghosts_rapidity_sorted; 0303 int _n_ghosts_phi; 0304 double _max_eta; 0305 bool _masses_to_zero; 0306 bool _use_nearby_hard; 0307 double _nearby_hard_radius; 0308 double _nearby_hard_factor; 0309 bool _fix_pseudorapidity; 0310 bool _scale_fourmomentum; 0311 std::vector<fastjet::PseudoJet> const *_hard_proxies; 0312 0313 std::vector<fastjet::PseudoJet> _ghosts; 0314 std::vector<double> _ghosts_area; 0315 std::vector<double> _ghosts_rapidities; // used for uniform grid 0316 double _grid_size_background_estimator; 0317 fastjet::Selector *_ghost_selector; 0318 fastjet::Selector *_particle_selector; 0319 fastjet::FunctionOfPseudoJet<double> *_rescaling; 0320 static LimitedWarning _warning_unused_rhom; 0321 }; 0322 0323 0324 0325 } // namespace contrib 0326 0327 FASTJET_END_NAMESPACE 0328 0329 #endif // __FASTJET_CONTRIB_CONSTITUENTSUBTRACTOR_HH__
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |