31#include "fastjet/tools/Subtractor.hh"
37FASTJET_BEGIN_NAMESPACE
39const double Subtractor::_invalid_rho = -numeric_limits<double>::infinity();
41LimitedWarning Subtractor::_unused_rho_m_warning;
45Subtractor::Subtractor(
double rho) : _bge(0), _rho(rho) {
46 if (
_rho<0.0)
throw Error(
"Subtractor(rho) was passed a negative rho value; rho should be >= 0");
53 if (
_rho<0.0)
throw Error(
"Subtractor(rho, rho_m) was passed a negative rho value; rho should be >= 0");
54 if (rho_m<0.0)
throw Error(
"Subtractor(rho, rho_m) was passed a negative rho_m value; rho_m should be >= 0");
74 throw Error(
"Subtractor::result(...): Trying to subtract a jet without area support");
84 vector<PseudoJet> constits_unknown, constits_known;
88 vector<PseudoJet> constits_known_lv, constits_known_pu;
96 known_lv = (constits_known_lv.size()!=0)
97 ? SelectorIdentity().
sum(constits_known_lv) : 0.0*jet;
98 known_pu = (constits_known_pu.size()!=0)
99 ? SelectorIdentity().
sum(constits_known_pu) : 0.0*jet;
100 if (constits_unknown.size()==0){
104 return subtracted_jet;
120 if (to_subtract.
pt2() < jet.
pt2() ) {
123 subtracted_jet -= to_subtract;
128 return subtracted_jet;
133 if (subtracted_jet.
pt2() < known_lv.
pt2()){
135 return subtracted_jet;
146 subtracted_jet.
phi(),
150 return subtracted_jet;
156 string desc =
"Subtractor that uses the following background estimator to determine rho: "+
_bge->
description();
157 if (
use_rho_m()) desc +=
"; including the rho_m correction";
158 if (
safe_mass()) desc +=
"; including mass safety tests";
165 ostr <<
"Subtractor that uses a fixed value of rho = " <<
_rho;
166 if (
use_rho_m()) ostr <<
" and rho_m = " << _rho_m;
169 return "Uninitialised subtractor";
182 rho = bg_estimate.
rho();
186 throw Error(
"Subtractor::_amount_to_subtract(...): default Subtractor does not have any information about the background, needed to perform the subtraction");
192 double const rho_m_warning_threshold = 1e-5;
200 if (!bg_estimate.
has_rho_m())
throw Error(
"Subtractor::_amount_to_subtract(...): requested subtraction with rho_m from a background estimator, but the estimator does not have rho_m support");
201 rho_m = bg_estimate.
rho_m();
205 throw Error(
"Subtractor::_amount_to_subtract(...): default Subtractor does not have any information about the background rho_m, needed to perform the rho_m subtraction");
207 to_subtract += rho_m *
PseudoJet(0.0, 0.0, area.pz(), area.E());
210 bg_estimate.
rho_m() > rho_m_warning_threshold * rho) {
211 _unused_rho_m_warning.
warn(
"Subtractor::_amount_to_subtract(...): Background estimator indicates non-zero rho_m, but use_rho_m()==false in subtractor; consider calling set_use_rho_m(true) to include the rho_m information");
/// a class that holds the result of the calculation
bool has_rho_m() const
true if this background estimate has a determination of rho_m.
double rho() const
background density per unit area
double rho_m() const
purely longitudinal (particle-mass-induced) component of the background density per unit area
virtual std::string description() const =0
returns a textual description of the background estimator
virtual BackgroundEstimate estimate() const =0
get the full set of background properties
base class corresponding to errors that can be thrown by FastJet
void warn(const char *warning)
outputs a warning to standard error (or the user's default warning stream if set)
Class to contain pseudojets, including minimal information of use to jet-clustering routines.
virtual std::vector< PseudoJet > constituents() const
retrieve the constituents.
double rap() const
returns the rapidity or some large value when the rapidity is infinite
void reset_momentum(double px, double py, double pz, double E)
reset the 4-momentum according to the supplied components but leave all other information (indices,...
double phi() const
returns phi (in the range 0..2pi)
virtual bool has_area() const
check if it has a defined area
double m2() const
returns the squared invariant mass // like CLHEP
double pt() const
returns the scalar transverse momentum
double m() const
returns the invariant mass (If m2() is negative then -sqrt(-m2()) is returned, as in CLHEP)
virtual PseudoJet area_4vector() const
return the jet 4-vector area.
double pt2() const
returns the squared transverse momentum
Class that encodes information about cuts and other selection criteria that can be applied to PseudoJ...
const SharedPtr< SelectorWorker > & worker() const
returns a (reference to) the underlying worker's shared pointer
PseudoJet sum(const std::vector< PseudoJet > &jets) const
Return the 4-vector sum of the objects that pass the selection.
std::string description() const
returns a textual description of the selector
void sift(const std::vector< PseudoJet > &jets, std::vector< PseudoJet > &jets_that_pass, std::vector< PseudoJet > &jets_that_fail) const
sift the input jets into two vectors – those that pass the selector and those that do not
bool use_rho_m() const
returns whether or not the rho_m component is used
PseudoJet _amount_to_subtract(const PseudoJet &jet) const
compute the 4-vector that should be subtracted from the given jet
Subtractor()
default constructor
static const double _invalid_rho
a value of rho that is used as a default to label that the stored rho is not valid for subtraction.
double _rho
the fixed value of rho and/or rho_m to use if the user has selected that option
virtual std::string description() const
class description
bool _safe_mass
ensures that the subtracted mass is +ve
bool safe_mass() const
returns whether or not safety tests on the mass are included
void set_use_rho_m(bool use_rho_m_in=true)
when 'use_rho_m' is true, include in the subtraction the correction from rho_m, the purely longitudin...
const BackgroundEstimatorBase * _bge
the tool used to estimate the background if has to be mutable in case its underlying selector takes a...
virtual PseudoJet result(const PseudoJet &jet) const
returns a jet that's subtracted
bool _use_rho_m
include the rho_m correction
Selector _sel_leading_vertex
amongst the particles with a known vertex origin, select those coming from the leading vertex
Selector _sel_known_vertex
selects the particles with a known vertex origin
void set_defaults()
reset all parameters to default values
PseudoJet PtYPhiM(double pt, double y, double phi, double m)
return a pseudojet with the given pt, y, phi and mass