31 #include "fastjet/tools/JetMedianBackgroundEstimator.hh"
32 #include <fastjet/ClusterSequenceArea.hh>
33 #include <fastjet/ClusterSequenceStructure.hh>
37 FASTJET_BEGIN_NAMESPACE
41 double BackgroundJetScalarPtDensity::result(
const PseudoJet & jet)
const {
46 for (
unsigned i = 0; i < constituents.size(); i++) {
47 scalar_pt += pow(constituents[i].perp(), _pt_power);
49 return scalar_pt / jet.
area();
53 std::string BackgroundJetScalarPtDensity::description()
const {
55 oss <<
"BackgroundScalarJetPtDensity";
56 if (_pt_power != 1.0) oss <<
" with pt_power = " << _pt_power;
62 double BackgroundRescalingYPolynomial::result(
const PseudoJet & jet)
const {
65 double rescaling = _a0 + _a1*y + _a2*y2 + _a3*y2*y + _a4*y2*y2;
72 LimitedWarning JetMedianBackgroundEstimator::_warnings_preliminary;
84 JetMedianBackgroundEstimator::JetMedianBackgroundEstimator(
const Selector &rho_range,
87 : _rho_range(rho_range), _jet_def(jet_def), _area_def(area_def){
93 _check_jet_alg_good_for_median();
120 throw Error(
"JetMedianBackgroundEstimator::set_particles can only be called if you set the jet (and area) definition explicitly through the class constructor");
166 _check_jet_alg_good_for_median();
170 throw Error(
"JetMedianBackgroundEstimator: either an area with explicit ghosts (recommended) or a Selector with finite area is needed (to allow for the computation of the empty area)");
187 throw Error(
"JetMedianBackgroundEstimator::JetMedianBackgroundEstimator: At least one jet is needed to compute the background properties");
192 if (! (jets[0].has_associated_cluster_sequence()) && (jets[0].has_area()))
193 throw Error(
"JetMedianBackgroundEstimator::JetMedianBackgroundEstimator: the jets used to estimate the background properties must be associated with a valid ClusterSequenceAreaBase");
195 _csi = jets[0].structure_shared_ptr();
199 for (
unsigned int i=1;i<jets.size(); i++){
200 if (! jets[i].has_associated_cluster_sequence())
201 throw Error(
"JetMedianBackgroundEstimator::set_jets(...): the jets used to estimate the background properties must be associated with a valid ClusterSequenceAreaBase");
203 if (jets[i].structure_shared_ptr().get() != _csi.get())
204 throw Error(
"JetMedianBackgroundEstimator::set_jets(...): all the jets used to estimate the background properties must share the same ClusterSequence");
208 _check_jet_alg_good_for_median();
212 throw Error(
"JetMedianBackgroundEstimator: either an area with explicit ghosts (recommended) or a Selector with finite area is needed (to allow for the computation of the empty area)");
217 _included_jets = jets;
231 throw Error(
"The background estimation is obtained from a selector that takes a reference jet. rho(PseudoJet) should be used in that case");
232 _recompute_if_needed();
239 throw Error(
"The background estimation is obtained from a selector that takes a reference jet. rho(PseudoJet) should be used in that case");
240 _recompute_if_needed();
251 _recompute_if_needed(jet);
252 double our_rho = _rho;
253 if (_rescaling_class != 0) {
254 our_rho *= (*_rescaling_class)(jet);
266 _recompute_if_needed(jet);
267 double our_sigma = _sigma;
268 if (_rescaling_class != 0) {
269 our_sigma *= (*_rescaling_class)(jet);
279 throw Error(
"JetMediamBackgroundEstimator: rho_m requested but rho_m calculation is disabled (either eplicitly or due to the presence of a jet density class).");
282 throw Error(
"The background estimation is obtained from a selector that takes a reference jet. rho(PseudoJet) should be used in that case");
283 _recompute_if_needed();
294 throw Error(
"JetMediamBackgroundEstimator: sigma_m requested but rho_m/sigma_m calculation is disabled (either explicitly or due to the presence of a jet density class).");
297 throw Error(
"The background estimation is obtained from a selector that takes a reference jet. rho(PseudoJet) should be used in that case");
298 _recompute_if_needed();
306 _recompute_if_needed(jet);
307 double our_rho = _rho_m;
308 if (_rescaling_class != 0) {
309 our_rho *= (*_rescaling_class)(jet);
319 _recompute_if_needed(jet);
320 double our_sigma = _sigma_m;
321 if (_rescaling_class != 0) {
322 our_sigma *= (*_rescaling_class)(jet);
338 _enable_rho_m =
true;
342 _rho_m = _sigma_m = 0.0;
343 _n_jets_used = _n_empty_jets = 0;
344 _empty_area = _mean_area = 0.0;
346 _included_jets.clear();
348 _jet_density_class = 0;
349 _rescaling_class = 0;
359 _warnings_preliminary.
warn(
"JetMedianBackgroundEstimator::set_jet_density_class: density classes are still preliminary in FastJet 3.1. Their interface may differ in future releases (without guaranteeing backward compatibility). Note that since FastJet 3.1, rho_m and sigma_m are accessible direclty in JetMedianBackgroundEstimator and GridMedianBackgroundEstimator(with no need for a density class).");
360 _jet_density_class = jet_density_class_in;
371 desc <<
"JetMedianBackgroundEstimator, using " << _jet_def.
description()
373 <<
" and selecting jets with " << _rho_range.
description();
385 void JetMedianBackgroundEstimator::_recompute_if_needed(
const PseudoJet &jet){
390 if (jet == _current_reference)
return;
398 _recompute_if_needed();
402 void JetMedianBackgroundEstimator::_compute()
const {
408 vector<double> vector_for_median_pt;
409 vector<double> vector_for_median_dt;
410 double total_area = 0.0;
414 vector<PseudoJet> selected_jets = _rho_range(_included_jets);
417 double median_input_pt, median_input_dt=0.0;
418 BackgroundJetPtMDensity m_density;
420 for (
unsigned i = 0; i < selected_jets.size(); i++) {
421 const PseudoJet & current_jet = selected_jets[i];
423 double this_area = (_use_area_4vector) ? current_jet.area_4vector().perp() : current_jet.area();
427 if (_jet_density_class == 0) {
428 median_input_pt = current_jet.perp()/this_area;
430 median_input_pt = (*_jet_density_class)(current_jet);
437 median_input_dt = m_density(current_jet);
440 if (_rescaling_class != 0) {
441 double resc = (*_rescaling_class)(current_jet);;
442 median_input_pt /= resc;
443 median_input_dt /= resc;
447 vector_for_median_pt.push_back(median_input_pt);
449 vector_for_median_dt.push_back(median_input_dt);
451 total_area += this_area;
454 _warnings_zero_area.
warn(
"JetMedianBackgroundEstimator::_compute(...): discarded jet with zero area. Zero-area jets may be due to (i) too large a ghost area (ii) a jet being outside the ghost range (iii) the computation not being done using an appropriate algorithm (kt;C/A).");
459 if (vector_for_median_pt.size() == 0) {
469 const ClusterSequenceAreaBase * csab = (
dynamic_cast<ClusterSequenceStructure*
>(_csi()))->validated_csab();
470 if (csab->has_explicit_ghosts()) {
474 _empty_area = csab->empty_area(_rho_range);
475 _n_empty_jets = csab->n_empty_jets(_rho_range);
478 double total_njets = _n_jets_used + _n_empty_jets;
479 total_area += _empty_area;
486 _mean_area = total_area / total_njets;
487 _sigma = stand_dev * sqrt(_mean_area);
493 _sigma_m = stand_dev * sqrt(_mean_area);
504 void JetMedianBackgroundEstimator::_check_csa_alive()
const{
505 ClusterSequenceStructure* csa =
dynamic_cast<ClusterSequenceStructure*
>(_csi());
507 throw Error(
"JetMedianBackgroundEstimator: there is no cluster sequence associated with the JetMedianBackgroundEstimator");
509 if (! dynamic_cast<ClusterSequenceStructure*>(_csi())->has_associated_cluster_sequence())
510 throw Error(
"JetMedianBackgroundEstimator: modifications are no longer possible as the underlying ClusterSequence has gone out of scope");
517 void JetMedianBackgroundEstimator::_check_jet_alg_good_for_median()
const{
518 const JetDefinition * jet_def = &_jet_def;
523 const ClusterSequence * cs =
dynamic_cast<ClusterSequenceStructure*
>(_csi())->validated_cs();
524 jet_def = &(cs->jet_def());
530 _warnings.
warn(
"JetMedianBackgroundEstimator: jet_def being used may not be suitable for estimating diffuse backgrounds (good alternatives are kt, cam)");
536 FASTJET_END_NAMESPACE
double rap() const
returns the rapidity or some large value when the rapidity is infinite
a version of cambridge with a special distance measure for particles whose pt is < extra_param(); thi...
General class for user to obtain ClusterSequence with additional area information.
the value for the jet algorithm in a JetDefinition for which no algorithm has yet been defined ...
JetAlgorithm jet_algorithm() const
return information about the definition...
bool takes_reference() const
returns true if this can be applied jet by jet
std::vector< PseudoJet > inclusive_jets(const double ptmin=0.0) const
return a vector of all jets (in the sense of the inclusive algorithm) with pt >= ptmin.
void _median_and_stddev(const std::vector< double > &quantity_vector, double n_empty_jets, double &median, double &stand_dev_if_gaussian, bool do_fj2_calculation=false) const
given a quantity in a vector (e.g.
virtual double area() const
return the jet (scalar) area.
the longitudinally invariant kt algorithm
class that holds a generic area definition
virtual bool has_explicit_ghosts() const
returns true if ghosts are explicitly included within jets for this ClusterSequence; ...
const SharedPtr< PseudoJetStructureBase > & structure_shared_ptr() const
retrieve a shared pointer to the wrapper to this ClusterSequence
Contains any information related to the clustering that should be directly accessible to PseudoJet...
void warn(const char *warning)
outputs a warning to standard error (or the user's default warning stream if set) ...
class to provide facilities for giving warnings up to some maximum number of times and to provide glo...
std::string description() const
return a textual description of the current jet definition
base class that sets interface for extensions of ClusterSequence that provide information about the a...
std::string description() const
return a description of the current area definition
base class corresponding to errors that can be thrown by FastJet
bool has_finite_area() const
returns true if it has a meaningful and finite area (i.e.
std::string description() const
returns a textual description of the selector
const Selector & set_reference(const PseudoJet &reference)
set the reference jet for this Selector
Class that encodes information about cuts and other selection criteria that can be applied to PseudoJ...
virtual const ClusterSequenceAreaBase * validated_csab() const
if the jet has valid area information then return a pointer to the associated ClusterSequenceAreaBase...
the longitudinally invariant variant of the cambridge algorithm (aka Aachen algoithm).
Selector SelectorIsPureGhost()
select objects that are (or are only made of) ghosts.
void delete_self_when_unused()
by calling this routine you tell the ClusterSequence to delete itself when all the Pseudojets associa...
virtual std::vector< PseudoJet > constituents() const
retrieve the constituents.
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