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]);
 
   86 double ClusterSequenceAreaBase::median_pt_per_unit_area(
const Selector & selector)
 const {
 
   87   return _median_pt_per_unit_area(selector);
 
   91 double ClusterSequenceAreaBase::_median_pt_per_unit_area(
const Selector & selector)
 const {
 
   92   return _median_pt_per_unit_something(selector,
false);
 
  102 double ClusterSequenceAreaBase::median_pt_per_unit_area_4vector(
const Selector & selector)
 const {
 
  103   return _median_pt_per_unit_area_4vector(selector);
 
  107 double ClusterSequenceAreaBase::_median_pt_per_unit_area_4vector(
const Selector & selector)
 const {
 
  108   return _median_pt_per_unit_something(selector,
true);
 
  118 double ClusterSequenceAreaBase::median_pt_per_unit_something(
 
  119                 const Selector & selector, 
bool use_area_4vector)
 const {
 
  120   return _median_pt_per_unit_something(selector, use_area_4vector);
 
  126 double ClusterSequenceAreaBase::_median_pt_per_unit_something(
 
  127                 const Selector & selector, 
bool use_area_4vector)
 const {
 
  128   double median, sigma, mean_area;
 
  129   _get_median_rho_and_sigma(selector, use_area_4vector, median, sigma, mean_area);
 
  138 void ClusterSequenceAreaBase::parabolic_pt_per_unit_area(
 
  139        double & a, 
double & b, 
const Selector & selector, 
 
  140        double exclude_above, 
bool use_area_4vector)
 const {
 
  141   return _parabolic_pt_per_unit_area(a, b, selector, exclude_above, use_area_4vector);
 
  144 void ClusterSequenceAreaBase::_parabolic_pt_per_unit_area(
 
  145        double & a, 
double & b, 
const Selector & selector, 
 
  146        double exclude_above, 
bool use_area_4vector)
 const {
 
  150   _check_selector_good_for_median(selector);
 
  154   double mean_f=0, mean_x2=0, mean_x4=0, mean_fx2=0; 
 
  156   vector<PseudoJet> incl_jets = inclusive_jets();
 
  158   for (
unsigned i = 0; i < incl_jets.size(); i++) {
 
  159     if (selector.
pass(incl_jets[i])) {
 
  161       if ( use_area_4vector ) {
 
  162           this_area = area_4vector(incl_jets[i]).perp();     
 
  164           this_area = area(incl_jets[i]);
 
  166       double f = incl_jets[i].perp()/this_area;
 
  167       if (exclude_above <= 0.0 || f < exclude_above) {
 
  168         double x = incl_jets[i].rap(); 
double x2 = x*x;
 
  192     b = (mean_f*mean_x2 - mean_fx2)/(mean_x2*mean_x2 - mean_x4);
 
  193     a = mean_f - b*mean_x2;
 
  200 void ClusterSequenceAreaBase::get_median_rho_and_sigma(
 
  201             const Selector & selector, 
bool use_area_4vector,
 
  202             double & median, 
double & sigma, 
double & mean_area)
 const {
 
  203   _get_median_rho_and_sigma(selector, use_area_4vector, median, sigma, mean_area);
 
  206 void ClusterSequenceAreaBase::_get_median_rho_and_sigma(
 
  207             const Selector & selector, 
bool use_area_4vector,
 
  208             double & median, 
double & sigma, 
double & mean_area)
 const {
 
  210   vector<PseudoJet> incl_jets = inclusive_jets();
 
  211   _get_median_rho_and_sigma(incl_jets, selector, use_area_4vector,
 
  212                             median, sigma, mean_area, 
true);
 
  215 void ClusterSequenceAreaBase::get_median_rho_and_sigma(
 
  216             const vector<PseudoJet> & all_jets,
 
  217             const Selector & selector, 
bool use_area_4vector,
 
  218             double & median, 
double & sigma, 
double & mean_area,
 
  219             bool all_are_incl)
 const {
 
  220   _get_median_rho_and_sigma(all_jets, selector, use_area_4vector,
 
  221                             median, sigma, mean_area, all_are_incl);
 
  224 void ClusterSequenceAreaBase::_get_median_rho_and_sigma(
 
  225             const vector<PseudoJet> & all_jets,
 
  226             const Selector & selector, 
bool use_area_4vector,
 
  227             double & median, 
double & sigma, 
double & mean_area,
 
  228             bool all_are_incl)
 const {
 
  230   _check_jet_alg_good_for_median();
 
  235   _check_selector_good_for_median(selector);
 
  237   vector<double> pt_over_areas;
 
  238   double total_area  = 0.0;
 
  239   double total_njets = 0;
 
  241   for (
unsigned i = 0; i < all_jets.size(); i++) {
 
  242     if (selector.pass(all_jets[i])) {
 
  244       if (use_area_4vector) {
 
  245           this_area = area_4vector(all_jets[i]).perp();
 
  247           this_area = area(all_jets[i]);
 
  251         pt_over_areas.push_back(all_jets[i].perp()/this_area);
 
  253         _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).");
 
  256       total_area  += this_area;
 
  262   if (pt_over_areas.size() == 0) {
 
  272   sort(pt_over_areas.begin(), pt_over_areas.end());
 
  276   double posn[2] = {0.5, (1.0-0.6827)/2.0};
 
  279   double n_empty, empty_a;
 
  280   if (has_explicit_ghosts()) {
 
  288   } 
else if (all_are_incl) {
 
  290     empty_a = empty_area(selector);
 
  291     n_empty = n_empty_jets(selector);
 
  296     empty_a = empty_area_from_jets(all_jets, selector);
 
  297     mean_area = total_area / total_njets; 
 
  298     n_empty   = empty_a / mean_area;
 
  302   total_njets += n_empty;
 
  303   total_area  += empty_a;
 
  308   int pt_over_areas_size = pt_over_areas.size();
 
  309   if (n_empty < -pt_over_areas_size/4.0)
 
  310     _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.");
 
  312   for (
int i = 0; i < 2; i++) {
 
  313     double nj_median_pos = 
 
  314       (pt_over_areas_size-1.0 + n_empty)*posn[i] - n_empty;
 
  315     double nj_median_ratio;
 
  316     if (nj_median_pos >= 0 && pt_over_areas_size > 1) {
 
  317       int int_nj_median = int(nj_median_pos);
 
  320       if (int_nj_median+1 > pt_over_areas_size-1){
 
  321         int_nj_median = pt_over_areas_size-2;
 
  322         nj_median_pos = pt_over_areas_size-1;
 
  326         pt_over_areas[int_nj_median] * (int_nj_median+1-nj_median_pos)
 
  327         + pt_over_areas[int_nj_median+1] * (nj_median_pos - int_nj_median);
 
  329       nj_median_ratio = 0.0;
 
  331     res[i] = nj_median_ratio;
 
  334   double error  = res[0] - res[1];
 
  335   mean_area = total_area / total_njets;
 
  336   sigma  = error * sqrt(mean_area);
 
  344 vector<PseudoJet> ClusterSequenceAreaBase::subtracted_jets(
const double rho,
 
  347   return _subtracted_jets(rho,ptmin);
 
  350 vector<PseudoJet> ClusterSequenceAreaBase::_subtracted_jets(
const double rho,
 
  353   vector<PseudoJet> sub_jets;
 
  354   vector<PseudoJet> jets_local = 
sorted_by_pt(inclusive_jets(ptmin));
 
  355   for (
unsigned i=0; i<jets_local.size(); i++) {
 
  356      PseudoJet sub_jet = _subtracted_jet(jets_local[i],rho);
 
  357      sub_jets.push_back(sub_jet);
 
  366 vector<PseudoJet> ClusterSequenceAreaBase::subtracted_jets(
 
  370   double rho = _median_pt_per_unit_area_4vector(selector);
 
  371   return _subtracted_jets(rho,ptmin);
 
  377                                                   const double rho)
 const {
 
  378   return _subtracted_jet(jet, rho);
 
  382                                                    const double rho)
 const {
 
  386   if (rho*area4vect.
perp() < jet.
perp() ) { 
 
  387     sub_jet = jet - rho*area4vect;
 
  388   } 
else { sub_jet = PseudoJet(0.0,0.0,0.0,0.0); }
 
  406   return _subtracted_jet(jet, selector);
 
  411   double rho = _median_pt_per_unit_area_4vector(selector);
 
  412   PseudoJet sub_jet = _subtracted_jet(jet, rho);
 
  418 double ClusterSequenceAreaBase::subtracted_pt(
const PseudoJet & jet,
 
  420                                               bool use_area_4vector)
 const {
 
  421   return _subtracted_pt(jet, rho, use_area_4vector);
 
  424 double ClusterSequenceAreaBase::_subtracted_pt(
const PseudoJet & jet,
 
  426                                               bool use_area_4vector)
 const {
 
  427   if ( use_area_4vector ) { 
 
  428      PseudoJet sub_jet = _subtracted_jet(jet,rho);
 
  429      return sub_jet.
perp();
 
  431      return jet.
perp() - rho*area(jet);
 
  439 double ClusterSequenceAreaBase::subtracted_pt(
const PseudoJet & jet,
 
  441                                               bool use_area_4vector)
 const {
 
  442   if ( use_area_4vector ) { 
 
  443      PseudoJet sub_jet = _subtracted_jet(jet,selector);
 
  444      return sub_jet.
perp();
 
  446      double rho = _median_pt_per_unit_area(selector);
 
  447      return _subtracted_pt(jet,rho,
false);
 
  453 void ClusterSequenceAreaBase::_check_selector_good_for_median(
const Selector &selector)
 const{
 
  456     throw Error(
"ClusterSequenceAreaBase: empty area can only be computed from selectors with a finite area");
 
  461     throw Error(
"ClusterSequenceAreaBase: empty area can only be computed from selectors that apply jet by jet");
 
  467 void ClusterSequenceAreaBase::_check_jet_alg_good_for_median()
 const {
 
  468   if (jet_def().jet_algorithm() != kt_algorithm
 
  469       && jet_def().jet_algorithm() != cambridge_algorithm
 
  470       && jet_def().jet_algorithm() !=  cambridge_for_passive_algorithm) {
 
  471     _warnings.warn(
"ClusterSequenceAreaBase: jet_def being used may not be suitable for estimating diffuse backgrounds (good options are kt, cam)");
 
  477 FASTJET_END_NAMESPACE