Page MenuHomeHEPForge

example.cc
No OneTemporary

example.cc

// To run this example, use the following command:
//
// ./example < ../data/single_event.dat
//
//----------------------------------------------------------------------
// $Id$
//
// Copyright (c) 2016,
// Jesse Thaler
//
//----------------------------------------------------------------------
// This file is part of FastJet contrib.
//
// It is free software; you can redistribute it and/or modify it under
// the terms of the GNU General Public License as published by the
// Free Software Foundation; either version 2 of the License, or (at
// your option) any later version.
//
// It is distributed in the hope that it will be useful, but WITHOUT
// ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
// or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public
// License for more details.
//
// You should have received a copy of the GNU General Public License
// along with this code. If not, see <http://www.gnu.org/licenses/>.
//----------------------------------------------------------------------
#include <iostream>
#include <sstream>
#include <iomanip>
#include "fastjet/PseudoJet.hh"
#include <sstream>
#include "CartesianJet.hh" // In external code, this should be fastjet/contrib/CartesianJet.hh
using namespace std;
using namespace fastjet;
using namespace fastjet::contrib;
// forward declaration to make things clearer
void read_event(vector<PseudoJet> &event);
//----------------------------------------------------------------------
int main(){
//----------------------------------------------------------
// read in input particles
vector<PseudoJet> event;
read_event(event);
cout << "# read an event with " << event.size() << " particles" << endl;
//----------------------------------------------------------
// illustrate how this CartesianJet contrib works
// Reading in CartesianCoordinates
vector<CartesianCoordinate> coords_input;
double random_scale_factor = 54.8; // Randomly scale up to check rescaling properties.
// For cross check, read in PseudoJets
vector<PseudoJet> pjs_input;
for (int i = 0; i < event.size(); i++) {
// From file, get x and y coordinates, rescale by 10.0 to avoid wrapping issue.
double x = event[i].rap()/10.0;
double y = event[i].phi()/10.0;
double intensity = event[i].pt();
int user_index = i;
coords_input.push_back(CartesianCoordinate(random_scale_factor*x,random_scale_factor*y,intensity,user_index));
// Make pseudojets for comparison
PseudoJet pj;
pj.reset_PtYPhiM(intensity,x,y);
pj.set_user_index(user_index);
pjs_input.push_back(pj);
}
double jet_radius = 0.1;
// Doing clustering with CartesianJetDefinition
CartesianJetDefinition jet_def(random_scale_factor*jet_radius);
vector<CartesianJet> cjs_output = jet_def(coords_input);
// Test that I get the same answer with antikt with boost invariant pt scheme
JetDefinition antikt(antikt_algorithm,jet_radius,BIpt_scheme);
vector<PseudoJet> pjs_output = antikt(pjs_input);
// Make sure there is the same number of jets
assert(cjs_output.size() == pjs_output.size());
cout << "Jet X Y Int. Const." << endl;
for (int i = 0; i < cjs_output.size() ; i++) {
CartesianJet my_jet = cjs_output[i];
cout << fixed << setprecision(2)
<< setw(3) << (i + 1) << " "
<< setw(7) << my_jet.x() << " "
<< setw(7) << my_jet.y() << " "
<< setw(7) << my_jet.intensity() << " "
<< setw(7) << my_jet.cartesian_constituents().size() << endl;
// Check for match in pTs
assert(std::abs(cjs_output[i].intensity() - pjs_output[i].pt()) < 10e-5);
}
return 0;
}
// read in input particles
void read_event(vector<PseudoJet> &event){
string line;
while (getline(cin, line)) {
istringstream linestream(line);
// take substrings to avoid problems when there is extra "pollution"
// characters (e.g. line-feed).
if (line.substr(0,4) == "#END") {return;}
if (line.substr(0,1) == "#") {continue;}
double px,py,pz,E;
linestream >> px >> py >> pz >> E;
PseudoJet particle(px,py,pz,E);
// push event onto back of full_event vector
event.push_back(particle);
}
}

File Metadata

Mime Type
text/x-c
Expires
Sat, May 3, 5:57 AM (1 d, 3 h)
Storage Engine
blob
Storage Format
Raw Data
Storage Handle
4982889
Default Alt Text
example.cc (4 KB)

Event Timeline