1 #ifndef __FASTJET_BACKGROUND_ESTIMATOR_BASE_HH__
2 #define __FASTJET_BACKGROUND_ESTIMATOR_BASE_HH__
34 #include <fastjet/ClusterSequenceAreaBase.hh>
35 #include <fastjet/FunctionOfPseudoJet.hh>
36 #include <fastjet/Selector.hh>
37 #include <fastjet/Error.hh>
40 FASTJET_BEGIN_NAMESPACE
68 virtual void set_particles(
const std::vector<PseudoJet> & particles) = 0;
77 virtual double rho()
const = 0;
83 throw Error(
"sigma() not supported for this Background Estimator");
91 virtual double rho(
const PseudoJet & jet) = 0;
96 throw Error(
"sigma(jet) not supported for this Background Estimator");
109 throw Error(
"rho_m() not supported for this Background Estimator");
117 throw Error(
"sigma_m() not supported for this Background Estimator");
122 throw Error(
"rho_m(jet) not supported for this Background Estimator");
127 throw Error(
"sigma_m(jet) not supported for this Background Estimator");
158 return _rescaling_class;
168 virtual std::string description()
const = 0;
188 void _median_and_stddev(
const std::vector<double> & quantity_vector,
191 double & stand_dev_if_gaussian,
192 bool do_fj2_calculation =
false
203 double _percentile(
const std::vector<double> & sorted_quantity_vector,
204 const double percentile,
205 const double nempty=0.0,
206 const bool do_fj2_calculation =
false)
const;
238 double a4=0) : _a0(a0), _a1(a1), _a2(a2), _a3(a3), _a4(a4) {}
241 virtual double result(
const PseudoJet & jet)
const;
243 double _a0, _a1, _a2, _a3, _a4;
250 FASTJET_END_NAMESPACE
252 #endif // __BACKGROUND_ESTIMATOR_BASE_HH__
virtual double sigma(const PseudoJet &)
get sigma, the background fluctuations per unit area, locally at the position of a given jet...
virtual double sigma_m() const
returns sigma_m, a measure of the fluctuations in the purely longitudinal, particle-mass-induced comp...
virtual double sigma() const
get sigma, the background fluctuations per unit area; must be multipled by sqrt(area) to get fluctuat...
virtual bool has_rho_m() const
Returns true if this background estimator has support for determination of rho_m. ...
BackgroundRescalingYPolynomial(double a0=1, double a1=0, double a2=0, double a3=0, double a4=0)
construct a background rescaling polynomial of the form a0 + a1*y + a2*y^2 + a3*y^3 + a4*y^4 ...
class to provide facilities for giving warnings up to some maximum number of times and to provide glo...
virtual void set_rescaling_class(const FunctionOfPseudoJet< double > *rescaling_class_in)
Set a pointer to a class that calculates the rescaling factor as a function of the jet (position)...
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 ~BackgroundEstimatorBase()
a default virtual destructor that does nothing
virtual double rho_m() const
returns rho_m, the purely longitudinal, particle-mass-induced component of the background density per...
A background rescaling that is a simple polynomial in y.
Class to contain pseudojets, including minimal information of use to jet-clustering routines...
virtual double rho_m(const PseudoJet &)
Returns rho_m locally at the jet position. As for rho(jet), it is non-const.
const FunctionOfPseudoJet< double > * rescaling_class() const
return the pointer to the jet density class
virtual double sigma_m(const PseudoJet &)
Returns sigma_m locally at the jet position. As for rho(jet), it is non-const.
virtual bool has_sigma()
returns true if this background estimator has support for determination of sigma