Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:57:16

0001 #ifndef __GRID_MEDIAN_BACKGROUND_ESTIMATOR_HH__
0002 #define __GRID_MEDIAN_BACKGROUND_ESTIMATOR_HH__
0003 
0004 //FJSTARTHEADER
0005 // $Id$
0006 //
0007 // Copyright (c) 2005-2021, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
0008 //
0009 //----------------------------------------------------------------------
0010 // This file is part of FastJet.
0011 //
0012 //  FastJet is free software; you can redistribute it and/or modify
0013 //  it under the terms of the GNU General Public License as published by
0014 //  the Free Software Foundation; either version 2 of the License, or
0015 //  (at your option) any later version.
0016 //
0017 //  The algorithms that underlie FastJet have required considerable
0018 //  development. They are described in the original FastJet paper,
0019 //  hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
0020 //  FastJet as part of work towards a scientific publication, please
0021 //  quote the version you use and include a citation to the manual and
0022 //  optionally also to hep-ph/0512210.
0023 //
0024 //  FastJet is distributed in the hope that it will be useful,
0025 //  but WITHOUT ANY WARRANTY; without even the implied warranty of
0026 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
0027 //  GNU General Public License for more details.
0028 //
0029 //  You should have received a copy of the GNU General Public License
0030 //  along with FastJet. If not, see <http://www.gnu.org/licenses/>.
0031 //----------------------------------------------------------------------
0032 //FJENDHEADER
0033 
0034 
0035 #include "fastjet/tools/BackgroundEstimatorBase.hh"
0036 #include "fastjet/RectangularGrid.hh"
0037 
0038 
0039 FASTJET_BEGIN_NAMESPACE      // defined in fastjet/internal/base.hh
0040 
0041 /// @ingroup tools_background
0042 /// \class GridMedianBackgroundEstimator
0043 /// 
0044 /// Background Estimator based on the median pt/area of a set of grid
0045 /// cells. 
0046 ///
0047 /// Description of the method:
0048 ///   This background estimator works by projecting the event onto a
0049 ///   grid in rapidity and azimuth. In each grid cell, the scalar pt
0050 ///   sum of the particles in the cell is computed. The background
0051 ///   density is then estimated by the median of (scalar pt sum/cell
0052 ///   area) for all cells.
0053 ///
0054 /// Parameters:
0055 ///   The class takes 2 arguments: the absolute rapidity extent of the 
0056 ///   cells and the size of the grid cells. Note that the size of the cell
0057 ///   will be adjusted in azimuth to satisfy the 2pi periodicity and
0058 ///   in rapidity to match the requested rapidity extent.
0059 ///
0060 /// Rescaling:
0061 ///   It is possible to use a rescaling profile. In this case, the
0062 ///   profile needs to be set before setting the particles and it will
0063 ///   be applied to each particle (i.e. not to each cell). 
0064 ///   Note also that in this case one needs to call rho(jet) instead of
0065 ///   rho() [Without rescaling, they are identical]
0066 ///
0067 class GridMedianBackgroundEstimator : public BackgroundEstimatorBase
0068                                     , public RectangularGrid
0069 {
0070 
0071 public:
0072   /// @name  constructors and destructors
0073   //\{
0074   //----------------------------------------------------------------
0075   ///   \param ymax   maximal absolute rapidity extent of the grid
0076   ///   \param requested_grid_spacing   size of the grid cell. The
0077   ///            "real" cell size could differ due e.g. to the 2pi
0078   ///             periodicity in azimuthal angle (size, not area)
0079   GridMedianBackgroundEstimator(double ymax, double requested_grid_spacing) :
0080     RectangularGrid(ymax, requested_grid_spacing),
0081     _enable_rho_m(true) {} 
0082 
0083   //----------------------------------------------------------------
0084   /// Constructor based on a user's fully specified RectangularGrid
0085   GridMedianBackgroundEstimator(const RectangularGrid & grid) :
0086     RectangularGrid(grid), _enable_rho_m(true) {
0087     if (!RectangularGrid::is_initialised()) 
0088       throw Error("attempt to construct GridMedianBackgroundEstimator with uninitialised RectangularGrid");
0089   }    
0090 
0091   //---------------------------------------------------------------- 
0092   /// Constructor with the explicit parameters for the underlying
0093   /// RectangularGrid
0094   ///
0095   ///  \param rapmin         the minimum rapidity extent of the grid
0096   ///  \param rapmax         the maximum rapidity extent of the grid
0097   ///  \param drap           the grid spacing in rapidity
0098   ///  \param dphi           the grid spacing in azimuth
0099   ///  \param tile_selector  optional (geometric) selector to specify 
0100   ///                        which tiles are good; a tile is good if
0101   ///                        a massless 4-vector at the center of the tile passes
0102   ///                        the selection
0103   GridMedianBackgroundEstimator(double rapmin_in, double rapmax_in, double drap_in, double dphi_in,
0104                                 Selector tile_selector = Selector()) :
0105     RectangularGrid(rapmin_in, rapmax_in, drap_in, dphi_in, tile_selector), _enable_rho_m(true) {}
0106 
0107   //\}
0108 
0109 
0110   /// @name setting a new event
0111   //\{
0112   //----------------------------------------------------------------
0113 
0114   /// tell the background estimator that it has a new event, composed
0115   /// of the specified particles.
0116   void set_particles(const std::vector<PseudoJet> & particles) FASTJET_OVERRIDE;
0117 
0118   /// determine whether the automatic calculation of rho_m and sigma_m
0119   /// is enabled (by default true)
0120   void set_compute_rho_m(bool enable){ _enable_rho_m = enable; }
0121 
0122   //\}
0123 
0124   /// return a pointer to a copy of this BGE; the user is responsible
0125   /// for eventually deleting the resulting object.
0126   BackgroundEstimatorBase * copy() const FASTJET_OVERRIDE {
0127     return new GridMedianBackgroundEstimator(*this);
0128   };
0129 
0130 
0131 
0132   /// @name  retrieving fundamental information
0133   //\{
0134   //----------------------------------------------------------------
0135   /// get the full set of background properties
0136   BackgroundEstimate estimate() const FASTJET_OVERRIDE;
0137   
0138   /// get the full set of background properties for a given reference jet
0139   BackgroundEstimate estimate(const PseudoJet &jet) const FASTJET_OVERRIDE;
0140 
0141   /// returns rho, the median background density per unit area
0142   double rho() const FASTJET_OVERRIDE;
0143 
0144   /// returns sigma, the background fluctuations per unit area; must be
0145   /// multipled by sqrt(area) to get fluctuations for a region of a
0146   /// given area.
0147   double sigma() const FASTJET_OVERRIDE;
0148 
0149   /// returns rho, the background density per unit area, locally at the
0150   /// position of a given jet. Note that this is not const, because a
0151   /// user may then wish to query other aspects of the background that
0152   /// could depend on the position of the jet last used for a rho(jet)
0153   /// determination.
0154   double rho(const PseudoJet & jet) FASTJET_OVERRIDE;
0155 
0156   /// returns sigma, the background fluctuations per unit area, locally at
0157   /// the position of a given jet. As for rho(jet), it is non-const.
0158   double sigma(const PseudoJet & jet) FASTJET_OVERRIDE;
0159 
0160   /// returns true if this background estimator has support for
0161   /// determination of sigma
0162   bool has_sigma() const FASTJET_OVERRIDE {return true;}
0163 
0164   //-----------------------------------------------------------------
0165   /// Returns rho_m, the purely longitudinal, particle-mass-induced
0166   /// component of the background density per unit area
0167   double rho_m() const FASTJET_OVERRIDE;
0168 
0169   /// returns sigma_m, a measure of the fluctuations in the purely
0170   /// longitudinal, particle-mass-induced component of the background
0171   /// density per unit area; must be multipled by sqrt(area) to get
0172   /// fluctuations for a region of a given area.
0173   double sigma_m() const FASTJET_OVERRIDE;
0174 
0175   /// Returns rho_m locally at the jet position. As for rho(jet), it is non-const.
0176   double rho_m(const PseudoJet & jet) FASTJET_OVERRIDE;
0177 
0178   /// Returns sigma_m locally at the jet position. As for rho(jet), it is non-const.
0179   double sigma_m(const PseudoJet & jet) FASTJET_OVERRIDE;
0180 
0181   /// Returns true if this background estimator has support for
0182   /// determination of rho_m.
0183   ///
0184   /// Note that support for sigma_m is automatic if one has sigma and
0185   /// rho_m support.
0186   bool has_rho_m() const FASTJET_OVERRIDE {return _enable_rho_m;}
0187 
0188 
0189   /// returns the area of the grid cells (all identical, but
0190   /// referred to as "mean" area for uniformity with JetMedianBGE).
0191   double mean_area() const {return mean_tile_area();}
0192   //\}
0193 
0194   /// @name configuring the behaviour
0195   //\{
0196   //----------------------------------------------------------------
0197 
0198   /// Set a pointer to a class that calculates the rescaling factor as
0199   /// a function of the jet (position). Note that the rescaling factor
0200   /// is used both in the determination of the "global" rho (the pt/A
0201   /// of each jet is divided by this factor) and when asking for a
0202   /// local rho (the result is multiplied by this factor).
0203   ///
0204   /// The BackgroundRescalingYPolynomial class can be used to get a
0205   /// rescaling that depends just on rapidity.
0206   ///
0207   /// Note that this has to be called BEFORE any attempt to do an
0208   /// actual computation
0209   ///
0210   /// The same profile will be used for both pt and mt (this is
0211   /// probabaly a good approximation since the particle density
0212   /// changes is what dominates the rapidity profile)
0213   virtual void set_rescaling_class(const FunctionOfPseudoJet<double> * rescaling_class) FASTJET_OVERRIDE;
0214 
0215   //\}
0216 
0217   /// @name description
0218   //\{
0219   //----------------------------------------------------------------
0220 
0221   /// returns a textual description of the background estimator
0222   std::string description() const FASTJET_OVERRIDE;
0223 
0224   //\}
0225 
0226 
0227 private:
0228 
0229   /// verify that particles have been set and throw an error if not
0230   void verify_particles_set() const;
0231 
0232   // information about the event
0233   //std::vector<double> _scalar_pt;
0234   //double _rho, _sigma, _rho_m, _sigma_m;
0235   BackgroundEstimate _cached_estimate;
0236   //bool _has_particles;
0237   bool _enable_rho_m;
0238 
0239   // various warnings to inform people of potential dangers
0240   LimitedWarning _warning_rho_of_jet;
0241   LimitedWarning _warning_rescaling;
0242 };
0243 
0244 FASTJET_END_NAMESPACE        // defined in fastjet/internal/base.hh
0245 
0246 #endif // __GRID_MEDIAN_BACKGROUND_ESTIMATOR_HH__