40 #include "fastjet/PseudoJet.hh"
41 #include "fastjet/ClusterSequenceArea.hh"
42 #include "fastjet/Selector.hh"
43 #include "fastjet/tools/JetMedianBackgroundEstimator.hh"
44 #include "fastjet/tools/Subtractor.hh"
46 #include <fastjet/config.h>
59 vector<PseudoJet> hard_event, full_event;
65 double particle_maxrap = 5.0;
69 while (getline(cin, line)) {
70 istringstream linestream(line);
73 if (line.substr(0,4) ==
"#END") {
break;}
74 if (line.substr(0,9) ==
"#SUBSTART") {
76 if (nsub == 1) hard_event = full_event;
79 if (line.substr(0,1) ==
"#") {
continue;}
81 linestream >> px >> py >> pz >> E;
86 if (abs(particle.rap()) <= particle_maxrap) full_event.push_back(particle);
90 if (nsub == 1) hard_event = full_event;
94 cerr <<
"Error: read empty event\n";
110 double ghost_maxrap = 6.0;
125 vector<PseudoJet> hard_jets =
sorted_by_pt(clust_seq_hard.inclusive_jets(ptmin));
126 vector<PseudoJet> full_jets =
sorted_by_pt(clust_seq_full.inclusive_jets(ptmin));
181 #if FASTJET_VERSION_NUMBER >= 30100
182 subtractor.set_use_rho_m(
true);
183 subtractor.set_safe_mass(
true);
192 bkgd_estimator.set_particles(full_event);
199 cout <<
"Main clustering:" << endl;
200 cout <<
" Ran: " << jet_def.description() << endl;
201 cout <<
" Area: " << area_def.description() << endl;
202 cout <<
" Particles up to |y|=" << particle_maxrap << endl;
205 cout <<
"Background estimation:" << endl;
206 cout <<
" " << bkgd_estimator.description() << endl << endl;;
207 cout <<
" Giving, for the full event" << endl;
209 cout <<
" rho = " << bkgd_estimate.
rho() << endl;
210 cout <<
" sigma = " << bkgd_estimate.
sigma() << endl;
211 cout <<
" rho_m = " << bkgd_estimate.
rho_m() << endl;
212 cout <<
" sigma_m = " << bkgd_estimate.
sigma_m() << endl;
215 cout <<
"Jets above " << ptmin <<
" GeV in the hard event (" << hard_event.size() <<
" particles)" << endl;
216 cout <<
"---------------------------------------\n";
217 printf(
"%5s %15s %15s %15s %15s %15s\n",
"jet #",
"rapidity",
"phi",
"pt",
"m",
"area");
218 for (
unsigned int i = 0; i < hard_jets.size(); i++) {
219 printf(
"%5u %15.8f %15.8f %15.8f %15.8f %15.8f\n", i,
220 hard_jets[i].rap(), hard_jets[i].phi(), hard_jets[i].pt(), hard_jets[i].m(),
221 hard_jets[i].area());
231 cout <<
"Jets above " << ptmin <<
" GeV in the full event (" << full_event.size() <<
" particles)" << endl;
232 cout <<
"---------------------------------------\n";
233 printf(
"%5s %15s %15s %15s %15s %15s %15s %15s %15s %15s\n",
"jet #",
"rapidity",
"phi",
"pt",
"m",
"area",
"rap_sub",
"phi_sub",
"pt_sub",
"m_sub");
237 vector<PseudoJet> subtracted_jets = subtractor(full_jets);
239 for (
unsigned int i=0; i<full_jets.size(); i++){
241 if (subtracted_jets[i].pt2() >= ptmin*ptmin){
242 printf(
"%5u %15.8f %15.8f %15.8f %15.8f %15.8f %15.8f %15.8f %15.8f %15.8f\n", idx,
243 full_jets[i].rap(), full_jets[i].phi(), full_jets[i].pt(), full_jets[i].m(),
245 subtracted_jets[i].rap(), subtracted_jets[i].phi(),
246 subtracted_jets[i].pt(),
247 subtracted_jets[i].m());
int main()
an example program showing how to use fastjet
class that holds a generic area definition
/// a class that holds the result of the calculation
double sigma_m() const
fluctuations in the purely longitudinal (particle-mass-induced) component of the background density p...
double rho() const
background density per unit area
double sigma() const
background fluctuations per unit square-root area must be multipled by sqrt(area) to get fluctuations...
double rho_m() const
purely longitudinal (particle-mass-induced) component of the background density per unit area
General class for user to obtain ClusterSequence with additional area information.
Parameters to configure the computation of jet areas using ghosts.
class that is intended to hold a full definition of the jet clusterer
Class to contain pseudojets, including minimal information of use to jet-clustering routines.
Class that encodes information about cuts and other selection criteria that can be applied to PseudoJ...
Class that helps perform jet background subtraction.
Selector SelectorNHardest(unsigned int n)
select the n hardest objects
Selector SelectorAbsRapMax(double absrapmax)
select objects with |rap| <= absrapmax
vector< PseudoJet > sorted_by_pt(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing kt2