FastJet  3.1.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Modules Pages
JetMedianBackgroundEstimator.cc
1 //FJSTARTHEADER
2 // $Id: JetMedianBackgroundEstimator.cc 3517 2014-08-01 14:23:13Z soyez $
3 //
4 // Copyright (c) 2005-2014, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
5 //
6 //----------------------------------------------------------------------
7 // This file is part of FastJet.
8 //
9 // FastJet is free software; you can redistribute it and/or modify
10 // it under the terms of the GNU General Public License as published by
11 // the Free Software Foundation; either version 2 of the License, or
12 // (at your option) any later version.
13 //
14 // The algorithms that underlie FastJet have required considerable
15 // development. They are described in the original FastJet paper,
16 // hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
17 // FastJet as part of work towards a scientific publication, please
18 // quote the version you use and include a citation to the manual and
19 // optionally also to hep-ph/0512210.
20 //
21 // FastJet is distributed in the hope that it will be useful,
22 // but WITHOUT ANY WARRANTY; without even the implied warranty of
23 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24 // GNU General Public License for more details.
25 //
26 // You should have received a copy of the GNU General Public License
27 // along with FastJet. If not, see <http://www.gnu.org/licenses/>.
28 //----------------------------------------------------------------------
29 //FJENDHEADER
30 
31 #include "fastjet/tools/JetMedianBackgroundEstimator.hh"
32 #include <fastjet/ClusterSequenceArea.hh>
33 #include <fastjet/ClusterSequenceStructure.hh>
34 #include <iostream>
35 #include <sstream>
36 
37 FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
38 
39 using namespace std;
40 
41 double BackgroundJetScalarPtDensity::result(const PseudoJet & jet) const {
42  // do not include the ghosts in the list of constituents to have a
43  // correct behaviour when _pt_power is <= 0
44  std::vector<PseudoJet> constituents = (!SelectorIsPureGhost())(jet.constituents());
45  double scalar_pt = 0;
46  for (unsigned i = 0; i < constituents.size(); i++) {
47  scalar_pt += pow(constituents[i].perp(), _pt_power);
48  }
49  return scalar_pt / jet.area();
50 }
51 
52 
53 std::string BackgroundJetScalarPtDensity::description() const {
54  ostringstream oss;
55  oss << "BackgroundScalarJetPtDensity";
56  if (_pt_power != 1.0) oss << " with pt_power = " << _pt_power;
57  return oss.str();
58 }
59 
60 
61 //----------------------------------------------------------------------
62 double BackgroundRescalingYPolynomial::result(const PseudoJet & jet) const {
63  double y = jet.rap();
64  double y2 = y*y;
65  double rescaling = _a0 + _a1*y + _a2*y2 + _a3*y2*y + _a4*y2*y2;
66  return rescaling;
67 }
68 
69 /// allow for warnings
70 LimitedWarning JetMedianBackgroundEstimator::_warnings;
71 LimitedWarning JetMedianBackgroundEstimator::_warnings_zero_area;
72 LimitedWarning JetMedianBackgroundEstimator::_warnings_preliminary;
73 
74 
75 //---------------------------------------------------------------------
76 // class JetMedianBackgroundEstimator
77 // Class to estimate the density of the background per unit area
78 //---------------------------------------------------------------------
79 
80 //----------------------------------------------------------------------
81 // ctors and dtors
82 //----------------------------------------------------------------------
83 // ctor that allows to set only the particles later on
84 JetMedianBackgroundEstimator::JetMedianBackgroundEstimator(const Selector &rho_range,
85  const JetDefinition &jet_def,
86  const AreaDefinition &area_def)
87  : _rho_range(rho_range), _jet_def(jet_def), _area_def(area_def){
88 
89  // initialise things decently
90  reset();
91 
92  // make a few checks
93  _check_jet_alg_good_for_median();
94 }
95 
96 
97 //----------------------------------------------------------------------
98 // ctor from a cluster sequence
99 // - csa the ClusterSequenceArea to use
100 // - rho_range the Selector specifying which jets will be considered
102  : _rho_range(rho_range), _jet_def(JetDefinition()){
103 
104  // initialise things properly
105  reset();
106 
107  // tell the BGE about the cluster sequence
109 }
110 
111 
112 //----------------------------------------------------------------------
113 // setting a new event
114 //----------------------------------------------------------------------
115 // tell the background estimator that it has a new event, composed
116 // of the specified particles.
117 void JetMedianBackgroundEstimator::set_particles(const vector<PseudoJet> & particles) {
118  // make sure that we have been provided a genuine jet definition
119  if (_jet_def.jet_algorithm() == undefined_jet_algorithm)
120  throw Error("JetMedianBackgroundEstimator::set_particles can only be called if you set the jet (and area) definition explicitly through the class constructor");
121 
122  // initialise things decently (including setting uptodate to false!)
123  //reset();
124  _uptodate=false;
125 
126  // cluster the particles
127  //
128  // One may argue that it is better to cache the particles and only
129  // do the clustering later but clustering the particles now has 2
130  // practical advantages:
131  // - it allows to une only '_included_jets' in all that follows
132  // - it avoids adding another flag to ensure particles are
133  // clustered only once
134  ClusterSequenceArea *csa = new ClusterSequenceArea(particles, _jet_def, _area_def);
135  _included_jets = csa->inclusive_jets();
136 
137  // store the CS for later on
138  _csi = csa->structure_shared_ptr();
140 }
141 
142 //----------------------------------------------------------------------
143 // (re)set the cluster sequence (with area support) to be used by
144 // future calls to rho() etc.
145 //
146 // \param csa the cluster sequence area
147 //
148 // Pre-conditions:
149 // - one should be able to estimate the "empty area" (i.e. the area
150 // not occupied by jets). This is feasible if at least one of the following
151 // conditions is satisfied:
152 // ( i) the ClusterSequence has explicit ghosts
153 // (ii) the range selected has a computable area.
154 // - the jet algorithm must be suited for median computation
155 // (otherwise a warning will be issues)
156 //
157 // Note that selectors with e.g. hardest-jets exclusion do not have
158 // a well-defined area. For this reasons, it is STRONGLY advised to
159 // use an area with explicit ghosts.
161  _csi = csa.structure_shared_ptr();
162 
163  // sanity checks
164  //---------------
165  // (i) check the alg is appropriate
166  _check_jet_alg_good_for_median();
167 
168  // (ii) check that, if there are no explicit ghosts, the selector has a finite area
169  if ((!csa.has_explicit_ghosts()) && (!_rho_range.has_finite_area())){
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)");
171  }
172 
173  // get the initial list of jets
174  _included_jets = csa.inclusive_jets();
175 
176  _uptodate = false;
177 }
178 
179 
180 //----------------------------------------------------------------------
181 // (re)set the jets (which must have area support) to be used by future
182 // calls to rho() etc.; for the conditions that must be satisfied
183 // by the jets, see the Constructor that takes jets.
184 void JetMedianBackgroundEstimator::set_jets(const vector<PseudoJet> &jets) {
185 
186  if (! jets.size())
187  throw Error("JetMedianBackgroundEstimator::JetMedianBackgroundEstimator: At least one jet is needed to compute the background properties");
188 
189  // sanity checks
190  //---------------
191  // (o) check that there is an underlying CS shared by all the jets
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");
194 
195  _csi = jets[0].structure_shared_ptr();
196  ClusterSequenceStructure * csi = dynamic_cast<ClusterSequenceStructure*>(_csi());
197  const ClusterSequenceAreaBase * csab = csi->validated_csab();
198 
199  for (unsigned int i=1;i<jets.size(); i++){
200  if (! jets[i].has_associated_cluster_sequence()) // area automatic if the next test succeeds
201  throw Error("JetMedianBackgroundEstimator::set_jets(...): the jets used to estimate the background properties must be associated with a valid ClusterSequenceAreaBase");
202 
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");
205  }
206 
207  // (i) check the alg is appropriate
208  _check_jet_alg_good_for_median();
209 
210  // (ii) check that, if there are no explicit ghosts, the selector has a finite area
211  if ((!csab->has_explicit_ghosts()) && (!_rho_range.has_finite_area())){
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)");
213  }
214 
215 
216  // get the initial list of jets
217  _included_jets = jets;
218 
219  // ensure recalculation of quantities that need it
220  _uptodate = false;
221 }
222 
223 
224 //----------------------------------------------------------------------
225 // retrieving fundamental information
226 //----------------------------------------------------------------
227 
228 // get rho, the median background density per unit area
230  if (_rho_range.takes_reference())
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();
233  return _rho;
234 }
235 
236 // get sigma, the background fluctuations per unit area
238  if (_rho_range.takes_reference())
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();
241  return _sigma;
242 }
243 
244 // get rho, the median background density per unit area, locally at
245 // the position of a given jet.
246 //
247 // If the Selector associated with the range takes a reference jet
248 // (i.e. is relocatable), then for subsequent operations the
249 // Selector has that jet set as its reference.
251  _recompute_if_needed(jet);
252  double our_rho = _rho;
253  if (_rescaling_class != 0) {
254  our_rho *= (*_rescaling_class)(jet);
255  }
256  return our_rho;
257 }
258 
259 // get sigma, the background fluctuations per unit area,
260 // locally at the position of a given jet.
261 //
262 // If the Selector associated with the range takes a reference jet
263 // (i.e. is relocatable), then for subsequent operations the
264 // Selector has that jet set as its reference.
266  _recompute_if_needed(jet);
267  double our_sigma = _sigma;
268  if (_rescaling_class != 0) {
269  our_sigma *= (*_rescaling_class)(jet);
270  }
271  return our_sigma;
272 }
273 
274 
275 //----------------------------------------------------------------------
276 // returns rho_m (particle-masses contribution to the 4-vector density)
278  if (! has_rho_m()){
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).");
280  }
281  if (_rho_range.takes_reference())
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();
284  return _rho_m;
285 }
286 
287 
288 //----------------------------------------------------------------------
289 // returns sigma_m (particle-masses contribution to the 4-vector
290 // density); must be multipled by sqrt(area) to get fluctuations
291 // for a region of a given area.
293  if (! has_rho_m()){
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).");
295  }
296  if (_rho_range.takes_reference())
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();
299  return _sigma_m;
300 }
301 
302 //----------------------------------------------------------------------
303 // returns rho_m locally at the position of a given jet. As for
304 // rho(jet), it is non-const.
306  _recompute_if_needed(jet);
307  double our_rho = _rho_m;
308  if (_rescaling_class != 0) {
309  our_rho *= (*_rescaling_class)(jet);
310  }
311  return our_rho;
312 }
313 
314 
315 //----------------------------------------------------------------------
316 // returns sigma_m locally at the position of a given jet. As for
317 // rho(jet), it is non-const.
319  _recompute_if_needed(jet);
320  double our_sigma = _sigma_m;
321  if (_rescaling_class != 0) {
322  our_sigma *= (*_rescaling_class)(jet);
323  }
324  return our_sigma;
325 }
326 
327 //----------------------------------------------------------------------
328 // configuring behaviour
329 //----------------------------------------------------------------------
330 // reset to default values
331 //
332 // set the variou options to their default values
334  // set the remaining default parameters
335  set_use_area_4vector(); // true by default
336  set_provide_fj2_sigma(false);
337 
338  _enable_rho_m = true;
339 
340  // reset the computed values
341  _rho = _sigma = 0.0;
342  _rho_m = _sigma_m = 0.0;
343  _n_jets_used = _n_empty_jets = 0;
344  _empty_area = _mean_area = 0.0;
345 
346  _included_jets.clear();
347 
348  _jet_density_class = 0; // null pointer
349  _rescaling_class = 0; // null pointer
350 
351  _uptodate = false;
352 }
353 
354 
355 // Set a pointer to a class that calculates the quantity whose
356 // median will be calculated; if the pointer is null then pt/area
357 // is used (as occurs also if this function is not called).
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;
361  _uptodate = false;
362 }
363 
364 
365 
366 //----------------------------------------------------------------------
367 // description
368 //----------------------------------------------------------------------
370  ostringstream desc;
371  desc << "JetMedianBackgroundEstimator, using " << _jet_def.description()
372  << " with " << _area_def.description()
373  << " and selecting jets with " << _rho_range.description();
374  return desc.str();
375 }
376 
377 
378 
379 //----------------------------------------------------------------------
380 // computation of the background properties
381 //----------------------------------------------------------------------
382 // for estimation using a relocatable selector (i.e. local range)
383 // this allows to set its position. Note that this HAS to be called
384 // before any attempt to compute the background properties
385 void JetMedianBackgroundEstimator::_recompute_if_needed(const PseudoJet &jet){
386  // if the range is relocatable, handles its relocation
387  if (_rho_range.takes_reference()){
388  // check that the reference is not the same as the previous one
389  // (would avoid an unnecessary recomputation)
390  if (jet == _current_reference) return;
391 
392  // relocate the range and make sure things get recomputed the next
393  // time one tries to get some information
394  _rho_range.set_reference(jet);
395  _uptodate=false;
396  }
397 
398  _recompute_if_needed();
399 }
400 
401 // do the actual job
402 void JetMedianBackgroundEstimator::_compute() const {
403  // check if the clustersequence is still valid
404  _check_csa_alive();
405 
406  // fill the vector of pt/area (or the quantity from the jet density class)
407  // - in the range
408  vector<double> vector_for_median_pt;
409  vector<double> vector_for_median_dt;
410  double total_area = 0.0;
411  _n_jets_used = 0;
412 
413  // apply the selector to the included jets
414  vector<PseudoJet> selected_jets = _rho_range(_included_jets);
415 
416  // compute the pt/area for the selected jets
417  double median_input_pt, median_input_dt=0.0;
418  BackgroundJetPtMDensity m_density;
419  bool do_rho_m = has_rho_m();
420  for (unsigned i = 0; i < selected_jets.size(); i++) {
421  const PseudoJet & current_jet = selected_jets[i];
422 
423  double this_area = (_use_area_4vector) ? current_jet.area_4vector().perp() : current_jet.area();
424  if (this_area>0){
425  // for the pt component, we either use pt or the user-provided
426  // density class
427  if (_jet_density_class == 0) {
428  median_input_pt = current_jet.perp()/this_area;
429  } else {
430  median_input_pt = (*_jet_density_class)(current_jet);
431  }
432 
433  // handle the rho_m part if requested
434  // note that we're using the scalar area as a normalisation inside the
435  // density class!
436  if (do_rho_m)
437  median_input_dt = m_density(current_jet);
438 
439  // perform rescaling if needed
440  if (_rescaling_class != 0) {
441  double resc = (*_rescaling_class)(current_jet);;
442  median_input_pt /= resc;
443  median_input_dt /= resc;
444  }
445 
446  // store the result for future computation of the median
447  vector_for_median_pt.push_back(median_input_pt);
448  if (do_rho_m)
449  vector_for_median_dt.push_back(median_input_dt);
450 
451  total_area += this_area;
452  _n_jets_used++;
453  } else {
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).");
455  }
456  }
457 
458  // there is nothing inside our region, so answer will always be zero
459  if (vector_for_median_pt.size() == 0) {
460  _rho = 0.0;
461  _sigma = 0.0;
462  _rho_m = 0.0;
463  _sigma_m = 0.0;
464  _mean_area = 0.0;
465  return;
466  }
467 
468  // determine the number of empty jets
469  const ClusterSequenceAreaBase * csab = (dynamic_cast<ClusterSequenceStructure*>(_csi()))->validated_csab();
470  if (csab->has_explicit_ghosts()) {
471  _empty_area = 0.0;
472  _n_empty_jets = 0;
473  } else {
474  _empty_area = csab->empty_area(_rho_range);
475  _n_empty_jets = csab->n_empty_jets(_rho_range);
476  }
477 
478  double total_njets = _n_jets_used + _n_empty_jets;
479  total_area += _empty_area;
480 
481  double stand_dev;
482  _median_and_stddev(vector_for_median_pt, _n_empty_jets, _rho, stand_dev,
483  _provide_fj2_sigma);
484 
485  // process and store the results (_rho was already stored above)
486  _mean_area = total_area / total_njets;
487  _sigma = stand_dev * sqrt(_mean_area);
488 
489  // compute the rho_m part now
490  if (do_rho_m){
491  _median_and_stddev(vector_for_median_dt, _n_empty_jets, _rho_m, stand_dev,
492  _provide_fj2_sigma);
493  _sigma_m = stand_dev * sqrt(_mean_area);
494  }
495 
496  // record that the computation has been performed
497  _uptodate = true;
498 }
499 
500 
501 
502 // check that the underlying structure is still alive;
503 // throw an error otherwise
504 void JetMedianBackgroundEstimator::_check_csa_alive() const{
505  ClusterSequenceStructure* csa = dynamic_cast<ClusterSequenceStructure*>(_csi());
506  if (csa == 0) {
507  throw Error("JetMedianBackgroundEstimator: there is no cluster sequence associated with the JetMedianBackgroundEstimator");
508  }
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");
511 }
512 
513 
514 // check that the algorithm used for the clustering is suitable for
515 // background estimation (i.e. either kt or C/A).
516 // Issue a warning otherwise
517 void JetMedianBackgroundEstimator::_check_jet_alg_good_for_median() const{
518  const JetDefinition * jet_def = &_jet_def;
519 
520  // if no explicit jet def has been provided, fall back on the
521  // cluster sequence
522  if (_jet_def.jet_algorithm() == undefined_jet_algorithm){
523  const ClusterSequence * cs = dynamic_cast<ClusterSequenceStructure*>(_csi())->validated_cs();
524  jet_def = &(cs->jet_def());
525  }
526 
527  if (jet_def->jet_algorithm() != kt_algorithm
528  && jet_def->jet_algorithm() != cambridge_algorithm
529  && jet_def->jet_algorithm() != cambridge_for_passive_algorithm) {
530  _warnings.warn("JetMedianBackgroundEstimator: jet_def being used may not be suitable for estimating diffuse backgrounds (good alternatives are kt, cam)");
531  }
532 }
533 
534 
535 
536 FASTJET_END_NAMESPACE
537 
538 
double rap() const
returns the rapidity or some large value when the rapidity is infinite
Definition: PseudoJet.hh:123
a version of cambridge with a special distance measure for particles whose pt is < extra_param(); thi...
std::string description() const
returns a textual description of the background estimator
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
Definition: Selector.hh:287
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.
Definition: PseudoJet.cc:717
the longitudinally invariant kt algorithm
JetMedianBackgroundEstimator(const Selector &rho_range, const JetDefinition &jet_def, const AreaDefinition &area_def)
Constructor that sets the rho range as well as the jet definition and area definition to be used to c...
class that holds a generic area definition
double sigma() const
get sigma, the background fluctuations per unit area
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...
void set_jet_density_class(const FunctionOfPseudoJet< double > *jet_density_class)
Set a pointer to a class that calculates the quantity whose median will be calculated; if the pointer...
void set_provide_fj2_sigma(bool provide_fj2_sigma=true)
The FastJet v2.X sigma calculation had a small spurious offset in the limit of a small number of jets...
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
virtual void set_particles(const std::vector< PseudoJet > &particles)
tell the background estimator that it has a new event, composed of the specified particles.
base class corresponding to errors that can be thrown by FastJet
Definition: Error.hh:47
void reset()
Resets the class to its default state, including the choice to use 4-vector areas.
void set_use_area_4vector(bool use_it=true)
By default when calculating pt/Area for a jet, it is the transverse component of the 4-vector area th...
void set_jets(const std::vector< PseudoJet > &jets)
(re)set the jets (which must have area support) to be used by future calls to rho() etc...
bool has_finite_area() const
returns true if it has a meaningful and finite area (i.e.
Definition: Selector.hh:250
std::string description() const
returns a textual description of the selector
Definition: Selector.hh:237
const Selector & set_reference(const PseudoJet &reference)
set the reference jet for this Selector
Definition: Selector.hh:292
Class that encodes information about cuts and other selection criteria that can be applied to PseudoJ...
Definition: Selector.hh:149
double rho() const
get rho, the median background density per unit area
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.
Definition: Selector.cc:1420
virtual double sigma_m() const
returns sigma_m, a measure of the fluctuations in the purely longitudinal, particle-mass-induced comp...
virtual double rho_m() const
returns rho_m, the purely longitudinal, particle-mass-induced component of the background density per...
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.
Definition: PseudoJet.cc:580
Class to contain pseudojets, including minimal information of use to jet-clustering routines...
Definition: PseudoJet.hh:67
void set_cluster_sequence(const ClusterSequenceAreaBase &csa)
(re)set the cluster sequence (with area support) to be used by future calls to rho() etc...
virtual bool has_rho_m() const
Returns true if this background estimator has support for determination of rho_m. ...
class that is intended to hold a full definition of the jet clusterer