FastJet 3.0beta1
JetMedianBackgroundEstimator.hh
00001 #ifndef __FASTJET_BACKGROUND_ESTIMATOR_HH__
00002 #define __FASTJET_BACKGROUND_ESTIMATOR_HH__
00003 
00004 //STARTHEADER
00005 // $Id: JetMedianBackgroundEstimator.hh 2472 2011-07-26 11:25:42Z cacciari $
00006 //
00007 // Copyright (c) 2005-2011, Matteo Cacciari, Gavin Salam and Gregory Soyez
00008 //
00009 //----------------------------------------------------------------------
00010 // This file is part of FastJet.
00011 //
00012 //  FastJet is free software; you can redistribute it and/or modify
00013 //  it under the terms of the GNU General Public License as published by
00014 //  the Free Software Foundation; either version 2 of the License, or
00015 //  (at your option) any later version.
00016 //
00017 //  The algorithms that underlie FastJet have required considerable
00018 //  development and are described in hep-ph/0512210. If you use
00019 //  FastJet as part of work towards a scientific publication, please
00020 //  include a citation to the FastJet paper.
00021 //
00022 //  FastJet is distributed in the hope that it will be useful,
00023 //  but WITHOUT ANY WARRANTY; without even the implied warranty of
00024 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00025 //  GNU General Public License for more details.
00026 //
00027 //  You should have received a copy of the GNU General Public License
00028 //  along with FastJet; if not, write to the Free Software
00029 //  Foundation, Inc.:
00030 //      59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00031 //----------------------------------------------------------------------
00032 //ENDHEADER
00033 
00034 #include <fastjet/ClusterSequenceAreaBase.hh>
00035 #include <fastjet/AreaDefinition.hh>
00036 #include <fastjet/FunctionOfPseudoJet.hh>
00037 #include <fastjet/Selector.hh>
00038 #include <fastjet/tools/BackgroundEstimatorBase.hh>
00039 #include <iostream>
00040 
00041 FASTJET_BEGIN_NAMESPACE     // defined in fastjet/internal/base.hh
00042 
00043 
00044 /// @ingroup tools_background
00045 /// \class JetMedianBackgroundEstimator
00046 ///
00047 /// Class to estimate the pt density of the background per unit area,
00048 /// using the median of the distribution of pt/area from jets that
00049 /// pass some selection criterion.
00050 ///
00051 /// Events are passed either in the form of the event particles (in
00052 /// which they're clustered by the class), a ClusterSequenceArea (in
00053 /// which case the jets used are those returned by "inclusive_jets()")
00054 /// or directly as a set of jets.
00055 ///
00056 /// The selection criterion is typically a geometrical one (e.g. all
00057 /// jets with |y|<2) sometimes supplemented with some kinematical
00058 /// restriction (e.g. exclusion of the two hardest jets). It is passed
00059 /// to the class through a Selector.
00060 ///
00061 /// Beware: 
00062 ///   by default, to correctly handle partially empty events, the
00063 ///   class attempts to calculate an "empty area", based
00064 ///   (schematically) on
00065 ///
00066 ///          range.total_area() - sum_{jets_in_range} jets.area()
00067 ///  
00068 ///   For ranges with small areas, this can be inaccurate (particularly 
00069 ///   relevant in dense events where empty_area should be zero and ends
00070 ///   up not being zero).
00071 ///
00072 ///   This calculation of empty area can be avoided if a
00073 ///   ClusterSequenceArea class with explicit ghosts
00074 ///   (ActiveAreaExplicitGhosts) is used.  This is _recommended_
00075 ///   unless speed requirements cause you to use Voronoi areas. For
00076 ///   speedy background estimation you could also consider using
00077 ///   GridMedianBackgroundEstimator.
00078 ///
00079 ///
00080 class JetMedianBackgroundEstimator : public BackgroundEstimatorBase {
00081 public:
00082   /// @name constructors and destructors
00083   //\{
00084   //----------------------------------------------------------------
00085   /// Constructor that sets the rho range as well as the jet
00086   /// definition and area definition to be used to cluster the
00087   /// particles. Prior to the estimation of rho, one has to provide
00088   /// the particles to cluster using set_particles(...)
00089   ///
00090   /// \param rho_range  the Selector specifying which jets will be considered
00091   /// \param jet_def    the jet definition to use for the clustering
00092   /// \param area_def   the area definition to use for the clustering
00093   JetMedianBackgroundEstimator(const Selector &rho_range,
00094                                const JetDefinition &jet_def,
00095                                const AreaDefinition &area_def);
00096 
00097   /// ctor from a ClusterSequenceAreaBase with area
00098   ///
00099   /// \param rho_range   the Selector specifying which jets will be considered
00100   /// \param csa         the ClusterSequenceArea to use
00101   ///
00102   /// Pre-conditions: 
00103   ///  - one should be able to estimate the "empty area" (i.e. the area
00104   ///    not occupied by jets). This is feasible if at least one of the following
00105   ///    conditions is satisfied:
00106   ///     ( i) the ClusterSequence has explicit ghosts
00107   ///     (ii) the range has a computable area.
00108   ///  - the jet algorithm must be suited for median computation
00109   ///    (otherwise a warning will be issues)
00110   ///
00111   /// Note that selectors with e.g. hardest-jets exclusion do not have
00112   /// a well-defined area. For this reasons, it is STRONGLY advised to
00113   /// use an area with explicit ghosts.
00114   JetMedianBackgroundEstimator(const Selector &rho_range,
00115                                const ClusterSequenceAreaBase &csa);
00116 
00117 
00118   /// Default constructor that optionally sets the rho range. The
00119   /// configuration must be done later calling
00120   /// set_cluster_sequence(...) or set_jets(...).
00121   ///
00122   /// \param rho_range   the Selector specifying which jets will be considered
00123   ///
00124   JetMedianBackgroundEstimator(const Selector &rho_range = SelectorIdentity())
00125     : _rho_range(rho_range), _jet_def(JetDefinition()) { reset(); }
00126   
00127 
00128   /// default dtor
00129   ~JetMedianBackgroundEstimator(){}
00130 
00131   //\}
00132 
00133 
00134   /// @name setting a new event
00135   //\{
00136   //----------------------------------------------------------------
00137 
00138   /// tell the background estimator that it has a new event, composed
00139   /// of the specified particles.
00140   virtual void set_particles(const std::vector<PseudoJet> & particles);
00141 
00142   /// (re)set the cluster sequence (with area support) to be used by
00143   /// future calls to rho() etc. 
00144   ///
00145   /// \param csa  the cluster sequence area
00146   ///
00147   /// Pre-conditions: 
00148   ///  - one should be able to estimate the "empty area" (i.e. the area
00149   ///    not occupied by jets). This is feasible if at least one of the following
00150   ///    conditions is satisfied:
00151   ///     ( i) the ClusterSequence has explicit ghosts
00152   ///     (ii) the range selected has a computable area.
00153   ///  - the jet algorithm must be suited for median computation
00154   ///    (otherwise a warning will be issues)
00155   ///
00156   /// Note that selectors with e.g. hardest-jets exclusion do not have
00157   /// a well-defined area. For this reasons, it is STRONGLY advised to
00158   /// use an area with explicit ghosts.
00159   void set_cluster_sequence(const ClusterSequenceAreaBase & csa);
00160 
00161   /// (re)set the jets (which must have area support) to be used by future
00162   /// calls to rho() etc.; for the conditions that must be satisfied
00163   /// by the jets, see the Constructor that takes jets.
00164   void set_jets(const std::vector<PseudoJet> &jets);
00165 
00166   /// (re)set the selector to be used for future calls to rho() etc.
00167   void set_selector(const Selector & rho_range_selector) {
00168     _rho_range = rho_range_selector;
00169     _uptodate = false;
00170   }
00171 
00172   //\}
00173 
00174 
00175   /// @ name  retrieving fundamental information
00176   //\{
00177   //----------------------------------------------------------------
00178 
00179   /// get rho, the median background density per unit area
00180   double rho() const;
00181 
00182   /// get sigma, the background fluctuations per unit area
00183   double sigma() const;
00184 
00185   /// get rho, the median background density per unit area, locally at
00186   /// the position of a given jet.
00187   ///
00188   /// If the Selector associated with the range takes a reference jet
00189   /// (i.e. is relocatable), then for subsequent operations the
00190   /// Selector has that jet set as its reference.
00191   double rho(const PseudoJet & jet);
00192 
00193   /// get sigma, the background fluctuations per unit area,
00194   /// locally at the position of a given jet.
00195   ///
00196   /// If the Selector associated with the range takes a reference jet
00197   /// (i.e. is relocatable), then for subsequent operations the
00198   /// Selector has that jet set as its reference.
00199   double sigma(const PseudoJet &jet);
00200 
00201   /// returns true if this background estimator has support for
00202   /// determination of sigma
00203   virtual bool has_sigma() {return true;}
00204 
00205   //\}
00206   
00207   /// @ name  retrieving additional useful information
00208   //\{
00209   //----------------------------------------------------------------
00210   /// Returns the mean area of the jets used to actually compute the
00211   /// background properties in the last call of rho() or sigma()
00212   double mean_area() const{
00213     _recompute_if_needed();
00214     return _mean_area;
00215   }
00216   
00217   /// returns the number of jets used to actually compute the
00218   /// background properties in the last call of rho() or sigma()
00219   unsigned int n_jets_used() const{
00220     _recompute_if_needed();
00221     return _n_jets_used;
00222   }
00223 
00224   /// Returns the estimate of the area (within the range defined by
00225   /// the selector) that is not occupied by jets. The value is that
00226   /// for the last call of rho() or sigma()
00227   ///
00228   /// The answer is defined to be zero if the area calculation
00229   /// involved explicit ghosts; if the area calculation was an active
00230   /// area, then use is made of the active area's internal list of
00231   /// pure ghost jets (taking those that pass the selector); otherwise
00232   /// it is based on the difference between the selector's total area
00233   /// and the area of the jets that pass the selector.
00234   ///
00235   /// The result here is just the cached result of the corresponding
00236   /// call to the ClusterSequenceAreaBase function.
00237   double empty_area() const{
00238     _recompute_if_needed();
00239     return _empty_area;
00240   }
00241 
00242   /// Returns the number of empty jets used when computing the
00243   /// background properties. The value is that for the last call of
00244   /// rho() or sigma().
00245   ///
00246   /// If the area has explicit ghosts the result is zero; for active
00247   /// areas it is the number of internal pure ghost jets that pass the
00248   /// selector; otherwise it is deduced from the empty area, divided by 
00249   /// \f$ 0.55 \pi R^2 \f$ (the average pure-ghost-jet area).
00250   ///
00251   /// The result here is just the cached result of the corresponding
00252   /// call to the ClusterSequenceAreaBase function.
00253   double n_empty_jets() const{
00254     _recompute_if_needed();
00255     return _n_empty_jets;
00256   }
00257 
00258   //}
00259 
00260 
00261   /// @name configuring behaviour
00262   //\{
00263   //----------------------------------------------------------------
00264 
00265   /// Resets the class to its default state, including the choice to
00266   /// use 4-vector areas.
00267   ///
00268   void reset();
00269 
00270   /// By default when calculating pt/Area for a jet, it is the
00271   /// transverse component of the 4-vector area that is used in the ratiof \f$p_t/A\f$. 
00272   /// Calling this function with a "false" argument causes the scalar area to
00273   /// be used instead. 
00274   ///
00275   /// While the difference between the two choices is usually small,
00276   /// for high-precision work it is usually the 4-vector area that is
00277   /// to be preferred.
00278   ///
00279   ///  \param use_it             whether one uses the 4-vector area or not (true by default)
00280   void set_use_area_4vector(bool use_it = true){
00281     _use_area_4vector = use_it;
00282     _uptodate = false;
00283   }  
00284 
00285   /// check if the estimator uses the 4-vector area or the scalar area
00286   bool use_area_4vector() const{ return _use_area_4vector;}
00287 
00288   /// The FastJet v2.X sigma calculation had a small spurious offset
00289   /// in the limit of a small number of jets. This is fixed by default
00290   /// in versions 3 upwards. The old behaviour can be obtained with a
00291   /// call to this function.
00292   void set_provide_fj2_sigma(bool provide_fj2_sigma = true) {
00293     _provide_fj2_sigma = provide_fj2_sigma;
00294     _uptodate = false;
00295   }
00296 
00297   /// Set a pointer to a class that calculates the quantity whose
00298   /// median will be calculated; if the pointer is null then pt/area
00299   /// is used (as occurs also if this function is not called).
00300   ///
00301   /// Note that this is still <i>preliminary</i> in FastJet 3.0 and
00302   /// that backward compatibility is not guaranteed in future releases
00303   /// of FastJet
00304   void set_jet_density_class(const FunctionOfPseudoJet<double> * jet_density_class);
00305 
00306   /// return the pointer to the jet density class
00307   const FunctionOfPseudoJet<double> *  jet_density_class() const{
00308     return _jet_density_class;
00309   }
00310 
00311   /// Set a pointer to a class that calculates the rescaling factor as
00312   /// a function of the jet (position). Note that the rescaling factor
00313   /// is used both in the determination of the "global" rho (the pt/A
00314   /// of each jet is divided by this factor) and when asking for a
00315   /// local rho (the result is multiplied by this factor).
00316   ///
00317   /// The BackgroundRescalingYPolynomial class can be used to get a
00318   /// rescaling that depends just on rapidity.
00319   virtual void set_rescaling_class(const FunctionOfPseudoJet<double> * rescaling_class) {
00320     BackgroundEstimatorBase::set_rescaling_class(rescaling_class);
00321     _uptodate = false;
00322   }
00323 
00324   //\}
00325 
00326   /// @name description
00327   //\{
00328   //----------------------------------------------------------------
00329 
00330   /// returns a textual description of the background estimator
00331   std::string description() const;
00332 
00333   //\}
00334 
00335 
00336 private:
00337 
00338   /// do the actual job
00339   void _compute() const;
00340    
00341   /// check if the properties need to be recomputed 
00342   /// and do so if needed
00343   void _recompute_if_needed() const {
00344     if (!_uptodate) _compute();
00345     _uptodate = true;
00346   }
00347 
00348   /// for estimation using a selector that takes a reference jet
00349   /// (i.e. a selector that can be relocated) this function allows one
00350   /// to set its position.
00351   ///
00352   /// Note that this HAS to be called before any attempt to compute
00353   /// the background properties. The call is, however, performed
00354   /// automatically by the functions rho(jet) and sigma(jet).
00355   void _recompute_if_needed(const PseudoJet &jet);
00356 
00357   /// check that the underlying structure is still alive
00358   /// throw an error otherwise
00359   void _check_csa_alive() const;
00360 
00361   /// check that the algorithm used for the clustering is adapted for
00362   /// background estimation (i.e. either kt or C/A)
00363   /// Issue a warning otherwise
00364   void _check_jet_alg_good_for_median() const;
00365   
00366   // the basic parameters of this class (passed through the variou ctors)
00367   Selector _rho_range;                   ///< range to compute the background in
00368   JetDefinition _jet_def;                ///< the jet def to use for teh clustering
00369   AreaDefinition _area_def;              ///< the area def to use for teh clustering
00370   std::vector<PseudoJet> _included_jets; ///< jets to be used
00371   
00372   // the tunable aprameters of the class
00373   bool _use_area_4vector;
00374   bool _provide_fj2_sigma;
00375   const FunctionOfPseudoJet<double> * _jet_density_class;
00376   //SharedPtr<BackgroundRescalingBase> _rescaling_class_sharedptr;
00377   
00378   // the actual results of the computation
00379   mutable double _rho;               ///< background estimated density per unit area
00380   mutable double _sigma;             ///< background estimated fluctuations
00381   mutable double _mean_area;         ///< mean area of the jets used to estimate the background
00382   mutable unsigned int _n_jets_used; ///< number of jets used to estimate the background
00383   mutable double _n_empty_jets;      ///< number of empty (pure-ghost) jets
00384   mutable double _empty_area;        ///< the empty (pure-ghost/unclustered) area!
00385 
00386   // internal variables
00387   SharedPtr<PseudoJetStructureBase> _csi; ///< allows to check if _csa is still valid
00388   PseudoJet _current_reference;           ///< current reference jet
00389   mutable bool _uptodate;                 ///< true when the background computation is up-to-date
00390 
00391   /// handle warning messages
00392   static LimitedWarning _warnings;
00393   static LimitedWarning _warnings_zero_area;
00394   static LimitedWarning _warnings_preliminary;
00395 };
00396 
00397 
00398 
00399 
00400 //----------------------------------------------------------------------
00401 /// @ingroup tools_background
00402 /// \class BackgroundJetPtDensity
00403 /// Class that implements pt/area_4vector.perp() for background estimation
00404 /// <i>(this is a preliminary class)</i>.
00405 class BackgroundJetPtDensity : public FunctionOfPseudoJet<double> {
00406 public:
00407   virtual double result(const PseudoJet & jet) const {
00408     return jet.perp() / jet.area_4vector().perp();
00409   }
00410   virtual std::string description() const {return "BackgroundJetPtDensity";}
00411 };
00412 
00413 
00414 //----------------------------------------------------------------------
00415 /// @ingroup tools_background
00416 /// \class BackgroundJetScalarPtDensity
00417 /// Class that implements (scalar pt sum of jet)/(scalar area of jet)
00418 /// for background estimation <i>(this is a preliminary class)</i>.
00419 ///
00420 /// Optionally it can return a quantity based on the sum of pt^n,
00421 /// e.g. for use in subtracting fragementation function moments.
00422 class BackgroundJetScalarPtDensity : public FunctionOfPseudoJet<double> {
00423 public:
00424   /// Default constructor provides background estimation with scalar pt sum
00425   BackgroundJetScalarPtDensity() : _pt_power(1) {}
00426 
00427   /// Constructor to provide background estimation based on 
00428   /// \f$ sum_{i\in jet} p_{ti}^{n} \f$
00429   BackgroundJetScalarPtDensity(double n) : _pt_power(n) {}
00430 
00431   virtual double result(const PseudoJet & jet) const;
00432 
00433   virtual std::string description() const {return "BackgroundScalarJetPtDensity";}
00434 
00435 private:
00436   double _pt_power;
00437 };
00438 
00439 //----------------------------------------------------------------------
00440 /// @ingroup tools_background
00441 /// \class BackgroundJetPtMDensity
00442 /// Class that implements
00443 /// \f$  \frac{1}{A} \sum_{i \in jet} (\sqrt{p_{ti}^2+m^2} - p_{ti}) \f$
00444 /// for background estimation <i>(this is a preliminary class)</i>.
00445 /// 
00446 ///
00447 /// This is useful for correcting jet masses in cases where the event
00448 /// involves massive particles.
00449 class BackgroundJetPtMDensity : public FunctionOfPseudoJet<double> {
00450 public:
00451   virtual double result(const PseudoJet & jet) const {
00452     std::vector<PseudoJet> constituents = jet.constituents();
00453     double scalar_ptm = 0;
00454     for (unsigned i = 0; i < constituents.size(); i++) {
00455       scalar_ptm += constituents[i].mperp() - constituents[i].perp();
00456     }
00457     return scalar_ptm / jet.area();
00458   }
00459 
00460   virtual std::string description() const {return "BackgroundPtMDensity";}
00461 };
00462 
00463 
00464 
00465 FASTJET_END_NAMESPACE
00466 
00467 #endif  // __BACKGROUND_ESTIMATOR_HH__
00468 
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends