Classes | Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Static Protected Attributes

fastjet::ClusterSequence Class Reference
[Fundamental FastJet classes]

deals with clustering More...

#include <ClusterSequence.hh>

Inheritance diagram for fastjet::ClusterSequence:
Inheritance graph
[legend]
Collaboration diagram for fastjet::ClusterSequence:
Collaboration graph
[legend]

List of all members.

Classes

class  Extras
 base class to store extra information that plugins may provide More...
struct  history_element
 a single element in the clustering history More...

Public Types

enum  JetType { Invalid = -3, InexistentParent = -2, BeamJet = -1 }
typedef ClusterSequenceStructure StructureType
 the structure type associated with a jet belonging to a ClusterSequence

Public Member Functions

 ClusterSequence ()
 default constructor
template<class L >
 ClusterSequence (const std::vector< L > &pseudojets, const double &R=1.0, const Strategy &strategy=Best, const bool &writeout_combinations=false)
 create a clustersequence starting from the supplied set of pseudojets and clustering them with the long-invariant kt algorithm (E-scheme recombination) with the supplied value for R.
template<class L >
 ClusterSequence (const std::vector< L > &pseudojets, const JetDefinition &jet_def, const bool &writeout_combinations=false)
 create a clustersequence starting from the supplied set of pseudojets and clustering them with jet definition specified by jet_def (which also specifies the clustering strategy)
std::vector< PseudoJetinclusive_jets (const double &ptmin=0.0) const
 return a vector of all jets (in the sense of the inclusive algorithm) with pt >= ptmin.
int n_exclusive_jets (const double &dcut) const
 return the number of jets (in the sense of the exclusive algorithm) that would be obtained when running the algorithm with the given dcut.
std::vector< PseudoJetexclusive_jets (const double &dcut) const
 return a vector of all jets (in the sense of the exclusive algorithm) that would be obtained when running the algorithm with the given dcut.
std::vector< PseudoJetexclusive_jets (const int &njets) const
 return a vector of all jets when the event is clustered (in the exclusive sense) to exactly njets.
double exclusive_dmerge (const int &njets) const
 return the dmin corresponding to the recombination that went from n+1 to n jets (sometimes known as d_{n n+1}).
double exclusive_dmerge_max (const int &njets) const
 return the maximum of the dmin encountered during all recombinations up to the one that led to an n-jet final state; identical to exclusive_dmerge, except in cases where the dmin do not increase monotonically.
double exclusive_ymerge (int njets) const
 return the ymin corresponding to the recombination that went from n+1 to n jets (sometimes known as y_{n n+1}).
double exclusive_ymerge_max (int njets) const
 same as exclusive_dmerge_max, but normalised to squared total energy
int n_exclusive_jets_ycut (double ycut) const
 the number of exclusive jets at the given ycut
std::vector< PseudoJetexclusive_jets_ycut (double ycut) const
 the exclusive jets obtained at the given ycut
std::vector< PseudoJetexclusive_subjets (const PseudoJet &jet, const double &dcut) const
 return a vector of all subjets of the current jet (in the sense of the exclusive algorithm) that would be obtained when running the algorithm with the given dcut.
int n_exclusive_subjets (const PseudoJet &jet, const double &dcut) const
 return the size of exclusive_subjets(...); still n ln n with same coefficient, but marginally more efficient than manually taking exclusive_subjets.size()
std::vector< PseudoJetexclusive_subjets (const PseudoJet &jet, int nsub) const
 return the list of subjets obtained by unclustering the supplied jet down to n subjets (or all constituents if there are fewer than n).
double exclusive_subdmerge (const PseudoJet &jet, int nsub) const
 return the dij that was present in the merging nsub+1 -> nsub subjets inside this jet.
double exclusive_subdmerge_max (const PseudoJet &jet, int nsub) const
 return the maximum dij that occurred in the whole event at the stage that the nsub+1 -> nsub merge of subjets occurred inside this jet.
double Q () const
 returns the sum of all energies in the event (relevant mainly for e+e-)
double Q2 () const
 return Q()^2
bool object_in_jet (const PseudoJet &object, const PseudoJet &jet) const
 returns true iff the object is included in the jet.
bool has_parents (const PseudoJet &jet, PseudoJet &parent1, PseudoJet &parent2) const
 if the jet has parents in the clustering, it returns true and sets parent1 and parent2 equal to them.
bool has_child (const PseudoJet &jet, PseudoJet &child) const
 if the jet has a child then return true and give the child jet otherwise return false and set the child to zero
bool has_child (const PseudoJet &jet, const PseudoJet *&childp) const
 Version of has_child that sets a pointer to the child if the child exists;.
bool has_partner (const PseudoJet &jet, PseudoJet &partner) const
 if this jet has a child (and so a partner) return true and give the partner, otherwise return false and set the partner to zero
std::vector< PseudoJetconstituents (const PseudoJet &jet) const
 return a vector of the particles that make up jet
void print_jets_for_root (const std::vector< PseudoJet > &jets, std::ostream &ostr=std::cout) const
 output the supplied vector of jets in a format that can be read by an appropriate root script; the format is: jet-n jet-px jet-py jet-pz jet-E particle-n particle-rap particle-phi particle-pt particle-n particle-rap particle-phi particle-pt ...
void print_jets_for_root (const std::vector< PseudoJet > &jets, const std::string &filename, const std::string &comment="") const
 print jets for root to the file labelled filename, with an optional comment at the beginning
void add_constituents (const PseudoJet &jet, std::vector< PseudoJet > &subjet_vector) const
 add on to subjet_vector the constituents of jet (for internal use mainly)
Strategy strategy_used () const
 return the enum value of the strategy used to cluster the event
std::string strategy_string () const
 return the name of the strategy used to cluster the event
std::string strategy_string (Strategy strategy_in) const
 return the name of the strategy associated with the enum strategy_in
const JetDefinitionjet_def () const
 return a reference to the jet definition
double jet_scale_for_algorithm (const PseudoJet &jet) const
 returns the scale associated with a jet as required for this clustering algorithm (kt^2 for the kt-algorithm, 1 for the Cambridge algorithm).
void plugin_record_ij_recombination (int jet_i, int jet_j, double dij, int &newjet_k)
 record the fact that there has been a recombination between jets()[jet_i] and jets()[jet_k], with the specified dij, and return the index (newjet_k) allocated to the new jet, whose momentum is assumed to be the 4-vector sum of that of jet_i and jet_j
void plugin_record_ij_recombination (int jet_i, int jet_j, double dij, const PseudoJet &newjet, int &newjet_k)
 as for the simpler variant of plugin_record_ij_recombination, except that the new jet is attributed the momentum and user_index of newjet
void plugin_record_iB_recombination (int jet_i, double diB)
 record the fact that there has been a recombination between jets()[jet_i] and the beam, with the specified diB; when looking for inclusive jets, any iB recombination will returned to the user as a jet.
void plugin_associate_extras (std::auto_ptr< Extras > extras_in)
 the plugin can associated some extra information with the ClusterSequence object by calling this function
bool plugin_activated () const
 returns true when the plugin is allowed to run the show.
const Extrasextras () const
 returns a pointer to the extras object (may be null)
template<class GBJ >
void plugin_simple_N2_cluster ()
 allows a plugin to run a templated clustering (nearest-neighbour heuristic)
const std::vector< PseudoJet > & jets () const
 allow the user to access the internally stored _jets() array, which contains both the initial particles and the various intermediate and final stages of recombination.
const std::vector
< history_element > & 
history () const
 allow the user to access the raw internal history.
unsigned int n_particles () const
 returns the number of particles that were provided to the clustering algorithm (helps the user find their way around the history and jets objects if they weren't paying attention beforehand).
std::vector< int > particle_jet_indices (const std::vector< PseudoJet > &) const
 returns a vector of size n_particles() which indicates, for each of the initial particles (in the order in which they were supplied), which of the supplied jets it belongs to; if it does not belong to any of the supplied jets, the index is set to -1;
std::vector< int > unique_history_order () const
 routine that returns an order in which to read the history such that clusterings that lead to identical jet compositions but different histories (because of degeneracies in the clustering order) will have matching constituents for each matching entry in the unique_history_order.
std::vector< PseudoJetunclustered_particles () const
 return the set of particles that have not been clustered.
bool contains (const PseudoJet &object) const
 returns true if the object (jet or particle) is contained by (ie belongs to) this cluster sequence.
void transfer_from_sequence (ClusterSequence &from_seq, bool transfer_ownership=true)
 transfer the sequence contained in other_seq into our own; any plugin "extras" contained in the from_seq will be lost from there.
const SharedPtr
< PseudoJetStructureBase > & 
structure_shared_ptr () const
 retrieve a shared pointer to the wrapper to this ClusterSequence
template<>
void _bj_set_jetinfo (EEBriefJet *const jetA, const int _jets_index) const
template<>
double _bj_dist (const EEBriefJet *const jeta, const EEBriefJet *const jetb) const

Static Public Member Functions

static void set_jet_algorithm (JetAlgorithm jet_algorithm)
 set the default (static) jet finder across all current and future ClusterSequence objects -- deprecated and obsolescent (i.e.
static void set_jet_finder (JetAlgorithm jet_algorithm)
 same as above for backward compatibility

Protected Member Functions

template<class L >
void _transfer_input_jets (const std::vector< L > &pseudojets)
 transfer the vector<L> of input jets into our own vector<PseudoJet> _jets (with some reserved space for future growth).
void _initialise_and_run (const JetDefinition &jet_def, const bool &writeout_combinations)
 This is the routine that will do all the initialisation and then run the clustering (may be called by various constructors).
void _initialise_and_run (const double &R, const Strategy &strategy, const bool &writeout_combinations)
 This is an alternative routine for initialising and running the clustering, provided for legacy purposes.
void _decant_options (const JetDefinition &jet_def, const bool &writeout_combinations)
 fills in the various member variables with "decanted" options from the jet_definition and writeout_combinations variables
void _fill_initial_history ()
 fill out the history (and jet cross refs) related to the initial set of jets (assumed already to have been "transferred"), without any clustering
void _do_ij_recombination_step (const int &jet_i, const int &jet_j, const double &dij, int &newjet_k)
 carry out the recombination between the jets numbered jet_i and jet_j, at distance scale dij; return the index newjet_k of the result of the recombination of i and j.
void _do_iB_recombination_step (const int &jet_i, const double &diB)
 carry out an recombination step in which _jets[jet_i] merges with the beam,
void get_subhist_set (std::set< const history_element * > &subhist, const PseudoJet &jet, double dcut, int maxjet) const
 set subhist to be a set pointers to history entries corresponding to the subjets of this jet; one stops going working down through the subjets either when

  • there is no further to go
  • one has found maxjet entries
  • max_dij_so_far <= dcut By setting maxjet=0 one can use just dcut; by setting dcut<0 one can use jet maxjet

Protected Attributes

JetDefinition _jet_def
std::vector< PseudoJet_jets
 This contains the physical PseudoJets; for each PseudoJet one can find the corresponding position in the _history by looking at _jets[i].cluster_hist_index().
std::vector< history_element_history
 this vector will contain the branching history; for each stage, _history[i].jetp_index indicates where to look in the _jets vector to get the physical PseudoJet.
bool _writeout_combinations
int _initial_n
double _Rparam
double _R2
double _invR2
double _Qtot
Strategy _strategy
JetAlgorithm _jet_algorithm
SharedPtr< PseudoJetStructureBase_structure_shared_ptr

Static Protected Attributes

static JetAlgorithm _default_jet_algorithm = kt_algorithm

Detailed Description

deals with clustering

Definition at line 111 of file ClusterSequence.hh.


Constructor & Destructor Documentation

template<class L >
fastjet::ClusterSequence::ClusterSequence ( const std::vector< L > &  pseudojets,
const double &  R = 1.0,
const Strategy strategy = Best,
const bool &  writeout_combinations = false 
)

create a clustersequence starting from the supplied set of pseudojets and clustering them with the long-invariant kt algorithm (E-scheme recombination) with the supplied value for R.

If strategy=DumbN3 a very stupid N^3 algorithm is used for the clustering; otherwise strategy = NlnN* uses cylinders algorithms with some number of pi coverage. If writeout_combinations=true a summary of the recombination sequence is written out

Definition at line 860 of file ClusterSequence.hh.

References _initialise_and_run(), and _transfer_input_jets().

                                                                      {

  // transfer the initial jets (type L) into our own array
  _transfer_input_jets(pseudojets);

  // run the clustering
  _initialise_and_run(R,strategy,writeout_combinations);
}

template<class L >
fastjet::ClusterSequence::ClusterSequence ( const std::vector< L > &  pseudojets,
const JetDefinition jet_def,
const bool &  writeout_combinations = false 
)

create a clustersequence starting from the supplied set of pseudojets and clustering them with jet definition specified by jet_def (which also specifies the clustering strategy)

constructor of a jet-clustering sequence from a vector of four-momenta, with the jet definition specified by jet_def

Definition at line 877 of file ClusterSequence.hh.

References _initialise_and_run(), and _transfer_input_jets().

                                                                      {

  // transfer the initial jets (type L) into our own array
  _transfer_input_jets(pseudojets);

  // run the clustering
  _initialise_and_run(jet_def,writeout_combinations);
}


Member Function Documentation

vector< PseudoJet > fastjet::ClusterSequence::inclusive_jets ( const double &  ptmin = 0.0  )  const

return a vector of all jets (in the sense of the inclusive algorithm) with pt >= ptmin.

Time taken should be of the order of the number of jets returned.

Definition at line 426 of file ClusterSequence.cc.

Referenced by fastjet::BackgroundEstimator::BackgroundEstimator(), fastjet::ClusterSequenceAreaBase::empty_area(), fastjet::ClusterSequenceAreaBase::get_median_rho_and_sigma(), main(), fastjet::ClusterSequenceAreaBase::parabolic_pt_per_unit_area(), fastjet::ClusterSequenceActiveArea::pt_per_unit_area(), and fastjet::ClusterSequenceAreaBase::subtracted_jets().

                                                                            {
  double dcut = ptmin*ptmin;
  int i = _history.size() - 1; // last jet
  vector<PseudoJet> jets;
  if (_jet_algorithm == kt_algorithm) {
    while (i >= 0) {
      // with our specific definition of dij and diB (i.e. R appears only in 
      // dij), then dij==diB is the same as the jet.perp2() and we can exploit
      // this in selecting the jets...
      if (_history[i].max_dij_so_far < dcut) {break;}
      if (_history[i].parent2 == BeamJet && _history[i].dij >= dcut) {
        // for beam jets
        int parent1 = _history[i].parent1;
        jets.push_back(_jets[_history[parent1].jetp_index]);}
      i--;
    }
  } else if (_jet_algorithm == cambridge_algorithm) {
    while (i >= 0) {
      // inclusive jets are all at end of clustering sequence in the
      // Cambridge algorithm -- so if we find a non-exclusive jet, then
      // we can exit
      if (_history[i].parent2 != BeamJet) {break;}
      int parent1 = _history[i].parent1;
      const PseudoJet & jet = _jets[_history[parent1].jetp_index];
      if (jet.perp2() >= dcut) {jets.push_back(jet);}
      i--;
    }
  } else if (_jet_algorithm == plugin_algorithm 
             || _jet_algorithm == ee_kt_algorithm
             || _jet_algorithm == antikt_algorithm
             || _jet_algorithm == genkt_algorithm
             || _jet_algorithm == ee_genkt_algorithm
             || _jet_algorithm == cambridge_for_passive_algorithm) {
    // for inclusive jets with a plugin algorithm, we make no
    // assumptions about anything (relation of dij to momenta,
    // ordering of the dij, etc.)
    while (i >= 0) {
      if (_history[i].parent2 == BeamJet) {

int fastjet::ClusterSequence::n_exclusive_jets ( const double &  dcut  )  const

return the number of jets (in the sense of the exclusive algorithm) that would be obtained when running the algorithm with the given dcut.

Definition at line 478 of file ClusterSequence.cc.

Referenced by exclusive_jets().

         {throw Error("cs::inclusive_jets(...): Unrecognized jet algorithm");}
  return jets;
}


//----------------------------------------------------------------------
// return the number of exclusive jets that would have been obtained
// running the algorithm in exclusive mode with the given dcut
int ClusterSequence::n_exclusive_jets (const double & dcut) const {

  // first locate the point where clustering would have stopped (i.e. the
  // first time max_dij_so_far > dcut)
  int i = _history.size() - 1; // last jet
  while (i >= 0) {

vector< PseudoJet > fastjet::ClusterSequence::exclusive_jets ( const double &  dcut  )  const

return a vector of all jets (in the sense of the exclusive algorithm) that would be obtained when running the algorithm with the given dcut.

Definition at line 497 of file ClusterSequence.cc.

Referenced by exclusive_jets(), and main().

vector< PseudoJet > fastjet::ClusterSequence::exclusive_jets ( const int &  njets  )  const

return a vector of all jets when the event is clustered (in the exclusive sense) to exactly njets.

If there are fewer than njets particles in the ClusterSequence then FastJet crashes since it is not able to return njets exclusive jets.

Definition at line 505 of file ClusterSequence.cc.

References exclusive_jets(), and n_exclusive_jets().

                                                                            {
  int njets = n_exclusive_jets(dcut);
  return exclusive_jets(njets);
}


//----------------------------------------------------------------------
// return the jets obtained by clustering the event to n jets.
vector<PseudoJet> ClusterSequence::exclusive_jets (const int & njets) const {

  // make sure the user does not ask for more than jets than there
  // were particles in the first place.
  assert (njets <= _initial_n);

  // provide a warning when extracting exclusive jets for algorithms 
  // that does not support it explicitly.
  // Native algorithm that support it are: kt, ee_kt, cambridge, 
  //   genkt and ee_genkt (both with p>=0)
  // For plugins, we check Plugin::exclusive_sequence_meaningful()
  if (( _jet_def.jet_algorithm() != kt_algorithm) &&
      ( _jet_def.jet_algorithm() != cambridge_algorithm) &&
      ( _jet_def.jet_algorithm() != ee_kt_algorithm) &&
      (((_jet_def.jet_algorithm() != genkt_algorithm) && 
        (_jet_def.jet_algorithm() != ee_genkt_algorithm)) || 
       (_jet_def.extra_param() <0)) &&
      ((_jet_def.jet_algorithm() != plugin_algorithm) ||
       (!_jet_def.plugin()->exclusive_sequence_meaningful())) &&
      (_n_exclusive_warnings < 5)) {
    _n_exclusive_warnings++;
    cerr << "FastJet WARNING: dcut and exclusive jets for jet-finders other than kt should be interpreted with care." << endl;
  }


  // calculate the point where we have to stop the clustering.
  // relation between stop_point, njets assumes one extra jet disappears
  // at each clustering.
  int stop_point = 2*_initial_n - njets;

  // some sanity checking to make sure that e+e- does not give us
  // surprises (should we ever implement e+e-)...
  if (2*_initial_n != static_cast<int>(_history.size())) {
    ostringstream err;
    err << "2*_initial_n != _history.size() -- this endangers internal assumptions!\n";
    throw Error(err.str());
    //assert(false);
  }

  // now go forwards and reconstitute the jets that we have --
  // basically for any history element, see if the parent jets to
  // which it refers were created before the stopping point -- if they
  // were then add them to the list, otherwise they are subsequent
  // recombinations of the jets that we are looking for.
  vector<PseudoJet> jets;
  for (unsigned int i = stop_point; i < _history.size(); i++) {
    int parent1 = _history[i].parent1;
    if (parent1 < stop_point) {
      jets.push_back(_jets[_history[parent1].jetp_index]);
    }
    int parent2 = _history[i].parent2;
    if (parent2 < stop_point && parent2 > 0) {
      jets.push_back(_jets[_history[parent2].jetp_index]);
    }
    
  }

  // sanity check...
  if (static_cast<int>(jets.size()) != njets) {

double fastjet::ClusterSequence::exclusive_dmerge ( const int &  njets  )  const

return the dmin corresponding to the recombination that went from n+1 to n jets (sometimes known as d_{n n+1}).

return the dmin corresponding to the recombination that went from n+1 to n jets

Definition at line 577 of file ClusterSequence.cc.

double fastjet::ClusterSequence::exclusive_dmerge_max ( const int &  njets  )  const

return the maximum of the dmin encountered during all recombinations up to the one that led to an n-jet final state; identical to exclusive_dmerge, except in cases where the dmin do not increase monotonically.

Definition at line 589 of file ClusterSequence.cc.

double fastjet::ClusterSequence::exclusive_ymerge ( int  njets  )  const [inline]

return the ymin corresponding to the recombination that went from n+1 to n jets (sometimes known as y_{n n+1}).

Definition at line 185 of file ClusterSequence.hh.

{return exclusive_dmerge(njets) / Q2();}

std::vector< PseudoJet > fastjet::ClusterSequence::exclusive_subjets ( const PseudoJet jet,
const double &  dcut 
) const

return a vector of all subjets of the current jet (in the sense of the exclusive algorithm) that would be obtained when running the algorithm with the given dcut.

Time taken is O(m ln m), where m is the number of subjets that are found. If m gets to be of order of the total number of constituents in the jet, this could be substantially slower than just getting that list of constituents.

Definition at line 601 of file ClusterSequence.cc.

Referenced by fastjet::ClusterSequenceStructure::exclusive_subjets().

                                                      {

  set<const history_element*> subhist;

  // get the set of history elements that correspond to subjets at
  // scale dcut
  get_subhist_set(subhist, jet, dcut, 0);

std::vector< PseudoJet > fastjet::ClusterSequence::exclusive_subjets ( const PseudoJet jet,
int  n 
) const

return the list of subjets obtained by unclustering the supplied jet down to n subjets (or all constituents if there are fewer than n).

This requires n ln n time

If the jet contains fewer than nsub particles (in which case it is not possible to return nsub subjets) then the vector of subjets that is returned is simply the list of particles (of size < nsub). Note that this behaviour differs from that of exclusive_jets().

Definition at line 637 of file ClusterSequence.cc.

                                        {

  set<const history_element*> subhist;

  // get the set of history elements that correspond to subjets at
  // scale dcut
  get_subhist_set(subhist, jet, -1.0, n);

double fastjet::ClusterSequence::exclusive_subdmerge ( const PseudoJet jet,
int  nsub 
) const

return the dij that was present in the merging nsub+1 -> nsub subjets inside this jet.

Returns 0 if there were nsub or fewer constituents in the jet.

If the jet has nsub or fewer constituents, it will return 0.

will be zero if nconst <= nsub, since highest will be an original particle have zero dij

Definition at line 661 of file ClusterSequence.cc.

Referenced by fastjet::ClusterSequenceStructure::exclusive_subdmerge().

                                                                                 {
  set<const history_element*> subhist;

  // get the set of history elements that correspond to subjets at

double fastjet::ClusterSequence::exclusive_subdmerge_max ( const PseudoJet jet,
int  nsub 
) const

return the maximum dij that occurred in the whole event at the stage that the nsub+1 -> nsub merge of subjets occurred inside this jet.

Returns 0 if there were nsub or fewer constituents in the jet.

If the jet has nsub or fewer constituents, it will return 0.

will be zero if nconst <= nsub, since highest will be an original particle have zero dij

Definition at line 682 of file ClusterSequence.cc.

Referenced by fastjet::ClusterSequenceStructure::exclusive_subdmerge_max().

                                                                                     {

  set<const history_element*> subhist;

  // get the set of history elements that correspond to subjets at

bool fastjet::ClusterSequence::object_in_jet ( const PseudoJet object,
const PseudoJet jet 
) const

returns true iff the object is included in the jet.

NB: this is only sensible if the object is already registered within the cluster sequence, so you cannot use it with an input particle to the CS (since the particle won't have the history index set properly).

For nice clustering structures it should run in O(ln(N)) time but in worst cases (certain cone plugins) it can take O(n) time, where n is the number of particles in the jet.

Definition at line 740 of file ClusterSequence.cc.

Referenced by fastjet::ClusterSequenceStructure::object_in_jet().

                                                                 {

  // make sure the object conceivably belongs to this clustering
  // sequence
  assert(contains(object) && contains(jet));

  const PseudoJet * this_object = &object;
  const PseudoJet * childp;
  while(true) {

bool fastjet::ClusterSequence::has_parents ( const PseudoJet jet,
PseudoJet parent1,
PseudoJet parent2 
) const

if the jet has parents in the clustering, it returns true and sets parent1 and parent2 equal to them.

if it has no parents it returns false and sets parent1 and parent2 to zero

Definition at line 766 of file ClusterSequence.cc.

Referenced by fastjet::ClusterSequenceStructure::has_parents().

                                                         {

  const history_element & hist = _history[jet.cluster_hist_index()];

  // make sure we do not run into any unexpected situations --
  // i.e. both parents valid, or neither
  assert ((hist.parent1 >= 0 && hist.parent2 >= 0) || 
          (hist.parent1 < 0 && hist.parent2 < 0));

  if (hist.parent1 < 0) {
    parent1 = PseudoJet(0.0,0.0,0.0,0.0);
    parent2 = parent1;

void fastjet::ClusterSequence::print_jets_for_root ( const std::vector< PseudoJet > &  jets,
std::ostream &  ostr = std::cout 
) const

output the supplied vector of jets in a format that can be read by an appropriate root script; the format is: jet-n jet-px jet-py jet-pz jet-E particle-n particle-rap particle-phi particle-pt particle-n particle-rap particle-phi particle-pt ...

END ... [i.e. above repeated]

double fastjet::ClusterSequence::jet_scale_for_algorithm ( const PseudoJet jet  )  const

returns the scale associated with a jet as required for this clustering algorithm (kt^2 for the kt-algorithm, 1 for the Cambridge algorithm).

[May become virtual at some point]

Definition at line 352 of file ClusterSequence.cc.

                      :
    strategy = "plugin strategy"; break;
  default:
    strategy = "Unrecognized";
  }
  return strategy;
}  


double ClusterSequence::jet_scale_for_algorithm(
                                  const PseudoJet & jet) const {
  if (_jet_algorithm == kt_algorithm)             {return jet.kt2();}
  else if (_jet_algorithm == cambridge_algorithm) {return 1.0;}
  else if (_jet_algorithm == antikt_algorithm) {
    double kt2=jet.kt2();
    return kt2 > 1e-300 ? 1.0/kt2 : 1e300;
  } else if (_jet_algorithm == genkt_algorithm) {
    double kt2 = jet.kt2();
    double p   = jet_def().extra_param();
    if (p <= 0 && kt2 < 1e-300) kt2 = 1e-300; // dodgy safety check

void fastjet::ClusterSequence::plugin_record_iB_recombination ( int  jet_i,
double  diB 
) [inline]

record the fact that there has been a recombination between jets()[jet_i] and the beam, with the specified diB; when looking for inclusive jets, any iB recombination will returned to the user as a jet.

Definition at line 361 of file ClusterSequence.hh.

Referenced by fastjet::PxConePlugin::run_clustering(), fastjet::JadePlugin::run_clustering(), fastjet::D0RunIIConePlugin::run_clustering(), and fastjet::CMSIterativeConePlugin::run_clustering().

                                                             {
    assert(plugin_activated());
    _do_iB_recombination_step(jet_i, diB);
  }

template<class GBJ >
void fastjet::ClusterSequence::plugin_simple_N2_cluster (  )  [inline]

allows a plugin to run a templated clustering (nearest-neighbour heuristic)

This has N^2 behaviour on "good" distance, but a worst case behaviour of N^3 (and many algs trigger the worst case behaviour)

For more details on how this works, see GenBriefJet below

Definition at line 398 of file ClusterSequence.hh.

                                                       {
    assert(plugin_activated());
    _simple_N2_cluster<GBJ>();
  }

static void fastjet::ClusterSequence::set_jet_algorithm ( JetAlgorithm  jet_algorithm  )  [inline, static]

set the default (static) jet finder across all current and future ClusterSequence objects -- deprecated and obsolescent (i.e.

may be suppressed in a future release).

Definition at line 408 of file ClusterSequence.hh.

{_default_jet_algorithm = jet_algorithm;}

const std::vector< PseudoJet > & fastjet::ClusterSequence::jets (  )  const [inline]

allow the user to access the internally stored _jets() array, which contains both the initial particles and the various intermediate and final stages of recombination.

The first n_particles() entries are the original particles, in the order in which they were supplied to the ClusterSequence constructor. It can be useful to access them for example when examining whether a given input object is part of a specific jet, via the objects_in_jet(...) member function (which only takes PseudoJets that are registered in the ClusterSequence).

One of the other (internal uses) is related to the fact because we don't seem to be able to access protected elements of the class for an object that is not "this" (at least in case where "this" is of a slightly different kind from the object, both derived from ClusterSequence).

Definition at line 890 of file ClusterSequence.hh.

References _jets.

Referenced by fastjet::PxConePlugin::run_clustering(), fastjet::JadePlugin::run_clustering(), fastjet::D0RunIIConePlugin::run_clustering(), fastjet::CMSIterativeConePlugin::run_clustering(), and fastjet::ClusterSequenceAreaBase::subtracted_jets().

                                                                 {
  return _jets;
}

const std::vector< ClusterSequence::history_element > & fastjet::ClusterSequence::history (  )  const [inline]

allow the user to access the raw internal history.

This is present (as for jets()) in part so that protected derived classes can access this information about other ClusterSequences.

A user who wishes to follow the details of the ClusterSequence can also make use of this information (and should consult the history_element documentation for more information), but should be aware that these internal structures may evolve in future FastJet versions.

Definition at line 894 of file ClusterSequence.hh.

References _history.

Referenced by fastjet::ClusterSequenceActiveArea::_transfer_ghost_free_history().

                                                                                         {
  return _history;
}

unsigned int fastjet::ClusterSequence::n_particles (  )  const [inline]

returns the number of particles that were provided to the clustering algorithm (helps the user find their way around the history and jets objects if they weren't paying attention beforehand).

Definition at line 898 of file ClusterSequence.hh.

Referenced by unique_history_order().

{return _initial_n;}

vector< int > fastjet::ClusterSequence::unique_history_order (  )  const

routine that returns an order in which to read the history such that clusterings that lead to identical jet compositions but different histories (because of degeneracies in the clustering order) will have matching constituents for each matching entry in the unique_history_order.

The order has the property that an entry's parents will always appear prior to that entry itself.

Roughly speaking the order is such that we first provide all steps that lead to the final jet containing particle 1; then we have the steps that lead to reconstruction of the jet containing the next-lowest-numbered unclustered particle, etc... [see GPS CCN28-12 for more info -- of course a full explanation here would be better...]

Definition at line 1028 of file ClusterSequence.cc.

References n_particles().

                                                        {

  // first construct an array that will tell us the lowest constituent
  // of a given jet -- this will always be one of the original
  // particles, whose order is well defined and so will help us to
  // follow the tree in a unique manner.
  valarray<int> lowest_constituent(_history.size());
  int hist_n = _history.size();
  lowest_constituent = hist_n; // give it a large number
  for (int i = 0; i < hist_n; i++) {
    // sets things up for the initial partons
    lowest_constituent[i] = min(lowest_constituent[i],i); 
    // propagates them through to the children of this parton
    if (_history[i].child > 0) lowest_constituent[_history[i].child] 
      = min(lowest_constituent[_history[i].child],lowest_constituent[i]);
  }

  // establish an array for what we have and have not extracted so far
  valarray<bool> extracted(_history.size()); extracted = false;
  vector<int> unique_tree;
  unique_tree.reserve(_history.size());

  // now work our way through the tree
  for (unsigned i = 0; i < n_particles(); i++) {

vector< PseudoJet > fastjet::ClusterSequence::unclustered_particles (  )  const

return the set of particles that have not been clustered.

For kt and cam/aachen algorithms this should always be null, but for cone type algorithms it can be non-null;

Definition at line 1083 of file ClusterSequence.cc.

Referenced by fastjet::ClusterSequenceActiveAreaExplicitGhosts::empty_area().

bool fastjet::ClusterSequence::contains ( const PseudoJet object  )  const

returns true if the object (jet or particle) is contained by (ie belongs to) this cluster sequence.

Tests performed: if thejet's interface is this cluster sequence and its cluster history index is in a consistent range.

Definition at line 1098 of file ClusterSequence.cc.

void fastjet::ClusterSequence::transfer_from_sequence ( ClusterSequence from_seq,
bool  transfer_ownership = true 
)

transfer the sequence contained in other_seq into our own; any plugin "extras" contained in the from_seq will be lost from there.

If transfer_ownership is true, it also sets the ClusterSequence pointers of the PseudoJets in the history to point to this ClusterSequence (true by default)

Definition at line 378 of file ClusterSequence.cc.

           {return 1.0;}
  } else {throw Error("Unrecognised jet algorithm");}
}


//----------------------------------------------------------------------
/// transfer the sequence contained in other_seq into our own;
/// any plugin "extras" contained in the from_seq will be lost
/// from there.
void ClusterSequence::transfer_from_sequence(ClusterSequence & from_seq, bool transfer_ownership) {

  // the metadata
  _jet_def                 = from_seq._jet_def                ;
  _writeout_combinations   = from_seq._writeout_combinations  ;
  _initial_n               = from_seq._initial_n              ;
  _Rparam                  = from_seq._Rparam                 ;
  _R2                      = from_seq._R2                     ;
  _invR2                   = from_seq._invR2                  ;
  _strategy                = from_seq._strategy               ;
  _jet_algorithm           = from_seq._jet_algorithm          ;
  _plugin_activated        = from_seq._plugin_activated       ;

  // the data
  _jets     = from_seq._jets;
  _history  = from_seq._history;
  // the following transferse ownership of the extras from the from_seq
  _extras   = from_seq._extras;

const SharedPtr<PseudoJetStructureBase>& fastjet::ClusterSequence::structure_shared_ptr (  )  const [inline]

retrieve a shared pointer to the wrapper to this ClusterSequence

this may turn useful if you want to track when this ClusterSequence goes out of scope

Definition at line 535 of file ClusterSequence.hh.

Referenced by fastjet::BackgroundEstimator::BackgroundEstimator().

                                                                        {
    return _structure_shared_ptr;
  }

template<class L >
void fastjet::ClusterSequence::_transfer_input_jets ( const std::vector< L > &  pseudojets  )  [protected]

transfer the vector<L> of input jets into our own vector<PseudoJet> _jets (with some reserved space for future growth).

Definition at line 841 of file ClusterSequence.hh.

References _jets.

Referenced by ClusterSequence().

                                                                        {

  // this will ensure that we can point to jets without difficulties
  // arising.
  _jets.reserve(pseudojets.size()*2);

  // insert initial jets this way so that any type L that can be
  // converted to a pseudojet will work fine (basically PseudoJet
  // and any type that has [] subscript access to the momentum
  // components, such as CLHEP HepLorentzVector).
  for (unsigned int i = 0; i < pseudojets.size(); i++) {
    _jets.push_back(pseudojets[i]);}
  
}

void fastjet::ClusterSequence::_initialise_and_run ( const JetDefinition jet_def,
const bool &  writeout_combinations 
) [protected]

This is the routine that will do all the initialisation and then run the clustering (may be called by various constructors).

It assumes _jets contains the momenta to be clustered.

Definition at line 72 of file ClusterSequence.cc.

References jet_def().

Referenced by ClusterSequence().

                                                                      {

  JetDefinition jet_def(_default_jet_algorithm, R, strategy);
  _initialise_and_run(jet_def, writeout_combinations);
}


//----------------------------------------------------------------------
void ClusterSequence::_initialise_and_run (
                                  const JetDefinition & jet_def,
                                  const bool & writeout_combinations) {

  // transfer all relevant info into internal variables
  _decant_options(jet_def, writeout_combinations);

  // set up the history entries for the initial particles (those
  // currently in _jets)
  _fill_initial_history();

  // don't run anything if the event is empty
  if (n_particles() == 0) return;

  // ----- deal with special cases: plugins & e+e- ------
  if (_jet_algorithm == plugin_algorithm) {
    // allows plugin_xyz() functions to modify cluster sequence
    _plugin_activated = true;
    // let the plugin do its work here
    _jet_def.plugin()->run_clustering( (*this) );
    _plugin_activated = false;
    return;
  } else if (_jet_algorithm == ee_kt_algorithm ||
             _jet_algorithm == ee_genkt_algorithm) {
    // ignore requested strategy
    _strategy = N2Plain;
    if (_jet_algorithm == ee_kt_algorithm) {
      // make sure that R is large enough so that "beam" recomb only
      // occurs when a single particle is left
      // Normally, this should be automatically set to 4 from JetDefinition
      assert(_Rparam > 2.0); 
      // this is used to renormalise the dij to get a "standard" form
      // and our convention in e+e- will be different from that
      // in long.inv case; NB: _invR2 name should be changed -> _renorm_dij?
      _invR2 = 1.0;
    } else {
      // as of 2009-01-09, choose R to be an angular distance, in
      // radians.  Since the algorithm uses 2(1-cos(theta)) as its
      // squared angular measure, make sure that the _R2 is defined
      // in a similar way.
      if (_Rparam > pi) {
        // choose a value that ensures that back-to-back particles will
        // always recombine 
        //_R2 = 4.0000000000001;
        _R2 = 2 * ( 3.0 + cos(_Rparam) );
      } else {
        _R2    = 2 * ( 1.0 - cos(_Rparam) );
      }
      _invR2 = 1.0/_R2;
    }
    //_simple_N2_cluster<EEBriefJet>();
    _simple_N2_cluster_EEBriefJet();
    return;
  }


  // automatically redefine the strategy according to N if that is
  // what the user requested -- transition points (and especially
  // their R-dependence) are based on empirical observations for a
  // R=0.4, 0.7 and 1.0, running on toth (3.4GHz, Pentium IV D [dual
  // core] with 2MB of cache).
  if (_strategy == Best) {
    int N = _jets.size();
    if (N <= 55*max(0.5,min(1.0,_Rparam))) {// empirical scaling with R
      _strategy = N2Plain;
    } else if (N > 6200/pow(_Rparam,2.0) && jet_def.jet_algorithm() == cambridge_algorithm) {
      _strategy = NlnNCam;
#ifndef DROP_CGAL
    } else if ((N > 16000/pow(_Rparam,1.15) && jet_def.jet_algorithm() != antikt_algorithm)
               || N > 35000/pow(_Rparam,1.15)) {
      _strategy = NlnN;
#endif  // DROP_CGAL
    } else if (N <= 450) {
      _strategy = N2Tiled;
    } else {                   
      _strategy = N2MinHeapTiled;
    }
  }

  // R >= 2pi is not supported by all clustering strategies owing to
  // periodicity issues (a particle might cluster with itself). When
  // R>=2pi, we therefore automatically switch to a strategy that is
  // known to work.
  if (_Rparam >= twopi) {
    if (   _strategy == NlnN
        || _strategy == NlnN3pi
        || _strategy == NlnNCam
        || _strategy == NlnNCam2pi2R
        || _strategy == NlnNCam4pi) {
#ifdef DROP_CGAL
      _strategy = N2MinHeapTiled;
#else
      _strategy = NlnN4pi;
#endif    
    }
    if (jet_def.strategy() != Best && _strategy != jet_def.strategy()) {
      ostringstream oss;
      oss << "Cluster strategy " << strategy_string(jet_def.strategy())
          << " automatically changed to " << strategy_string()
          << " because the former is not supported for R = " << _Rparam
          << " >= 2pi";
      _changed_strategy_warning.warn(oss.str());
    }
  }


  // run the code containing the selected strategy
  // 
  // We order the strategies stqrting from the ones used by the Best
  // strategy in the order of increasing N, then the remaining ones
  // again in the order of increasing N.
  if (_strategy == N2Plain) {
    // BriefJet provides standard long.invariant kt alg.
    this->_simple_N2_cluster_BriefJet();
  } else if (_strategy == N2Tiled) {
    this->_faster_tiled_N2_cluster();
  } else if (_strategy == N2MinHeapTiled) {
    this->_minheap_faster_tiled_N2_cluster();
  } else if (_strategy == NlnN) {
    this->_delaunay_cluster();
  } else if (_strategy == NlnNCam) {
    this->_CP2DChan_cluster_2piMultD();
  } else if (_strategy == NlnN3pi || _strategy == NlnN4pi ) {
    this->_delaunay_cluster();
  } else if (_strategy ==  N3Dumb ) {
    this->_really_dumb_cluster();
  } else if (_strategy == N2PoorTiled) {
    this->_tiled_N2_cluster();
  } else if (_strategy == NlnNCam4pi) {

void fastjet::ClusterSequence::_initialise_and_run ( const double &  R,
const Strategy strategy,
const bool &  writeout_combinations 
) [protected]

This is an alternative routine for initialising and running the clustering, provided for legacy purposes.

The jet finder is that specified in the static member _default_jet_algorithm.

Definition at line 61 of file ClusterSequence.cc.

void fastjet::ClusterSequence::_do_ij_recombination_step ( const int &  jet_i,
const int &  jet_j,
const double &  dij,
int &  newjet_k 
) [protected]

carry out the recombination between the jets numbered jet_i and jet_j, at distance scale dij; return the index newjet_k of the result of the recombination of i and j.

carries out the bookkeeping associated with the step of recombining jet_i and jet_j (assuming a distance dij) and returns the index of the recombined jet, newjet_k.

Definition at line 1140 of file ClusterSequence.cc.

Referenced by fastjet::ClusterSequenceActiveArea::_transfer_ghost_free_history().

                                               {

  // create the new jet by recombining the first two
  PseudoJet newjet;
  _jet_def.recombiner()->recombine(_jets[jet_i], _jets[jet_j], newjet);
  _jets.push_back(newjet);
  // original version...
  //_jets.push_back(_jets[jet_i] + _jets[jet_j]);

  // get its index
  newjet_k = _jets.size()-1;

  // get history index
  int newstep_k = _history.size();
  // and provide jet with the info
  _jets[newjet_k].set_cluster_hist_index(newstep_k);

void fastjet::ClusterSequence::_do_iB_recombination_step ( const int &  jet_i,
const double &  diB 
) [protected]

carry out an recombination step in which _jets[jet_i] merges with the beam,

carries out the bookkeeping associated with the step of recombining jet_i with the beam

Definition at line 1173 of file ClusterSequence.cc.

Referenced by fastjet::ClusterSequenceActiveArea::_transfer_ghost_free_history().

                                                                         {


Member Data Documentation

std::vector<PseudoJet> fastjet::ClusterSequence::_jets [protected]

This contains the physical PseudoJets; for each PseudoJet one can find the corresponding position in the _history by looking at _jets[i].cluster_hist_index().

Definition at line 589 of file ClusterSequence.hh.

Referenced by _transfer_input_jets(), and jets().

this vector will contain the branching history; for each stage, _history[i].jetp_index indicates where to look in the _jets vector to get the physical PseudoJet.

Definition at line 595 of file ClusterSequence.hh.

Referenced by fastjet::ClusterSequenceActiveArea::_transfer_ghost_free_history(), and history().


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