ClusterSequenceAreaBase Class Reference

base class that sets interface for extensions of ClusterSequence that provide information about the area of each jet; More...

#include <ClusterSequenceAreaBase.hh>

Inheritance diagram for ClusterSequenceAreaBase:
Inheritance graph
[legend]
Collaboration diagram for ClusterSequenceAreaBase:
Collaboration graph
[legend]

List of all members.

Public Member Functions

template<class L >
 ClusterSequenceAreaBase (const std::vector< L > &pseudojets, const JetDefinition &jet_def, const bool &writeout_combinations=false)
 a constructor which just carries out the construction of the parent class
 ClusterSequenceAreaBase ()
 default constructor
virtual ~ClusterSequenceAreaBase ()
 destructor
virtual double area (const PseudoJet &) const
 return the area associated with the given jet; this base class returns 0.
virtual double area_error (const PseudoJet &) const
 return the error (uncertainty) associated with the determination of the area of this jet; this base class returns 0.
virtual PseudoJet area_4vector (const PseudoJet &) const
 return a PseudoJet whose 4-vector is defined by the following integral
virtual bool is_pure_ghost (const PseudoJet &) const
 true if a jet is made exclusively of ghosts
virtual bool has_explicit_ghosts () const
 returns true if ghosts are explicitly included within jets for this ClusterSequence;
virtual double empty_area (const RangeDefinition &range) const
 return the total area, within range, that is free of jets, in general based on the inclusive jets
double empty_area_from_jets (const std::vector< PseudoJet > &all_jets, const RangeDefinition &range) const
 return the total area, within range, that is free of jets, based on the supplied all_jets
virtual double n_empty_jets (const RangeDefinition &range) const
 return something similar to the number of pure ghost jets in the given range in an active area case.
double median_pt_per_unit_area (const RangeDefinition &range) const
 the median of (pt/area) for jets contained within range, making use also of the info on n_empty_jets
double median_pt_per_unit_area_4vector (const RangeDefinition &range) const
 the median of (pt/area_4vector) for jets contained within making use also of the info on n_empty_jets
double median_pt_per_unit_something (const RangeDefinition &range, bool use_area_4vector) const
 the function that does the work for median_pt_per_unit_area and median_pt_per_unit_area_4vector:

  • something_is_area_4vect = false -> use plain area
  • something_is_area_4vect = true -> use 4-vector area

virtual void get_median_rho_and_sigma (const RangeDefinition &range, bool use_area_4vector, double &median, double &sigma, double &mean_area) const
 using jets withing range (and with 4-vector areas if use_area_4vector), calculate the median pt/area, as well as an "error" (uncertainty), which is defined as the 1-sigma half-width of the distribution of pt/A, obtained by looking for the point below which we have (1-0.6827)/2 of the jets (including empty jets).
virtual void get_median_rho_and_sigma (const std::vector< PseudoJet > &all_jets, const RangeDefinition &range, bool use_area_4vector, double &median, double &sigma, double &mean_area, bool all_are_inclusive=false) const
 a more advanced version of get_median_rho_and_sigma, which allows one to use any "view" of the event containing all jets (so that, e.g.
virtual void get_median_rho_and_sigma (const RangeDefinition &range, bool use_area_4vector, double &median, double &sigma) const
 same as the full version of get_median_rho_and_error, but without access to the mean_area
virtual void parabolic_pt_per_unit_area (double &a, double &b, const RangeDefinition &range, double exclude_above=-1.0, bool use_area_4vector=false) const
 fits a form pt_per_unit_area(y) = a + b*y^2 in the range "range".
std::vector< PseudoJetsubtracted_jets (const double rho, const double ptmin=0.0) const
 return a vector of all subtracted jets, using area_4vector, given rho.
std::vector< PseudoJetsubtracted_jets (const RangeDefinition &range, const double ptmin=0.0) const
 return a vector of subtracted jets, using area_4vector.
PseudoJet subtracted_jet (const PseudoJet &jet, const double rho) const
 return a subtracted jet, using area_4vector, given rho
PseudoJet subtracted_jet (const PseudoJet &jet, const RangeDefinition &range) const
 return a subtracted jet, using area_4vector; note that this is potentially inefficient if repeatedly used for many different jets, because rho will be recalculated each time around.
double subtracted_pt (const PseudoJet &jet, const double rho, bool use_area_4vector=false) const
 return the subtracted pt, given rho
double subtracted_pt (const PseudoJet &jet, const RangeDefinition &range, bool use_area_4vector=false) const
 return the subtracted pt; note that this is potentially inefficient if repeatedly used for many different jets, because rho will be recalculated each time around.

Private Member Functions

void _check_jet_alg_good_for_median () const
 check the jet algorithm is suitable (and if not issue a warning)

Static Private Attributes

static LimitedWarning _warnings
 handle warning messages
static LimitedWarning _warnings_zero_area

Detailed Description

base class that sets interface for extensions of ClusterSequence that provide information about the area of each jet;

the virtual functions here all return 0, since no area determination is implemented.

Definition at line 45 of file ClusterSequenceAreaBase.hh.


Constructor & Destructor Documentation

template<class L >
ClusterSequenceAreaBase::ClusterSequenceAreaBase ( const std::vector< L > &  pseudojets,
const JetDefinition jet_def,
const bool &  writeout_combinations = false 
) [inline]

a constructor which just carries out the construction of the parent class

Definition at line 51 of file ClusterSequenceAreaBase.hh.

00053                                                       :
00054            ClusterSequence(pseudojets, jet_def, writeout_combinations) {}

ClusterSequenceAreaBase::ClusterSequenceAreaBase (  )  [inline]

default constructor

Definition at line 58 of file ClusterSequenceAreaBase.hh.

00058 {}

virtual ClusterSequenceAreaBase::~ClusterSequenceAreaBase (  )  [inline, virtual]

destructor

Definition at line 62 of file ClusterSequenceAreaBase.hh.

00062 {}


Member Function Documentation

void ClusterSequenceAreaBase::_check_jet_alg_good_for_median (  )  const [private]

check the jet algorithm is suitable (and if not issue a warning)

Definition at line 353 of file ClusterSequenceAreaBase.cc.

References _warnings, cambridge_algorithm, cambridge_for_passive_algorithm, ClusterSequence::jet_def(), kt_algorithm, and LimitedWarning::warn().

00353                                                                    {
00354   if (jet_def().jet_algorithm() != kt_algorithm
00355       && jet_def().jet_algorithm() != cambridge_algorithm
00356       && jet_def().jet_algorithm() !=  cambridge_for_passive_algorithm) {
00357     _warnings.warn("ClusterSequenceAreaBase: jet_def being used may not be suitable for estimating diffuse backgrounds (good options are kt, cam)");
00358   }
00359 }

virtual double ClusterSequenceAreaBase::area ( const PseudoJet  )  const [inline, virtual]

return the area associated with the given jet; this base class returns 0.

Reimplemented in ClusterSequenceActiveArea, ClusterSequenceActiveAreaExplicitGhosts, ClusterSequenceArea, and ClusterSequenceVoronoiArea.

Definition at line 67 of file ClusterSequenceAreaBase.hh.

Referenced by empty_area_from_jets(), parabolic_pt_per_unit_area(), and subtracted_pt().

00067 {return 0.0;}

virtual PseudoJet ClusterSequenceAreaBase::area_4vector ( const PseudoJet  )  const [inline, virtual]

return a PseudoJet whose 4-vector is defined by the following integral

drap d PseudoJet("rap,phi,pt=one") * Theta("rap,phi inside jet boundary")

where PseudoJet("rap,phi,pt=one") is a 4-vector with the given rapidity (rap), azimuth (phi) and pt=1, while Theta("rap,phi inside jet boundary") is a function that is 1 when rap,phi define a direction inside the jet boundary and 0 otherwise.

This base class returns a null 4-vector.

Reimplemented in ClusterSequenceActiveArea, ClusterSequenceActiveAreaExplicitGhosts, ClusterSequenceArea, and ClusterSequenceVoronoiArea.

Definition at line 84 of file ClusterSequenceAreaBase.hh.

Referenced by parabolic_pt_per_unit_area(), and subtracted_jet().

00084                                                            {
00085     return PseudoJet(0.0,0.0,0.0,0.0);}

virtual double ClusterSequenceAreaBase::area_error ( const PseudoJet  )  const [inline, virtual]

return the error (uncertainty) associated with the determination of the area of this jet; this base class returns 0.

Reimplemented in ClusterSequenceActiveArea, ClusterSequenceArea, and ClusterSequenceVoronoiArea.

Definition at line 71 of file ClusterSequenceAreaBase.hh.

00071 {return 0.0;}

double ClusterSequenceAreaBase::empty_area ( const RangeDefinition range  )  const [virtual]

return the total area, within range, that is free of jets, in general based on the inclusive jets

return the total area, within range, that is free of jets.

Calculate this as (range area) - {i in range} A_i

for ClusterSequences with explicit ghosts, assume that there will never be any empty area, i.e. it is always filled in by pure ghosts jets. This holds for seq.rec. algorithms

Reimplemented in ClusterSequenceActiveArea, ClusterSequenceActiveAreaExplicitGhosts, ClusterSequenceArea, and ClusterSequencePassiveArea.

Definition at line 55 of file ClusterSequenceAreaBase.cc.

References empty_area_from_jets(), has_explicit_ghosts(), and ClusterSequence::inclusive_jets().

Referenced by n_empty_jets().

00055                                                                               {
00056 
00057   if (has_explicit_ghosts()) {return 0.0;}
00058   else { return empty_area_from_jets(inclusive_jets(0.0), range);}
00059 
00060 }

double ClusterSequenceAreaBase::empty_area_from_jets ( const std::vector< PseudoJet > &  all_jets,
const RangeDefinition range 
) const

return the total area, within range, that is free of jets, based on the supplied all_jets

return the total area, within range, that is free of jets.

Calculate this as (range area) - {i in range} A_i

Definition at line 67 of file ClusterSequenceAreaBase.cc.

References area(), RangeDefinition::area(), and RangeDefinition::is_in_range().

Referenced by empty_area().

00069                                                            {
00070 
00071   double empty = range.area();
00072   for (unsigned i = 0; i < all_jets.size(); i++) {
00073     if (range.is_in_range(all_jets[i])) empty -= area(all_jets[i]);
00074   }
00075   return empty;
00076 }

virtual void ClusterSequenceAreaBase::get_median_rho_and_sigma ( const RangeDefinition range,
bool  use_area_4vector,
double &  median,
double &  sigma 
) const [inline, virtual]

same as the full version of get_median_rho_and_error, but without access to the mean_area

Reimplemented in ClusterSequenceArea.

Definition at line 188 of file ClusterSequenceAreaBase.hh.

References get_median_rho_and_sigma().

00190                                                                        {
00191     double mean_area;
00192     get_median_rho_and_sigma(range,  use_area_4vector,
00193                              median,  sigma, mean_area);
00194   }

virtual void ClusterSequenceAreaBase::get_median_rho_and_sigma ( const std::vector< PseudoJet > &  all_jets,
const RangeDefinition range,
bool  use_area_4vector,
double &  median,
double &  sigma,
double &  mean_area,
bool  all_are_inclusive = false 
) const [virtual]

a more advanced version of get_median_rho_and_sigma, which allows one to use any "view" of the event containing all jets (so that, e.g.

one might use Cam on a different resolution scale without have to rerun the algorithm).

By default it will assume that "all" are not inclusive jets, so that in dealing with empty area it has to calculate the number of empty jets based on the empty area and the the observed <area> of jets rather than a surmised area

Note that for small effective radii, this can cause problems because the harder jets get an area >> <ghost-jet-area> and so the estimate comes out all wrong. In these situations it is highly advisable to use an area with explicit ghosts, since then the "empty" jets are actually visible.

Reimplemented in ClusterSequenceArea.

void ClusterSequenceAreaBase::get_median_rho_and_sigma ( const RangeDefinition range,
bool  use_area_4vector,
double &  median,
double &  sigma,
double &  mean_area 
) const [virtual]

using jets withing range (and with 4-vector areas if use_area_4vector), calculate the median pt/area, as well as an "error" (uncertainty), which is defined as the 1-sigma half-width of the distribution of pt/A, obtained by looking for the point below which we have (1-0.6827)/2 of the jets (including empty jets).

The subtraction for a jet with uncorrected pt pt^U and area A is

pt^S = pt^U - median*A +- sigma*sqrt(A)

where the error is only that associated with the fluctuations in the noise and not that associated with the noise having caused changes in the hard-particle content of the jet.

NB: subtraction may also be done with 4-vector area of course, and this is recommended for jets with larger values of R, as long as rho has also been determined with a 4-vector area; using a scalar area causes one to neglect terms of relative order $R^2/8$ in the jet $p_t$.

Reimplemented in ClusterSequenceArea.

Definition at line 157 of file ClusterSequenceAreaBase.cc.

References get_median_rho_and_sigma(), and ClusterSequence::inclusive_jets().

Referenced by get_median_rho_and_sigma(), and median_pt_per_unit_something().

00159                                                                        {
00160 
00161   vector<PseudoJet> incl_jets = inclusive_jets();
00162   get_median_rho_and_sigma(incl_jets, range, use_area_4vector,
00163                            median, sigma, mean_area, true);
00164 }

virtual bool ClusterSequenceAreaBase::has_explicit_ghosts (  )  const [inline, virtual]

returns true if ghosts are explicitly included within jets for this ClusterSequence;

Derived classes that do include explicit ghosts should provide an alternative version of this routine and set it properly.

Reimplemented in ClusterSequenceActiveAreaExplicitGhosts, and ClusterSequenceArea.

Definition at line 101 of file ClusterSequenceAreaBase.hh.

Referenced by empty_area().

00101                                            {
00102     return false;
00103   }

virtual bool ClusterSequenceAreaBase::is_pure_ghost ( const PseudoJet  )  const [inline, virtual]

true if a jet is made exclusively of ghosts

NB: most area classes do not give any explicit ghost jets, but some do, and they should replace this function with their own version.

Reimplemented in ClusterSequenceActiveAreaExplicitGhosts, and ClusterSequenceArea.

Definition at line 92 of file ClusterSequenceAreaBase.hh.

00092                                                        {
00093     return false;
00094   }

double ClusterSequenceAreaBase::median_pt_per_unit_area ( const RangeDefinition range  )  const

the median of (pt/area) for jets contained within range, making use also of the info on n_empty_jets

Definition at line 78 of file ClusterSequenceAreaBase.cc.

References median_pt_per_unit_something().

Referenced by subtracted_pt().

00078                                                                                            {
00079   return median_pt_per_unit_something(range,false);
00080 }

double ClusterSequenceAreaBase::median_pt_per_unit_area_4vector ( const RangeDefinition range  )  const

the median of (pt/area_4vector) for jets contained within making use also of the info on n_empty_jets

Definition at line 82 of file ClusterSequenceAreaBase.cc.

References median_pt_per_unit_something().

Referenced by subtracted_jet(), and subtracted_jets().

00082                                                                                                    {
00083   return median_pt_per_unit_something(range,true);
00084 }

double ClusterSequenceAreaBase::median_pt_per_unit_something ( const RangeDefinition range,
bool  use_area_4vector 
) const

the function that does the work for median_pt_per_unit_area and median_pt_per_unit_area_4vector:

  • something_is_area_4vect = false -> use plain area
  • something_is_area_4vect = true -> use 4-vector area

the median of (pt/area) for jets contained within range, counting the empty area as if it were made up of a collection of empty jets each of area (0.55 * pi R^2).

Definition at line 91 of file ClusterSequenceAreaBase.cc.

References get_median_rho_and_sigma().

Referenced by median_pt_per_unit_area(), and median_pt_per_unit_area_4vector().

00092                                                                             {
00093 
00094   double median, sigma, mean_area;
00095   get_median_rho_and_sigma(range, use_area_4vector, median, sigma, mean_area);
00096   return median;
00097 
00098 }

virtual double ClusterSequenceAreaBase::n_empty_jets ( const RangeDefinition range  )  const [inline, virtual]

return something similar to the number of pure ghost jets in the given range in an active area case.

For the local implementation we return empty_area/(0.55 pi R^2), based on measured properties of ghost jets with kt and cam. Note that the number returned is a double.

Reimplemented in ClusterSequence1GhostPassiveArea, ClusterSequenceActiveArea, and ClusterSequenceArea.

Definition at line 119 of file ClusterSequenceAreaBase.hh.

References empty_area(), ClusterSequence::jet_def(), fastjet::pi, and JetDefinition::R().

00119                                                                    {
00120     double R = jet_def().R();
00121     return empty_area(range)/(0.55*pi*R*R);
00122   }

void ClusterSequenceAreaBase::parabolic_pt_per_unit_area ( double &  a,
double &  b,
const RangeDefinition range,
double  exclude_above = -1.0,
bool  use_area_4vector = false 
) const [virtual]

fits a form pt_per_unit_area(y) = a + b*y^2 in the range "range".

fits a form pt_per_unit_area(y) = a + b*y^2 for jets in range.

exclude_above allows one to exclude large values of pt/area from fit. (if negative, the cut is discarded) use_area_4vector = true uses the 4vector areas.

exclude_above allows one to exclude large values of pt/area from fit. use_area_4vector = true uses the 4vector areas.

Reimplemented in ClusterSequenceArea.

Definition at line 105 of file ClusterSequenceAreaBase.cc.

References area(), area_4vector(), ClusterSequence::inclusive_jets(), RangeDefinition::is_in_range(), and PseudoJet::perp().

00107                                                           {
00108   
00109   int n=0;
00110   int n_excluded = 0;
00111   double mean_f=0, mean_x2=0, mean_x4=0, mean_fx2=0; 
00112 
00113   vector<PseudoJet> incl_jets = inclusive_jets();
00114 
00115   for (unsigned i = 0; i < incl_jets.size(); i++) {
00116     if (range.is_in_range(incl_jets[i])) {
00117       double this_area;
00118       if ( use_area_4vector ) {
00119           this_area = area_4vector(incl_jets[i]).perp();     
00120       } else {
00121           this_area = area(incl_jets[i]);
00122       }
00123       double f = incl_jets[i].perp()/this_area;
00124       if (exclude_above <= 0.0 || f < exclude_above) {
00125         double x = incl_jets[i].rap(); double x2 = x*x;
00126         mean_f   += f;
00127         mean_x2  += x2;
00128         mean_x4  += x2*x2;
00129         mean_fx2 += f*x2;
00130         n++;
00131       } else {
00132         n_excluded++;
00133       }
00134     }
00135   }
00136 
00137   if (n <= 1) {
00138     // meaningful results require at least two jets inside the
00139     // area -- mind you if there are empty jets we should be in 
00140     // any case doing something special...
00141     a = 0.0;
00142     b = 0.0;
00143   } else {
00144     mean_f   /= n;
00145     mean_x2  /= n;
00146     mean_x4  /= n;
00147     mean_fx2 /= n;
00148     
00149     b = (mean_f*mean_x2 - mean_fx2)/(mean_x2*mean_x2 - mean_x4);
00150     a = mean_f - b*mean_x2;
00151   }
00152   //cerr << "n_excluded = "<< n_excluded << endl;
00153 }

PseudoJet ClusterSequenceAreaBase::subtracted_jet ( const PseudoJet jet,
const RangeDefinition range 
) const

return a subtracted jet, using area_4vector; note that this is potentially inefficient if repeatedly used for many different jets, because rho will be recalculated each time around.

Definition at line 315 of file ClusterSequenceAreaBase.cc.

References median_pt_per_unit_area_4vector(), and subtracted_jet().

00316                                                                             {
00317   double rho = median_pt_per_unit_area_4vector(range);
00318   PseudoJet sub_jet = subtracted_jet(jet, rho);
00319   return sub_jet;
00320 }

PseudoJet ClusterSequenceAreaBase::subtracted_jet ( const PseudoJet jet,
const double  rho 
) const

return a subtracted jet, using area_4vector, given rho

Definition at line 294 of file ClusterSequenceAreaBase.cc.

References area_4vector(), PseudoJet::cluster_hist_index(), PseudoJet::perp(), PseudoJet::set_cluster_hist_index(), PseudoJet::set_user_index(), and PseudoJet::user_index().

Referenced by subtracted_jet(), subtracted_jets(), and subtracted_pt().

00295                                                                           {
00296   PseudoJet area4vect = area_4vector(jet);
00297   PseudoJet sub_jet;
00298   // sanity check
00299   if (rho*area4vect.perp() < jet.perp() ) { 
00300     sub_jet = jet - rho*area4vect;
00301   } else { sub_jet = PseudoJet(0.0,0.0,0.0,0.0); }
00302   
00303   // make sure the subtracted jet has the same index (cluster and user)
00304   // (i.e. "looks like") the original jet
00305   sub_jet.set_cluster_hist_index(jet.cluster_hist_index());
00306   sub_jet.set_user_index(jet.user_index());
00307   
00308   return sub_jet;
00309 }

vector< PseudoJet > ClusterSequenceAreaBase::subtracted_jets ( const RangeDefinition range,
const double  ptmin = 0.0 
) const

return a vector of subtracted jets, using area_4vector.

Only inclusive_jets above ptmin are subtracted and returned. the ordering is the same as that of sorted_by_pt(cs.inclusive_jets()), i.e. not necessarily ordered in pt once subtracted

Definition at line 284 of file ClusterSequenceAreaBase.cc.

References median_pt_per_unit_area_4vector(), and subtracted_jets().

00287                                                        {
00288   double rho = median_pt_per_unit_area_4vector(range);
00289   return subtracted_jets(rho,ptmin);
00290 }

vector< PseudoJet > ClusterSequenceAreaBase::subtracted_jets ( const double  rho,
const double  ptmin = 0.0 
) const

return a vector of all subtracted jets, using area_4vector, given rho.

Only inclusive_jets above ptmin are subtracted and returned. the ordering is the same as that of sorted_by_pt(cs.inclusive_jets()), i.e. not necessarily ordered in pt once subtracted

Definition at line 268 of file ClusterSequenceAreaBase.cc.

References ClusterSequence::inclusive_jets(), ClusterSequence::jets(), sorted_by_pt(), and subtracted_jet().

Referenced by subtracted_jets().

00270                                                                  {
00271   vector<PseudoJet> sub_jets;
00272   vector<PseudoJet> jets = sorted_by_pt(inclusive_jets(ptmin));
00273   for (unsigned i=0; i<jets.size(); i++) {
00274      PseudoJet sub_jet = subtracted_jet(jets[i],rho);
00275      sub_jets.push_back(sub_jet);
00276   }
00277   return sub_jets;
00278 }

double ClusterSequenceAreaBase::subtracted_pt ( const PseudoJet jet,
const RangeDefinition range,
bool  use_area_4vector = false 
) const

return the subtracted pt; note that this is potentially inefficient if repeatedly used for many different jets, because rho will be recalculated each time around.

Definition at line 339 of file ClusterSequenceAreaBase.cc.

References median_pt_per_unit_area(), PseudoJet::perp(), subtracted_jet(), and subtracted_pt().

00341                                                                            {
00342   if ( use_area_4vector ) { 
00343      PseudoJet sub_jet = subtracted_jet(jet,range);
00344      return sub_jet.perp();
00345   } else {
00346      double rho = median_pt_per_unit_area(range);
00347      return subtracted_pt(jet,rho,false);
00348   }
00349 }  

double ClusterSequenceAreaBase::subtracted_pt ( const PseudoJet jet,
const double  rho,
bool  use_area_4vector = false 
) const

return the subtracted pt, given rho

Definition at line 324 of file ClusterSequenceAreaBase.cc.

References area(), PseudoJet::perp(), and subtracted_jet().

Referenced by subtracted_pt().

00326                                                                            {
00327   if ( use_area_4vector ) { 
00328      PseudoJet sub_jet = subtracted_jet(jet,rho);
00329      return sub_jet.perp();
00330   } else {
00331      return jet.perp() - rho*area(jet);
00332   }
00333 }  


Member Data Documentation

handle warning messages

allow for warnings

Reimplemented in ClusterSequenceActiveAreaExplicitGhosts.

Definition at line 246 of file ClusterSequenceAreaBase.hh.

Referenced by _check_jet_alg_good_for_median().

Definition at line 247 of file ClusterSequenceAreaBase.hh.


The documentation for this class was generated from the following files:

Generated on 26 Feb 2010 for fastjet by  doxygen 1.6.1