1 #ifndef __FASTJET_BACKGROUND_ESTIMATOR_HH__     2 #define __FASTJET_BACKGROUND_ESTIMATOR_HH__    34 #include <fastjet/ClusterSequenceAreaBase.hh>    35 #include <fastjet/AreaDefinition.hh>    36 #include <fastjet/FunctionOfPseudoJet.hh>    37 #include <fastjet/Selector.hh>    38 #include <fastjet/tools/BackgroundEstimatorBase.hh>    41 FASTJET_BEGIN_NAMESPACE     
   126       _enable_rho_m(true){ reset(); }
   141   virtual void set_particles(
const std::vector<PseudoJet> & particles);
   165   void set_jets(
const std::vector<PseudoJet> &jets);
   169     _rho_range = rho_range_selector;
   188   double sigma() 
const;
   215   virtual double rho_m() 
const;
   221   virtual double sigma_m() 
const;
   224   virtual double rho_m(
const PseudoJet & );
   227   virtual double sigma_m(
const PseudoJet & );
   237   virtual bool has_rho_m()
 const {
return _enable_rho_m && (_jet_density_class == 0);}
   248       throw Error(
"JetMedianBackgroundEstimator::mean_area(): one may not retrieve information about the last call to rho() or sigma() when the configuration has changed in the meantime.");
   258       throw Error(
"JetMedianBackgroundEstimator::n_jets_used(): one may not retrieve information about the last call to rho() or sigma() when the configuration has changed in the meantime.");
   266     if (!_uptodate) 
throw Error(
"JetMedianBackgroundEstimator::n_jets_used(): one may not retrieve information about the last call to rho() or sigma() when the configuration has changed in the meantime.");
   268     std::vector<PseudoJet> tmp_jets = _rho_range(_included_jets);
   269     std::vector<PseudoJet> used_jets;
   270     for (
unsigned int i=0; i<tmp_jets.size(); i++){
   271       if (tmp_jets[i].area()>0) used_jets.push_back(tmp_jets[i]);
   292       throw Error(
"JetMedianBackgroundEstimator::empty_area(): one may not retrieve information about the last call to rho() or sigma() when the configuration has changed in the meantime.");
   311       throw Error(
"JetMedianBackgroundEstimator::n_empty_jets(): one may not retrieve information about the last call to rho() or sigma() when the configuration has changed in the meantime.");
   313     return _n_empty_jets;
   339     _use_area_4vector = use_it;
   351     _provide_fj2_sigma = provide_fj2_sigma;
   366     return _jet_density_class;
   378     BackgroundEstimatorBase::set_rescaling_class(rescaling_class_in);
   389   std::string description() 
const;
   397   void _compute() 
const;
   401   void _recompute_if_needed()
 const {
   402     if (!_uptodate) _compute();
   413   void _recompute_if_needed(
const PseudoJet &jet);
   417   void _check_csa_alive() 
const;
   422   void _check_jet_alg_good_for_median() 
const;
   428   std::vector<PseudoJet> _included_jets; 
   431   bool _use_area_4vector;
   432   bool _provide_fj2_sigma;
   439   mutable double _sigma;             
   440   mutable double _rho_m;             
   441   mutable double _sigma_m;           
   442   mutable double _mean_area;         
   443   mutable unsigned int _n_jets_used; 
   444   mutable double _n_empty_jets;      
   445   mutable double _empty_area;        
   450   mutable bool _uptodate;                 
   471   virtual std::string 
description()
 const {
return "BackgroundJetPtDensity";}
   492   virtual double result(
const PseudoJet & jet) 
const;
   494   virtual std::string description() 
const;
   513     std::vector<PseudoJet> constituents = jet.
constituents();
   514     double scalar_ptm = 0;
   515     for (
unsigned i = 0; i < constituents.size(); i++) {
   516       scalar_ptm += constituents[i].mperp() - constituents[i].perp();
   518     return scalar_ptm / jet.
area();
   521   virtual std::string 
description()
 const {
return "BackgroundPtMDensity";}
   526 FASTJET_END_NAMESPACE
   528 #endif  // __BACKGROUND_ESTIMATOR_HH__ 
BackgroundJetScalarPtDensity(double n)
Constructor to provide background estimation based on . 
virtual PseudoJet area_4vector() const
return the jet 4-vector area. 
BackgroundJetScalarPtDensity()
Default constructor provides background estimation with scalar pt sum. 
Class that implements pt/area_4vector.perp() for background estimation (this is a preliminary class)...
virtual std::string description() const
returns a description of the function (an empty string by default) 
virtual double result(const PseudoJet &jet) const
the action of the function this has to be overloaded in derived classes 
Class that implements  for background estimation (this is a preliminary class). 
virtual double area() const
return the jet (scalar) area. 
virtual double result(const PseudoJet &jet) const
the action of the function this has to be overloaded in derived classes 
class that holds a generic area definition 
class to provide facilities for giving warnings up to some maximum number of times and to provide glo...
virtual std::vector< PseudoJet > constituents() const
retrieve the constituents. 
double perp() const
returns the scalar transverse momentum 
base class that sets interface for extensions of ClusterSequence that provide information about the a...
base class corresponding to errors that can be thrown by FastJet 
Abstract base class that provides the basic interface for classes that estimate levels of background ...
virtual std::string description() const
returns a description of the function (an empty string by default) 
Class that implements (scalar pt sum of jet)/(scalar area of jet) for background estimation (this is ...
an implementation of C++0x shared pointers (or boost's) 
Class that encodes information about cuts and other selection criteria that can be applied to PseudoJ...
Class to contain pseudojets, including minimal information of use to jet-clustering routines...
class that is intended to hold a full definition of the jet clusterer