Skip to content

Commit e681162

Browse files
Aleksandra MulewiczMyFavoriteGitHub
authored andcommitted
Update MC and jet histograms
1 parent afe2378 commit e681162

1 file changed

Lines changed: 56 additions & 22 deletions

File tree

PWGJE/Tasks/jetHadronsPID.cxx

Lines changed: 56 additions & 22 deletions
Original file line numberDiff line numberDiff line change
@@ -98,6 +98,16 @@ struct PIDHadronsInJets {
9898
if (setMCDefaultItsParams) {
9999
itsResponse.setMCDefaultParameters();
100100
}
101+
registryData.add("n_events", "Event counter", HistType::kTH1F, {{1, 0.5, 1.5, "N_{events}"}});
102+
registryData.add("n_events_raw", "All events", HistType::kTH1F, {{1, 0.5, 1.5, ""}});
103+
104+
registryData.add("jet_pt_raw", "Jet pT raw", HistType::kTH1F, {{200, 0.0, 200.0, "#it{p}_{T}^{raw} (GeV/#it{c})"}});
105+
registryData.add("jet_pt_subtracted", "Jet pT subtracted", HistType::kTH1F, {{200, 0.0, 200.0, "#it{p}_{T}^{sub} (GeV/#it{c})"}});
106+
registryData.add("jet_pt_raw_vs_sub", "Raw vs sub jet pT", HistType::kTH2F, {{200, 0, 200}, {200, 0, 200}});
107+
registryData.add("jet_eta", "Jet eta", HistType::kTH1F, {{100, -1.0, 1.0, "#eta_{jet}"}});
108+
registryData.add("jet_phi", "Jet phi", HistType::kTH1F, {{100, 0.0, TwoPI, "#phi_{jet}"}});
109+
registryData.add("jet_area", "Jet area", HistType::kTH1F, {{100, 0.0, 1.5, "Area"}});
110+
registryData.add("jet_n_constituents", "Jet multiplicity", HistType::kTH1I, {{100, 0, 100, "N_{constituents}"}});
101111

102112
registryData.add("pion_pure_tpc", "TPC Pion PID", HistType::kTH2F, {{120, 0.0, 4.0, "#it{p}_{T} (GeV/#it{c})"}, {200, -3.0, 3.0, "n#sigma_{TPC}"}});
103113
registryData.add("pion_pure_tof", "TOF Pion PID", HistType::kTH2F, {{120, 0.0, 4.0, "#it{p}_{T} (GeV/#it{c})"}, {200, -3.0, 3.0, "n#sigma_{TOF}"}});
@@ -120,7 +130,6 @@ struct PIDHadronsInJets {
120130
registryData.add("proton_pure_dcaz", "Proton DCAz", HistType::kTH1F, {{100, -0.1, 0.1, "DCA_{z} (cm)" }});
121131
registryData.add("proton_pure_dcaxy", "Proton DCAxy", HistType::kTH1F, {{100, -0.1, 0.1, "DCA_{xy} (cm)" }});
122132

123-
124133
registryData.add("z_vtx", "Z-Vertex Distribution", HistType::kTH1F, {{200, -20.0, 20.0, "Z-Vertex (cm)"}});
125134
registryData.add("tracks_in_jets", "Number of Tracks Inside Jets", HistType::kTH1F, {{100, 0, 100, "N_{tracks}"}});
126135
registryData.add("tracks_outside_jets", "Number of Tracks Outside Jets", HistType::kTH1F, {{500, 0, 500, "N_{tracks}"}});
@@ -163,7 +172,17 @@ struct PIDHadronsInJets {
163172

164173
registryData.add("mc_gen_proton_pt", "Generated Proton pT", HistType::kTH1F, {{120, 0.0, 4.0, "#it{p}_{T} (GeV/#it{c})" }});
165174
registryData.add("mc_rec_proton_pt", "Reconstructed Proton pT", HistType::kTH1F, {{120, 0.0, 4.0, "#it{p}_{T} (GeV/#it{c})" }});
166-
175+
176+
registryData.add("rec_pion_all", "All identified pions", HistType::kTH1F, {{120, 0.0, 4.0, "#it{p}_{T}"}});
177+
registryData.add("rec_kaon_all", "All identified kaons", HistType::kTH1F, {{120, 0.0, 4.0, "#it{p}_{T}"}});
178+
registryData.add("rec_proton_all", "All identified protons", HistType::kTH1F, {{120, 0.0, 4.0, "#it{p}_{T}"}});
179+
registryData.add("contamination_matrix_pion", "Pion Contamination Matrix;True absolute PDG;Rec #it{p}_{T} (GeV/#it{c})", HistType::kTH2F, {{3000, -0.5, 2999.5}, {120, 0.0, 4.0}});
180+
registryData.add("contamination_matrix_kaon", "Kaon Contamination Matrix;True absolute PDG;Rec #it{p}_{T} (GeV/#it{c})", HistType::kTH2F, {{3000, -0.5, 2999.5}, {120, 0.0, 4.0}});
181+
registryData.add("contamination_matrix_proton", "Proton Contamination Matrix;True absolute PDG;Rec #it{p}_{T} (GeV/#it{c})", HistType::kTH2F, {{3000, -0.5, 2999.5}, {120, 0.0, 4.0}});
182+
183+
registryData.add("mc_sec_pion_pt", "Reconstructed Pion Secondaries", HistType::kTH1F, {{120, 0.0, 4.0, "#it{p}_{T} (GeV/#it{c})" }});
184+
registryData.add("mc_sec_kaon_pt", "Reconstructed Kaon Secondaries", HistType::kTH1F, {{120, 0.0, 4.0, "#it{p}_{T} (GeV/#it{c})" }});
185+
registryData.add("mc_sec_proton_pt", "Reconstructed Proton Secondaries", HistType::kTH1F, {{120, 0.0, 4.0, "#it{p}_{T} (GeV/#it{c})" }});
167186
}
168187

169188
void getPerpendicularDirections(const TVector3& p, TVector3& u1, TVector3& u2)
@@ -256,6 +275,8 @@ struct PIDHadronsInJets {
256275

257276
void processForJets(SelectedCollisions::iterator const& collision, HadronTracks const& tracks)
258277
{
278+
registryData.fill(HIST("n_events_raw"), 1);
279+
259280
if (!collision.sel8() || std::fabs(collision.posZ()) > zVtx) return;
260281
if (rejectITSROFBorder && !collision.selection_bit(o2::aod::evsel::kNoITSROFrameBorder)) return;
261282
if (rejectTFBorder && !collision.selection_bit(o2::aod::evsel::kNoTimeFrameBorder)) return;
@@ -264,9 +285,8 @@ struct PIDHadronsInJets {
264285
if (requireIsGoodZvtxFT0VsPV && !collision.selection_bit(o2::aod::evsel::kIsGoodZvtxFT0vsPV)) return;
265286
if (requireIsVertexTOFmatched && !collision.selection_bit(o2::aod::evsel::kIsVertexTOFmatched)) return;
266287

267-
268-
269288
//pure histograms
289+
registryData.fill(HIST("n_events"), 1);
270290
registryData.fill(HIST("z_vtx"), collision.posZ());
271291

272292
for (auto const& track : tracks) {
@@ -278,7 +298,6 @@ struct PIDHadronsInJets {
278298
double dcaxy = track.dcaXY();
279299
double dcaz = track.dcaZ();
280300

281-
282301
bool passTpcPi = (std::abs(track.tpcNSigmaPi()) <= 3.0);
283302
bool passTpcKa = (std::abs(track.tpcNSigmaKa()) <= 3.0);
284303
bool passTpcPr = (std::abs(track.tpcNSigmaPr()) <= 3.0);
@@ -336,7 +355,6 @@ struct PIDHadronsInJets {
336355

337356
}
338357

339-
// creating jets
340358
int id(-1);
341359
std::vector<fastjet::PseudoJet> fjParticles;
342360
for (auto const& track : tracks) {
@@ -362,15 +380,20 @@ struct PIDHadronsInJets {
362380

363381
std::set<int> tracksInJetsSet;
364382

365-
// loop over jets and fill histograms
366-
367383
for (const auto& jet : jets) {
384+
registryData.fill(HIST("jet_pt_raw"), jet.pt());
385+
registryData.fill(HIST("jet_eta"), jet.eta());
386+
registryData.fill(HIST("jet_phi"), jet.phi());
387+
registryData.fill(HIST("jet_area"), jet.area());
388+
registryData.fill(HIST("jet_n_constituents"), jet.constituents().size());
368389

369390
if (!isppRefAnalysis && ((std::fabs(jet.eta()) + rJet) > (maxEta - deltaEtaEdge))) continue;
370391
if (isppRefAnalysis && std::fabs(jet.eta()) > cfgEtaJetMax) continue;
371392

372393
auto jetForSub = jet;
373394
fastjet::PseudoJet jetMinusBkg = backgroundSub.doRhoAreaSub(jetForSub, rhoPerp, rhoMPerp);
395+
registryData.fill(HIST("jet_pt_subtracted"), jetMinusBkg.pt());
396+
registryData.fill(HIST("jet_pt_raw_vs_sub"), jet.pt(), jetMinusBkg.pt());
374397

375398
if (isppRefAnalysis && (jet.pt() < minJetPt || jet.pt() > maxJetPt)) continue;
376399
if (!isppRefAnalysis && (jetMinusBkg.pt() < minJetPt || jetMinusBkg.pt() > maxJetPt)) continue;
@@ -388,8 +411,6 @@ struct PIDHadronsInJets {
388411

389412
std::vector<fastjet::PseudoJet> jetConstituents = jet.constituents();
390413

391-
392-
//loop over particles in jets
393414
for (const auto& particle : jetConstituents) {
394415

395416
int trackIdx = particle.user_index();
@@ -462,8 +483,6 @@ struct PIDHadronsInJets {
462483

463484
}
464485

465-
466-
//loop for UE
467486
for (auto const& track : tracks) {
468487

469488
if (!passedTrackSelection(track)) continue;
@@ -482,12 +501,10 @@ struct PIDHadronsInJets {
482501
maxConeRadius = std::sqrt(maxNormalizedJetArea) * rJet;
483502
}
484503

485-
// Only process tracks inside the UE cone
486504
if (deltaRUe1 > maxConeRadius && deltaRUe2 > maxConeRadius) continue;
487505

488506
double pt = track.pt();
489507

490-
491508
bool passTpcPi = (std::abs(track.tpcNSigmaPi()) <= 3.0);
492509
bool passTpcKa = (std::abs(track.tpcNSigmaKa()) <= 3.0);
493510
bool passTpcPr = (std::abs(track.tpcNSigmaPr()) <= 3.0);
@@ -549,7 +566,7 @@ struct PIDHadronsInJets {
549566
}
550567

551568
}
552-
} // end jet loop
569+
}
553570

554571
registryData.fill(HIST("tracks_in_jets"), tracksInJetsSet.size());
555572

@@ -566,8 +583,6 @@ struct PIDHadronsInJets {
566583

567584
PROCESS_SWITCH(PIDHadronsInJets, processForJets, "Process jets", true);
568585

569-
570-
571586
void processMC(SelectedCollisions::iterator const& collision, aod::McParticles const& mcParticles, HadronTracksMC const& tracks)
572587
{
573588
if (!collision.sel8() || std::fabs(collision.posZ()) > zVtx) return;
@@ -613,20 +628,39 @@ struct PIDHadronsInJets {
613628

614629
if (!track.has_mcParticle()) continue;
615630
auto const& trueParticle = track.mcParticle();
616-
if (!trueParticle.isPhysicalPrimary()) continue;
617-
631+
618632
int pdg = std::abs(trueParticle.pdgCode());
633+
bool isPrimary = trueParticle.isPhysicalPrimary();
619634

620635
if (passTpcPi && (pt < ptThreshold || passTofPi)) {
621-
if (pdg == 211) registryData.fill(HIST("mc_rec_pion_pt"), pt);
636+
registryData.fill(HIST("rec_pion_all"), pt); //reconstructed after cuts
637+
registryData.fill(HIST("contamination_matrix_pion"), pdg, pt);
638+
639+
if (isPrimary) {
640+
if (pdg == 211) registryData.fill(HIST("mc_rec_pion_pt"), pt); // Efficiency
641+
} else {
642+
registryData.fill(HIST("mc_sec_pion_pt"), pt);
643+
}
622644
}
623645

624646
if (passTpcKa && (pt < ptThreshold || passTofKa)) {
625-
if (pdg == 321) registryData.fill(HIST("mc_rec_kaon_pt"), pt);
647+
registryData.fill(HIST("rec_kaon_all"), pt);
648+
registryData.fill(HIST("contamination_matrix_kaon"), pdg, pt);
649+
if (isPrimary) {
650+
if (pdg == 321) registryData.fill(HIST("mc_rec_kaon_pt"), pt);
651+
} else {
652+
registryData.fill(HIST("mc_sec_kaon_pt"), pt);
653+
}
626654
}
627655

628656
if (passTpcPr && (pt < ptThreshold || passTofPr)) {
629-
if (pdg == 2212) registryData.fill(HIST("mc_rec_proton_pt"), pt);
657+
registryData.fill(HIST("rec_proton_all"), pt);
658+
registryData.fill(HIST("contamination_matrix_proton"), pdg, pt);
659+
if (isPrimary) {
660+
if (pdg == 2212) registryData.fill(HIST("mc_rec_proton_pt"), pt);
661+
} else {
662+
registryData.fill(HIST("mc_sec_proton_pt"), pt);
663+
}
630664
}
631665
}
632666
}

0 commit comments

Comments
 (0)