|
|
|||
File indexing completed on 2026-06-21 08:23:37
0001 // 0002 // SignalFreeBackgroundEstimator package 0003 // Questions/comments: peter.berta@cern.ch 0004 // 0005 // Copyright (c) 2024-, Peter Berta 0006 // 0007 //---------------------------------------------------------------------- 0008 // This file is part of FastJet contrib. 0009 // 0010 // It is free software; you can redistribute it and/or modify it under 0011 // the terms of the GNU General Public License as published by the 0012 // Free Software Foundation; either version 2 of the License, or (at 0013 // your option) any later version. 0014 // 0015 // It is distributed in the hope that it will be useful, but WITHOUT 0016 // ANY WARRANTY; without even the implied warranty of MERCHANTABILITY 0017 // or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public 0018 // License for more details. 0019 // 0020 // You should have received a copy of the GNU General Public License 0021 // along with this code. If not, see <http://www.gnu.org/licenses/>. 0022 //---------------------------------------------------------------------- 0023 0024 0025 #ifndef __SIGNAL_FREE_BACKGROUND_ESTIMATOR_H__ 0026 #define __SIGNAL_FREE_BACKGROUND_ESTIMATOR_H__ 0027 0028 #include "fastjet/tools/BackgroundEstimatorBase.hh" 0029 #include "fastjet/RectangularGrid.hh" 0030 #include "fastjet/tools/GridMedianBackgroundEstimator.hh" 0031 #include "fastjet/ClusterSequenceArea.hh" 0032 0033 FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh 0034 0035 0036 namespace contrib{ 0037 0038 /// \class SignalFreeBackgroundEstimator 0039 /// Background estimator class based on the method described in arXiv:2304.08383. The design of this class is based on the GridMedianBackgroundEstimator class located in the official FastJet software: ./tools/fastjet/tools/GridMedianBackgroundEstimator.hh 0040 0041 0042 class SignalFreeBackgroundEstimator : public fastjet::BackgroundEstimatorBase, 0043 public fastjet::RectangularGrid { 0044 public: 0045 SignalFreeBackgroundEstimator(double rapidity_max, 0046 double requested_grid_spacing): fastjet::RectangularGrid(rapidity_max, requested_grid_spacing), _rapidity_max(rapidity_max) {} 0047 0048 0049 /// Return a pointer to a copy of this estimator 0050 SignalFreeBackgroundEstimator * copy() const override { 0051 return new SignalFreeBackgroundEstimator(*this); 0052 }; 0053 0054 /// 0055 /// default dtor 0056 virtual ~SignalFreeBackgroundEstimator(){} 0057 0058 /// Tell the background estimator that it has a new event, composed of the specified particles. 0059 void set_particles(const std::vector<fastjet::PseudoJet>& particles) override { 0060 std::vector<fastjet::PseudoJet> empty; 0061 this->set_particles(particles,empty); 0062 } 0063 0064 /// Tell the background estimator that it has a new event, composed of the specified particles and also include estimated signal particles. If the variable measureOfPileup is set to -1, then the minimal jet rho threshold is constant and set to the parameter _jet_rho_min_intercept. If it is >0, then it is not constant and varies linearly with measureOfPileup. If signalChargedParticles are also provided, then the signal jets are obtained also from clustering of these particles (separate jet clustering is used). 0065 void set_particles(const std::vector<fastjet::PseudoJet>& particles, const std::vector<fastjet::PseudoJet>& signalParticles, const double& measureOfPileup=-1 /* can be nPU, nPV, mu, rho*/, const std::vector<fastjet::PseudoJet>& signalChargedParticles=std::vector<fastjet::PseudoJet>()); 0066 0067 /// Define signal seeds from the user. This function must be called before "set_particles" 0068 void add_seeds_from_user(const std::vector<fastjet::PseudoJet> &additional_seeds); 0069 0070 /// Get the full set of background properties 0071 fastjet::BackgroundEstimate estimate() const override; 0072 0073 /// Get the full set of background properties for a given reference jet 0074 fastjet::BackgroundEstimate estimate(const fastjet::PseudoJet& jet) const override; 0075 0076 /// Returns rho, the average background density per unit area 0077 double rho() const override; 0078 0079 /// Returns rho, the average background density per unit area, locally at the position of a given jet 0080 double rho(const fastjet::PseudoJet & jet) override; 0081 0082 /// Returns rho_m, the purely longitudinal, particle-mass-induced component of the background density per unit area 0083 double rho_m() const FASTJET_OVERRIDE; 0084 0085 /// Returns rho_m locally at the jet position. As for rho(jet), it is non-const. 0086 double rho_m(const PseudoJet & jet) FASTJET_OVERRIDE; 0087 0088 /// determine whether the automatic calculation of rho_m (by default true) 0089 void set_compute_rho_m(bool enable){ _enable_rho_m = enable; } 0090 0091 /// Returns true if this background estimator has support for determination of rho_m. 0092 bool has_rho_m() const FASTJET_OVERRIDE {return _enable_rho_m;} 0093 0094 /// Returns a textual description of the background estimator. 0095 std::string description() const override; 0096 0097 /// Set a pointer to a rescaling function 0098 virtual void set_rescaling_class(const fastjet::FunctionOfPseudoJet<double>* rescaling_class) override; 0099 0100 /// Set the parameters for signal seeds (distance parameter for clustering of signal particles into jets, and deltaR distance around signal jets to exclude areas with signal) 0101 void set_signal_seed_parameters(const double &signal_jets_R, const double &exclusion_deltaR){ _signal_jets_R = signal_jets_R; _exclusion_deltaR = exclusion_deltaR; } 0102 0103 /// Set the parameters for the minimal jet rho when selecting jets from the signal jets 0104 void set_jet_rho_min(const double &jet_rho_min_intercept, const double &jet_rho_min_slope){ _jet_rho_min_intercept = jet_rho_min_intercept; _jet_rho_min_slope=jet_rho_min_slope;} 0105 0106 /// Set the parameter for the minimal jet rho when selecting jets obtained from charged signal particles 0107 void set_jet_rho_min_charged(const double &jet_rho_min_charged){ _jet_rho_min_charged = jet_rho_min_charged;} 0108 0109 /// Set the parameter for the spacing of ghosts 0110 void set_ghost_size(const double &ghost_size){ _ghost_size = ghost_size; } 0111 0112 /// Set the parameter for the minimal tile area 0113 void set_tile_area_min(const double &tile_area_min){ _tile_area_min = tile_area_min; } 0114 0115 /// Set parameters for the window to be used to compute the weighted median 0116 void set_window_parameters(const double ¢er, const double &half_window, const double ®ulator_for_floating_center=-1, const double &signal_fraction_min=0.005, const double &shift_max=0.15){ _center = center; _half_window = half_window; _regulator_for_floating_center=regulator_for_floating_center; _signal_fraction_min=signal_fraction_min; _shift_max=shift_max; } 0117 0118 /// Set to use weighted tiles for rho computation 0119 void set_use_weighted_tiles(const bool &use_weighted_tiles){ _use_weighted_tiles = use_weighted_tiles; } 0120 0121 private: 0122 // function to compute the weighted median 0123 double _compute_weighted_median(std::vector<std::pair<double,double> > &tileRho) const; 0124 0125 /// Store eta extend. 0126 double _rapidity_max=-1; 0127 0128 /// Distance parameter for jet clustering with anti-kt algorithm to obtain the signal seeds 0129 double _signal_jets_R=0.3; 0130 0131 /// Parameter for the exclusion delta R around the signal seeds 0132 double _exclusion_deltaR=0.4; 0133 0134 /// Parameters to compute the minimal jet rho threshold: intercept and slope of the linear function as a function of some measure of pileup specified in the set_particles function: 0135 double _jet_rho_min_intercept=5; 0136 double _jet_rho_min_slope=8; 0137 0138 /// Parameter for the minimal jet rho threshold for the jets clustered from charged particles. The default value is 20 GeV/(unit area), if the input charged particles are provided in GeV. 0139 double _jet_rho_min_charged=20; 0140 0141 /// Parameter for the minimal tile area after excluding signal-contaminated areas. Tiles with smaller areas are excluded from the weighted median computation. Can be modified with function set_tile_area_min 0142 double _tile_area_min=0.00001; 0143 0144 /// Parameter for the spacing of ghosts to evaluate tile areas after excluding signal-contaminated areas. The smaller ghostSize, the more precise is the computation, but also slower. The default is 0.04. Can be modified with function set_ghost_size 0145 double _ghost_size=0.04; 0146 0147 /// Window parameters: 0148 /// Parameters for the window center and the size of half window used to compute rho 0149 double _center=0.5; 0150 double _half_window=0.1; 0151 /// Parameters for the floating center of the window to compute rho 0152 double _regulator_for_floating_center=-1; 0153 double _signal_fraction_min=0.005; 0154 double _shift_max=0.15; 0155 /// Fraction of scalar pT from signal particles used to move the window to get the median. It is computed based on the provided signal particles. 0156 double _signal_fraction; 0157 0158 /// Sum of all tile areas which are used as weights to compute the weighted median 0159 double _weights_sum; 0160 0161 /// bool to determine whether the rho should be computed using weigted tiles 0162 bool _use_weighted_tiles=true; 0163 0164 /// verify that particles have been set and throw an error if not 0165 void verify_particles_set() const; 0166 0167 /// Rescaling warning 0168 fastjet::LimitedWarning _warning_rescaling; 0169 0170 bool _enable_rho_m=true; 0171 std::vector<fastjet::PseudoJet> _signal_seeds_from_user; 0172 0173 }; 0174 0175 } 0176 0177 FASTJET_END_NAMESPACE 0178 0179 #endif // __SIGNAL_FREE_BACKGROUND_ESTIMATOR_H__ 0180
| [ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
|
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |
|