35 #include "fastjet/ClusterSequenceAreaBase.hh"
38 FASTJET_BEGIN_NAMESPACE
44 LimitedWarning ClusterSequenceAreaBase::_warnings;
45 LimitedWarning ClusterSequenceAreaBase::_warnings_zero_area;
46 LimitedWarning ClusterSequenceAreaBase::_warnings_empty_area;
57 double ClusterSequenceAreaBase::empty_area(
const Selector & selector)
const {
59 if (has_explicit_ghosts()) {
return 0.0;}
60 else {
return empty_area_from_jets(inclusive_jets(0.0), selector);}
69 double ClusterSequenceAreaBase::empty_area_from_jets(
70 const std::vector<PseudoJet> & all_jets,
72 _check_selector_good_for_median(selector);
74 double empty = selector.
area();
75 for (
unsigned i = 0; i < all_jets.size(); i++) {
76 if (selector.
pass(all_jets[i])) empty -= area(all_jets[i]);
81 double ClusterSequenceAreaBase::median_pt_per_unit_area(
const Selector & selector)
const {
82 return median_pt_per_unit_something(selector,
false);
85 double ClusterSequenceAreaBase::median_pt_per_unit_area_4vector(
const Selector & selector)
const {
86 return median_pt_per_unit_something(selector,
true);
94 double ClusterSequenceAreaBase::median_pt_per_unit_something(
95 const Selector & selector,
bool use_area_4vector)
const {
97 double median, sigma, mean_area;
98 get_median_rho_and_sigma(selector, use_area_4vector, median, sigma, mean_area);
108 void ClusterSequenceAreaBase::parabolic_pt_per_unit_area(
109 double & a,
double & b,
const Selector & selector,
110 double exclude_above,
bool use_area_4vector)
const {
114 _check_selector_good_for_median(selector);
118 double mean_f=0, mean_x2=0, mean_x4=0, mean_fx2=0;
120 vector<PseudoJet> incl_jets = inclusive_jets();
122 for (
unsigned i = 0; i < incl_jets.size(); i++) {
123 if (selector.
pass(incl_jets[i])) {
125 if ( use_area_4vector ) {
126 this_area = area_4vector(incl_jets[i]).perp();
128 this_area = area(incl_jets[i]);
130 double f = incl_jets[i].perp()/this_area;
131 if (exclude_above <= 0.0 || f < exclude_above) {
132 double x = incl_jets[i].rap();
double x2 = x*x;
156 b = (mean_f*mean_x2 - mean_fx2)/(mean_x2*mean_x2 - mean_x4);
157 a = mean_f - b*mean_x2;
164 void ClusterSequenceAreaBase::get_median_rho_and_sigma(
165 const Selector & selector,
bool use_area_4vector,
166 double & median,
double & sigma,
double & mean_area)
const {
168 vector<PseudoJet> incl_jets = inclusive_jets();
169 get_median_rho_and_sigma(incl_jets, selector, use_area_4vector,
170 median, sigma, mean_area,
true);
174 void ClusterSequenceAreaBase::get_median_rho_and_sigma(
175 const vector<PseudoJet> & all_jets,
176 const Selector & selector,
bool use_area_4vector,
177 double & median,
double & sigma,
double & mean_area,
178 bool all_are_incl)
const {
180 _check_jet_alg_good_for_median();
185 _check_selector_good_for_median(selector);
187 vector<double> pt_over_areas;
188 double total_area = 0.0;
189 double total_njets = 0;
191 for (
unsigned i = 0; i < all_jets.size(); i++) {
192 if (selector.
pass(all_jets[i])) {
194 if (use_area_4vector) {
195 this_area = area_4vector(all_jets[i]).perp();
197 this_area = area(all_jets[i]);
201 pt_over_areas.push_back(all_jets[i].perp()/this_area);
203 _warnings_zero_area.warn(
"ClusterSequenceAreaBase::get_median_rho_and_sigma(...): 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).");
206 total_area += this_area;
212 if (pt_over_areas.size() == 0) {
222 sort(pt_over_areas.begin(), pt_over_areas.end());
226 double posn[2] = {0.5, (1.0-0.6827)/2.0};
229 double n_empty, empty_a;
230 if (has_explicit_ghosts()) {
238 }
else if (all_are_incl) {
240 empty_a = empty_area(selector);
241 n_empty = n_empty_jets(selector);
246 empty_a = empty_area_from_jets(all_jets, selector);
247 mean_area = total_area / total_njets;
248 n_empty = empty_a / mean_area;
252 total_njets += n_empty;
253 total_area += empty_a;
258 int pt_over_areas_size = pt_over_areas.size();
259 if (n_empty < -pt_over_areas_size/4.0)
260 _warnings_empty_area.warn(
"ClusterSequenceAreaBase::get_median_rho_and_sigma(...): the estimated empty area is suspiciously large and negative and may lead to an over-estimation of rho. This may be due to (i) a rare statistical fluctuation or (ii) too small a range used to estimate the background properties.");
262 for (
int i = 0; i < 2; i++) {
263 double nj_median_pos =
264 (pt_over_areas_size-1.0 + n_empty)*posn[i] - n_empty;
265 double nj_median_ratio;
266 if (nj_median_pos >= 0 && pt_over_areas_size > 1) {
267 int int_nj_median = int(nj_median_pos);
270 if (int_nj_median+1 > pt_over_areas_size-1){
271 int_nj_median = pt_over_areas_size-2;
272 nj_median_pos = pt_over_areas_size-1;
276 pt_over_areas[int_nj_median] * (int_nj_median+1-nj_median_pos)
277 + pt_over_areas[int_nj_median+1] * (nj_median_pos - int_nj_median);
279 nj_median_ratio = 0.0;
281 res[i] = nj_median_ratio;
284 double error = res[0] - res[1];
285 mean_area = total_area / total_njets;
286 sigma = error * sqrt(mean_area);
294 vector<PseudoJet> ClusterSequenceAreaBase::subtracted_jets(
const double rho,
297 vector<PseudoJet> sub_jets;
298 vector<PseudoJet> jets_local =
sorted_by_pt(inclusive_jets(ptmin));
299 for (
unsigned i=0; i<jets_local.size(); i++) {
300 PseudoJet sub_jet = subtracted_jet(jets_local[i],rho);
301 sub_jets.push_back(sub_jet);
310 vector<PseudoJet> ClusterSequenceAreaBase::subtracted_jets(
314 double rho = median_pt_per_unit_area_4vector(selector);
315 return subtracted_jets(rho,ptmin);
321 const double rho)
const {
325 if (rho*area4vect.
perp() < jet.
perp() ) {
326 sub_jet = jet - rho*area4vect;
327 }
else { sub_jet =
PseudoJet(0.0,0.0,0.0,0.0); }
345 double rho = median_pt_per_unit_area_4vector(selector);
346 PseudoJet sub_jet = subtracted_jet(jet, rho);
352 double ClusterSequenceAreaBase::subtracted_pt(
const PseudoJet & jet,
354 bool use_area_4vector)
const {
355 if ( use_area_4vector ) {
356 PseudoJet sub_jet = subtracted_jet(jet,rho);
357 return sub_jet.
perp();
359 return jet.
perp() - rho*area(jet);
367 double ClusterSequenceAreaBase::subtracted_pt(
const PseudoJet & jet,
369 bool use_area_4vector)
const {
370 if ( use_area_4vector ) {
371 PseudoJet sub_jet = subtracted_jet(jet,selector);
372 return sub_jet.
perp();
374 double rho = median_pt_per_unit_area(selector);
375 return subtracted_pt(jet,rho,
false);
381 void ClusterSequenceAreaBase::_check_selector_good_for_median(
const Selector &selector)
const{
384 throw Error(
"ClusterSequenceAreaBase: empty area can only be computed from selectors with a finite area");
389 throw Error(
"ClusterSequenceAreaBase: empty area can only be computed from selectors that apply jet by jet");
395 void ClusterSequenceAreaBase::_check_jet_alg_good_for_median()
const {
396 if (jet_def().jet_algorithm() != kt_algorithm
397 && jet_def().jet_algorithm() != cambridge_algorithm
398 && jet_def().jet_algorithm() != cambridge_for_passive_algorithm) {
399 _warnings.warn(
"ClusterSequenceAreaBase: jet_def being used may not be suitable for estimating diffuse backgrounds (good options are kt, cam)");
405 FASTJET_END_NAMESPACE
vector< PseudoJet > sorted_by_pt(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing kt2
bool applies_jet_by_jet() const
returns true if this can be applied jet by jet
void set_user_index(const int index)
set the user_index, intended to allow the user to add simple identifying information to a particle/je...
const SharedPtr< PseudoJetStructureBase > & structure_shared_ptr() const
return a reference to the shared pointer to the PseudoJetStructureBase associated with this PseudoJet...
int cluster_hist_index() const
return the cluster_hist_index, intended to be used by clustering routines.
void set_cluster_hist_index(const int index)
set the cluster_hist_index, intended to be used by clustering routines.
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.
Class that encodes information about cuts and other selection criteria that can be applied to PseudoJ...
double area() const
returns the rapidity-phi area associated with the Selector (throws InvalidArea if the area does not m...
int user_index() const
return the user_index,
void set_structure_shared_ptr(const SharedPtr< PseudoJetStructureBase > &structure)
set the associated structure
double perp() const
returns the scalar transverse momentum
bool pass(const PseudoJet &jet) const
return true if the jet passes the selection
Class to contain pseudojets, including minimal information of use to jet-clustering routines...