Skip to content

Commit 2a7f691

Browse files
zovargaalibuild
andauthored
[PWGJE] Change MiniCollision index column declaration for jetLundPlane.cxx (#15633)
Co-authored-by: ALICE Action Bot <alibuild@cern.ch>
1 parent 2a0b11a commit 2a7f691

File tree

1 file changed

+10
-5
lines changed

1 file changed

+10
-5
lines changed

PWGJE/Tasks/jetLundPlane.cxx

Lines changed: 10 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -53,7 +53,7 @@ DECLARE_SOA_TABLE(MiniCollisions, "AOD", "MINICOLL",
5353
MiniCollTag);
5454

5555
// MiniJets -> MiniCollisions
56-
DECLARE_SOA_INDEX_COLUMN(MiniCollision, miniCollision);
56+
DECLARE_SOA_INDEX_COLUMN_CUSTOM(MiniCollision, miniCollision, "MINICOLLS");
5757

5858
// Jet payload
5959
DECLARE_SOA_COLUMN(Level, level, uint8_t); // JetLevel::Det=reco(det), JetLevel::Part=truth(part)
@@ -96,6 +96,9 @@ DECLARE_SOA_TABLE(MiniJetMatches, "AOD", "MINIMCH",
9696
namespace
9797
{
9898
constexpr float kTiny = 1e-12f;
99+
constexpr uint64_t collisionKeyShift = 1ULL;
100+
constexpr uint64_t partCollisionKeyTag = 1ULL;
101+
constexpr size_t MinConstituentsForJet = 2;
99102

100103
struct JetLevel {
101104
enum Type : uint8_t {
@@ -416,7 +419,7 @@ struct JetLundPlaneUnfolding {
416419
std::vector<SplittingObs> getPrimarySplittings(JetRowT const& jet, ConstituentTableT const&)
417420
{
418421
auto fjInputs = buildFastJetInputs(jet.template tracks_as<ConstituentTableT>(), trackPtMin.value);
419-
if (fjInputs.size() < 2) {
422+
if (fjInputs.size() < MinConstituentsForJet) {
420423
return {};
421424
}
422425

@@ -583,7 +586,8 @@ struct JetLundPlaneUnfolding {
583586
fillSplittingQAHists(spl, /*isTruth*/ true, partJet.pt());
584587

585588
if (writeMiniAOD.value) {
586-
const uint64_t partCollKey = (static_cast<uint64_t>(partJet.mcCollisionId()) << 1U) | 1ULL;
589+
const uint64_t partCollKey =
590+
(static_cast<uint64_t>(partJet.mcCollisionId()) << collisionKeyShift) | partCollisionKeyTag;
587591
int partMiniCollIdx = -1;
588592
auto collIt = partMiniCollByKey.find(partCollKey);
589593
if (collIt == partMiniCollByKey.end()) {
@@ -650,7 +654,8 @@ struct JetLundPlaneUnfolding {
650654
fillSplittingQAHists(detSpl, /*isTruth*/ false, detJet.pt());
651655

652656
if (writeMiniAOD.value) {
653-
const uint64_t detCollKey = (static_cast<uint64_t>(detJet.collisionId()) << 1U);
657+
const uint64_t detCollKey =
658+
(static_cast<uint64_t>(detJet.collisionId()) << collisionKeyShift);
654659
int detMiniCollIdx = -1;
655660
auto collIt = detMiniCollByKey.find(detCollKey);
656661
if (collIt == detMiniCollByKey.end()) {
@@ -795,5 +800,5 @@ struct JetLundPlaneUnfolding {
795800
WorkflowSpec defineDataProcessing(ConfigContext const& cfgc)
796801
{
797802
return WorkflowSpec{
798-
adaptAnalysisTask<JetLundPlaneUnfolding>(cfgc, TaskName{"jet-lund-plane"})};
803+
adaptAnalysisTask<JetLundPlaneUnfolding>(cfgc)};
799804
}

0 commit comments

Comments
 (0)