|
| 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 | +/// |
| 13 | +/// \file qaPIDTPCSignal.cxx |
| 14 | +/// \author Nicolò Jacazio nicolo.jacazio@cern.ch |
| 15 | +/// \since 24/03/2023 |
| 16 | +/// \brief Implementation for QA tasks of the TPC signal (lite task, refer to qaPIDTPC.cxx for full TPC QA for PID) |
| 17 | +/// |
| 18 | + |
| 19 | +#include "Framework/AnalysisTask.h" |
| 20 | +#include "Framework/runDataProcessing.h" |
| 21 | +#include "Framework/HistogramRegistry.h" |
| 22 | +#include "Framework/StaticFor.h" |
| 23 | +#include "Common/DataModel/TrackSelectionTables.h" |
| 24 | +#include "Common/DataModel/EventSelection.h" |
| 25 | + |
| 26 | +using namespace o2; |
| 27 | +using namespace o2::framework; |
| 28 | +using namespace o2::framework::expressions; |
| 29 | + |
| 30 | +/// Task to produce the TPC QA plots |
| 31 | +struct tpcPidQaSignal { |
| 32 | + HistogramRegistry histos{"Histos", {}, OutputObjHandlingPolicy::AnalysisObject}; |
| 33 | + |
| 34 | + Configurable<int> logAxis{"logAxis", 1, "Flag to use a log momentum axis"}; |
| 35 | + Configurable<int> nBinsP{"nBinsP", 3000, "Number of bins for the momentum"}; |
| 36 | + Configurable<float> minP{"minP", 0.01, "Minimum momentum in range"}; |
| 37 | + Configurable<float> maxP{"maxP", 20, "Maximum momentum in range"}; |
| 38 | + ConfigurableAxis dEdxBins{"dEdxBins", {5000, 0.f, 5000.f}, "Binning in dE/dx"}; |
| 39 | + Configurable<int> applyEvSel{"applyEvSel", 2, "Flag to apply event selection cut: 0 -> no event selection, 1 -> Run 2 event selection, 2 -> Run 3 event selection"}; |
| 40 | + Configurable<float> minTPCNcls{"minTPCNcls", 0.f, "Minimum number or TPC Clusters for tracks"}; |
| 41 | + Configurable<float> minNClsCrossedRows{"minNClsCrossedRows", 70.f, "Minimum number or TPC crossed rows for tracks"}; |
| 42 | + |
| 43 | + void init(o2::framework::InitContext&) |
| 44 | + { |
| 45 | + const AxisSpec vtxZAxis{100, -20, 20, "Vtx_{z} (cm)"}; |
| 46 | + AxisSpec pAxis{nBinsP, minP, maxP, "#it{p}/|Z| (GeV/#it{c})"}; |
| 47 | + if (logAxis) { |
| 48 | + pAxis.makeLogarithmic(); |
| 49 | + } |
| 50 | + const AxisSpec dedxAxis{dEdxBins, "d#it{E}/d#it{x} Arb. units"}; |
| 51 | + const AxisSpec chargeAxis{2, -2.f, 2.f, "Charge"}; |
| 52 | + |
| 53 | + // Event properties |
| 54 | + auto h = histos.add<TH1>("event/evsel", "", kTH1D, {{10, 0.5, 10.5, "Ev. Sel."}}); |
| 55 | + h->GetXaxis()->SetBinLabel(1, "Events read"); |
| 56 | + h->GetXaxis()->SetBinLabel(2, "Passed ev. sel."); |
| 57 | + h->GetXaxis()->SetBinLabel(3, "Passed vtx Z"); |
| 58 | + |
| 59 | + histos.add("event/vertexz", "", kTH1D, {vtxZAxis}); |
| 60 | + histos.add("event/tpcsignal", "", kTH3F, {pAxis, dedxAxis, chargeAxis}); |
| 61 | + LOG(info) << "QA PID TPC histograms:"; |
| 62 | + histos.print(); |
| 63 | + } |
| 64 | + |
| 65 | + Filter trackFilterEta = (nabs(aod::track::eta) < 0.8f); |
| 66 | + Filter trackFilterITS = (aod::track::itsChi2NCl < 36.f); |
| 67 | + Filter trackFilterTPC = ((aod::track::tpcChi2NCl < 4.f)); |
| 68 | + using TrackCandidates = soa::Join<aod::TracksIU, aod::TracksExtra>; |
| 69 | + void processEvSel(soa::Join<aod::Collisions, aod::EvSels>::iterator const& collision, |
| 70 | + soa::Filtered<TrackCandidates> const& tracks) |
| 71 | + { |
| 72 | + |
| 73 | + histos.fill(HIST("event/evsel"), 1); |
| 74 | + if (!collision.sel8()) { |
| 75 | + return; |
| 76 | + } |
| 77 | + |
| 78 | + histos.fill(HIST("event/evsel"), 2); |
| 79 | + |
| 80 | + if (abs(collision.posZ()) > 10.f) { |
| 81 | + return; |
| 82 | + } |
| 83 | + histos.fill(HIST("event/evsel"), 3); |
| 84 | + histos.fill(HIST("event/vertexz"), collision.posZ()); |
| 85 | + |
| 86 | + for (const auto& t : tracks) { |
| 87 | + if (t.itsNCls() < 6) { |
| 88 | + continue; |
| 89 | + } |
| 90 | + if (t.tpcNClsFindable() <= 0) { |
| 91 | + continue; |
| 92 | + } |
| 93 | + if (t.tpcNClsFound() < minTPCNcls) { |
| 94 | + continue; |
| 95 | + } |
| 96 | + if (t.tpcCrossedRowsOverFindableCls() < 0.8f) { |
| 97 | + continue; |
| 98 | + } |
| 99 | + if (t.tpcNClsCrossedRows() < minNClsCrossedRows) { |
| 100 | + continue; |
| 101 | + } |
| 102 | + histos.fill(HIST("event/tpcsignal"), t.tpcInnerParam(), t.tpcSignal(), t.sign()); |
| 103 | + } |
| 104 | + } |
| 105 | + PROCESS_SWITCH(tpcPidQaSignal, processEvSel, "Process with event selection", false); |
| 106 | + |
| 107 | + void processNoEvSel(aod::Collision const& collision, |
| 108 | + soa::Filtered<TrackCandidates> const& tracks) |
| 109 | + { |
| 110 | + |
| 111 | + histos.fill(HIST("event/evsel"), 1); |
| 112 | + histos.fill(HIST("event/evsel"), 2); |
| 113 | + |
| 114 | + if (abs(collision.posZ()) > 10.f) { |
| 115 | + return; |
| 116 | + } |
| 117 | + histos.fill(HIST("event/evsel"), 3); |
| 118 | + histos.fill(HIST("event/vertexz"), collision.posZ()); |
| 119 | + |
| 120 | + for (const auto& t : tracks) { |
| 121 | + if (t.itsNCls() < 6) { |
| 122 | + continue; |
| 123 | + } |
| 124 | + if (t.tpcNClsFindable() <= 0) { |
| 125 | + continue; |
| 126 | + } |
| 127 | + if (t.tpcNClsFound() < minTPCNcls) { |
| 128 | + continue; |
| 129 | + } |
| 130 | + if (t.tpcCrossedRowsOverFindableCls() < 0.8f) { |
| 131 | + continue; |
| 132 | + } |
| 133 | + if (t.tpcNClsCrossedRows() < minNClsCrossedRows) { |
| 134 | + continue; |
| 135 | + } |
| 136 | + histos.fill(HIST("event/tpcsignal"), t.tpcInnerParam(), t.tpcSignal(), t.sign()); |
| 137 | + } |
| 138 | + } |
| 139 | + PROCESS_SWITCH(tpcPidQaSignal, processNoEvSel, "Process without event selection", true); |
| 140 | +}; |
| 141 | + |
| 142 | +WorkflowSpec defineDataProcessing(ConfigContext const& cfgc) { return WorkflowSpec{adaptAnalysisTask<tpcPidQaSignal>(cfgc)}; } |
0 commit comments