Index: contrib/contribs/LundPlane/trunk/example_dpsi_slice.cc =================================================================== --- contrib/contribs/LundPlane/trunk/example_dpsi_slice.cc (revision 1329) +++ contrib/contribs/LundPlane/trunk/example_dpsi_slice.cc (revision 1330) @@ -1,173 +1,174 @@ //---------------------------------------------------------------------- /// \file example_dpsi_slice.cc /// /// This example program is meant to illustrate how the /// fastjet::contrib::RecursiveLundEEGenerator class is used. /// /// Run this example with /// /// \verbatim /// ./example_dpsi_slice < ../data/single-ee-event.dat /// \endverbatim //---------------------------------------------------------------------- // $Id$ // // Copyright (c) 2018-, Frederic A. Dreyer, Keith Hamilton, Alexander Karlberg, // Gavin P. Salam, Ludovic Scyboz, Gregory Soyez, Rob Verheyen // //---------------------------------------------------------------------- // 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 . //---------------------------------------------------------------------- #include #include #include #include "fastjet/PseudoJet.hh" #include "fastjet/EECambridgePlugin.hh" #include #include "RecursiveLundEEGenerator.hh" // In external code, this should be fastjet/contrib/RecursiveLundEEGenerator.hh using namespace std; using namespace fastjet; // Definitions for the slice observable (ymax = 1, zcut = 0.1) double abs_rap_slice = 1; double z2_cut = 0.1; // forward declaration to make things clearer void read_event(vector &event); // returns true if PseudoJet p is in the central slice of half-width abs_rap_slice // w.r.t. the event axis pref bool in_slice(const PseudoJet & p, const PseudoJet & pref) { // Rotate p, and use the rapidity to determine if p is in the slice contrib::lund_plane::Matrix3 rotmat = contrib::lund_plane::Matrix3::from_direction(pref).transpose(); const PseudoJet p_rot = rotmat*p; return fabs(p_rot.rap()) < abs_rap_slice; } //---------------------------------------------------------------------- int main(){ //---------------------------------------------------------- // read in input particles vector event; read_event(event); cout << "# read an event with " << event.size() << " particles" << endl; //---------------------------------------------------------- // create an instance of RecursiveLundEEGenerator, with default options int depth = -1; bool dynamic_psi_reference = true; fastjet::contrib::RecursiveLundEEGenerator lund(depth, dynamic_psi_reference); // first get some C/A jets double y3_cut = 1.0; JetDefinition::Plugin* ee_plugin = new EECambridgePlugin(y3_cut); JetDefinition jet_def(ee_plugin); ClusterSequence cs(event, jet_def); // Get the event axis (to be used for the slice) std::vector excl = cs.exclusive_jets(2); assert(excl.size() == 2); // order the two jets according to momentum along z axis if (excl[0].pz() < excl[1].pz()) { std::swap(excl[0],excl[1]); } const PseudoJet ev_axis = excl[0]-excl[1]; // Get the list of primary declusterings const vector declusts = lund.result(cs); // find the declustering that throws something into a fixed slice // (from a parent that was not in the slice) and that throws the // largest pt (relative to the z axis) into that slice int index_of_max_pt_in_slice = -1; double max_pt_in_slice = 0.0; - double psi_1; + double psi_1 = 0.0; //< init just to avoid a compiler warning (the + //< test in L115 makes it safe already) for (unsigned i = 0; i < declusts.size(); i++) { const auto & decl = declusts[i]; if (!in_slice(decl.harder(), ev_axis) && in_slice(decl.softer(), ev_axis)) { if (decl.softer().pt() > max_pt_in_slice) { index_of_max_pt_in_slice = i; max_pt_in_slice = decl.softer().pt(); psi_1 = decl.psibar(); } } } if (index_of_max_pt_in_slice < 0) return 0; // establish what we need to follow and the reference psi_1 int iplane_to_follow = declusts[index_of_max_pt_in_slice].leaf_iplane(); vector secondaries; for (const auto & declust: declusts){ if (declust.iplane() == iplane_to_follow) secondaries.push_back(&declust); } int index_of_max_kt_secondary = -1; double dpsi; for (uint i_secondary=0; i_secondaryz() > z2_cut) { index_of_max_kt_secondary = i_secondary; double psi_2 = secondaries[i_secondary]->psibar(); dpsi = contrib::lund_plane::map_to_pi(psi_2 - psi_1); break; } } if (index_of_max_kt_secondary < 0) return 0; cout << "Primary in the central slice, with Lund coordinates ( ln 1/Delta, ln kt, psibar ):" << endl; pair coords = declusts[index_of_max_pt_in_slice].lund_coordinates(); cout << "index [" << index_of_max_pt_in_slice << "](" << coords.first << ", " << coords.second << ", " << declusts[index_of_max_pt_in_slice].psibar() << ")" << endl; cout << endl << "with Lund coordinates for the (highest-kT) secondary plane that passes the zcut of " << z2_cut << endl; coords = secondaries[index_of_max_kt_secondary]->lund_coordinates(); cout << "index [" << index_of_max_kt_secondary << "](" << coords.first << ", " << coords.second << ", " << secondaries[index_of_max_kt_secondary]->psibar() << ")"; cout << " --> delta_psi,slice = " << dpsi << endl; // Delete the EECambridge plugin delete ee_plugin; return 0; } // read in input particles void read_event(vector &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); } } Index: contrib/contribs/LundPlane/trunk/example.cc =================================================================== --- contrib/contribs/LundPlane/trunk/example.cc (revision 1329) +++ contrib/contribs/LundPlane/trunk/example.cc (revision 1330) @@ -1,119 +1,119 @@ //---------------------------------------------------------------------- /// \file example.cc /// /// This example program is meant to illustrate how the /// fastjet::contrib::LundGenerator class is used. /// /// Run this example with /// /// \verbatim /// ./example < ../data/single-event.dat /// \endverbatim //---------------------------------------------------------------------- // $Id$ // // Copyright (c) 2018-, Frederic A. Dreyer, Keith Hamilton, Alexander Karlberg, // Gavin P. Salam, Ludovic Scyboz, Gregory Soyez, Rob Verheyen // //---------------------------------------------------------------------- // 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 . //---------------------------------------------------------------------- #include #include #include #include "fastjet/PseudoJet.hh" #include #include "LundGenerator.hh" // In external code, this should be fastjet/contrib/LundGenerator.hh #include "LundJSON.hh" // In external code, this should be fastjet/contrib/LundJSON.hh using namespace std; using namespace fastjet; // forward declaration to make things clearer void read_event(vector &event); //---------------------------------------------------------------------- int main(){ //---------------------------------------------------------- // read in input particles vector event; read_event(event); string filename = "jets.json"; cout << "# read an event with " << event.size() << " particles" << endl; cout << "# writing declusterings of primary and secondary plane to file " << filename << endl; ofstream outfile; outfile.open(filename.c_str()); // first get some anti-kt jets double R = 1.0, ptmin = 100.0; JetDefinition jet_def(antikt_algorithm, R); ClusterSequence cs(event, jet_def); vector jets = sorted_by_pt(cs.inclusive_jets(ptmin)); //---------------------------------------------------------- // create an instance of LundGenerator, with default options contrib::LundGenerator lund; cout << lund.description() << endl; for (unsigned ijet = 0; ijet < jets.size(); ijet++) { cout << endl << "Lund coordinates ( ln 1/Delta, ln kt ) of declusterings of jet " << ijet << " are:" << endl; vector declusts = lund(jets[ijet]); - for (int idecl = 0; idecl < declusts.size(); idecl++) { + for (unsigned int idecl = 0; idecl < declusts.size(); ++idecl) { pair coords = declusts[idecl].lund_coordinates(); cout << "(" << coords.first << ", " << coords.second << ")"; if (idecl < declusts.size() - 1) cout << "; "; } cout << endl; // outputs the primary Lund plane lund_to_json(outfile, declusts); outfile << endl; // outputs the full Lund tree //to_json(cout, lund_gen, jets[ijet]); cout << endl; } cout << endl << "File " << filename << " written." << endl; return 0; } // read in input particles void read_event(vector &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); } } Index: contrib/contribs/LundPlane/trunk/example_dpsi_collinear.cc =================================================================== --- contrib/contribs/LundPlane/trunk/example_dpsi_collinear.cc (revision 1329) +++ contrib/contribs/LundPlane/trunk/example_dpsi_collinear.cc (revision 1330) @@ -1,145 +1,145 @@ //---------------------------------------------------------------------- /// \file example_dpsi_collinear.cc /// /// This example program is meant to illustrate how the /// fastjet::contrib::RecursiveLundEEGenerator class is used. /// /// Run this example with /// /// \verbatim /// ./example_dpsi_collinear < ../data/single-ee-event.dat /// \endverbatim //---------------------------------------------------------------------- // $Id$ // // Copyright (c) 2018-, Frederic A. Dreyer, Keith Hamilton, Alexander Karlberg, // Gavin P. Salam, Ludovic Scyboz, Gregory Soyez, Rob Verheyen // //---------------------------------------------------------------------- // 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 . //---------------------------------------------------------------------- #include #include #include #include "fastjet/PseudoJet.hh" #include "fastjet/EECambridgePlugin.hh" #include #include "RecursiveLundEEGenerator.hh" // In external code, this should be fastjet/contrib/RecursiveLundEEGenerator.hh using namespace std; using namespace fastjet; // Definitions for the collinear observable double z1_cut = 0.1; double z2_cut = 0.1; // forward declaration to make things clearer void read_event(vector &event); //---------------------------------------------------------------------- int main(){ //---------------------------------------------------------- // read in input particles vector event; read_event(event); cout << "# read an event with " << event.size() << " particles" << endl; //---------------------------------------------------------- // create an instance of RecursiveLundEEGenerator, with default options int depth = -1; bool dynamic_psi_reference = true; fastjet::contrib::RecursiveLundEEGenerator lund(depth, dynamic_psi_reference); // first get some C/A jets double y3_cut = 1.0; JetDefinition::Plugin* ee_plugin = new EECambridgePlugin(y3_cut); JetDefinition jet_def(ee_plugin); ClusterSequence cs(event, jet_def); // Get the list of primary declusterings const vector declusts = lund.result(cs); // Find the highest-kt primary declustering (ordered in kt by default) int i_primary = -1; double psi_1; - for (int i=0; i z1_cut) { - i_primary = i; - psi_1 = declusts[i].psibar(); - break; + i_primary = i; + psi_1 = declusts[i].psibar(); + break; } } if (i_primary < 0) return 0; // Find the highest-kt secondary associated to that Lund leaf int iplane_to_follow = declusts[i_primary].leaf_iplane(); vector secondaries; for (const auto & declust: declusts){ - if (declust.iplane() == iplane_to_follow) secondaries.push_back(&declust); + if (declust.iplane() == iplane_to_follow) secondaries.push_back(&declust); } if(secondaries.size() < 1) return 0; int i_secondary = -1; double dpsi; - for (int i=0; iz() > z2_cut) { - i_secondary = i; - double psi_2 = secondaries[i]->psibar(); - dpsi = contrib::lund_plane::map_to_pi(psi_2-psi_1); - break; - } + for (unsigned int i=0; iz() > z2_cut) { + i_secondary = i; + double psi_2 = secondaries[i]->psibar(); + dpsi = contrib::lund_plane::map_to_pi(psi_2-psi_1); + break; + } } if (i_secondary < 0) return 0; cout << "Primary that passes the zcut of " << z1_cut << ", with Lund coordinates ( ln 1/Delta, ln kt, psibar ):" << endl; pair coords = declusts[i_primary].lund_coordinates(); cout << "index [" << i_primary << "](" << coords.first << ", " << coords.second << ", " << declusts[i_primary].psibar() << ")" << endl; cout << endl << "with Lund coordinates for the secondary plane that passes the zcut of " << z2_cut << endl; coords = secondaries[i_secondary]->lund_coordinates(); cout << "index [" << i_secondary << "](" << coords.first << ", " << coords.second << ", " << secondaries[i_secondary]->psibar() << ")"; cout << " --> delta_psi = " << dpsi << endl; // Delete the EECambridge plugin delete ee_plugin; return 0; } // read in input particles void read_event(vector &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); } } Index: contrib/contribs/LundPlane/trunk/example_secondary.cc =================================================================== --- contrib/contribs/LundPlane/trunk/example_secondary.cc (revision 1329) +++ contrib/contribs/LundPlane/trunk/example_secondary.cc (revision 1330) @@ -1,117 +1,117 @@ //---------------------------------------------------------------------- /// \file example_secondary.cc /// /// This example program is meant to illustrate how the /// fastjet::contrib::LundWithSecondary class is used. /// /// Run this example with /// /// \verbatim /// ./example_secondary < ../data/single-event.dat /// \endverbatim //---------------------------------------------------------------------- // $Id$ // // Copyright (c) 2018-, Frederic A. Dreyer, Keith Hamilton, Alexander Karlberg, // Gavin P. Salam, Ludovic Scyboz, Gregory Soyez, Rob Verheyen // //---------------------------------------------------------------------- // 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 . //---------------------------------------------------------------------- #include #include #include #include "fastjet/PseudoJet.hh" #include #include "LundWithSecondary.hh" // In external code, this should be fastjet/contrib/LundWithSecondary.hh using namespace std; using namespace fastjet; // forward declaration to make things clearer void read_event(vector &event); //---------------------------------------------------------------------- int main(){ //---------------------------------------------------------- // read in input particles vector event; read_event(event); cout << "# read an event with " << event.size() << " particles" << endl; // first get some anti-kt jets double R = 1.0, ptmin = 100.0; JetDefinition jet_def(antikt_algorithm, R); ClusterSequence cs(event, jet_def); vector jets = sorted_by_pt(cs.inclusive_jets(ptmin)); //---------------------------------------------------------- // create an instance of LundWithSecondary, with default options contrib::SecondaryLund_mMDT secondary; contrib::LundWithSecondary lund(&secondary); cout << lund.description() << endl; for (unsigned ijet = 0; ijet < jets.size(); ijet++) { cout << endl << "Lund coordinates ( ln 1/Delta, ln kt ) of declusterings of jet " << ijet << " are:" << endl; vector declusts = lund.primary(jets[ijet]); - for (int idecl = 0; idecl < declusts.size(); idecl++) { + for (unsigned int idecl = 0; idecl < declusts.size(); ++idecl) { pair coords = declusts[idecl].lund_coordinates(); cout << "["<< idecl<< "](" << coords.first << ", " << coords.second << ")"; if (idecl < declusts.size() - 1) cout << "; "; } cout << endl; vector sec_declusts = lund.secondary(declusts); cout << endl << "with Lund coordinates for the secondary plane (from primary declustering [" << lund.secondary_index(declusts) <<"]):" << endl; - for (int idecl = 0; idecl < sec_declusts.size(); ++idecl) { + for (unsigned int idecl = 0; idecl < sec_declusts.size(); ++idecl) { pair coords = sec_declusts[idecl].lund_coordinates(); cout << "["<< idecl<< "](" << coords.first << ", " << coords.second << ")"; if (idecl < sec_declusts.size() - 1) cout << "; "; } cout << endl; } return 0; } // read in input particles void read_event(vector &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); } } Index: contrib/contribs/LundPlane/trunk/ChangeLog =================================================================== --- contrib/contribs/LundPlane/trunk/ChangeLog (revision 1329) +++ contrib/contribs/LundPlane/trunk/ChangeLog (revision 1330) @@ -1,115 +1,124 @@ +2022-10-04 Gregory Soyez + + * example.cc: + * example_secondary.cc: + * example_dpsi_collinear.cc: + * example_dpsi_slice.cc: + removed (harmless) compilation warnings (+ minor uniformisation + of indentation) + 2022-10-04 Tue + Ludo Scyboz * VERSION: * NEWS: prepared for release 2.0.3 * EEHelpers.hh -> LundEEHelpers.hh: EEHelpers.hh was not being installed; given that it's a potentially common name, renamed it before including among the installation targets. Also placed whole contents in a new contrib::lund_plane namespace, because of certain potentially common names like cross_product * Makefile: replaced EEHelpers.hh -> LundEEHelpers.hh and made sure that LundEEHelpers.hh gets installed (without which RecursiveLundEEGenerator.hh cannot be used) * example_dpsi_collinear.cc: * example_dpsi_slice.cc: use of things from LundEEHelpers.hh updated to use of new namespace * RecursiveLundEEGenerator.hh: * RecursiveLundEEGenerator.cc: updated to use LundEEHelpers.hh (and associated namespace); also added a (protected) default constructor to LundEEDeclustering 2022-08-20 Gavin Salam * VERSION: * NEWS: prepared for release 2.0.2 * Makefile: updated some dependencies * EEHelpers.hh: added #include , as per request from Andy Buckley for compilation with g++-12 on some systems 2021-12-06 Gavin Salam * NEWS: * VERSION: prepared for release 2.0.1 * AUTHORS: fixed missing names and publication info 2021-12-06 Gregory Soyez * SecondaryLund.cc: fixed int v. unsigned int in loop over vector indices 2021-12-06 Gavin Salam * example.py: fixed name of executable in comments about how to execute this (thanks to Matteo Cacciari) 2021-11-09 Ludovic Scyboz * VERSION: preparing for release of 2.0.0 * RecursiveLundEEGenerator.hh: * RecursiveLundEEGenerator.cc: class for recursive Lund declustering in e+e- * example_dpsi_collinear.cc: spin-sensitive collinear observable from 2103.16526 * example_dpsi_slice.cc: spin-sensitive non-global observable from 2111.01161 2020-02-23 Gavin Salam * NEWS: * VERSION: preparing for release of 1.0.3 * example.cc: changed outfile open(filename) to outfile.open(filename.c_str()); to attempt to solve issue reported by Steven Schramm. 2018-10-26 Gavin Salam * read_lund_json.py: removed extraneous normalisation of zeroth bin in the LundImage class. Added documentation. 2018-08-30 Gavin Salam * VERSION: * NEWS: Release of version 1.0.1 2018-08-23 Gavin Salam * LundWithSecondary.hh: * LundWithSecondary.cc: added secondary_index(...), removed virtual qualifier from various functions * example_secondary.cc: * example_secondary.ref: example now prints out index of the primary declustering being used for the secondary. Referemce file updated accordingly. 2018-08-09 Frédéric Dreyer First version of LundPlane.