fastjet 2.4.5
Functions
fastjet_areas.cc File Reference
#include "fastjet/PseudoJet.hh"
#include "fastjet/ClusterSequenceArea.hh"
#include "fastjet/ClusterSequencePassiveArea.hh"
#include "fastjet/config.h"
#include "fastjet/SISConePlugin.hh"
#include <iostream>
#include <sstream>
#include <vector>
Include dependency graph for fastjet_areas.cc:

Go to the source code of this file.

Functions

void print_jets (const fastjet::ClusterSequenceAreaBase &clust_seq, const vector< fastjet::PseudoJet > &unsorted_jets)
 a function that pretty prints a list of jets
int main (int argc, char **argv)
 an example program showing how to use fastjet

Function Documentation

int main ( int  argc,
char **  argv 
)

an example program showing how to use fastjet

Definition at line 67 of file fastjet_areas.cc.

References Best, fastjet::AreaDefinition::description(), fastjet::JetDefinition::description(), fastjet::ClusterSequence::inclusive_jets(), fastjet::kt_algorithm, fastjet::passive_area, print_jets(), fastjet::ClusterSequence::strategy_string(), and fastjet::ClusterSequence::unclustered_particles().

                                  {
  
  vector<fastjet::PseudoJet> input_particles;
  
  // read in input particles
  double px, py , pz, E;
  while (cin >> px >> py >> pz >> E) {
    // create a fastjet::PseudoJet with these components and put it onto
    // back of the input_particles vector
    input_particles.push_back(fastjet::PseudoJet(px,py,pz,E)); 
  }

  // create an object that represents your choice of jet algorithm, and 
  // the associated parameters
  double Rparam = 1.0;
  fastjet::Strategy strategy = fastjet::Best;
  fastjet::JetDefinition jet_def(fastjet::kt_algorithm, Rparam, strategy);
  //fastjet::JetDefinition jet_def(fastjet::cambridge_algorithm, Rparam, strategy);
  //fastjet::JetDefinition jet_def(fastjet::antikt_algorithm, Rparam, strategy);
  //fastjet::JetDefinition jet_def(new fastjet::SISConePlugin(Rparam,0.75));

  // create an object that specifies how we to define the area
  fastjet::AreaDefinition area_def;
  bool use_voronoi = false;
  if (!use_voronoi) {
    double ghost_etamax = 7.0;
    double ghost_area    = 0.05;
    int    active_area_repeats = 1;
    // alternative settings for more precision:
    // reducing ghost area gives better sensitivity to the exact edges of the jet
    //double ghost_area    = 0.01;
    // increasing the repeats is useful in sparse events
    //int    active_area_repeats = 100; 

    // now create the object that holds info about ghosts, and from that
    // get an area definition
    fastjet::GhostedAreaSpec ghost_spec(ghost_etamax, active_area_repeats, 
                                        ghost_area);
    area_def = fastjet::AreaDefinition(fastjet::passive_area,ghost_spec);
    //area_def = fastjet::AreaDefinition(fastjet::active_area,ghost_spec);
  } else {
    double effective_Rfact = 1.0;
    area_def = fastjet::VoronoiAreaSpec(effective_Rfact);
  }

  // run the jet clustering with the above jet definition
  fastjet::ClusterSequenceArea clust_seq(input_particles, 
                                             jet_def, area_def);
  // you can also run the individual area classes directly
  //fastjet::ClusterSequencePassiveArea clust_seq(input_particles, jet_def, 
  //                                              area_def.ghost_spec());

  // you may want to find out how much area in a given range (|y|<range)
  // is empty of real jets (or corresponds to pure "ghost" jets).
  //double range = 4.0;
  //cout << clust_seq.empty_area(range) << endl;
  //cout << clust_seq.n_empty_jets(range) << endl;

  // tell the user what was done
  cout << "Jet definition was: " << jet_def.description() << endl;
  cout << "Area definition was: " << area_def.description() << endl;
  cout << "Strategy adopted by FastJet was "<<
       clust_seq.strategy_string()<<endl<<endl;

  // extract the inclusive jets with pt > 5 GeV, sorted by pt
  double ptmin = 5.0;
  vector<fastjet::PseudoJet> inclusive_jets = clust_seq.inclusive_jets(ptmin);

  // print them out
  cout << "Printing inclusive jets with pt > "<< ptmin<<" GeV\n";
  cout << "---------------------------------------\n";
  print_jets(clust_seq, inclusive_jets);
  cout << endl;

  
  cout << "Number of unclustered particles: " 
       << clust_seq.unclustered_particles().size() << endl;


}
void print_jets ( const fastjet::ClusterSequenceAreaBase clust_seq,
const vector< fastjet::PseudoJet > &  unsorted_jets 
)

a function that pretty prints a list of jets

a function that pretty prints a list of jets, and performs the subtraction in two different ways, using a generic ClusterSequenceAreaBase type object.

Definition at line 151 of file fastjet_areas.cc.

References fastjet::ClusterSequenceAreaBase::area(), fastjet::ClusterSequenceAreaBase::area_error(), fastjet::d0::inline_maths::phi(), and fastjet::sorted_by_pt().

Referenced by main(), and run_jet_finder().

                                                                 {

  // sort jets into increasing pt
  vector<fastjet::PseudoJet> jets = sorted_by_pt(unsorted_jets);  

  printf(" ijet   rap      phi        Pt         area  +-   err\n");
  for (size_t j = 0; j < jets.size(); j++) {

    double area     = clust_seq.area(jets[j]);
    double area_error = clust_seq.area_error(jets[j]);

    printf("%5u %9.5f %8.5f %10.3f %8.3f +- %6.3f\n",j,jets[j].rap(),
           jets[j].phi(),jets[j].perp(), area, area_error);
  }


}
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines