fastjet 2.4.3
Functions

PseudoJet.cc File Reference

#include "fastjet/Error.hh"
#include "fastjet/PseudoJet.hh"
#include <valarray>
#include <iostream>
#include <sstream>
#include <cmath>
#include <algorithm>
Include dependency graph for PseudoJet.cc:

Go to the source code of this file.

Functions

PseudoJet operator+ (const PseudoJet &jet1, const PseudoJet &jet2)
PseudoJet operator- (const PseudoJet &jet1, const PseudoJet &jet2)
PseudoJet operator* (double coeff, const PseudoJet &jet)
PseudoJet operator* (const PseudoJet &jet, double coeff)
PseudoJet operator/ (const PseudoJet &jet, double coeff)
bool have_same_momentum (const PseudoJet &jeta, const PseudoJet &jetb)
 returns true if the momenta of the two input jets are identical
PseudoJet PtYPhiM (double pt, double y, double phi, double m)
 return a pseudojet with the given pt, y, phi and mass
void sort_indices (vector< int > &indices, const vector< double > &values)
template<class T >
vector< T > objects_sorted_by_values (const vector< T > &objects, const vector< double > &values)
 given a vector of values with a one-to-one correspondence with the vector of objects, sort objects into an order such that the associated values would be in increasing order
vector< PseudoJet > sorted_by_pt (const vector< PseudoJet > &jets)
 return a vector of jets sorted into decreasing kt2
vector< PseudoJet > sorted_by_rapidity (const vector< PseudoJet > &jets)
 return a vector of jets sorted into increasing rapidity
vector< PseudoJet > sorted_by_E (const vector< PseudoJet > &jets)
 return a vector of jets sorted into decreasing energy
vector< PseudoJet > sorted_by_pz (const vector< PseudoJet > &jets)
 return a vector of jets sorted into increasing pz

Function Documentation

bool have_same_momentum ( const PseudoJet &  jeta,
const PseudoJet &  jetb 
)

returns true if the momenta of the two input jets are identical

Definition at line 284 of file PseudoJet.cc.

Referenced by SISConeSphericalPlugin::run_clustering(), and SISConePlugin::run_clustering().

                                                                        {
  return jeta.px() == jetb.px()
    &&   jeta.py() == jetb.py()
    &&   jeta.pz() == jetb.pz()
    &&   jeta.E()  == jetb.E();
}
template<class T >
vector<T> objects_sorted_by_values ( const vector< T > &  objects,
const vector< double > &  values 
)

given a vector of values with a one-to-one correspondence with the vector of objects, sort objects into an order such that the associated values would be in increasing order

Definition at line 345 of file PseudoJet.cc.

References sort_indices().

Referenced by sorted_by_E(), sorted_by_pt(), sorted_by_pz(), and sorted_by_rapidity().

                                                      {

  assert(objects.size() == values.size());

  // get a vector of indices
  vector<int> indices(values.size());
  for (size_t i = 0; i < indices.size(); i++) {indices[i] = i;}
  
  // sort the indices
  sort_indices(indices, values);
  
  // copy the objects 
  vector<T> objects_sorted(objects.size());
  
  // place the objects in the correct order
  for (size_t i = 0; i < indices.size(); i++) {
    objects_sorted[i] = objects[indices[i]];
  }

  return objects_sorted;
}
PseudoJet operator* ( const PseudoJet &  jet,
double  coeff 
)

Definition at line 178 of file PseudoJet.cc.

                                                          {
  return coeff*jet;
} 
PseudoJet operator* ( double  coeff,
const PseudoJet &  jet 
)

Definition at line 168 of file PseudoJet.cc.

                                                          {
  //return PseudoJet(coeff*jet.four_mom());
  // the following code is hopefully more efficient
  PseudoJet coeff_times_jet(jet);
  coeff_times_jet *= coeff;
  return coeff_times_jet;
} 
PseudoJet operator+ ( const PseudoJet &  jet1,
const PseudoJet &  jet2 
)

Definition at line 148 of file PseudoJet.cc.

                                                                     {
  //return PseudoJet(jet1.four_mom()+jet2.four_mom());
  return PseudoJet(jet1.px()+jet2.px(),
                   jet1.py()+jet2.py(),
                   jet1.pz()+jet2.pz(),
                   jet1.E() +jet2.E()  );
} 
PseudoJet operator- ( const PseudoJet &  jet1,
const PseudoJet &  jet2 
)

Definition at line 158 of file PseudoJet.cc.

                                                                     {
  //return PseudoJet(jet1.four_mom()-jet2.four_mom());
  return PseudoJet(jet1.px()-jet2.px(),
                   jet1.py()-jet2.py(),
                   jet1.pz()-jet2.pz(),
                   jet1.E() -jet2.E()  );
} 
PseudoJet operator/ ( const PseudoJet &  jet,
double  coeff 
)

Definition at line 184 of file PseudoJet.cc.

                                                          {
  return (1.0/coeff)*jet;
} 
PseudoJet PtYPhiM ( double  pt,
double  y,
double  phi,
double  m 
)

return a pseudojet with the given pt, y, phi and mass

Definition at line 294 of file PseudoJet.cc.

                                                             {
  double ptm = sqrt(pt*pt+m*m);
  return PseudoJet(pt*cos(phi), pt*sin(phi), ptm*sinh(y), ptm*cosh(y));
}
void sort_indices ( vector< int > &  indices,
const vector< double > &  values 
)

Definition at line 335 of file PseudoJet.cc.

Referenced by objects_sorted_by_values(), and CDFJetCluPlugin::run_clustering().

                                                        {
  IndexedSortHelper index_sort_helper(&values);
  sort(indices.begin(), indices.end(), index_sort_helper);
}
vector<PseudoJet> sorted_by_E ( const vector< PseudoJet > &  jets)

return a vector of jets sorted into decreasing energy

Definition at line 387 of file PseudoJet.cc.

References objects_sorted_by_values().

Referenced by main(), and print_jets().

                                                              {
  vector<double> energies(jets.size());
  for (size_t i = 0; i < jets.size(); i++) {energies[i] = -jets[i].E();}
  return objects_sorted_by_values(jets, energies);
}
vector<PseudoJet> sorted_by_pt ( const vector< PseudoJet > &  jets)

return a vector of jets sorted into decreasing kt2

Definition at line 371 of file PseudoJet.cc.

References objects_sorted_by_values().

Referenced by main(), print_jets(), print_jets_and_sub(), and ClusterSequenceAreaBase::subtracted_jets().

                                                               {
  vector<double> minus_kt2(jets.size());
  for (size_t i = 0; i < jets.size(); i++) {minus_kt2[i] = -jets[i].kt2();}
  return objects_sorted_by_values(jets, minus_kt2);
}
vector<PseudoJet> sorted_by_pz ( const vector< PseudoJet > &  jets)

return a vector of jets sorted into increasing pz

Definition at line 395 of file PseudoJet.cc.

References objects_sorted_by_values().

                                                               {
  vector<double> pz(jets.size());
  for (size_t i = 0; i < jets.size(); i++) {pz[i] = jets[i].pz();}
  return objects_sorted_by_values(jets, pz);
}
vector<PseudoJet> sorted_by_rapidity ( const vector< PseudoJet > &  jets)

return a vector of jets sorted into increasing rapidity

Definition at line 379 of file PseudoJet.cc.

References objects_sorted_by_values().

Referenced by main().

                                                                     {
  vector<double> rapidities(jets.size());
  for (size_t i = 0; i < jets.size(); i++) {rapidities[i] = jets[i].rap();}
  return objects_sorted_by_values(jets, rapidities);
}