Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

feat: adding groomed jet radius and momentum balance variables for softdrop groomed jet #312

Draft
wants to merge 6 commits into
base: main
Choose a base branch
from
Draft
Changes from 2 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
10 changes: 9 additions & 1 deletion src/_ext.cpp
Original file line number Diff line number Diff line change
@@ -1658,6 +1658,8 @@ PYBIND11_MODULE(_ext, m) {
std::vector<double> jet_groomed_m;
std::vector<double> jet_groomed_E;
std::vector<double> jet_groomed_pz;
std::vector<double> jet_groomed_delta_R;
std::vector<double> jet_groomed_symmetry;

fastjet::contrib::RecursiveSymmetryCutBase::SymmetryMeasure sym_meas = fastjet::contrib::RecursiveSymmetryCutBase::SymmetryMeasure::scalar_z;
if (symmetry_measure == "scalar_z") {
@@ -1703,6 +1705,8 @@ PYBIND11_MODULE(_ext, m) {
jet_groomed_m.push_back(soft.m());
jet_groomed_E.push_back(soft.E());
jet_groomed_pz.push_back(soft.pz());
jet_groomed_delta_R.push_back(soft.structure_of<fastjet::contrib::SoftDrop>().delta_R());
jet_groomed_symmetry.push_back(soft.structure_of<fastjet::contrib::SoftDrop>().symmetry());
for (unsigned int k = 0; k < soft.constituents().size(); k++){
consts_groomed_px.push_back(soft.constituents()[k].px());
consts_groomed_py.push_back(soft.constituents()[k].py());
@@ -1723,6 +1727,8 @@ PYBIND11_MODULE(_ext, m) {
auto jet_m = py::array(jet_groomed_m.size(), jet_groomed_m.data());
auto jet_E = py::array(jet_groomed_E.size(), jet_groomed_E.data());
auto jet_pz = py::array(jet_groomed_pz.size(), jet_groomed_pz.data());
auto jet_delta_R = py::array(jet_groomed_delta_R.size(), jet_groomed_delta_R.data());
auto jet_symmetry = py::array(jet_groomed_symmetry.size(), jet_groomed_symmetry.data());

return std::make_tuple(
consts_px,
@@ -1735,7 +1741,9 @@ PYBIND11_MODULE(_ext, m) {
jet_phi,
jet_m,
jet_E,
jet_pz
jet_pz,
jet_delta_R,
jet_symmetry
);
}, R"pbdoc(
Performs softdrop pruning on jets.
4 changes: 4 additions & 0 deletions src/fastjet/_generalevent.py
Original file line number Diff line number Diff line change
@@ -491,6 +491,8 @@ def exclusive_jets_softdrop_grooming(
jetmass = ak.Array(ak.contents.NumpyArray(np_results[8]))
jetE = ak.Array(ak.contents.NumpyArray(np_results[9]))
jetpz = ak.Array(ak.contents.NumpyArray(np_results[10]))
jetdeltaR = ak.Array(ak.contents.NumpyArray(np_results[11]))
jetsymmetry = ak.Array(ak.contents.NumpyArray(np_results[12]))

self._out.append(
ak.zip(
@@ -504,6 +506,8 @@ def exclusive_jets_softdrop_grooming(
"phisoftdrop": jetphi,
"Esoftdrop": jetE,
"pzsoftdrop": jetpz,
"deltaRsoftdrop": jetdeltaR,
"symmetrysoftdrop": jetsymmetry,
},
depth_limit=1,
)
4 changes: 4 additions & 0 deletions src/fastjet/_multievent.py
Original file line number Diff line number Diff line change
@@ -264,6 +264,8 @@ def exclusive_jets_softdrop_grooming(
jetmass = ak.Array(ak.contents.NumpyArray(np_results[8]))
jetE = ak.Array(ak.contents.NumpyArray(np_results[9]))
jetpz = ak.Array(ak.contents.NumpyArray(np_results[10]))
jetdeltaR = ak.Array(ak.contents.NumpyArray(np_results[11]))
jetsymmetry = ak.Array(ak.contents.NumpyArray(np_results[12]))

out = ak.zip(
{
@@ -276,6 +278,8 @@ def exclusive_jets_softdrop_grooming(
"phisoftdrop": jetphi,
"Esoftdrop": jetE,
"pzsoftdrop": jetpz,
"deltaRsoftdrop": jetdeltaR,
"symmetrysoftdrop": jetsymmetry,
},
depth_limit=1,
)
4 changes: 4 additions & 0 deletions src/fastjet/_singleevent.py
Original file line number Diff line number Diff line change
@@ -248,6 +248,8 @@ def exclusive_jets_softdrop_grooming(
jetmass = ak.Array(ak.contents.NumpyArray(np_results[8]))
jetE = ak.Array(ak.contents.NumpyArray(np_results[9]))
jetpz = ak.Array(ak.contents.NumpyArray(np_results[10]))
jetdeltaR = ak.Array(ak.contents.NumpyArray(np_results[11]))
jetsymmetry = ak.Array(ak.contents.NumpyArray(np_results[12]))

out = ak.zip(
{
@@ -260,6 +262,8 @@ def exclusive_jets_softdrop_grooming(
"phisoftdrop": jetphi,
"Esoftdrop": jetE,
"pzsoftdrop": jetpz,
"deltaRsoftdrop": jetdeltaR,
"symmetrysoftdrop": jetsymmetry,
},
depth_limit=1,
)
14 changes: 14 additions & 0 deletions tests/test_002-exclusive_jets.py
Original file line number Diff line number Diff line change
@@ -352,6 +352,8 @@ def test_exclusive_jets_softdrop_grooming_multi():
],
),
"msoftdrop": ak.Array([488.2395243115817, 488.2395243115817]),
"deltaRsoftdrop": ak.Array([0.009873899817126915, 0.009873899817126915]),
"symmetrysoftdrop": ak.Array([0.49727522889673303, 0.49727522889673303]),
"ptsoftdrop": ak.Array([142.88274528437645, 142.88274528437645]),
"etasoftdrop": ak.Array([2.726117171791057, 2.726117171791057]),
"phisoftdrop": ak.Array([1.1012644074821902, 1.1012644074821902]),
@@ -393,6 +395,18 @@ def test_exclusive_jets_softdrop_grooming_multi():
rtol=1e-12,
atol=0,
),
ak.isclose(
ak.Array([softdrop_output.deltaRsoftdrop]),
ak.Array([softdrop.deltaRsoftdrop]),
rtol=1e-12,
atol=0,
),
ak.isclose(
ak.Array([softdrop_output.symmetrysoftdrop]),
ak.Array([softdrop.symmetrysoftdrop]),
rtol=1e-12,
atol=0,
),
ak.isclose(
ak.Array([softdrop_output.ptsoftdrop]),
ak.Array([softdrop.ptsoftdrop]),
18 changes: 18 additions & 0 deletions tests/test_008-dask.py
Original file line number Diff line number Diff line change
@@ -179,6 +179,12 @@ def dask_multi_test_exclusive_jets_softdrop_grooming():
{"px": 32.45, "py": 63.21, "pz": 543.14, "E": 599.56},
],
"msoftdrop": ak.Array([488.2395243115817, 488.2395243115817]),
"deltaRsoftdrop": ak.Array(
[0.009873899817126915, 0.009873899817126915]
),
"symmetrysoftdrop": ak.Array(
[0.49727522889673303, 0.49727522889673303]
),
"ptsoftdrop": ak.Array([142.88274528437645, 142.88274528437645]),
"etasoftdrop": ak.Array([2.726117171791057, 2.726117171791057]),
"phisoftdrop": ak.Array([1.1012644074821902, 1.1012644074821902]),
@@ -220,6 +226,18 @@ def dask_multi_test_exclusive_jets_softdrop_grooming():
rtol=1e-12,
atol=0,
),
ak.isclose(
ak.Array([softdrop_output.deltaRsoftdrop]),
ak.Array([softdrop.deltaRsoftdrop]),
rtol=1e-12,
atol=0,
),
ak.isclose(
ak.Array([softdrop_output.symmetrysoftdrop]),
ak.Array([softdrop.symmetrysoftdrop]),
rtol=1e-12,
atol=0,
),
ak.isclose(
ak.Array([softdrop_output.ptsoftdrop]),
ak.Array([softdrop.ptsoftdrop]),