//---------------------------------------------------------------------- // Example on how to use this contribution // // run it with // ./example_areasub < ../data/Pythia-Zp2jets-lhc-pileup-1ev.dat //---------------------------------------------------------------------- // $Id: example.cc 3001 2013-01-29 10:41:40Z soyez $ // // Copyright (c) 2012-, Matteo Cacciari, Jihun Kim, Gavin P. Salam and Gregory Soyez // //---------------------------------------------------------------------- // 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 . //---------------------------------------------------------------------- //---------------------------------------------------------- // GS notes for CHS events //---------------------------------------------------------- // Is we want to show an example of what happens for CHS events. we // need to introduce a UserInfo class similar to what is in example 09 // in FastJet. For area wsubtraction, we'd just discard charged-PU // tracks and that's mostly it (watching out that the Selector needs // to be passed to the subtractor to properly handle the safety // tests). But for NpC, we'd need a transformer (unless we use a loop) // that rescales aprticles. // // Question: do we do that at all? Only Area? In this example or in // another? #include #include #include "fastjet/PseudoJet.hh" #include "fastjet/ClusterSequenceArea.hh" #include "fastjet/Selector.hh" #include "fastjet/tools/JetMedianBackgroundEstimator.hh" //#include "fastjet/tools/Subtractor.hh" #include "SafeSubtractor.hh" // we'll use the implementation provided // in this contrib #include "GenericSubtractor.hh" using namespace std; using namespace fastjet; void read_event(vector &hard_event, vector &full_event); void print_jet(const PseudoJet &jet); //---------------------------------------------------------------------- int main(){ // read in input particles // // since we use here simulated data we can split the hard event // from the full (i.e. with pileup added) one // // (see also example 07 in FastJet) //---------------------------------------------------------- vector hard_event, full_event; read_event(hard_event, full_event); // keep the particles up to 4 units in rapidity hard_event = SelectorAbsRapMax(4.0)(hard_event); full_event = SelectorAbsRapMax(4.0)(full_event); // do the clustering and get the jets //---------------------------------------------------------- JetDefinition jet_def(antikt_algorithm, 0.7); AreaDefinition area_def(active_area_explicit_ghosts, GhostedAreaSpec(SelectorAbsRapMax(4.0))); ClusterSequenceArea clust_seq_hard(hard_event, jet_def, area_def); ClusterSequenceArea clust_seq_full(full_event, jet_def, area_def); Selector sel_jets = SelectorNHardest(2) * SelectorAbsRapMax(3.0); vector hard_jets = sel_jets(clust_seq_hard.inclusive_jets()); vector full_jets = sel_jets(clust_seq_full.inclusive_jets()); // print the results without subtraction //---------------------------------------------------------- cout << setprecision(4); cout << "# original hard jets" << endl; for (unsigned int i=0; i &hard_event, vector &full_event){ string line; int nsub = 0; // counter to keep track of which sub-event we're reading while (getline(cin, line)) { istringstream linestream(line); // take substrings to avoid problems when there are extra "pollution" // characters (e.g. line-feed). if (line.substr(0,4) == "#END") {break;} if (line.substr(0,9) == "#SUBSTART") { // if more sub events follow, make copy of first one (the hard one) here if (nsub == 1) hard_event = full_event; nsub += 1; } 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 full_event.push_back(particle); } // if we have read in only one event, copy it across here... if (nsub == 1) hard_event = full_event; // if there was nothing in the event if (nsub == 0) { cerr << "Error: read empty event\n"; exit(-1); } cout << "# " << nsub-1 << " pileup events on top of the hard event" << endl; } //------------------------------------------------------------------------ // print iut some info about a giben jet void print_jet(const PseudoJet &jet){ cout << "pt = " << jet.pt() << ", rap = " << jet.rap() << ", mass = " << jet.m() << endl; }