Skip to content

Commit d749fc8

Browse files
pdhankhePreeti DhankheralibuildPreeti Dhankher
authored
[PWGJE] Add new D0 substructure task (#15311)
Co-authored-by: Preeti Dhankher <preetidhanker@Preetis-MacBook-Air-8.local> Co-authored-by: ALICE Action Bot <alibuild@cern.ch> Co-authored-by: Preeti Dhankher <preetidhanker@dhcp-135-171.nikhef.nl>
1 parent bf0341e commit d749fc8

File tree

2 files changed

+212
-0
lines changed

2 files changed

+212
-0
lines changed

PWGJE/Tasks/CMakeLists.txt

Lines changed: 4 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -395,6 +395,10 @@ if(FastJet_FOUND)
395395
SOURCES jetDsSpectrumAndSubstructure.cxx
396396
PUBLIC_LINK_LIBRARIES O2::Framework O2Physics::PWGJECore O2Physics::AnalysisCore
397397
COMPONENT_NAME Analysis)
398+
o2physics_add_dpl_workflow(jet-d0-ang-substructure
399+
SOURCES jetD0AngSubstructure.cxx
400+
PUBLIC_LINK_LIBRARIES O2::Framework O2Physics::PWGJECore O2Physics::AnalysisCore
401+
COMPONENT_NAME Analysis)
398402
o2physics_add_dpl_workflow(bjet-cent-mult
399403
SOURCES bjetCentMult.cxx
400404
PUBLIC_LINK_LIBRARIES O2::Framework O2Physics::PWGJECore O2Physics::AnalysisCore
Lines changed: 208 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,208 @@
1+
// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2+
// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3+
// All rights not expressly granted are reserved.
4+
//
5+
// This software is distributed under the terms of the GNU General Public
6+
// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7+
//
8+
// In applying this license CERN does not waive the privileges and immunities
9+
// granted to it by virtue of its status as an Intergovernmental Organization
10+
// or submit itself to any jurisdiction.
11+
//
12+
// \file JetD0AngSubstructure.cxx
13+
//
14+
// \brief Analysis task for the reconstruction and study of charged jets
15+
// containing D_0 mesons in pp collisions.
16+
// \inherited from D0 fragmentation and Ds
17+
// \P. Dhankher
18+
19+
#include "PWGHF/Core/DecayChannels.h"
20+
#include "PWGJE/Core/JetDerivedDataUtilities.h"
21+
#include "PWGJE/Core/JetUtilities.h"
22+
#include "PWGJE/DataModel/Jet.h"
23+
#include "PWGJE/DataModel/JetReducedData.h"
24+
#include "PWGJE/DataModel/JetSubtraction.h"
25+
26+
#include "Common/Core/RecoDecay.h"
27+
#include "Common/DataModel/TrackSelectionTables.h"
28+
29+
#include "Framework/ASoA.h"
30+
#include "Framework/AnalysisDataModel.h"
31+
#include "Framework/HistogramRegistry.h"
32+
#include <CommonConstants/MathConstants.h>
33+
#include <Framework/AnalysisTask.h>
34+
#include <Framework/ConfigContext.h>
35+
#include <Framework/Configurable.h>
36+
#include <Framework/DataProcessorSpec.h>
37+
#include <Framework/HistogramSpec.h>
38+
#include <Framework/InitContext.h>
39+
#include <Framework/runDataProcessing.h>
40+
41+
#include "TVector3.h"
42+
43+
#include <cmath>
44+
#include <cstdlib>
45+
#include <string>
46+
#include <vector>
47+
48+
using namespace o2;
49+
using namespace o2::framework;
50+
using namespace o2::framework::expressions;
51+
52+
// Definition of a custom AOD table to store jet–D0 quantities
53+
namespace o2::aod
54+
{
55+
namespace jet_obj
56+
{
57+
// Jet-related quantities
58+
DECLARE_SOA_COLUMN(JetHfDist, jetHfDist, float);
59+
DECLARE_SOA_COLUMN(JetPt, jetPt, float);
60+
DECLARE_SOA_COLUMN(JetEta, jetEta, float);
61+
DECLARE_SOA_COLUMN(JetPhi, jetPhi, float);
62+
DECLARE_SOA_COLUMN(JetNConst, jetNConst, int);
63+
DECLARE_SOA_COLUMN(JetAng, jetAng, float);
64+
// D0 candidate quantities
65+
DECLARE_SOA_COLUMN(HfPt, hfPt, float);
66+
DECLARE_SOA_COLUMN(HfEta, hfEta, float);
67+
DECLARE_SOA_COLUMN(HfPhi, hfPhi, float);
68+
DECLARE_SOA_COLUMN(HfMass, hfMass, float);
69+
DECLARE_SOA_COLUMN(HfY, hfY, float);
70+
// ML scores
71+
DECLARE_SOA_COLUMN(HfMlScore0, hfMlScore0, float);
72+
DECLARE_SOA_COLUMN(HfMlScore1, hfMlScore1, float);
73+
DECLARE_SOA_COLUMN(HfMlScore2, hfMlScore2, float);
74+
} // namespace jet_obj
75+
// AOD table definition
76+
DECLARE_SOA_TABLE(JetObjTable, "AOD", "JETOBJTABLE",
77+
jet_obj::JetHfDist,
78+
jet_obj::JetPt,
79+
jet_obj::JetEta,
80+
jet_obj::JetPhi,
81+
jet_obj::JetNConst,
82+
jet_obj::JetAng,
83+
jet_obj::HfPt,
84+
jet_obj::HfEta,
85+
jet_obj::HfPhi,
86+
jet_obj::HfMass,
87+
jet_obj::HfY,
88+
jet_obj::HfMlScore0,
89+
jet_obj::HfMlScore1,
90+
jet_obj::HfMlScore2);
91+
} // namespace o2::aod
92+
struct JetD0AngSubstructure {
93+
/**
94+
* Histogram registry
95+
*
96+
* Contains:
97+
* - Event and track histograms
98+
* - Jet kinematic distributions
99+
* - D0–jet substructure observables
100+
*/
101+
HistogramRegistry registry{"registry",
102+
{{"h_collision_counter", "# of collisions;", {HistType::kTH1F, {{2, 0., 2.}}}},
103+
{"h_jet_counter", ";# of D^{0} jets;", {HistType::kTH1F, {{6, 0., 3.0}}}},
104+
{"h_d0_jet_projection", ";z^{D^{0},jet}_{||};dN/dz^{D^{0},jet}_{||}", {HistType::kTH1F, {{1000, 0., 10.}}}},
105+
{"h_d0_jet_distance_vs_projection", ";#DeltaR_{D^{0},jet};z^{D^{0},jet}_{||}", {HistType::kTH2F, {{1000, 0., 10.}, {1000, 0., 10.}}}},
106+
{"h_d0_jet_distance", ";#DeltaR_{D^{0},jet};dN/d(#DeltaR)", {HistType::kTH1F, {{1000, 0., 10.}}}},
107+
{"h_d0_jet_pt", ";p_{T,D^{0} jet};dN/dp_{T,D^{0} jet}", {HistType::kTH1F, {{200, 0., 10.}}}},
108+
{"h_d0_jet_eta", ";#eta_{T,D^{0} jet};dN/d#eta_{D^{0} jet}", {HistType::kTH1F, {{250, -5., 5.}}}},
109+
{"h_d0_jet_phi", ";#phi_{T,D^{0} jet};dN/d#phi_{D^{0} jet}", {HistType::kTH1F, {{250, -10., 10.}}}},
110+
{"h_d0_mass", ";m_{D^{0}});dN/dm_{D^{0}}", {HistType::kTH1F, {{1000, 0., 10.}}}},
111+
{"h_d0_eta", ";#eta_{D^{0}});dN/d#eta_{D_{}}", {HistType::kTH1F, {{250, -5., 5.}}}},
112+
{"h_d0_phi", ";#phi_{D^{0}});dN/d#phi_{D^{0}}", {HistType::kTH1F, {{250, -10., 10.}}}},
113+
{"h_d0_ang", ";#lambda_{#kappa}^{#alpha};counts", {HistType::kTH1F, {{100, 0., 1.}}}}}};
114+
115+
// Configurables
116+
Configurable<float> vertexZCut{"vertexZCut", 10.0f, "Accepted z-vertex range"};
117+
118+
Configurable<float> jetPtMin{"jetPtMin", 5.0, "minimum jet pT cut"};
119+
Configurable<float> jetR{"jetR", 0.4, "jet resolution parameter"};
120+
121+
Configurable<std::string> eventSelections{"eventSelections", "sel8", "choose event selection"};
122+
Configurable<std::string> trackSelections{"trackSelections", "globalTracks", "set track selections"};
123+
Configurable<float> kappa{"kappa", 1.0, "angularity kappa"}; // to do: configurable from json
124+
Configurable<float> alpha{"alpha", 1.0, "angularity alpha"};
125+
126+
std::vector<int> eventSelectionBits;
127+
int trackSelection = -1;
128+
129+
// Output table producer
130+
Produces<aod::JetObjTable> ObjJetTable;
131+
132+
float angularity;
133+
134+
void init(o2::framework::InitContext&)
135+
{
136+
eventSelectionBits = jetderiveddatautilities::initialiseEventSelectionBits(static_cast<std::string>(eventSelections));
137+
trackSelection = jetderiveddatautilities::initialiseTrackSelection(static_cast<std::string>(trackSelections));
138+
}
139+
140+
template <typename T, typename U>
141+
void jetCalculateAngularity(T const& jet, U const& /*tracks*/)
142+
{
143+
angularity = 0.0;
144+
for (auto& constituent : jet.template tracks_as<U>()) {
145+
angularity += std::pow(constituent.pt(), kappa) * std::pow(jetutilities::deltaR(jet, constituent) / (jet.r() / 100.f), alpha);
146+
}
147+
angularity /= std::pow(jet.pt(), kappa);
148+
}
149+
150+
void processDataChargedSubstructure(aod::JetCollision const& collision,
151+
soa::Join<aod::D0ChargedJets, aod::D0ChargedJetConstituents> const& jets,
152+
aod::CandidatesD0Data const&, aod::JetTracks const& tracks)
153+
{
154+
// apply event selection and fill histograms for sanity check
155+
registry.fill(HIST("h_collision_counter"), 0.5);
156+
157+
if (!jetderiveddatautilities::selectCollision(collision, eventSelectionBits) || !(std::abs(collision.posZ()) < vertexZCut)) {
158+
return;
159+
}
160+
registry.fill(HIST("h_collision_counter"), 1.5);
161+
162+
// Loop over jets containing D0 candidates
163+
for (const auto& jet : jets) {
164+
// number of charged jets with D0
165+
registry.fill(HIST("h_jet_counter"), 0.5);
166+
// obtaining jet 3-vector
167+
TVector3 jetVector(jet.px(), jet.py(), jet.pz());
168+
169+
// Loop over D0 candidates associated to the jet
170+
for (const auto& d0Candidate : jet.candidates_as<aod::CandidatesD0Data>()) {
171+
// obtaining jet 3-vector
172+
TVector3 d0Vector(d0Candidate.px(), d0Candidate.py(), d0Candidate.pz());
173+
174+
// calculating fraction of the jet momentum carried by the D0 along the direction of the jet axis
175+
double zParallel = (jetVector * d0Vector) / (jetVector * jetVector);
176+
177+
// calculating angular distance in eta-phi plane
178+
double axisDistance = jetutilities::deltaR(jet, d0Candidate);
179+
180+
jetCalculateAngularity(jet, tracks);
181+
182+
// filling histograms
183+
registry.fill(HIST("h_d0_jet_projection"), zParallel);
184+
registry.fill(HIST("h_d0_jet_distance_vs_projection"), axisDistance, zParallel);
185+
registry.fill(HIST("h_d0_jet_distance"), axisDistance);
186+
registry.fill(HIST("h_d0_jet_pt"), jet.pt());
187+
registry.fill(HIST("h_d0_jet_eta"), jet.eta());
188+
registry.fill(HIST("h_d0_jet_phi"), jet.phi());
189+
registry.fill(HIST("h_d0_mass"), d0Candidate.m());
190+
registry.fill(HIST("h_d0_eta"), d0Candidate.eta());
191+
registry.fill(HIST("h_d0_phi"), d0Candidate.phi());
192+
registry.fill(HIST("h_d0_ang"), angularity); // add more axis
193+
194+
// filling table
195+
ObjJetTable(axisDistance,
196+
jet.pt(), jet.eta(), jet.phi(), jet.tracks_as<aod::JetTracks>().size(), angularity,
197+
d0Candidate.pt(), d0Candidate.eta(), d0Candidate.phi(), d0Candidate.m(), d0Candidate.y(), d0Candidate.mlScores()[0], d0Candidate.mlScores()[1], d0Candidate.mlScores()[2]);
198+
199+
break; // get out of candidates' loop after first HF particle is found in jet
200+
} // end of D0 candidates loop
201+
202+
} // end of jets loop
203+
204+
} // end of process function
205+
PROCESS_SWITCH(JetD0AngSubstructure, processDataChargedSubstructure, "charged HF jet substructure", false);
206+
};
207+
// Workflow definition
208+
WorkflowSpec defineDataProcessing(ConfigContext const& cfgc) { return WorkflowSpec{adaptAnalysisTask<JetD0AngSubstructure>(cfgc, TaskName{"jet-d0-ang-substructure"})}; }

0 commit comments

Comments
 (0)