Page MenuHomeHEPForge

No OneTemporary

Index: contrib/contribs/LundPlane/trunk/RecursiveLundEEGenerator.hh
===================================================================
--- contrib/contribs/LundPlane/trunk/RecursiveLundEEGenerator.hh (revision 1364)
+++ contrib/contribs/LundPlane/trunk/RecursiveLundEEGenerator.hh (revision 1365)
@@ -1,347 +1,349 @@
// $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 <http://www.gnu.org/licenses/>.
//----------------------------------------------------------------------
#ifndef __FASTJET_CONTRIB_RECURSIVELUNDEEGENERATOR_HH__
#define __FASTJET_CONTRIB_RECURSIVELUNDEEGENERATOR_HH__
#include "LundEEHelpers.hh"
#include <fastjet/internal/base.hh>
#include "fastjet/tools/Recluster.hh"
#include "fastjet/JetDefinition.hh"
#include "fastjet/PseudoJet.hh"
#include <string>
#include <vector>
#include <utility>
#include <queue>
using namespace std;
FASTJET_BEGIN_NAMESPACE
namespace contrib{
//----------------------------------------------------------------------
/// \class LundEEDeclustering
/// Contains the declustering variables associated with a single node
/// on the LundEE plane
class LundEEDeclustering {
public:
/// return the pair PseudoJet, i.e. sum of the two subjets
const PseudoJet & pair() const {return pair_;}
/// returns the subjet with larger transverse momentum
const PseudoJet & harder() const {return harder_;}
/// returns the subjet with smaller transverse momentum
const PseudoJet & softer() const {return softer_;}
/// returns pair().m() [cached]
double m() const {return m_;}
/// returns the effective pseudorapidity of the emission [cached]
double eta() const {return eta_;}
/// returns sin(theta) of the branching [cached]
double sin_theta() const {return sin_theta_;}
/// returns softer().modp() / (softer().modp() + harder().modp()) [cached]
double z() const {return z_;}
/// returns softer().modp() * sin(theta()) [cached]
double kt() const {return kt_;}
/// returns ln(softer().modp() * sin(theta())) [cached]
double lnkt() const {return lnkt_;}
/// returns z() * Delta() [cached]
double kappa() const {return kappa_;}
/// returns the index of the plane to which this branching belongs
int iplane() const {return iplane_;}
/// returns the depth of the plane on which this declustering
/// occurred. 0 is the primary plane, 1 is the first set of leaves, etc.
int depth() const {return depth_;}
/// returns iplane (plane index) of the leaf associated with the
/// potential further declustering of the softer of the objects in
/// this splitting
int leaf_iplane() const {return leaf_iplane_;}
/// Returns sign_s, indicating the initial parent jet index of this splitting
int sign_s() const {return sign_s_;}
/// returns an azimuthal angle psibar associated with this declustering.
/// The actual value of psibar is arbitrary, but differences in psibar
/// values between different clusterings are meaningful.
double psibar() const {return psibar_;}
/// (DEPRECATED)
/// returns an azimuthal type angle between this declustering plane and the previous one
/// Note: one should use psibar() instead, since we found that this definition of psi is
/// not invariant under rotations of the event
double psi() const {return psi_;}
/// update the azimuthal angle (deprecated)
void set_psi(double psi) {psi_ = psi;}
/// returns the coordinates in the Lund plane
std::pair<double,double> const lund_coordinates() const {
return std::pair<double,double>(eta_,lnkt_);
}
virtual ~LundEEDeclustering() {}
private:
int iplane_;
double psi_, psibar_, lnkt_, eta_;
double m_, z_, kt_, kappa_, sin_theta_;
PseudoJet pair_, harder_, softer_;
int depth_ = -1, leaf_iplane_ = -1;
int sign_s_;
protected:
/// default ctor (protected, should not normally be needed by users,
/// but can be useful for derived classes)
LundEEDeclustering() {}
/// the constructor is protected, because users will not generally be
/// constructing a LundEEDeclustering element themselves.
LundEEDeclustering(const PseudoJet& pair,
const PseudoJet& j1, const PseudoJet& j2,
int iplane = -1, double psi = 0.0, double psibar = 0.0, int depth = -1, int leaf_iplane = -1, int sign_s = 1);
friend class RecursiveLundEEGenerator;
};
/// Default comparison operator for LundEEDeclustering, using kt as the ordering.
/// Useful when including declusterings in structures like priority queues
inline bool operator<(const LundEEDeclustering& d1, const LundEEDeclustering& d2) {
return d1.kt() < d2.kt();
}
//----------------------------------------------------------------------
/// Class to carry out Lund declustering to get anything from the
/// primary Lund plane declusterings to the full Lund diagram with all
/// its leaves, etc.
class RecursiveLundEEGenerator {
public:
/// constructs a RecursiveLundEEGenerator with the specified depth.
/// - depth = 0 means only primary declusterings are registered
/// - depth = 1 means the first set of leaves are declustered
/// - ...
/// - depth < 0 means no limit, i.e. recurse through all leaves
///
/// The psibar values that are set in the result Lund tree have the
/// following property:
///
/// - if the jet with the larger pz has splittings, then its
/// first splitting has psibar = 0
/// - otherwise the first splitting of the other jet has psibar = 0
///
/// Note that this makes psibar IR unsafe (because an arbitrarily soft
/// splitting can be the one that gets the reference psibar=0 value),
/// but differences between psibar values are IR safe.
///
/// NB: The dynamical_psi_ref option relates to the deprecated definition of psi
/// New code should use the psibar() function and dynamical_psi_ref
/// is irrelevant.
RecursiveLundEEGenerator(int max_depth = 0, bool dynamical_psi_ref = false) :
max_depth_(max_depth), nx_(1,0,0,0), ny_(0,1,0,0), dynamical_psi_ref_(dynamical_psi_ref)
{}
/// destructor
virtual ~RecursiveLundEEGenerator() {}
/// This takes a cluster sequence with an e+e- C/A style algorithm, e.g.
/// EECambridgePlugin(ycut=1.0).
///
/// The output is a vector of LundEEDeclustering objects, ordered
/// according to kt
virtual std::vector<LundEEDeclustering> result(const ClusterSequence & cs) const {
std::vector<PseudoJet> exclusive_jets = cs.exclusive_jets(2);
assert(exclusive_jets.size() == 2);
// order the two jets according to momentum along z axis
if (exclusive_jets[0].pz() < exclusive_jets[1].pz()) {
std::swap(exclusive_jets[0],exclusive_jets[1]);
}
PseudoJet d_ev = exclusive_jets[0] - exclusive_jets[1];
lund_plane::Matrix3 rotmat = lund_plane::Matrix3::from_direction(d_ev);
std::vector<LundEEDeclustering> declusterings;
int depth = 0;
int max_iplane_sofar = 1;
// 2024-01: new code, that fixes up issue of psibar differences
// between hemispheres. If RLEEG_NEWPSIBAR is false, answers
// will come out wrong.
#define RLEEG_NEWPSIBAR
#ifdef RLEEG_NEWPSIBAR
// 2024-01 -- attempt at new definition of psibar
PseudoJet ref_plane;
double last_psibar = 0.;
bool first_time = true;
for (unsigned ijet = 0; ijet < exclusive_jets.size(); ijet++) {
int sign_s = ijet == 0? +1 : -1;
append_to_vector(declusterings, exclusive_jets[ijet], depth, ijet, max_iplane_sofar,
rotmat, sign_s, ref_plane, last_psibar, first_time);
}
#else
for (unsigned ijet = 0; ijet < exclusive_jets.size(); ijet++) {
// reference direction for psibar calculation
PseudoJet axis = d_ev/sqrt(d_ev.modp2());
PseudoJet ref_plane = axis;
int sign_s = ijet == 0? +1 : -1;
bool
// We can pass a vector normal to a plane of reference for phi definitions
append_to_vector(declusterings, exclusive_jets[ijet], depth, ijet, max_iplane_sofar,
rotmat, sign_s, ref_plane, 0., true);
}
#endif
// a typedef to save typing below
typedef LundEEDeclustering LD;
// sort so that declusterings are returned in order of decreasing
// kt (if result of the lambda is true, then first object appears
// before the second one in the final sorted list)
sort(declusterings.begin(), declusterings.end(),
[](const LD & d1, const LD & d2){return d1.kt() > d2.kt();});
return declusterings;
}
private:
/// internal routine to recursively carry out the declusterings,
/// adding each one to the declusterings vector; the primary
/// ones are dealt with first (from large to small angle),
/// and then secondary ones take place.
void append_to_vector(std::vector<LundEEDeclustering> & declusterings,
const PseudoJet & jet, int depth,
int iplane, int & max_iplane_sofar,
const lund_plane::Matrix3 & rotmat, int sign_s,
PseudoJet & psibar_ref_plane,
const double & last_psibar, bool & first_time) const {
PseudoJet j1, j2;
if (!jet.has_parents(j1, j2)) return;
if (j1.modp2() < j2.modp2()) std::swap(j1,j2);
// calculation of azimuth psi
lund_plane::Matrix3 new_rotmat;
if (dynamical_psi_ref_) {
new_rotmat = lund_plane::Matrix3::from_direction(rotmat.transpose()*(sign_s*jet)) * rotmat;
} else {
new_rotmat = rotmat;
}
PseudoJet rx = new_rotmat * nx_;
PseudoJet ry = new_rotmat * ny_;
PseudoJet u1 = j1/j1.modp(), u2 = j2/j2.modp();
PseudoJet du = u2 - u1;
double x = du.px() * rx.px() + du.py() * rx.py() + du.pz() * rx.pz();
double y = du.px() * ry.px() + du.py() * ry.py() + du.pz() * ry.pz();
double psi = atan2(y,x);
// calculation of psibar
double psibar = 0.;
PseudoJet n1, n2;
// First psibar for this jet
if (first_time) {
#ifdef RLEEG_NEWPSIBAR
// 2024-01: new code, that fixes up issue of psibar differences
// between hemispheres
assert(last_psibar == 0.0);
psibar = 0.0;
n2 = lund_plane::cross_product(j1,j2);
n2 /= n2.modp();
psibar_ref_plane = n2;
first_time = false;
#else
// Compute the angle between the planes spanned by (some axis,j1) and by (j1,j2)
n1 = lund_plane::cross_product(psibar_ref_plane,j1);
n2 = lund_plane::cross_product(j1,j2);
double signed_angle = 0.;
n2 /= n2.modp();
if (n1.modp() != 0) {
n1 /= n1.modp();
signed_angle = lund_plane::signed_angle_between_planes(n1,n2,j1);
}
psibar = lund_plane::map_to_pi(j1.phi() + signed_angle);
#endif
}
// Else take the value of psibar_i and the plane from the last splitting to define psibar_{i+1}
+ // The signed angle is multiplied by sign_s (+1 for the +z hemisphere, -1 for the -z hemisphere)
+ // to take the correct orientation into account.
else {
n2 = lund_plane::cross_product(j1,j2);
n2 /= n2.modp();
- psibar = lund_plane::map_to_pi(last_psibar + lund_plane::signed_angle_between_planes(psibar_ref_plane, n2, j1));
+ psibar = lund_plane::map_to_pi(last_psibar + sign_s*lund_plane::signed_angle_between_planes(psibar_ref_plane, n2, j1));
}
int leaf_iplane = -1;
// we will recurse into the softer "parent" only if the depth is
// not yet at its limit or if there is no limit on the depth (max_depth<0)
bool recurse_into_softer = (depth < max_depth_ || max_depth_ < 0);
if (recurse_into_softer) {
max_iplane_sofar += 1;
leaf_iplane = max_iplane_sofar;
}
LundEEDeclustering declust(jet, j1, j2, iplane, psi, psibar, depth, leaf_iplane, sign_s);
declusterings.push_back(declust);
// now recurse
// for the definition of psibar, we recursively pass the last splitting plane (normal to n2) and the last value
// of psibar
bool lcl_first_time = false;
append_to_vector(declusterings, j1, depth, iplane, max_iplane_sofar, new_rotmat, sign_s, n2, psibar, lcl_first_time);
if (recurse_into_softer) {
append_to_vector(declusterings, j2, depth+1, leaf_iplane, max_iplane_sofar, new_rotmat, sign_s, n2, psibar, lcl_first_time);
}
}
int max_depth_ = 0;
/// vectors used to help define psi
PseudoJet nx_;
PseudoJet ny_;
bool dynamical_psi_ref_;
};
} // namespace contrib
FASTJET_END_NAMESPACE
/// for output of declustering information
std::ostream & operator<<(std::ostream & ostr, const fastjet::contrib::LundEEDeclustering & d);
#endif // __FASTJET_CONTRIB_RECURSIVELUNDEEGENERATOR_HH__
Index: contrib/contribs/LundPlane/trunk/NEWS
===================================================================
--- contrib/contribs/LundPlane/trunk/NEWS (revision 1364)
+++ contrib/contribs/LundPlane/trunk/NEWS (revision 1365)
@@ -1,36 +1,36 @@
2024/01/XX: release of version 2.1.0
- changed psibar calculation in RecursiveLundEEGenerator to
get correct answers for delta psibar between opposite hemisphere
(changes all psibar values, but leaves in-jet delta psibar results
- unchanged).
+ unchanged up to a sign).
- fixed uint -> unsigned int to avoid compilation failures on some systems.
2023/09/22: release of version 2.0.4
- fixed const-correctness issue
- added LundPlane.hh, which includes the main other headers
2022/10/04: release of version 2.0.3
- renamed EEHelpers.hh -> LundEEHelpers.hh, put its contents in
a dedicated lund_plane namespace and ensure it gets installed
2022/08/20: release of version 2.0.2
- fixed missing header for g++-12
- fixed missing Makefile dependencies
2021/12/06: release of version 2.0.1
- fixed some compilation warnings (signed/unsigned)
- fixed up missing info in AUTHORS
- fixed commented python example lines
2021/11/08: release of version 2.0.0
- Recursive e+e- generator with spin-correlation observables
with the RecursiveLundEEGenerator
- See example_dpsi_collinear.cc (implements analysis of arXiv:2103.16526)
- and example_dpsi_slice.cc (implements analysis of arXiv:2111.01161)
2020/02/23: release of version 1.0.3
- fixes issue with C++98 compatibility
2018/10/29: release of version 1.0.2
2018/08/30: release of version 1.0.1
Index: contrib/contribs/LundPlane/trunk/example_dpsi_collinear.ref
===================================================================
--- contrib/contribs/LundPlane/trunk/example_dpsi_collinear.ref (revision 1364)
+++ contrib/contribs/LundPlane/trunk/example_dpsi_collinear.ref (revision 1365)
@@ -1,22 +1,22 @@
# read an event with 70 particles
#--------------------------------------------------------------------------
# FastJet release 3.4.2
# M. Cacciari, G.P. Salam and G. Soyez
# A software package for jet finding and analysis at colliders
# http://fastjet.fr
#
# Please cite EPJC72(2012)1896 [arXiv:1111.6097] if you use this package
# for scientific work and optionally PLB641(2006)57 [hep-ph/0512210].
#
# FastJet is provided without warranty under the GNU GPL v2 or higher.
# It uses T. Chan's closest pair algorithm, S. Fortune's Voronoi code
# and 3rd party plugin jet algorithms. See COPYING file for details.
#--------------------------------------------------------------------------
Highest-kt primary that passes the zcut of 0.1, with Lund coordinates ( ln 1/Delta, ln kt, psibar ):
-index [0](0.461518, 2.74826, -2.25293)
+index [0](0.461518, 2.74826, 2.25293)
with Lund coordinates for the secondary plane that passes the zcut of 0.1
-index [0](0.860316, 1.51421, -0.949104) --> delta_psi(primary,secondary) = 1.30383
+index [0](0.860316, 1.51421, 0.949104) --> delta_psi(primary,secondary) = -1.30383
with Lund coordinates for the second primary plane that passes the zcut of 0.1
-index [2](1.99207, 1.49759, -1.91726) --> delta_psi(primary_1, primary_2) = 0.335672
+index [2](1.99207, 1.49759, -1.91726) --> delta_psi(primary_1, primary_2) = 2.11299
Index: contrib/contribs/LundPlane/trunk/python/gen-test-event.py
===================================================================
--- contrib/contribs/LundPlane/trunk/python/gen-test-event.py (revision 1364)
+++ contrib/contribs/LundPlane/trunk/python/gen-test-event.py (revision 1365)
@@ -1,38 +1,56 @@
#!/usr/bin/env python3
#
# Script to help produce simple small events for testing the
# azimuthal angle aspects of the LundPlane contrib
#
# Usage (from the LundPlane directory):
#
# python3 python/gen-test-event.py | ./example_dpsi_collinear
#
import fastjet as fj
from math import *
-# produce an event with a collinear branching on each side
-# and some angle phi between the two branches.
+# produce an event with a collinear branching on each side,
+# some angle phi between the two branches, and an additional
+# particle in the z>0 hemisphere at phi!=0 (to test whether
+# signs are set properly across the two hemispheres).
+#
+# (That particle generates the highest-kt primary on that side,
+# but that primary does not pass the z-cut. Hence the two
+# opposite-side primaries will have psibar != 0)
#
# - smaller value of scale make the event more collinear
# - phi = 0 causes the soft particles to be on opposite sides
-scale = 0.01
+# - mirror_event rotates the event by pi (to test whether
+# seeding the declustering with the other jet makes a difference)
+scale = 0.001
phi = 0.1
-particles = [
+mirror_event = False
+
+particles_in = [
fj.PseudoJet(0,0,+1,1),
fj.PseudoJet(0,0,-1,1),
- fj.PtYPhiM(0.12*scale, 1.0*(1 - log(scale)), 0.0),
- fj.PtYPhiM(0.10*scale, -1.0*(1 - log(scale)), -pi - phi),
+ fj.PtYPhiM(0.12*scale, 1.0*(1 - log(scale)), -phi),
+ fj.PtYPhiM(0.10*scale, -1.0*(1 - log(scale)), pi),
+ fj.PtYPhiM(0.15*scale, 0.92*(1 - log(scale)), 1.15),
]
# dumb way to get sum of all momenta -- clustering with e+e- Cambridge/Aachen
# with a large radius
-psum = fj.JetDefinition(fj.ee_genkt_algorithm, 2*pi, 0.0)(particles)[0]
+psum = fj.JetDefinition(fj.ee_genkt_algorithm, 2*pi, 0.0)(particles_in)[0]
#print(psum)
# then boost to balance the event
-for p in particles:
+for p in particles_in:
p.unboost(psum)
+# and mirror the event if requested
+particles = []
+for p in particles_in:
+ if mirror_event:
+ p = fj.PseudoJet(-p.px(), -p.py(), -p.pz(), p.E())
+ particles.append(p)
+
for p in particles:
- print(p.px(), p.py(), p.pz(), p.E())
\ No newline at end of file
+ print(p.px(), p.py(), p.pz(), p.E())
Index: contrib/contribs/LundPlane/trunk/example_dpsi_slice.ref
===================================================================
--- contrib/contribs/LundPlane/trunk/example_dpsi_slice.ref (revision 1364)
+++ contrib/contribs/LundPlane/trunk/example_dpsi_slice.ref (revision 1365)
@@ -1,19 +1,19 @@
# read an event with 70 particles
#--------------------------------------------------------------------------
# FastJet release 3.4.2
# M. Cacciari, G.P. Salam and G. Soyez
# A software package for jet finding and analysis at colliders
# http://fastjet.fr
#
# Please cite EPJC72(2012)1896 [arXiv:1111.6097] if you use this package
# for scientific work and optionally PLB641(2006)57 [hep-ph/0512210].
#
# FastJet is provided without warranty under the GNU GPL v2 or higher.
# It uses T. Chan's closest pair algorithm, S. Fortune's Voronoi code
# and 3rd party plugin jet algorithms. See COPYING file for details.
#--------------------------------------------------------------------------
Primary in the central slice, with Lund coordinates ( ln 1/Delta, ln kt, psibar ):
-index [0](0.461518, 2.74826, -2.25293)
+index [0](0.461518, 2.74826, 2.25293)
with Lund coordinates for the (highest-kT) secondary plane that passes the zcut of 0.1
-index [0](0.860316, 1.51421, -0.949104) --> delta_psi,slice = 1.30383
+index [0](0.860316, 1.51421, 0.949104) --> delta_psi,slice = -1.30383

File Metadata

Mime Type
text/x-diff
Expires
Sat, May 3, 5:52 AM (17 h, 56 m)
Storage Engine
blob
Storage Format
Raw Data
Storage Handle
4982844
Default Alt Text
(21 KB)

Event Timeline