Skip to content

Commit 92931f4

Browse files
Preeti DhankherPreeti Dhankher
authored andcommitted
Add new D0 substructure task
1 parent 77b3b50 commit 92931f4

File tree

2 files changed

+270
-0
lines changed

2 files changed

+270
-0
lines changed

PWGJE/Tasks/CMakeLists.txt

Lines changed: 4 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -389,6 +389,10 @@ if(FastJet_FOUND)
389389
SOURCES jetDsSpectrumAndSubstructure.cxx
390390
PUBLIC_LINK_LIBRARIES O2::Framework O2Physics::PWGJECore O2Physics::AnalysisCore
391391
COMPONENT_NAME Analysis)
392+
o2physics_add_dpl_workflow(jet-d0-substructure
393+
SOURCES jetD0Substructure.cxx
394+
PUBLIC_LINK_LIBRARIES O2::Framework O2Physics::PWGJECore O2Physics::AnalysisCore
395+
COMPONENT_NAME Analysis)
392396
o2physics_add_dpl_workflow(bjet-cent-mult
393397
SOURCES bjetCentMult.cxx
394398
PUBLIC_LINK_LIBRARIES O2::Framework O2Physics::PWGJECore O2Physics::AnalysisCore

PWGJE/Tasks/jetD0Substructure.cxx

Lines changed: 266 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,266 @@
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 jetD0Substructure.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+
20+
#include "PWGJE/DataModel/Jet.h"
21+
#include "PWGJE/DataModel/JetSubtraction.h"
22+
#include "PWGJE/DataModel/JetReducedData.h"
23+
24+
#include "PWGJE/Core/JetDerivedDataUtilities.h"
25+
#include "PWGJE/Core/JetUtilities.h"
26+
#include "PWGHF/Core/DecayChannels.h"
27+
28+
#include "Common/Core/RecoDecay.h"
29+
#include "Common/DataModel/TrackSelectionTables.h"
30+
31+
#include "Framework/ASoA.h"
32+
#include "Framework/AnalysisDataModel.h"
33+
#include <Framework/AnalysisTask.h>
34+
#include "Framework/HistogramRegistry.h"
35+
#include <CommonConstants/MathConstants.h>
36+
#include <Framework/Configurable.h>
37+
#include <Framework/HistogramSpec.h>
38+
#include <Framework/InitContext.h>
39+
#include <Framework/ConfigContext.h>
40+
#include <Framework/DataProcessorSpec.h>
41+
#include <Framework/runDataProcessing.h>
42+
43+
44+
#include "TVector3.h"
45+
#include <cmath>
46+
#include <cstdlib>
47+
#include <string>
48+
#include <vector>
49+
50+
using namespace o2;
51+
using namespace o2::framework;
52+
using namespace o2::framework::expressions;
53+
54+
// Definition of a custom AOD table to store jet–D0 quantities
55+
namespace o2::aod
56+
{
57+
namespace jet_obj
58+
{
59+
// Jet-related quantities
60+
DECLARE_SOA_COLUMN(JetHfDist, jetHfDist, float);
61+
DECLARE_SOA_COLUMN(JetPt, jetPt, float);
62+
DECLARE_SOA_COLUMN(JetEta, jetEta, float);
63+
DECLARE_SOA_COLUMN(JetPhi, jetPhi, float);
64+
DECLARE_SOA_COLUMN(JetNConst, jetNConst, int);
65+
DECLARE_SOA_COLUMN(JetAng, jetAng, float);
66+
// D0 candidate quantities
67+
DECLARE_SOA_COLUMN(HfPt, hfPt, float);
68+
DECLARE_SOA_COLUMN(HfEta, hfEta, float);
69+
DECLARE_SOA_COLUMN(HfPhi, hfPhi, float);
70+
DECLARE_SOA_COLUMN(HfMass, hfMass, float);
71+
DECLARE_SOA_COLUMN(HfY, hfY, float);
72+
// ML scores
73+
DECLARE_SOA_COLUMN(HfMlScore0, hfMlScore0, float);
74+
DECLARE_SOA_COLUMN(HfMlScore1, hfMlScore1, float);
75+
DECLARE_SOA_COLUMN(HfMlScore2, hfMlScore2, float);
76+
} // namespace jet_obj
77+
// AOD table definition
78+
DECLARE_SOA_TABLE(JetObjTable, "AOD", "JETOBJTABLE",
79+
jet_obj::JetHfDist,
80+
jet_obj::JetPt,
81+
jet_obj::JetEta,
82+
jet_obj::JetPhi,
83+
jet_obj::JetNConst,
84+
jet_obj::JetAng,
85+
jet_obj::HfPt,
86+
jet_obj::HfEta,
87+
jet_obj::HfPhi,
88+
jet_obj::HfMass,
89+
jet_obj::HfY,
90+
jet_obj::HfMlScore0,
91+
jet_obj::HfMlScore1,
92+
jet_obj::HfMlScore2);
93+
}
94+
struct JetD0Substructure {
95+
/**
96+
* Histogram registry
97+
*
98+
* Contains:
99+
* - Event and track histograms
100+
* - Jet kinematic distributions
101+
* - D0–jet substructure observables
102+
*/
103+
HistogramRegistry registry{"registry",
104+
{{"h_collisions", "event status;event status;entries", {HistType::kTH1F, {{4, 0.0, 4.0}}}},
105+
{"h_track_pt", "track pT;#it{p}_{T,track} (GeV/#it{c});entries", {HistType::kTH1F, {{200, 0., 200.}}}},
106+
{"h_track_eta", "track #eta;#eta_{track};entries", {HistType::kTH1F, {{100, -1.0, 1.0}}}},
107+
{"h_track_phi", "track #varphi;#varphi_{track};entries", {HistType::kTH1F, {{80, -1.0, 7.}}}},
108+
{"h_jet_pt", "jet pT;#it{p}_{T,jet} (GeV/#it{c});entries", {HistType::kTH1F, {{200, 0., 200.}}}},
109+
{"h_jet_eta", "jet #eta;#eta_{jet};entries", {HistType::kTH1F, {{100, -1.0, 1.0}}}},
110+
{"h_jet_phi", "jet #phi;#phi_{jet};entries", {HistType::kTH1F, {{80, -1.0, 7.}}}},
111+
{"h_collision_counter", "# of collisions;", {HistType::kTH1F, {{2, 0., 2.}}}},
112+
{"h_jet_counter", ";# of D^{0} jets;", {HistType::kTH1F, {{6, 0., 3.0}}}},
113+
{"h_d0_jet_projection", ";z^{D^{0},jet}_{||};dN/dz^{D^{0},jet}_{||}", {HistType::kTH1F, {{1000, 0., 10.}}}},
114+
{"h_d0_jet_distance_vs_projection", ";#DeltaR_{D^{0},jet};z^{D^{0},jet}_{||}", {HistType::kTH2F, {{1000, 0., 10.}, {1000, 0., 10.}}}},
115+
{"h_d0_jet_distance", ";#DeltaR_{D^{0},jet};dN/d(#DeltaR)", {HistType::kTH1F, {{1000, 0., 10.}}}},
116+
{"h_d0_jet_pt", ";p_{T,D^{0} jet};dN/dp_{T,D^{0} jet}", {HistType::kTH1F, {{200, 0., 10.}}}},
117+
{"h_d0_jet_eta", ";#eta_{T,D^{0} jet};dN/d#eta_{D^{0} jet}", {HistType::kTH1F, {{250, -5., 5.}}}},
118+
{"h_d0_jet_phi", ";#phi_{T,D^{0} jet};dN/d#phi_{D^{0} jet}", {HistType::kTH1F, {{250, -10., 10.}}}},
119+
{"h_d0_mass", ";m_{D^{0}} (GeV/c^{2});dN/dm_{D^{0}}", {HistType::kTH1F, {{1000, 0., 10.}}}},
120+
{"h_d0_eta", ";#eta_{D^{0}} (GeV/c^{2});dN/d#eta_{D_{}}", {HistType::kTH1F, {{250, -5., 5.}}}},
121+
{"h_d0_phi", ";#phi_{D^{0}} (GeV/c^{2});dN/d#phi_{D^{0}}", {HistType::kTH1F, {{250, -10., 10.}}}},
122+
{"h_d0_ang", ";#lambda_{#kappa}^{#alpha};counts", {HistType::kTH1F, {{100, 0., 1.}}}}}
123+
};
124+
125+
// Configurables
126+
Configurable<float> vertexZCut{"vertexZCut", 10.0f, "Accepted z-vertex range"};
127+
128+
Configurable<float> jetPtMin{"jetPtMin", 5.0, "minimum jet pT cut"};
129+
Configurable<float> jetR{"jetR", 0.4, "jet resolution parameter"};
130+
131+
Configurable<std::string> eventSelections{"eventSelections", "sel8", "choose event selection"};
132+
Configurable<std::string> trackSelections{"trackSelections", "globalTracks", "set track selections"};
133+
Configurable<float> kappa{"kappa", 1.0, "angularity kappa"}; //to do: configurable from json
134+
Configurable<float> alpha{"alpha", 1.0, "angularity alpha"};
135+
136+
std::vector<int> eventSelectionBits;
137+
int trackSelection = -1;
138+
139+
// Output table producer
140+
Produces<aod::JetObjTable> ObjJetTable;
141+
142+
float angularity;
143+
float leadingConstituentPt;
144+
145+
void init(o2::framework::InitContext&)
146+
{
147+
eventSelectionBits = jetderiveddatautilities::initialiseEventSelectionBits(static_cast<std::string>(eventSelections));
148+
trackSelection = jetderiveddatautilities::initialiseTrackSelection(static_cast<std::string>(trackSelections));
149+
}
150+
151+
Filter jetCuts = aod::jet::pt > jetPtMin&& aod::jet::r == nround(jetR.node() * 100.0f);
152+
Filter collisionFilter = nabs(aod::jcollision::posZ) < vertexZCut;
153+
154+
template <typename T, typename U>
155+
void jetCalculateAngularity(T const& jet, U const& /*tracks*/)
156+
{
157+
angularity = 0.0;
158+
leadingConstituentPt = 0.0;
159+
for (auto& constituent : jet.template tracks_as<U>()) {
160+
if (constituent.pt() >= leadingConstituentPt) {
161+
leadingConstituentPt = constituent.pt();
162+
}
163+
angularity += std::pow(constituent.pt(), kappa) * std::pow(jetutilities::deltaR(jet, constituent), alpha);
164+
}
165+
angularity /= (jet.pt() * (jet.r() / 100.f));
166+
}
167+
// Collision-level and for Tracks QA
168+
void processCollisions(aod::JetCollision const& collision, aod::JetTracks const& tracks)
169+
{
170+
171+
registry.fill(HIST("h_collisions"), 0.5);
172+
if (!jetderiveddatautilities::selectCollision(collision, eventSelectionBits)) {
173+
return;
174+
}
175+
registry.fill(HIST("h_collisions"), 1.5);
176+
177+
// Loop over tracks and apply track selection
178+
for (auto const& track : tracks) {
179+
if (!jetderiveddatautilities::selectTrack(track, trackSelection)) {
180+
continue;
181+
}
182+
registry.fill(HIST("h_track_pt"), track.pt());
183+
registry.fill(HIST("h_track_eta"), track.eta());
184+
registry.fill(HIST("h_track_phi"), track.phi());
185+
}
186+
}
187+
PROCESS_SWITCH(JetD0Substructure, processCollisions, "process JE collisions", false);
188+
189+
// Charged jet processing (no HF requirement)
190+
void processDataCharged(soa::Filtered<aod::JetCollisions>::iterator const& collision, soa::Filtered<aod::ChargedJets> const& jets)
191+
{
192+
if (!jetderiveddatautilities::selectCollision(collision, eventSelectionBits)) {
193+
return;
194+
}
195+
//jets -> charged jet
196+
for (auto& jet : jets) {
197+
registry.fill(HIST("h_jet_pt"), jet.pt());
198+
registry.fill(HIST("h_jet_eta"), jet.eta());
199+
registry.fill(HIST("h_jet_phi"), jet.phi());
200+
}
201+
}
202+
PROCESS_SWITCH(JetD0Substructure, processDataCharged, "charged jets in data", false);
203+
204+
void processDataChargedSubstructure(aod::JetCollision const& collision,
205+
soa::Join<aod::D0ChargedJets, aod::D0ChargedJetConstituents> const& jets,
206+
aod::CandidatesD0Data const&, aod::JetTracks const& tracks)
207+
{
208+
// apply event selection and fill histograms for sanity check
209+
registry.fill(HIST("h_collision_counter"), 0.5);
210+
211+
if (!jetderiveddatautilities::selectCollision(collision, eventSelectionBits) || !(std::abs(collision.posZ()) < vertexZCut)) {
212+
return;
213+
}
214+
registry.fill(HIST("h_collision_counter"), 1.5);
215+
216+
// Loop over jets containing D0 candidates
217+
for (const auto& jet : jets) {
218+
//number of charged jets with D0
219+
registry.fill(HIST("h_jet_counter"), 0.5);
220+
// obtaining jet 3-vector
221+
TVector3 jetVector(jet.px(), jet.py(), jet.pz());
222+
223+
//Loop over D0 candidates associated to the jet
224+
for (const auto& d0Candidate : jet.candidates_as<aod::CandidatesD0Data>()) {
225+
// obtaining jet 3-vector
226+
TVector3 d0Vector(d0Candidate.px(), d0Candidate.py(), d0Candidate.pz());
227+
228+
// calculating fraction of the jet momentum carried by the D0 along the direction of the jet axis
229+
double zParallel = (jetVector * d0Vector) / (jetVector * jetVector);
230+
231+
// calculating angular distance in eta-phi plane
232+
double axisDistance = jetutilities::deltaR(jet, d0Candidate);
233+
234+
jetCalculateAngularity(jet,tracks);
235+
236+
// filling histograms
237+
registry.fill(HIST("h_d0_jet_projection"), zParallel);
238+
registry.fill(HIST("h_d0_jet_distance_vs_projection"), axisDistance, zParallel);
239+
registry.fill(HIST("h_d0_jet_distance"), axisDistance);
240+
registry.fill(HIST("h_d0_jet_pt"), jet.pt());
241+
registry.fill(HIST("h_d0_jet_eta"), jet.eta());
242+
registry.fill(HIST("h_d0_jet_phi"), jet.phi());
243+
registry.fill(HIST("h_d0_mass"), d0Candidate.m());
244+
registry.fill(HIST("h_d0_eta"), d0Candidate.eta());
245+
registry.fill(HIST("h_d0_phi"), d0Candidate.phi());
246+
registry.fill(HIST("h_d0_ang"), angularity); // add more axis
247+
248+
249+
// filling table
250+
ObjJetTable(axisDistance,
251+
jet.pt(), jet.eta(), jet.phi(), jet.tracks_as<aod::JetTracks>().size(), angularity,
252+
d0Candidate.pt(), d0Candidate.eta(), d0Candidate.phi(), d0Candidate.m(), d0Candidate.y(), d0Candidate.mlScores()[0], d0Candidate.mlScores()[1], d0Candidate.mlScores()[2]);
253+
254+
break; // get out of candidates' loop after first HF particle is found in jet
255+
} // end of D0 candidates loop
256+
257+
} // end of jets loop
258+
259+
} // end of process function
260+
PROCESS_SWITCH(JetD0Substructure, processDataChargedSubstructure, "charged HF jet substructure", false);
261+
262+
263+
264+
};
265+
// Workflow definition
266+
WorkflowSpec defineDataProcessing(ConfigContext const& cfgc) { return WorkflowSpec{adaptAnalysisTask<JetD0Substructure>(cfgc, TaskName{"jet-d0-substructure"})}; }

0 commit comments

Comments
 (0)