Project
Loading...
Searching...
No Matches
TPCTimeSeriesSpec.cxx
Go to the documentation of this file.
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
16
17#include "Framework/Task.h"
18#include "Framework/Logger.h"
23#include "TPCBase/Mapper.h"
31#include "MathUtils/Tsallis.h"
39#include <random>
40#include <chrono>
43#include "TROOT.h"
48
49using namespace o2::globaltracking;
53using namespace o2::framework;
54
55namespace o2
56{
57namespace tpc
58{
59
60class TPCTimeSeries : public Task
61{
62 public:
64 TPCTimeSeries(std::shared_ptr<o2::base::GRPGeomRequest> req, const bool disableWriter, const o2::base::Propagator::MatCorrType matType, const bool enableUnbinnedWriter, const bool tpcOnly, std::shared_ptr<o2::globaltracking::DataRequest> dr) : mCCDBRequest(req), mDisableWriter(disableWriter), mMatType(matType), mUnbinnedWriter(enableUnbinnedWriter), mTPCOnly(tpcOnly), mDataRequest(dr) {};
65
67 {
69 mNMaxTracks = ic.options().get<int>("max-tracks");
70 mMinMom = ic.options().get<float>("min-momentum");
71 mMinNCl = ic.options().get<int>("min-cluster");
72 mMaxTgl = ic.options().get<float>("max-tgl");
73 mMaxQPt = ic.options().get<float>("max-qPt");
74 mCoarseStep = ic.options().get<float>("coarse-step");
75 mFineStep = ic.options().get<float>("fine-step");
76 mCutDCA = ic.options().get<float>("cut-DCA-median");
77 mCutRMS = ic.options().get<float>("cut-DCA-RMS");
78 mRefXSec = ic.options().get<float>("refX-for-sector");
79 mTglBins = ic.options().get<int>("tgl-bins");
80 mPhiBins = ic.options().get<int>("phi-bins");
81 mQPtBins = ic.options().get<int>("qPt-bins");
82 mNThreads = ic.options().get<int>("threads");
83 maxITSTPCDCAr = ic.options().get<float>("max-ITS-TPC-DCAr");
84 maxITSTPCDCAz = ic.options().get<float>("max-ITS-TPC-DCAz");
85 maxITSTPCDCAr_comb = ic.options().get<float>("max-ITS-TPC-DCAr_comb");
86 maxITSTPCDCAz_comb = ic.options().get<float>("max-ITS-TPC-DCAz_comb");
87 mTimeWindowMUS = ic.options().get<float>("time-window-mult-mus");
88 mMIPdEdx = ic.options().get<float>("MIP-dedx");
89 mMaxSnp = ic.options().get<float>("max-snp");
90 mXCoarse = ic.options().get<float>("mX-coarse");
91 mSqrt = ic.options().get<float>("sqrts");
92 mMultBins = ic.options().get<int>("mult-bins");
93 mMultMax = ic.options().get<int>("mult-max");
94 mMinTracksPerVertex = ic.options().get<int>("min-tracks-per-vertex");
95 mMaxdEdxRatio = ic.options().get<float>("max-dedx-ratio");
96 mMaxdEdxRegionRatio = ic.options().get<float>("max-dedx-region-ratio");
97 mSamplingFactor = ic.options().get<float>("sampling-factor");
98 mSampleTsallis = ic.options().get<bool>("sample-unbinned-tsallis");
99 mXOuterMatching = ic.options().get<float>("refX-for-outer-ITS");
100 mUseMinBiasTrigger = !ic.options().get<bool>("disable-min-bias-trigger");
101 mMaxOccupancyHistBins = ic.options().get<int>("max-occupancy-bins");
102
103 if (mUnbinnedWriter) {
104 for (int iThread = 0; iThread < mNThreads; ++iThread) {
105 mGenerator.emplace_back(std::mt19937(std::random_device{}()));
106 }
107 }
108 mBufferVals.resize(mNThreads);
109 mBufferDCA.setBinning(mPhiBins, mTglBins, mQPtBins, mMultBins, mMaxTgl, mMaxQPt, mMultMax);
110 if (mUnbinnedWriter) {
111 std::string outfile = ic.options().get<std::string>("out-file-unbinned");
112 if (mNThreads > 1) {
113 ROOT::EnableThreadSafety();
114 }
115 mStreamer.resize(mNThreads);
116 for (int iThread = 0; iThread < mNThreads; ++iThread) {
117 std::string outfileThr = outfile;
118 outfileThr.replace(outfileThr.length() - 5, outfileThr.length(), fmt::format("_{}.root", iThread));
119 LOGP(info, "Writing unbinned data to: {}", outfileThr);
120 mStreamer[iThread] = std::make_unique<o2::utils::TreeStreamRedirector>(outfileThr.data(), "recreate");
121 }
122 }
123 }
124
125 void run(ProcessingContext& pc) final
126 {
128 mTPCVDriftHelper.extractCCDBInputs(pc);
129 mPTHelper.extractCCDBInputs(pc);
130 if (mTPCVDriftHelper.isUpdated()) {
131 mTPCVDriftHelper.acknowledgeUpdate();
132 mVDrift = mTPCVDriftHelper.getVDriftObject().getVDrift();
133 LOGP(info, "Updated reference drift velocity to: {}", mVDrift);
134 }
135 mBufferDCA.mVDrift = mVDrift;
136
137 const int nBins = getNBins();
138
141 mBufferDCA.mTemperature = mPTHelper.getMeanTemperature(mTimeMS);
142 mBufferDCA.mPressure = mPTHelper.getPressure(mTimeMS);
143
144 // init only once
145 if (mAvgADCAr.size() != nBins) {
146 mBufferDCA.resize(nBins);
147 mAvgADCAr.resize(nBins);
148 mAvgCDCAr.resize(nBins);
149 mAvgADCAz.resize(nBins);
150 mAvgCDCAz.resize(nBins);
151 mAvgMeffA.resize(nBins);
152 mAvgMeffC.resize(nBins);
153 mAvgChi2MatchA.resize(nBins);
154 mAvgChi2MatchC.resize(nBins);
155 mMIPdEdxRatioQMaxA.resize(nBins);
156 mMIPdEdxRatioQMaxC.resize(nBins);
157 mMIPdEdxRatioQTotA.resize(nBins);
158 mMIPdEdxRatioQTotC.resize(nBins);
159 mTPCChi2A.resize(nBins);
160 mTPCChi2C.resize(nBins);
161 mTPCNClA.resize(nBins);
162 mTPCNClC.resize(nBins);
163 mLogdEdxQTotA.resize(nBins);
164 mLogdEdxQTotC.resize(nBins);
165 mLogdEdxQMaxA.resize(nBins);
166 mLogdEdxQMaxC.resize(nBins);
167 mITSPropertiesA.resize(nBins);
168 mITSPropertiesC.resize(nBins);
169 mITSTPCDeltaPA.resize(nBins);
170 mITSTPCDeltaPC.resize(nBins);
171 mSigmaYZA.resize(nBins);
172 mSigmaYZC.resize(nBins);
173 }
174
175 RecoContainer recoData;
176 recoData.collectData(pc, *mDataRequest.get());
177
178 // getting tracks
179 auto tracksTPC = recoData.getTPCTracks();
180 auto tracksITSTPC = mTPCOnly ? gsl::span<o2::dataformats::TrackTPCITS>() : recoData.getTPCITSTracks();
181 auto tracksITS = mTPCOnly ? gsl::span<o2::its::TrackITS>() : recoData.getITSTracks();
182
183 // getting the vertices
184 auto vertices = mTPCOnly ? gsl::span<o2::dataformats::PrimaryVertex>() : recoData.getPrimaryVertices();
185 auto primMatchedTracks = mTPCOnly ? gsl::span<o2::dataformats::VtxTrackIndex>() : recoData.getPrimaryVertexMatchedTracks(); // Global ID's for associated tracks
186 auto primMatchedTracksRef = mTPCOnly ? gsl::span<o2::dataformats::VtxTrackRef>() : recoData.getPrimaryVertexMatchedTrackRefs(); // references from vertex to these track IDs
187
188 // get occupancy map
189 mBufferDCA.mOccupancyMapTPC = std::vector<unsigned int>(recoData.occupancyMapTPC.begin(), recoData.occupancyMapTPC.end());
190 if (mBufferDCA.mOccupancyMapTPC.size() > mMaxOccupancyHistBins) {
191 mBufferDCA.mOccupancyMapTPC.resize(mMaxOccupancyHistBins);
192 }
193
194 // TOF clusters
195 const auto& tofClusters = mTPCOnly ? gsl::span<o2::tof::Cluster>() : recoData.getTOFClusters();
196
197 LOGP(info, "Processing {} vertices, {} primary matched vertices, {} TPC tracks, {} ITS tracks, {} ITS-TPC tracks, {} TOF clusters", vertices.size(), primMatchedTracks.size(), tracksTPC.size(), tracksITS.size(), tracksITSTPC.size(), tofClusters.size());
198
199 // calculate mean vertex, RMS and count vertices
200 auto indicesITSTPC_vtx = processVertices(vertices, primMatchedTracks, primMatchedTracksRef, recoData);
201
202 // storing indices to ITS-TPC tracks and vertex ID for tpc track
203 std::unordered_map<unsigned int, std::array<int, 2>> indicesITSTPC; // TPC track index -> ITS-TPC track index, vertex ID
204 // loop over all ITS-TPC tracks
205 for (int i = 0; i < tracksITSTPC.size(); ++i) {
206 auto it = indicesITSTPC_vtx.find(i);
207 // check if ITS-TPC track has attached vertex
208 const auto idxVtx = (it != indicesITSTPC_vtx.end()) ? (it->second) : -1;
209 // store TPC index and ITS-TPC+vertex index
210 indicesITSTPC[tracksITSTPC[i].getRefTPC().getIndex()] = {i, idxVtx};
211 }
212
213 std::vector<std::tuple<int, float, float, o2::track::TrackLTIntegral, double, float, unsigned int>> idxTPCTrackToTOFCluster; // store for each tpc track index the index to the TOF cluster
214
215 // get matches to TOF in case skimmed data is produced
216 if (mUnbinnedWriter) {
217 // getLTIntegralOut(), ///< L,TOF integral calculated during the propagation
218 // getSignal() mSignal = 0.0; ///< TOF time in ps
220 idxTPCTrackToTOFCluster = std::vector<std::tuple<int, float, float, o2::track::TrackLTIntegral, double, float, unsigned int>>(tracksTPC.size(), {-1, -999, -999, defLT, 0, 0, 0});
221 const std::vector<gsl::span<const o2::dataformats::MatchInfoTOF>> tofMatches{recoData.getTPCTOFMatches(), recoData.getTPCTRDTOFMatches(), recoData.getITSTPCTOFMatches(), recoData.getITSTPCTRDTOFMatches()};
222
223 const auto& ft0rec = recoData.getFT0RecPoints();
224 // fill available FT0-AC event times vs BClong
225 std::map<ULong64_t, short> t0array;
226 for (const auto& t0 : ft0rec) {
227 if (!(t0.isValidTime(1) && t0.isValidTime(2))) { // skip if !(A & C)
228 continue;
229 }
230
231 auto bclong = t0.mIntRecord.differenceInBC(recoData.startIR);
232 if (t0array.find(bclong) == t0array.end()) { // add if it doesn't exist
233 t0array.emplace(std::make_pair(bclong, t0.getCollisionTime(0)));
234 }
235 }
236
237 static const double BC_TIME_INPS_INV = 1E-3 / o2::constants::lhc::LHCBunchSpacingNS;
238
239 // loop over ITS-TPC-TRD-TOF and ITS-TPC-TOF tracks an store for each ITS-TPC track the TOF track index
240 for (const auto& tofMatch : tofMatches) {
241 for (const auto& tpctofmatch : tofMatch) {
242 auto refTPC = recoData.getTPCContributorGID(tpctofmatch.getTrackRef());
243 if (refTPC.isIndexSet()) {
244 o2::track::TrackLTIntegral ltIntegral = tpctofmatch.getLTIntegralOut();
245 ULong64_t bclongtof = (tpctofmatch.getSignal() - 10000) * BC_TIME_INPS_INV;
246 double t0 = 0; // bclongtof * o2::constants::lhc::LHCBunchSpacingNS * 1E3; // if you want to subtract also the BC uncomment this part (-> tofsignal can be a float)
247 unsigned int mask = 0;
248 if (!(t0array.find(bclongtof) == t0array.end())) { // subtract FT0-AC if it exists in the same BC
249 t0 += t0array.find(bclongtof)->second;
250 mask |= o2::dataformats::MatchInfoTOF::QualityFlags::hasT0sameBC; // 8th bit if FT0-AC in same BC
251 }
252
253 double signal = tpctofmatch.getSignal() - t0;
254 float deltaT = tpctofmatch.getDeltaT();
255
256 float dy = tpctofmatch.getDYatTOF(); // residual orthogonal to the strip (it should be close to zero)
257 bool isMultiHitZ = tpctofmatch.getHitPatternUpDown();
258 bool isMultiHitX = tpctofmatch.getHitPatternLeftRight();
259 bool isMultiStripMatch = tpctofmatch.getChi2() < 1E-9;
260 float chi2 = tpctofmatch.getChi2();
261 bool hasT0_1BCbefore = (t0array.find(bclongtof - 1) != t0array.end());
262 bool hasT0_2BCbefore = (t0array.find(bclongtof - 2) != t0array.end());
263
264 if (isMultiHitX) { // 1nd bit on if multiple hits along X
266 }
267 if (isMultiHitZ) { // 2nd bit on if multiple hits along Z
269 }
270 if (fabs(dy) > 0.5) { // 3rd bit on if Y-residual too large
272 }
273 if (isMultiStripMatch) { // 4th bit on if two strips fired
275 }
276 if (chi2 > 1E-4) { // 5th bit on if chi2 > 1E-4 -> not inside the pad
278 }
279 if (chi2 > 3) { // 6th bit on if chi2 > 3
281 }
282 if (chi2 > 5) { // 7th bit on if chi2 > 5
284 }
285 if (hasT0_1BCbefore) { // 9th bit if FT0-AC also BC before
287 }
288 if (hasT0_2BCbefore) { // 10th bit if FT0-AC also 2BCs before
290 }
291
292 idxTPCTrackToTOFCluster[refTPC] = {tpctofmatch.getIdxTOFCl(), tpctofmatch.getDXatTOF(), tpctofmatch.getDZatTOF(), ltIntegral, signal, deltaT, mask};
293 }
294 }
295 }
296 }
297
298 // find nearest vertex of tracks which have no vertex assigned
299 findNearesVertex(tracksTPC, vertices);
300
301 // getting cluster references for cluster bitmask
302 if (mUnbinnedWriter) {
303 mTPCTrackClIdx = pc.inputs().get<gsl::span<o2::tpc::TPCClRefElem>>("trackTPCClRefs");
304 mFirstTFOrbit = processing_helpers::getFirstTForbit(pc);
305 }
306
307 // get local multiplicity - count neighbouring tracks
308 findNNeighbourTracks(tracksTPC);
309
310 // reset buffers
311 for (int i = 0; i < nBins; ++i) {
312 for (int type = 0; type < mAvgADCAr[i].size(); ++type) {
313 mAvgADCAr[i][type].clear();
314 mAvgCDCAr[i][type].clear();
315 mAvgADCAz[i][type].clear();
316 mAvgCDCAz[i][type].clear();
317 }
318 for (int type = 0; type < mMIPdEdxRatioQMaxA[i].size(); ++type) {
319 mMIPdEdxRatioQMaxA[i][type].clear();
320 mMIPdEdxRatioQMaxC[i][type].clear();
321 mMIPdEdxRatioQTotA[i][type].clear();
322 mMIPdEdxRatioQTotC[i][type].clear();
323 mTPCChi2A[i][type].clear();
324 mTPCChi2C[i][type].clear();
325 mTPCNClA[i][type].clear();
326 mTPCNClC[i][type].clear();
327 mMIPdEdxRatioQMaxA[i][type].setUseWeights(false);
328 mMIPdEdxRatioQMaxC[i][type].setUseWeights(false);
329 mMIPdEdxRatioQTotA[i][type].setUseWeights(false);
330 mMIPdEdxRatioQTotC[i][type].setUseWeights(false);
331 mTPCChi2A[i][type].setUseWeights(false);
332 mTPCChi2C[i][type].setUseWeights(false);
333 mTPCNClA[i][type].setUseWeights(false);
334 mTPCNClC[i][type].setUseWeights(false);
335 }
336 for (int type = 0; type < mLogdEdxQTotA[i].size(); ++type) {
337 mLogdEdxQTotA[i][type].clear();
338 mLogdEdxQTotC[i][type].clear();
339 mLogdEdxQMaxA[i][type].clear();
340 mLogdEdxQMaxC[i][type].clear();
341 mLogdEdxQTotA[i][type].setUseWeights(false);
342 mLogdEdxQTotC[i][type].setUseWeights(false);
343 mLogdEdxQMaxA[i][type].setUseWeights(false);
344 mLogdEdxQMaxC[i][type].setUseWeights(false);
345 }
346 for (int j = 0; j < mITSPropertiesA[i].size(); ++j) {
347 mITSPropertiesA[i][j].clear();
348 mITSPropertiesC[i][j].clear();
349 mITSPropertiesA[i][j].setUseWeights(false);
350 mITSPropertiesC[i][j].setUseWeights(false);
351 }
352
353 for (int j = 0; j < mITSTPCDeltaPA[i].size(); ++j) {
354 mITSTPCDeltaPA[i][j].clear();
355 mITSTPCDeltaPC[i][j].clear();
356 mITSTPCDeltaPA[i][j].setUseWeights(false);
357 mITSTPCDeltaPC[i][j].setUseWeights(false);
358 }
359
360 for (int j = 0; j < mSigmaYZA[i].size(); ++j) {
361 mSigmaYZA[i][j].clear();
362 mSigmaYZC[i][j].clear();
363 mSigmaYZA[i][j].setUseWeights(false);
364 mSigmaYZC[i][j].setUseWeights(false);
365 }
366
367 for (int j = 0; j < mAvgMeffA[i].size(); ++j) {
368 mAvgMeffA[i][j].clear();
369 mAvgMeffC[i][j].clear();
370 mAvgChi2MatchA[i][j].clear();
371 mAvgChi2MatchC[i][j].clear();
372 mAvgMeffA[i][j].setUseWeights(false);
373 mAvgMeffC[i][j].setUseWeights(false);
374 mAvgChi2MatchA[i][j].setUseWeights(false);
375 mAvgChi2MatchC[i][j].setUseWeights(false);
376 }
377 }
378
379 for (int i = 0; i < mNThreads; ++i) {
380 mBufferVals[i].front().clear();
381 mBufferVals[i].back().clear();
382 }
383
384 // define number of tracks which are used
385 const auto nTracks = tracksTPC.size();
386 const size_t loopEnd = (mNMaxTracks < 0) ? nTracks : ((mNMaxTracks > nTracks) ? nTracks : size_t(mNMaxTracks));
387
388 // reserve memory
389 for (int i = 0; i < nBins; ++i) {
390 const int lastIdxPhi = mBufferDCA.mTSTPC.getIndexPhi(mPhiBins);
391 const int lastIdxTgl = mBufferDCA.mTSTPC.getIndexTgl(mTglBins);
392 const int lastIdxQPt = mBufferDCA.mTSTPC.getIndexqPt(mQPtBins);
393 const int lastIdxMult = mBufferDCA.mTSTPC.getIndexMult(mMultBins);
394 const int firstIdx = mBufferDCA.mTSTPC.getIndexInt();
395 int resMem = 0;
396 if (i < lastIdxPhi) {
397 resMem = loopEnd / mPhiBins;
398 } else if (i < lastIdxTgl) {
399 resMem = loopEnd / mTglBins;
400 } else if (i < lastIdxQPt) {
401 resMem = loopEnd / mQPtBins;
402 } else if (i < lastIdxMult) {
403 resMem = loopEnd / mMultBins;
404 } else {
405 resMem = loopEnd;
406 }
407 // Divide by 2 for A-C-side
408 resMem /= 2;
409 for (int type = 0; type < mAvgADCAr[i].size(); ++type) {
410 mAvgADCAr[i][type].reserve(resMem);
411 mAvgCDCAr[i][type].reserve(resMem);
412 mAvgADCAz[i][type].reserve(resMem);
413 mAvgCDCAz[i][type].reserve(resMem);
414 }
415 for (int type = 0; type < mMIPdEdxRatioQMaxA[i].size(); ++type) {
416 mMIPdEdxRatioQMaxA[i][type].reserve(resMem);
417 mMIPdEdxRatioQMaxC[i][type].reserve(resMem);
418 mMIPdEdxRatioQTotA[i][type].reserve(resMem);
419 mMIPdEdxRatioQTotC[i][type].reserve(resMem);
420 mTPCChi2A[i][type].reserve(resMem);
421 mTPCChi2C[i][type].reserve(resMem);
422 mTPCNClA[i][type].reserve(resMem);
423 mTPCNClC[i][type].reserve(resMem);
424 }
425 for (int j = 0; j < mAvgMeffA[i].size(); ++j) {
426 mAvgMeffA[i][j].reserve(resMem);
427 mAvgMeffC[i][j].reserve(resMem);
428 mAvgChi2MatchA[i][j].reserve(resMem);
429 mAvgChi2MatchC[i][j].reserve(resMem);
430 }
431 for (int j = 0; j < mITSPropertiesA[i].size(); ++j) {
432 mITSPropertiesA[i][j].reserve(resMem);
433 mITSPropertiesC[i][j].reserve(resMem);
434 }
435 for (int j = 0; j < mITSTPCDeltaPA[i].size(); ++j) {
436 mITSTPCDeltaPA[i][j].reserve(resMem);
437 mITSTPCDeltaPC[i][j].reserve(resMem);
438 }
439 for (int j = 0; j < mSigmaYZA[i].size(); ++j) {
440 mSigmaYZA[i][j].reserve(resMem);
441 mSigmaYZC[i][j].reserve(resMem);
442 }
443 for (int j = 0; j < mAvgMeffA[i].size(); ++j) {
444 mLogdEdxQTotA[i][j].reserve(resMem);
445 mLogdEdxQTotC[i][j].reserve(resMem);
446 mLogdEdxQMaxA[i][j].reserve(resMem);
447 mLogdEdxQMaxC[i][j].reserve(resMem);
448 }
449 }
450 for (int iThread = 0; iThread < mNThreads; ++iThread) {
451 const int resMem = (mNThreads > 0) ? loopEnd / mNThreads : loopEnd;
452 mBufferVals[iThread].front().reserve(loopEnd, 1);
453 mBufferVals[iThread].back().reserve(loopEnd, 0);
454 }
455
456 using timer = std::chrono::high_resolution_clock;
457 auto startTotal = timer::now();
458
459 // loop over tracks and calculate DCAs
460 if (loopEnd < nTracks) {
461 // draw random tracks
462 std::vector<size_t> ind(nTracks);
463 std::iota(ind.begin(), ind.end(), 0);
464 std::minstd_rand rng(std::time(nullptr));
465 std::shuffle(ind.begin(), ind.end(), rng);
466
467 auto myThread = [&](int iThread) {
468 for (size_t i = iThread; i < loopEnd; i += mNThreads) {
469 if (acceptTrack(tracksTPC[i])) {
470 fillDCA(tracksTPC, tracksITSTPC, vertices, i, iThread, indicesITSTPC, tracksITS, idxTPCTrackToTOFCluster, tofClusters);
471 }
472 }
473 };
474
475 std::vector<std::thread> threads(mNThreads);
476 for (int i = 0; i < mNThreads; i++) {
477 threads[i] = std::thread(myThread, i);
478 }
479
480 for (auto& th : threads) {
481 th.join();
482 }
483 } else {
484 auto myThread = [&](int iThread) {
485 for (size_t i = iThread; i < loopEnd; i += mNThreads) {
486 if (acceptTrack(tracksTPC[i])) {
487 fillDCA(tracksTPC, tracksITSTPC, vertices, i, iThread, indicesITSTPC, tracksITS, idxTPCTrackToTOFCluster, tofClusters);
488 }
489 }
490 };
491
492 std::vector<std::thread> threads(mNThreads);
493 for (int i = 0; i < mNThreads; i++) {
494 threads[i] = std::thread(myThread, i);
495 }
496
497 for (auto& th : threads) {
498 th.join();
499 }
500 }
501
502 // fill DCA values from buffer
503 for (const auto& vals : mBufferVals) {
504 for (int type = 0; type < vals.size(); ++type) {
505 const auto& val = vals[type];
506 const auto nPoints = val.side.size();
507 for (int i = 0; i < nPoints; ++i) {
508 const auto tglBin = val.tglBin[i];
509 const auto phiBin = val.phiBin[i];
510 const auto qPtBin = val.qPtBin[i];
511 const auto multBin = val.multBin[i];
512 const auto dcar = val.dcar[i];
513 const auto dcaz = val.dcaz[i];
514 const auto dcarW = val.dcarW[i];
515 const int binInt = nBins - 1;
516 const bool fillCombDCA = ((type == 1) && (val.dcarcomb[i] != -1) && (val.dcazcomb[i] != -1));
517 const bool fillDCAR = (type == 1) ? (dcar != -999) : true;
518 const std::array<int, 5> bins{tglBin, phiBin, qPtBin, multBin, binInt};
519 // fill bins
520 for (auto bin : bins) {
521 if (val.side[i] == Side::C) {
522 if (fillDCAR) {
523 mAvgCDCAr[bin][type].addValue(dcar, dcarW);
524 }
525 if (fillCombDCA) {
526 mAvgCDCAr[bin][2].addValue(val.dcarcomb[i], dcarW);
527 mAvgCDCAz[bin][2].addValue(val.dcazcomb[i], dcarW);
528 }
529 // fill only in case of valid value
530 if (dcaz != 0) {
531 mAvgCDCAz[bin][type].addValue(dcaz, dcarW);
532 }
533 } else {
534 if (fillDCAR) {
535 mAvgADCAr[bin][type].addValue(dcar, dcarW);
536 }
537 if (fillCombDCA) {
538 mAvgADCAr[bin][2].addValue(val.dcarcomb[i], dcarW);
539 mAvgADCAz[bin][2].addValue(val.dcazcomb[i], dcarW);
540 }
541 // fill only in case of valid value
542 if (dcaz != 0) {
543 mAvgADCAz[bin][type].addValue(dcaz, dcarW);
544 }
545 }
546 }
547 }
548 }
549 }
550
551 // calculate statistics and store values
552 // loop over TPC sides
553 for (int type = 0; type < 2; ++type) {
554 // loop over phi and tgl bins
555 for (int slice = 0; slice < nBins; ++slice) {
556 auto& bufferDCA = (type == 0) ? mBufferDCA.mTSTPC : mBufferDCA.mTSITSTPC;
557
558 const auto dcaAr = mAvgADCAr[slice][type].filterPointsMedian(mCutDCA, mCutRMS);
559 bufferDCA.mDCAr_A_Median[slice] = std::get<0>(dcaAr);
560 bufferDCA.mDCAr_A_WeightedMean[slice] = std::get<1>(dcaAr);
561 bufferDCA.mDCAr_A_RMS[slice] = std::get<2>(dcaAr);
562 bufferDCA.mDCAr_A_NTracks[slice] = std::get<3>(dcaAr);
563
564 const auto dcaAz = mAvgADCAz[slice][type].filterPointsMedian(mCutDCA, mCutRMS);
565 bufferDCA.mDCAz_A_Median[slice] = std::get<0>(dcaAz);
566 bufferDCA.mDCAz_A_WeightedMean[slice] = std::get<1>(dcaAz);
567 bufferDCA.mDCAz_A_RMS[slice] = std::get<2>(dcaAz);
568 bufferDCA.mDCAz_A_NTracks[slice] = std::get<3>(dcaAz);
569
570 const auto dcaCr = mAvgCDCAr[slice][type].filterPointsMedian(mCutDCA, mCutRMS);
571 bufferDCA.mDCAr_C_Median[slice] = std::get<0>(dcaCr);
572 bufferDCA.mDCAr_C_WeightedMean[slice] = std::get<1>(dcaCr);
573 bufferDCA.mDCAr_C_RMS[slice] = std::get<2>(dcaCr);
574 bufferDCA.mDCAr_C_NTracks[slice] = std::get<3>(dcaCr);
575
576 const auto dcaCz = mAvgCDCAz[slice][type].filterPointsMedian(mCutDCA, mCutRMS);
577 bufferDCA.mDCAz_C_Median[slice] = std::get<0>(dcaCz);
578 bufferDCA.mDCAz_C_WeightedMean[slice] = std::get<1>(dcaCz);
579 bufferDCA.mDCAz_C_RMS[slice] = std::get<2>(dcaCz);
580 bufferDCA.mDCAz_C_NTracks[slice] = std::get<3>(dcaCz);
581 // store combined ITS-TPC DCAs
582 if (type == 1) {
583 const auto dcaArComb = mAvgADCAr[slice][2].filterPointsMedian(mCutDCA, mCutRMS);
584 mBufferDCA.mDCAr_comb_A_Median[slice] = std::get<0>(dcaArComb);
585 mBufferDCA.mDCAr_comb_A_RMS[slice] = std::get<2>(dcaArComb);
586
587 const auto dcaAzCom = mAvgADCAz[slice][2].filterPointsMedian(mCutDCA, mCutRMS);
588 mBufferDCA.mDCAz_comb_A_Median[slice] = std::get<0>(dcaAzCom);
589 mBufferDCA.mDCAz_comb_A_RMS[slice] = std::get<2>(dcaAzCom);
590
591 const auto dcaCrComb = mAvgCDCAr[slice][2].filterPointsMedian(mCutDCA, mCutRMS);
592 mBufferDCA.mDCAr_comb_C_Median[slice] = std::get<0>(dcaCrComb);
593 mBufferDCA.mDCAr_comb_C_RMS[slice] = std::get<2>(dcaCrComb);
594
595 const auto dcaCzComb = mAvgCDCAz[slice][2].filterPointsMedian(mCutDCA, mCutRMS);
596 mBufferDCA.mDCAz_comb_C_Median[slice] = std::get<0>(dcaCzComb);
597 mBufferDCA.mDCAz_comb_C_RMS[slice] = std::get<2>(dcaCzComb);
598 }
599 }
600 }
601
602 // calculate matching eff
603 for (const auto& vals : mBufferVals) {
604 const auto& val = vals.front();
605 const auto nPoints = val.side.size();
606 for (int i = 0; i < nPoints; ++i) {
607 const auto tglBin = val.tglBin[i];
608 const auto phiBin = val.phiBin[i];
609 const auto qPtBin = val.qPtBin[i];
610 const auto multBin = val.multBin[i];
611 const auto dcar = val.dcar[i];
612 const auto dcaz = val.dcaz[i];
613 const auto hasITS = val.hasITS[i];
614 const auto chi2Match = val.chi2Match[i];
615 const auto dedxRatioqMax = val.dedxRatioqMax[i];
616 const auto dedxRatioqTot = val.dedxRatioqTot[i];
617 const auto sqrtChi2TPC = val.sqrtChi2TPC[i];
618 const auto nClTPC = val.nClTPC[i];
619 const int binInt = nBins - 1;
620 const Side side = val.side[i];
621 const bool isCSide = (side == Side::C);
622 const auto& bufferDCARMSR = isCSide ? mBufferDCA.mTSTPC.mDCAr_C_RMS : mBufferDCA.mTSTPC.mDCAr_A_RMS;
623 const auto& bufferDCARMSZ = isCSide ? mBufferDCA.mTSTPC.mDCAz_C_RMS : mBufferDCA.mTSTPC.mDCAz_A_RMS;
624 const auto& bufferDCAMedR = isCSide ? mBufferDCA.mTSTPC.mDCAr_C_Median : mBufferDCA.mTSTPC.mDCAr_A_Median;
625 const auto& bufferDCAMedZ = isCSide ? mBufferDCA.mTSTPC.mDCAz_C_Median : mBufferDCA.mTSTPC.mDCAz_A_Median;
626 auto& mAvgEff = isCSide ? mAvgMeffC : mAvgMeffA;
627 auto& mAvgChi2Match = isCSide ? mAvgChi2MatchC : mAvgChi2MatchA;
628 auto& mAvgmMIPdEdxRatioqMax = isCSide ? mMIPdEdxRatioQMaxC : mMIPdEdxRatioQMaxA;
629 auto& mAvgmMIPdEdxRatioqTot = isCSide ? mMIPdEdxRatioQTotC : mMIPdEdxRatioQTotA;
630 auto& mAvgmTPCChi2 = isCSide ? mTPCChi2C : mTPCChi2A;
631 auto& mAvgmTPCNCl = isCSide ? mTPCNClC : mTPCNClA;
632 auto& mAvgmdEdxRatioQMax = isCSide ? mLogdEdxQMaxC : mLogdEdxQMaxA;
633 auto& mAvgmdEdxRatioQTot = isCSide ? mLogdEdxQTotC : mLogdEdxQTotA;
634 auto& mITSProperties = isCSide ? mITSPropertiesC : mITSPropertiesA;
635 auto& mSigmaYZ = isCSide ? mSigmaYZC : mSigmaYZA;
636 auto& mITSTPCDeltaP = isCSide ? mITSTPCDeltaPC : mITSTPCDeltaPA;
637
638 const std::array<int, 5> bins{tglBin, phiBin, qPtBin, multBin, binInt};
639 // fill bins
640 for (auto bin : bins) {
641 // make DCA cut - select only good tracks
642 if ((std::abs(dcar - bufferDCAMedR[bin]) < (bufferDCARMSR[bin] * mCutRMS)) && (std::abs(dcaz - bufferDCAMedZ[bin]) < (bufferDCARMSZ[bin] * mCutRMS))) {
643 const auto gID = val.gID[i];
644 mAvgEff[bin][0].addValue(hasITS);
645 // count tpc only tracks not matched
646 if (!hasITS) {
647 mAvgEff[bin][1].addValue(hasITS);
648 mAvgEff[bin][2].addValue(hasITS);
649 }
650 // count tracks from ITS standalone and afterburner
652 mAvgEff[bin][1].addValue(hasITS);
654 mAvgEff[bin][2].addValue(hasITS);
655 }
656 if (chi2Match > 0) {
657 mAvgChi2Match[bin][0].addValue(chi2Match);
659 mAvgChi2Match[bin][1].addValue(chi2Match);
661 mAvgChi2Match[bin][2].addValue(chi2Match);
662 }
663 }
664 if (dedxRatioqMax > 0) {
665 mAvgmMIPdEdxRatioqMax[bin][0].addValue(dedxRatioqMax);
666 }
667 if (dedxRatioqTot > 0) {
668 mAvgmMIPdEdxRatioqTot[bin][0].addValue(dedxRatioqTot);
669 }
670 mAvgmTPCChi2[bin][0].addValue(sqrtChi2TPC);
671 mAvgmTPCNCl[bin][0].addValue(nClTPC);
672 if (hasITS) {
673 if (dedxRatioqMax > 0) {
674 mAvgmMIPdEdxRatioqMax[bin][1].addValue(dedxRatioqMax);
675 }
676 if (dedxRatioqTot > 0) {
677 mAvgmMIPdEdxRatioqTot[bin][1].addValue(dedxRatioqTot);
678 }
679 mAvgmTPCChi2[bin][1].addValue(sqrtChi2TPC);
680 mAvgmTPCNCl[bin][1].addValue(nClTPC);
681 }
682
683 float dedxNormQMax = val.dedxValsqMax[i].dedxNorm;
684 if (dedxNormQMax > 0) {
685 mAvgmdEdxRatioQMax[bin][0].addValue(dedxNormQMax);
686 }
687
688 float dedxNormQTot = val.dedxValsqTot[i].dedxNorm;
689 if (dedxNormQTot > 0) {
690 mAvgmdEdxRatioQTot[bin][0].addValue(dedxNormQTot);
691 }
692
693 float dedxIROCQMax = val.dedxValsqMax[i].dedxIROC;
694 if (dedxIROCQMax > 0) {
695 mAvgmdEdxRatioQMax[bin][1].addValue(dedxIROCQMax);
696 }
697
698 float dedxIROCQTot = val.dedxValsqTot[i].dedxIROC;
699 if (dedxIROCQTot > 0) {
700 mAvgmdEdxRatioQTot[bin][1].addValue(dedxIROCQTot);
701 }
702
703 float dedxOROC1QMax = val.dedxValsqMax[i].dedxOROC1;
704 if (dedxOROC1QMax > 0) {
705 mAvgmdEdxRatioQMax[bin][2].addValue(dedxOROC1QMax);
706 }
707
708 float dedxOROC1QTot = val.dedxValsqTot[i].dedxOROC1;
709 if (dedxOROC1QTot > 0) {
710 mAvgmdEdxRatioQTot[bin][2].addValue(dedxOROC1QTot);
711 }
712
713 float dedxOROC2QMax = val.dedxValsqMax[i].dedxOROC2;
714 if (dedxOROC2QMax > 0) {
715 mAvgmdEdxRatioQMax[bin][3].addValue(dedxOROC2QMax);
716 }
717
718 float dedxOROC2QTot = val.dedxValsqTot[i].dedxOROC2;
719 if (dedxOROC2QTot > 0) {
720 mAvgmdEdxRatioQTot[bin][3].addValue(dedxOROC2QTot);
721 }
722
723 float dedxOROC3QMax = val.dedxValsqMax[i].dedxOROC3;
724 if (dedxOROC3QMax > 0) {
725 mAvgmdEdxRatioQMax[bin][4].addValue(dedxOROC3QMax);
726 }
727
728 float dedxOROC3QTot = val.dedxValsqTot[i].dedxOROC3;
729 if (dedxOROC3QTot > 0) {
730 mAvgmdEdxRatioQTot[bin][4].addValue(dedxOROC3QTot);
731 }
732
733 float nClITS = val.nClITS[i];
734 if (nClITS > 0) {
735 mITSProperties[bin][0].addValue(nClITS);
736 }
737 float chi2ITS = val.chi2ITS[i];
738 if (chi2ITS > 0) {
739 mITSProperties[bin][1].addValue(chi2ITS);
740 }
741
742 float sigmay2 = val.sigmaY2[i];
743 if (sigmay2 > 0) {
744 mSigmaYZ[bin][0].addValue(sigmay2);
745 }
746 float sigmaz2 = val.sigmaZ2[i];
747 if (sigmaz2 > 0) {
748 mSigmaYZ[bin][1].addValue(sigmaz2);
749 }
750
751 float deltaP2 = val.deltaP2[i];
752 if (deltaP2 != -999) {
753 mITSTPCDeltaP[bin][0].addValue(deltaP2);
754 }
755
756 float deltaP3 = val.deltaP3[i];
757 if (deltaP3 != -999) {
758 mITSTPCDeltaP[bin][1].addValue(deltaP3);
759 }
760
761 float deltaP4 = val.deltaP4[i];
762 if (deltaP4 != -999) {
763 mITSTPCDeltaP[bin][2].addValue(deltaP4);
764 }
765 }
766 }
767 }
768 }
769
770 // store matching eff
771 for (int slice = 0; slice < nBins; ++slice) {
772 for (int i = 0; i < mAvgMeffA[slice].size(); ++i) {
773 auto& itsBuf = (i == 0) ? mBufferDCA.mITSTPCAll : ((i == 1) ? mBufferDCA.mITSTPCStandalone : mBufferDCA.mITSTPCAfterburner);
774 itsBuf.mITSTPC_A_MatchEff[slice] = mAvgMeffA[slice][i].getMean();
775 itsBuf.mITSTPC_C_MatchEff[slice] = mAvgMeffC[slice][i].getMean();
776 itsBuf.mITSTPC_A_Chi2Match[slice] = mAvgChi2MatchA[slice][i].getMean();
777 itsBuf.mITSTPC_C_Chi2Match[slice] = mAvgChi2MatchC[slice][i].getMean();
778 }
779
780 // loop over TPC and ITS-TPC tracks
781 for (int i = 0; i < mMIPdEdxRatioQMaxC[slice].size(); ++i) {
782 auto& buff = (i == 0) ? mBufferDCA.mTSTPC : mBufferDCA.mTSITSTPC;
783 buff.mMIPdEdxRatioQMaxA[slice] = mMIPdEdxRatioQMaxA[slice][i].getMean();
784 buff.mMIPdEdxRatioQMaxC[slice] = mMIPdEdxRatioQMaxC[slice][i].getMean();
785 buff.mMIPdEdxRatioQTotA[slice] = mMIPdEdxRatioQTotA[slice][i].getMean();
786 buff.mMIPdEdxRatioQTotC[slice] = mMIPdEdxRatioQTotC[slice][i].getMean();
787 buff.mTPCChi2C[slice] = mTPCChi2C[slice][i].getMean();
788 buff.mTPCChi2A[slice] = mTPCChi2A[slice][i].getMean();
789 buff.mTPCNClC[slice] = mTPCNClC[slice][i].getMean();
790 buff.mTPCNClA[slice] = mTPCNClA[slice][i].getMean();
791 }
792
793 // loop over qMax and qTot
794 for (int type = 0; type < 2; ++type) {
795 auto& logdEdxA = (type == 0) ? mLogdEdxQMaxA : mLogdEdxQTotA;
796 auto& buffer = (type == 0) ? mBufferDCA.mdEdxQMax : mBufferDCA.mdEdxQTot;
797 // fill A-side
798 buffer.mLogdEdx_A_Median[slice] = logdEdxA[slice][0].getMedian();
799 buffer.mLogdEdx_A_RMS[slice] = logdEdxA[slice][0].getStdDev();
800 buffer.mLogdEdx_A_IROC_Median[slice] = logdEdxA[slice][1].getMedian();
801 buffer.mLogdEdx_A_IROC_RMS[slice] = logdEdxA[slice][1].getStdDev();
802 buffer.mLogdEdx_A_OROC1_Median[slice] = logdEdxA[slice][2].getMedian();
803 buffer.mLogdEdx_A_OROC1_RMS[slice] = logdEdxA[slice][2].getStdDev();
804 buffer.mLogdEdx_A_OROC2_Median[slice] = logdEdxA[slice][3].getMedian();
805 buffer.mLogdEdx_A_OROC2_RMS[slice] = logdEdxA[slice][3].getStdDev();
806 buffer.mLogdEdx_A_OROC3_Median[slice] = logdEdxA[slice][4].getMedian();
807 buffer.mLogdEdx_A_OROC3_RMS[slice] = logdEdxA[slice][4].getStdDev();
808 // fill C-side
809 auto& logdEdxC = (type == 0) ? mLogdEdxQMaxC : mLogdEdxQTotC;
810 buffer.mLogdEdx_C_Median[slice] = logdEdxC[slice][0].getMedian();
811 buffer.mLogdEdx_C_RMS[slice] = logdEdxC[slice][0].getStdDev();
812 buffer.mLogdEdx_C_IROC_Median[slice] = logdEdxC[slice][1].getMedian();
813 buffer.mLogdEdx_C_IROC_RMS[slice] = logdEdxC[slice][1].getStdDev();
814 buffer.mLogdEdx_C_OROC1_Median[slice] = logdEdxC[slice][2].getMedian();
815 buffer.mLogdEdx_C_OROC1_RMS[slice] = logdEdxC[slice][2].getStdDev();
816 buffer.mLogdEdx_C_OROC2_Median[slice] = logdEdxC[slice][3].getMedian();
817 buffer.mLogdEdx_C_OROC2_RMS[slice] = logdEdxC[slice][3].getStdDev();
818 buffer.mLogdEdx_C_OROC3_Median[slice] = logdEdxC[slice][4].getMedian();
819 buffer.mLogdEdx_C_OROC3_RMS[slice] = logdEdxC[slice][4].getStdDev();
820 }
821
822 // ITS properties
823 // A-side
824 mBufferDCA.mITS_A_NCl_Median[slice] = mITSPropertiesA[slice][0].getMedian();
825 mBufferDCA.mITS_A_NCl_RMS[slice] = mITSPropertiesA[slice][0].getStdDev();
826 mBufferDCA.mSqrtITSChi2_Ncl_A_Median[slice] = mITSPropertiesA[slice][1].getMedian();
827 mBufferDCA.mSqrtITSChi2_Ncl_A_RMS[slice] = mITSPropertiesA[slice][1].getStdDev();
828 // C-side
829 mBufferDCA.mITS_C_NCl_Median[slice] = mITSPropertiesC[slice][0].getMedian();
830 mBufferDCA.mITS_C_NCl_RMS[slice] = mITSPropertiesC[slice][0].getStdDev();
831 mBufferDCA.mSqrtITSChi2_Ncl_C_Median[slice] = mITSPropertiesC[slice][1].getMedian();
832 mBufferDCA.mSqrtITSChi2_Ncl_C_RMS[slice] = mITSPropertiesC[slice][1].getStdDev();
833
834 //...
835 mBufferDCA.mITSTPCDeltaP2_A_Median[slice] = mITSTPCDeltaPA[slice][0].getMedian();
836 mBufferDCA.mITSTPCDeltaP3_A_Median[slice] = mITSTPCDeltaPA[slice][1].getMedian();
837 mBufferDCA.mITSTPCDeltaP4_A_Median[slice] = mITSTPCDeltaPA[slice][2].getMedian();
838 mBufferDCA.mITSTPCDeltaP2_C_Median[slice] = mITSTPCDeltaPC[slice][0].getMedian();
839 mBufferDCA.mITSTPCDeltaP3_C_Median[slice] = mITSTPCDeltaPC[slice][1].getMedian();
840 mBufferDCA.mITSTPCDeltaP4_C_Median[slice] = mITSTPCDeltaPC[slice][2].getMedian();
841 mBufferDCA.mITSTPCDeltaP2_A_RMS[slice] = mITSTPCDeltaPA[slice][0].getStdDev();
842 mBufferDCA.mITSTPCDeltaP3_A_RMS[slice] = mITSTPCDeltaPA[slice][1].getStdDev();
843 mBufferDCA.mITSTPCDeltaP4_A_RMS[slice] = mITSTPCDeltaPA[slice][2].getStdDev();
844 mBufferDCA.mITSTPCDeltaP2_C_RMS[slice] = mITSTPCDeltaPC[slice][0].getStdDev();
845 mBufferDCA.mITSTPCDeltaP3_C_RMS[slice] = mITSTPCDeltaPC[slice][1].getStdDev();
846 mBufferDCA.mITSTPCDeltaP4_C_RMS[slice] = mITSTPCDeltaPC[slice][2].getStdDev();
847 mBufferDCA.mTPCSigmaY2A_Median[slice] = mSigmaYZA[slice][0].getMedian();
848 mBufferDCA.mTPCSigmaZ2A_Median[slice] = mSigmaYZA[slice][1].getMedian();
849 mBufferDCA.mTPCSigmaY2C_Median[slice] = mSigmaYZC[slice][0].getMedian();
850 mBufferDCA.mTPCSigmaZ2C_Median[slice] = mSigmaYZC[slice][1].getMedian();
851 mBufferDCA.mTPCSigmaY2A_RMS[slice] = mSigmaYZA[slice][0].getStdDev();
852 mBufferDCA.mTPCSigmaZ2A_RMS[slice] = mSigmaYZA[slice][1].getStdDev();
853 mBufferDCA.mTPCSigmaY2C_RMS[slice] = mSigmaYZC[slice][0].getStdDev();
854 mBufferDCA.mTPCSigmaZ2C_RMS[slice] = mSigmaYZC[slice][1].getStdDev();
855 }
856
857 auto stop = timer::now();
858 std::chrono::duration<float> time = stop - startTotal;
859 LOGP(info, "Time series creation took {}", time.count());
860
861 // send data
862 sendOutput(pc);
863 }
864
866 {
867 for (auto& streamer : mStreamer) {
868 streamer->Close();
869 }
870 eos.services().get<ControlService>().readyToQuit(QuitRequest::Me);
871 }
872
873 void finaliseCCDB(o2::framework::ConcreteDataMatcher& matcher, void* obj) final
874 {
875 mTPCVDriftHelper.accountCCDBInputs(matcher, obj);
876 mPTHelper.accountCCDBInputs(matcher, obj);
878 }
879
880 private:
882 struct ValsdEdx {
883 float dedxNorm = 0;
884 float dedxIROC = 0;
885 float dedxOROC1 = 0;
886 float dedxOROC2 = 0;
887 float dedxOROC3 = 0;
888 };
889
890 struct FillVals {
891 void reserve(int n, int type)
892 {
893 side.reserve(n);
894 tglBin.reserve(n);
895 phiBin.reserve(n);
896 qPtBin.reserve(n);
897 multBin.reserve(n);
898 dcar.reserve(n);
899 dcaz.reserve(n);
900 dcarW.reserve(n);
901 dedxRatioqTot.reserve(n);
902 dedxRatioqMax.reserve(n);
903 sqrtChi2TPC.reserve(n);
904 nClTPC.reserve(n);
905 if (type == 1) {
906 hasITS.reserve(n);
907 chi2Match.reserve(n);
908 gID.reserve(n);
909 dedxValsqTot.reserve(n);
910 dedxValsqMax.reserve(n);
911 nClITS.reserve(n);
912 chi2ITS.reserve(n);
913 deltaP2.reserve(n);
914 deltaP3.reserve(n);
915 deltaP4.reserve(n);
916 sigmaY2.reserve(n);
917 sigmaZ2.reserve(n);
918 } else if (type == 0) {
919 dcarcomb.reserve(n);
920 dcazcomb.reserve(n);
921 }
922 }
923
924 void clear()
925 {
926 side.clear();
927 tglBin.clear();
928 phiBin.clear();
929 qPtBin.clear();
930 multBin.clear();
931 dcar.clear();
932 dcaz.clear();
933 dcarW.clear();
934 hasITS.clear();
935 chi2Match.clear();
936 dedxRatioqTot.clear();
937 dedxRatioqMax.clear();
938 sqrtChi2TPC.clear();
939 nClTPC.clear();
940 gID.clear();
941 dcarcomb.clear();
942 dcazcomb.clear();
943 nClITS.clear();
944 chi2ITS.clear();
945 dedxValsqTot.clear();
946 dedxValsqMax.clear();
947 deltaP2.clear();
948 deltaP3.clear();
949 deltaP4.clear();
950 sigmaY2.clear();
951 sigmaZ2.clear();
952 }
953
954 void emplace_back(Side sideTmp, int tglBinTmp, int phiBinTmp, int qPtBinTmp, int multBinTmp, float dcarTmp, float dcazTmp, float dcarWTmp, float dedxRatioqTotTmp, float dedxRatioqMaxTmp, float sqrtChi2TPCTmp, float nClTPCTmp, o2::dataformats::GlobalTrackID::Source gIDTmp, float chi2MatchTmp, int hasITSTmp, int nClITSTmp, float chi2ITSTmp, const ValsdEdx& dedxValsqTotTmp, const ValsdEdx& dedxValsqMaxTmp, float sigmaY2Tmp, float sigmaZ2Tmp)
955 {
956 side.emplace_back(sideTmp);
957 tglBin.emplace_back(tglBinTmp);
958 phiBin.emplace_back(phiBinTmp);
959 qPtBin.emplace_back(qPtBinTmp);
960 multBin.emplace_back(multBinTmp);
961 dcar.emplace_back(dcarTmp);
962 dcaz.emplace_back(dcazTmp);
963 dcarW.emplace_back(dcarWTmp);
964 dedxRatioqTot.emplace_back(dedxRatioqTotTmp);
965 dedxRatioqMax.emplace_back(dedxRatioqMaxTmp);
966 sqrtChi2TPC.emplace_back(sqrtChi2TPCTmp);
967 nClTPC.emplace_back(nClTPCTmp);
968 chi2Match.emplace_back(chi2MatchTmp);
969 hasITS.emplace_back(hasITSTmp);
970 gID.emplace_back(gIDTmp);
971 dedxValsqMax.emplace_back(dedxValsqMaxTmp);
972 dedxValsqTot.emplace_back(dedxValsqTotTmp);
973 nClITS.emplace_back(nClITSTmp);
974 chi2ITS.emplace_back(chi2ITSTmp);
975 sigmaY2.emplace_back(sigmaY2Tmp);
976 sigmaZ2.emplace_back(sigmaZ2Tmp);
977 deltaP2.emplace_back(-999);
978 deltaP3.emplace_back(-999);
979 deltaP4.emplace_back(-999);
980 }
981
982 void setDeltaParam(float deltaP2Tmp, float deltaP3Tmp, float deltaP4Tmp)
983 {
984 if (!deltaP2.empty()) {
985 deltaP2.back() = deltaP2Tmp;
986 deltaP3.back() = deltaP3Tmp;
987 deltaP4.back() = deltaP4Tmp;
988 }
989 }
990
991 void emplace_back_ITSTPC(Side sideTmp, int tglBinTmp, int phiBinTmp, int qPtBinTmp, int multBinTmp, float dcarTmp, float dcazTmp, float dcarWTmp, float dedxRatioqTotTmp, float dedxRatioqMaxTmp, float sqrtChi2TPCTmp, float nClTPCTmp, float dcarCombTmp, float dcazCombTmp)
992 {
993 side.emplace_back(sideTmp);
994 tglBin.emplace_back(tglBinTmp);
995 phiBin.emplace_back(phiBinTmp);
996 qPtBin.emplace_back(qPtBinTmp);
997 multBin.emplace_back(multBinTmp);
998 dcar.emplace_back(dcarTmp);
999 dcaz.emplace_back(dcazTmp);
1000 dcarW.emplace_back(dcarWTmp);
1001 dedxRatioqTot.emplace_back(dedxRatioqTotTmp);
1002 dedxRatioqMax.emplace_back(dedxRatioqMaxTmp);
1003 sqrtChi2TPC.emplace_back(sqrtChi2TPCTmp);
1004 nClTPC.emplace_back(nClTPCTmp);
1005 dcarcomb.emplace_back(dcarCombTmp);
1006 dcazcomb.emplace_back(dcazCombTmp);
1007 }
1008
1009 std::vector<Side> side;
1010 std::vector<int> tglBin;
1011 std::vector<int> phiBin;
1012 std::vector<int> qPtBin;
1013 std::vector<int> multBin;
1014 std::vector<float> dcar;
1015 std::vector<float> dcaz;
1016 std::vector<float> dcarW;
1017 std::vector<bool> hasITS;
1018 std::vector<float> chi2Match;
1019 std::vector<float> dedxRatioqTot;
1020 std::vector<float> dedxRatioqMax;
1021 std::vector<float> sqrtChi2TPC;
1022 std::vector<float> nClTPC;
1023 std::vector<float> dcarcomb;
1024 std::vector<float> dcazcomb;
1025 std::vector<ValsdEdx> dedxValsqTot;
1026 std::vector<ValsdEdx> dedxValsqMax;
1027 std::vector<int> nClITS;
1028 std::vector<float> chi2ITS;
1029 std::vector<o2::dataformats::GlobalTrackID::Source> gID;
1030 std::vector<float> deltaP2;
1031 std::vector<float> deltaP3;
1032 std::vector<float> deltaP4;
1033 std::vector<float> sigmaY2;
1034 std::vector<float> sigmaZ2;
1035 };
1036 std::shared_ptr<o2::base::GRPGeomRequest> mCCDBRequest;
1037 const bool mDisableWriter{false};
1039 const bool mUnbinnedWriter{false};
1040 const bool mTPCOnly{false};
1041 std::shared_ptr<o2::globaltracking::DataRequest> mDataRequest;
1042 int mPhiBins = SECTORSPERSIDE;
1043 int mTglBins{3};
1044 int mQPtBins{20};
1045 TimeSeriesITSTPC mBufferDCA;
1046 std::vector<std::array<RobustAverage, 3>> mAvgADCAr;
1047 std::vector<std::array<RobustAverage, 3>> mAvgCDCAr;
1048 std::vector<std::array<RobustAverage, 3>> mAvgADCAz;
1049 std::vector<std::array<RobustAverage, 3>> mAvgCDCAz;
1050 std::vector<std::array<RobustAverage, 2>> mMIPdEdxRatioQMaxA;
1051 std::vector<std::array<RobustAverage, 2>> mMIPdEdxRatioQMaxC;
1052 std::vector<std::array<RobustAverage, 2>> mMIPdEdxRatioQTotA;
1053 std::vector<std::array<RobustAverage, 2>> mMIPdEdxRatioQTotC;
1054 std::vector<std::array<RobustAverage, 2>> mTPCChi2A;
1055 std::vector<std::array<RobustAverage, 2>> mTPCChi2C;
1056 std::vector<std::array<RobustAverage, 2>> mTPCNClA;
1057 std::vector<std::array<RobustAverage, 2>> mTPCNClC;
1058 std::vector<std::array<RobustAverage, 3>> mAvgMeffA;
1059 std::vector<std::array<RobustAverage, 3>> mAvgMeffC;
1060 std::vector<std::array<RobustAverage, 3>> mAvgChi2MatchA;
1061 std::vector<std::array<RobustAverage, 3>> mAvgChi2MatchC;
1062 std::vector<std::array<RobustAverage, 10>> mLogdEdxQTotA;
1063 std::vector<std::array<RobustAverage, 10>> mLogdEdxQTotC;
1064 std::vector<std::array<RobustAverage, 10>> mLogdEdxQMaxA;
1065 std::vector<std::array<RobustAverage, 10>> mLogdEdxQMaxC;
1066 std::vector<std::array<RobustAverage, 2>> mITSPropertiesA;
1067 std::vector<std::array<RobustAverage, 2>> mITSPropertiesC;
1068 std::vector<std::array<RobustAverage, 3>> mITSTPCDeltaPA;
1069 std::vector<std::array<RobustAverage, 3>> mITSTPCDeltaPC;
1070 std::vector<std::array<RobustAverage, 2>> mSigmaYZA;
1071 std::vector<std::array<RobustAverage, 2>> mSigmaYZC;
1072 int mNMaxTracks{-1};
1073 float mMinMom{1};
1074 int mMinNCl{80};
1075 float mMaxTgl{1};
1076 float mMaxQPt{5};
1077 float mCoarseStep{1};
1078 float mFineStep{0.005};
1079 float mCutDCA{5};
1080 float mCutRMS{5};
1081 float mRefXSec{108.475};
1082 int mNThreads{1};
1083 float maxITSTPCDCAr{0.2};
1084 float maxITSTPCDCAz{10};
1085 float maxITSTPCDCAr_comb{0.2};
1086 float maxITSTPCDCAz_comb{0.2};
1087 gsl::span<const TPCClRefElem> mTPCTrackClIdx{};
1088 std::vector<std::array<FillVals, 2>> mBufferVals;
1089 uint32_t mFirstTFOrbit{0};
1090 float mTimeWindowMUS{50};
1091 float mMIPdEdx{50};
1092 std::vector<int> mNTracksWindow;
1093 std::vector<int> mNearestVtxTPC;
1094 o2::tpc::VDriftHelper mTPCVDriftHelper{};
1095 float mVDrift{2.64};
1096 float mMaxSnp{0.85};
1097 float mXCoarse{40};
1098 float mSqrt{13600};
1099 int mMultBins{20};
1100 int mMultMax{80000};
1101 PIDResponse mPID;
1102 int mMinTracksPerVertex{5};
1103 float mMaxdEdxRatio{0.3};
1104 float mMaxdEdxRegionRatio{0.5};
1105 float mSamplingFactor{0.1};
1106 bool mSampleTsallis{false};
1107 std::vector<std::mt19937> mGenerator;
1108 std::vector<std::unique_ptr<o2::utils::TreeStreamRedirector>> mStreamer;
1109 float mXOuterMatching{60};
1110 bool mUseMinBiasTrigger{false};
1111 long mTimeMS{};
1112 int mRun{};
1113 int mMaxOccupancyHistBins{912};
1114 PressureTemperatureHelper mPTHelper;
1115
1117 bool acceptTrack(const TrackTPC& track) const { return std::abs(track.getTgl()) < mMaxTgl; }
1118
1119 bool checkTrack(const TrackTPC& track) const
1120 {
1121 const bool isGoodTrack = ((track.getNClusters() < mMinNCl) || (track.getP() < mMinMom)) ? false : true;
1122 return isGoodTrack;
1123 }
1124
1125 void fillDCA(const gsl::span<const TrackTPC> tracksTPC, const gsl::span<const o2::dataformats::TrackTPCITS> tracksITSTPC, const gsl::span<const o2::dataformats::PrimaryVertex> vertices, const int iTrk, const int iThread, const std::unordered_map<unsigned int, std::array<int, 2>>& indicesITSTPC, const gsl::span<const o2::its::TrackITS> tracksITS, const std::vector<std::tuple<int, float, float, o2::track::TrackLTIntegral, double, float, unsigned int>>& idxTPCTrackToTOFCluster, const gsl::span<const o2::tof::Cluster> tofClusters)
1126 {
1127 const auto& trackFull = tracksTPC[iTrk];
1128 const bool isGoodTrack = checkTrack(trackFull);
1129
1130 // check for min bias trigger - sample flat -
1131 bool minBiasOk = false;
1132 const float factorMinBias = 0.1 * mSamplingFactor;
1133 if (mUnbinnedWriter && mUseMinBiasTrigger) {
1134 std::uniform_real_distribution<> distr(0., 1.);
1135 if (distr(mGenerator[iThread]) < factorMinBias) {
1136 minBiasOk = true;
1137 }
1138 }
1139
1140 // check if at least one check passed
1141 if (!isGoodTrack && !minBiasOk) {
1142 return;
1143 }
1144
1145 o2::track::TrackParCov track = tracksTPC[iTrk];
1146
1147 // propagate track to the DCA and fill in slice
1148 auto propagator = o2::base::Propagator::Instance();
1149
1150 // propagate track to DCA
1151 std::array<float, 2> dca;
1152 const o2::math_utils::Point3D<float> refPoint{0, 0, 0};
1153
1154 // coarse propagation
1155 if (!propagator->PropagateToXBxByBz(track, mXCoarse, mMaxSnp, mCoarseStep, mMatType)) {
1156 return;
1157 }
1158
1159 // fine propagation with Bz only
1160 if (!propagator->propagateToDCA(refPoint, track, propagator->getNominalBz(), mFineStep, mMatType, &dca)) {
1161 return;
1162 }
1163
1164 o2::track::TrackPar trackTmp(tracksTPC[iTrk]);
1165
1166 // coarse propagation to centre of IROC for phi bin
1167 if (!propagator->propagateTo(trackTmp, mRefXSec, false, mMaxSnp, mCoarseStep, mMatType)) {
1168 return;
1169 }
1170
1171 const int tglBin = mTglBins * std::abs(trackTmp.getTgl()) / mMaxTgl + mPhiBins;
1172 const int phiBin = mPhiBins * trackTmp.getPhi() / o2::constants::math::TwoPI;
1173
1174 const int offsQPtBin = mPhiBins + mTglBins;
1175 const int qPtBin = offsQPtBin + mQPtBins * (trackTmp.getQ2Pt() + mMaxQPt) / (2 * mMaxQPt);
1176 const int localMult = mNTracksWindow[iTrk];
1177
1178 const int offsMult = offsQPtBin + mQPtBins;
1179 const int multBin = offsMult + mMultBins * localMult / mMultMax;
1180 const int nBins = getNBins();
1181
1182 if ((phiBin < 0) || (phiBin > mPhiBins) || (tglBin < mPhiBins) || (tglBin > offsQPtBin) || (qPtBin < offsQPtBin) || (qPtBin > offsMult) || (multBin < offsMult) || (multBin > offsMult + mMultBins)) {
1183 return;
1184 }
1185
1186 float sigmaY2 = 0;
1187 float sigmaZ2 = 0;
1188 const int sector = o2::math_utils::angle2Sector(trackTmp.getPhiPos());
1189 // find possible ITS-TPC track and vertex index
1190 auto it = indicesITSTPC.find(iTrk);
1191 const auto idxITSTPC = (it != indicesITSTPC.end()) ? (it->second) : std::array<int, 2>{-1, -1};
1192
1193 // get vertex (check if vertex ID is valid). In case no vertex is assigned return nearest vertex or else default vertex
1194 const auto vertex = (idxITSTPC.back() != -1) ? vertices[idxITSTPC.back()] : ((mNearestVtxTPC[iTrk] != -1) ? vertices[mNearestVtxTPC[iTrk]] : o2::dataformats::PrimaryVertex{});
1195
1196 // calculate DCAz: (time TPC track - time vertex) * vDrift + sign_side * vertexZ
1197 const float signSide = trackFull.hasCSideClustersOnly() ? -1 : 1; // invert sign for C-side
1198 const float dcaZFromDeltaTime = (vertex.getTimeStamp().getTimeStamp() == 0) ? 0 : (o2::tpc::ParameterElectronics::Instance().ZbinWidth * trackFull.getTime0() - vertex.getTimeStamp().getTimeStamp()) * mVDrift + signSide * vertex.getZ();
1199
1200 // for weight of DCA
1201 const float resCl = std::min(trackFull.getNClusters(), static_cast<int>(Mapper::PADROWS)) / static_cast<float>(Mapper::PADROWS);
1202
1203 const float div = (resCl * track.getPt());
1204 if (div == 0) {
1205 return;
1206 }
1207
1208 const float fB = 0.2 / div;
1209 const float fA = 0.15 + 0.15; // = 0.15 with additional misalignment error
1210 const float dcarW = 1. / std::sqrt(fA * fA + fB * fB); // Weight of DCA: Rms2 ~ 0.15^2 + k/(L^2*pt) → 0.15**2 + (0.2/((NCl/152)*pt)^2);
1211
1212 // store values for TPC DCA only for A- or C-side only tracks
1213 const bool hasITSTPC = idxITSTPC.front() != -1;
1214
1215 // get ratio of chi2 in case ITS-TPC track has been found
1216 const float chi2 = hasITSTPC ? tracksITSTPC[idxITSTPC.front()].getChi2Match() : -1;
1218 // check for source in case of ITS-TPC
1219 if (hasITSTPC) {
1220 const auto src = tracksITSTPC[idxITSTPC.front()].getRefITS().getSource();
1225 }
1226 }
1227
1228 const float chi2Match = (chi2 > 0) ? std::sqrt(chi2) : -1;
1229 const float sqrtChi2TPC = (trackFull.getChi2() > 0) ? std::sqrt(trackFull.getChi2()) : 0;
1230 const float nClTPC = trackFull.getNClusters();
1231
1232 // const float dedx = mUseQMax ? track.getdEdx().dEdxMaxTPC : track.getdEdx().dEdxTotTPC;
1233 const float dedxRatioqTot = (trackFull.getdEdx().dEdxTotTPC > 0) ? (mMIPdEdx / trackFull.getdEdx().dEdxTotTPC) : -1;
1234 const float dedxRatioqMax = (trackFull.getdEdx().dEdxMaxTPC > 0) ? (mMIPdEdx / trackFull.getdEdx().dEdxMaxTPC) : -1;
1235
1236 const auto dedxQTotVars = getdEdxVars(0, trackFull);
1237 const auto dedxQMaxVars = getdEdxVars(1, trackFull);
1238
1239 // make check to avoid crash in case no or less ITS tracks have been found!
1240 const int idxITSTrack = (hasITSTPC && (gID == o2::dataformats::GlobalTrackID::Source::ITS)) ? tracksITSTPC[idxITSTPC.front()].getRefITS().getIndex() : -1;
1241 const bool idxITSCheck = (idxITSTrack != -1);
1242
1243 const int nClITS = idxITSCheck ? tracksITS[idxITSTrack].getNClusters() : -1;
1244 float chi2ITS = idxITSCheck ? tracksITS[idxITSTrack].getChi2() : -1;
1245 if ((chi2ITS > 0) && (nClITS > 0)) {
1246 chi2ITS = std::sqrt(chi2ITS / nClITS);
1247 }
1248 sigmaY2 = track.getSigmaY2();
1249 sigmaZ2 = track.getSigmaZ2();
1250 if (isGoodTrack) {
1251 if (trackFull.hasCSideClustersOnly()) {
1252 mBufferVals[iThread].front().emplace_back(Side::C, tglBin, phiBin, qPtBin, multBin, dca[0], dcaZFromDeltaTime, dcarW, dedxRatioqTot, dedxRatioqMax, sqrtChi2TPC, nClTPC, gID, chi2Match, hasITSTPC, nClITS, chi2ITS, dedxQTotVars, dedxQMaxVars, sigmaY2, sigmaZ2);
1253 } else if (trackFull.hasASideClustersOnly()) {
1254 mBufferVals[iThread].front().emplace_back(Side::A, tglBin, phiBin, qPtBin, multBin, dca[0], dcaZFromDeltaTime, dcarW, dedxRatioqTot, dedxRatioqMax, sqrtChi2TPC, nClTPC, gID, chi2Match, hasITSTPC, nClITS, chi2ITS, dedxQTotVars, dedxQMaxVars, sigmaY2, sigmaZ2);
1255 }
1256 }
1257
1258 // make propagation for ITS-TPC Track
1259 // check if the track was assigned to ITS track
1260 std::array<float, 2> dcaITSTPC{0, 0};
1261 float deltaP0 = -999;
1262 float deltaP1 = -999;
1263 float deltaP2 = -999;
1264 float deltaP3 = -999;
1265 float deltaP4 = -999;
1266 float phiITSTPCAtVertex = -999; // phi of ITS-TPC track at vertex
1267 float dcaTPCAtVertex = -999;
1268 if (hasITSTPC) {
1269 // propagate ITS-TPC track to (0,0)
1270 auto trackITSTPCTmp = tracksITSTPC[idxITSTPC.front()];
1271 // fine propagation with Bz only
1272 if (propagator->propagateToDCA(refPoint, trackITSTPCTmp, propagator->getNominalBz(), mFineStep, mMatType, &dcaITSTPC)) {
1273 // make cut on abs(DCA)
1274 if ((std::abs(dcaITSTPC[0]) < maxITSTPCDCAr) && (std::abs(dcaITSTPC[1]) < maxITSTPCDCAz)) {
1275 // store TPC only DCAs
1276 // propagate to vertex in case the track belongs to vertex
1277 const bool contributeToVertex = (idxITSTPC.back() != -1);
1278 std::array<float, 2> dcaITSTPCTmp{-1, -1};
1279
1280 if (contributeToVertex) {
1281 if (propagator->propagateToDCA(vertex.getXYZ(), trackITSTPCTmp, propagator->getNominalBz(), mFineStep, mMatType, &dcaITSTPCTmp)) {
1282 phiITSTPCAtVertex = trackITSTPCTmp.getPhi();
1283 dcaITSTPC = dcaITSTPCTmp;
1284 }
1285
1286 // propagate TPC track to vertex
1287 std::array<float, 2> dcaTPCTmp{-1, -1};
1288 if (propagator->propagateToDCA(vertex.getXYZ(), track, propagator->getNominalBz(), mFineStep, mMatType, &dcaTPCTmp)) {
1289 dcaTPCAtVertex = dcaTPCTmp[0];
1290 }
1291 }
1292
1293 // make cut around DCA to vertex due to gammas
1294 if ((std::abs(dcaITSTPCTmp[0]) < maxITSTPCDCAr_comb) && (std::abs(dcaITSTPCTmp[1]) < maxITSTPCDCAz_comb)) {
1295 // propagate TPC track to ITS track and store delta track parameters
1296 if (idxITSTrack >= 0 && track.rotate(tracksITS[idxITSTrack].getAlpha()) && propagator->propagateTo(track, trackITSTPCTmp.getX(), false, mMaxSnp, mFineStep, mMatType)) {
1297 o2::track::TrackPar trackITS(tracksITS[idxITSTrack]);
1298 const bool propITSOk = propagator->propagateTo(trackITS, trackITSTPCTmp.getX(), false, mMaxSnp, mFineStep, mMatType);
1299 if (propITSOk) {
1300 deltaP0 = track.getParam(0) - trackITS.getParam(0);
1301 deltaP1 = track.getParam(1) - trackITS.getParam(1);
1302 deltaP2 = track.getParam(2) - trackITS.getParam(2);
1303 deltaP3 = track.getParam(3) - trackITS.getParam(3);
1304 deltaP4 = track.getParam(4) - trackITS.getParam(4);
1305 mBufferVals[iThread].front().setDeltaParam(deltaP2, deltaP3, deltaP4);
1306 }
1307 }
1308 } else {
1309 dcaITSTPCTmp[0] = -1;
1310 dcaITSTPCTmp[1] = -1;
1311 }
1312
1313 if (isGoodTrack) {
1314 if (trackFull.hasCSideClustersOnly()) {
1315 mBufferVals[iThread].back().emplace_back_ITSTPC(Side::C, tglBin, phiBin, qPtBin, multBin, dcaTPCAtVertex, dcaZFromDeltaTime, dcarW, dedxRatioqTot, dedxRatioqMax, sqrtChi2TPC, nClTPC, dcaITSTPCTmp[0], dcaITSTPCTmp[1]);
1316 } else if (trackFull.hasASideClustersOnly()) {
1317 mBufferVals[iThread].back().emplace_back_ITSTPC(Side::A, tglBin, phiBin, qPtBin, multBin, dcaTPCAtVertex, dcaZFromDeltaTime, dcarW, dedxRatioqTot, dedxRatioqMax, sqrtChi2TPC, nClTPC, dcaITSTPCTmp[0], dcaITSTPCTmp[1]);
1318 }
1319 }
1320 }
1321 }
1322 }
1323
1324 if (mUnbinnedWriter && mStreamer[iThread]) {
1325 const float factorPt = mSamplingFactor;
1326 bool writeData = true;
1327 float weight = 0;
1328 if (mSampleTsallis) {
1329 std::uniform_real_distribution<> distr(0., 1.);
1330 writeData = o2::math_utils::Tsallis::downsampleTsallisCharged(tracksTPC[iTrk].getPt(), factorPt, mSqrt, weight, distr(mGenerator[iThread]));
1331 }
1332 if (writeData || minBiasOk) {
1333 auto clusterMask = makeClusterBitMask(trackFull);
1334 const auto& trkOrig = tracksTPC[iTrk];
1335 const bool isNearestVtx = (idxITSTPC.back() == -1); // is nearest vertex in case no vertex was found
1336 const float mx_ITS = hasITSTPC ? tracksITSTPC[idxITSTPC.front()].getX() : -1;
1337 const float pt_ITS = hasITSTPC ? tracksITSTPC[idxITSTPC.front()].getQ2Pt() : -1;
1338 const float chi2match_ITSTPC = hasITSTPC ? tracksITSTPC[idxITSTPC.front()].getChi2Match() : -1;
1339 const int nClITS = idxITSCheck ? tracksITS[idxITSTrack].getNClusters() : -1;
1340 const int chi2ITS = idxITSCheck ? tracksITS[idxITSTrack].getChi2() : -1;
1341 int typeSide = 2; // A- and C-Side cluster
1342 if (trackFull.hasASideClustersOnly()) {
1343 typeSide = 0;
1344 } else if (trackFull.hasCSideClustersOnly()) {
1345 typeSide = 1;
1346 }
1347
1348 // check for TOF and propagate TPC track to TOF cluster
1349 bool hasTOFCluster = (std::get<0>(idxTPCTrackToTOFCluster[iTrk]) != -1);
1350 auto tofCl = hasTOFCluster ? tofClusters[std::get<0>(idxTPCTrackToTOFCluster[iTrk])] : o2::tof::Cluster();
1351
1352 float tpcYDeltaAtTOF = -999;
1353 float tpcZDeltaAtTOF = -999;
1354 if (hasTOFCluster) {
1355 o2::track::TrackPar trackTmpOut(tracksTPC[iTrk].getParamOut());
1356 if (trackTmpOut.rotate(o2::math_utils::sector2Angle(tofCl.getSector())) && propagator->propagateTo(trackTmpOut, tofCl.getX(), false, mMaxSnp, mFineStep, mMatType)) {
1357 tpcYDeltaAtTOF = trackTmpOut.getY() - tofCl.getY();
1358 tpcZDeltaAtTOF = signSide * (o2::tpc::ParameterElectronics::Instance().ZbinWidth * trackFull.getTime0() - vertex.getTimeStamp().getTimeStamp()) * mVDrift - trackTmpOut.getZ() + tofCl.getZ();
1359 }
1360 }
1361
1362 // get delta parameter between inner and outer
1363 float deltaTPCParamInOutTgl = trackFull.getTgl() - trackFull.getParamOut().getTgl();
1364 float deltaTPCParamInOutQPt = trackFull.getQ2Pt() - trackFull.getParamOut().getQ2Pt();
1365
1366 // propagate TPC track and ITS-TPC track to outer matching at 60cm
1367 float deltaP0OuterITS = -999;
1368 float deltaP1OuterITS = -999;
1369 float deltaP2OuterITS = -999;
1370 float deltaP3OuterITS = -999;
1371 float deltaP4OuterITS = -999;
1372 if (idxITSCheck) {
1373 o2::track::TrackPar trackTmpOut(tracksITS[idxITSTrack].getParamOut());
1374 const bool propITSOk = propagator->propagateTo(trackTmpOut, mXOuterMatching, false, mMaxSnp, mCoarseStep, mMatType);
1375 if (propITSOk && trackTmp.rotate(trackTmpOut.getAlpha())) {
1376 const bool propTPCOk = propagator->propagateTo(trackTmp, mXOuterMatching, false, mMaxSnp, mCoarseStep, mMatType);
1377 if (propTPCOk) {
1378 // store delta parameters
1379 deltaP0OuterITS = trackTmp.getParam(0) - trackTmpOut.getParam(0);
1380 deltaP1OuterITS = trackTmp.getParam(1) - trackTmpOut.getParam(2);
1381 deltaP2OuterITS = trackTmp.getParam(2) - trackTmpOut.getParam(2);
1382 deltaP3OuterITS = trackTmp.getParam(3) - trackTmpOut.getParam(3);
1383 deltaP4OuterITS = trackTmp.getParam(4) - trackTmpOut.getParam(4);
1384 }
1385 }
1386 }
1387 const int triggerMask = 0x1 * minBiasOk + 0x2 * writeData;
1388
1389 float deltaP2ConstrVtx = -999;
1390 float deltaP3ConstrVtx = -999;
1391 float deltaP4ConstrVtx = -999;
1392
1393 // cov of TPC track constrained at vertex
1394 float covTPCConstrVtxP2 = -999;
1395 float covTPCConstrVtxP3 = -999;
1396 float covTPCConstrVtxP4 = -999;
1397
1398 // cov of ITS-TPC track at vertex
1399 float covITSTPCConstrVtxP2 = -999;
1400 float covITSTPCConstrVtxP3 = -999;
1401 float covITSTPCConstrVtxP4 = -999;
1402
1403 float covTPCAtVertex0 = -999;
1404 float covTPCAtVertex1 = -999;
1405
1406 const bool contributeToVertex = (idxITSTPC.back() != -1);
1407 if (hasITSTPC && contributeToVertex) {
1408 o2::track::TrackParCov trackITSTPCTmp = tracksITSTPC[idxITSTPC.front()];
1409 std::array<float, 2> dcaITSTPCTmp{-1, -1};
1410 if (propagator->propagateToDCA(vertex.getXYZ(), trackITSTPCTmp, propagator->getNominalBz(), mFineStep, mMatType, &dcaITSTPCTmp)) {
1411 o2::track::TrackParCov trackTPC = tracksTPC[iTrk];
1412 if (trackTPC.rotate(trackITSTPCTmp.getAlpha()) && propagator->propagateTo(trackTPC, trackITSTPCTmp.getX(), false, mMaxSnp, mFineStep, mMatType)) {
1413 // store covariance of TPC track at vertex
1414 covTPCAtVertex0 = trackTPC.getCovarElem(0, 0);
1415 covTPCAtVertex1 = trackTPC.getCovarElem(1, 1);
1416
1417 trackTPC.update(vertex);
1418 deltaP2ConstrVtx = trackTPC.getParam(2) - trackITSTPCTmp.getParam(2);
1419 deltaP3ConstrVtx = trackTPC.getParam(3) - trackITSTPCTmp.getParam(3);
1420 deltaP4ConstrVtx = trackTPC.getParam(4) - trackITSTPCTmp.getParam(4);
1421 covTPCConstrVtxP2 = trackTPC.getCovarElem(2, 2);
1422 covTPCConstrVtxP3 = trackTPC.getCovarElem(3, 3);
1423 covTPCConstrVtxP4 = trackTPC.getCovarElem(4, 4);
1424 covITSTPCConstrVtxP2 = trackITSTPCTmp.getCovarElem(2, 2);
1425 covITSTPCConstrVtxP3 = trackITSTPCTmp.getCovarElem(3, 3);
1426 covITSTPCConstrVtxP4 = trackITSTPCTmp.getCovarElem(4, 4);
1427 }
1428 }
1429 }
1430 double vertexTime = vertex.getTimeStamp().getTimeStamp();
1431 double trackTime0 = trackFull.getTime0();
1432 *mStreamer[iThread] << "treeTimeSeries"
1433 // DCAs
1434 << "triggerMask=" << triggerMask
1435 << "factorMinBias=" << factorMinBias
1436 << "factorPt=" << factorPt
1437 << "weight=" << weight
1438 << "dcar_tpc_vertex=" << dcaTPCAtVertex
1439 << "dcar_tpc=" << dca[0]
1440 << "dcaz_tpc=" << dca[1]
1441 << "dcar_itstpc=" << dcaITSTPC[0]
1442 << "dcaz_itstpc=" << dcaITSTPC[1]
1443 << "dcarW=" << dcarW
1444 << "dcaZFromDeltaTime=" << dcaZFromDeltaTime
1445 << "hasITSTPC=" << hasITSTPC
1446 // vertex
1447 << "vertex_x=" << vertex.getX()
1448 << "vertex_y=" << vertex.getY()
1449 << "vertex_z=" << vertex.getZ()
1450 << "vertex_time=" << vertex.getTimeStamp().getTimeStamp()
1451 << "vertex_nContributors=" << vertex.getNContributors()
1452 << "isNearestVertex=" << isNearestVtx
1453 // tpc track properties
1454 << "pt=" << trkOrig.getPt()
1455 << "qpt_ITSTPC=" << pt_ITS
1456 << "tpc_timebin=" << trkOrig.getTime0()
1457 << "qpt=" << trkOrig.getParam(4)
1458 << "ncl=" << trkOrig.getNClusters()
1459 << "ncl_shared=" << trkOrig.getNClusters()
1460 << "tgl=" << trkOrig.getTgl()
1461 << "side_type=" << typeSide
1462 << "phi=" << trkOrig.getPhi()
1463 << "clusterMask=" << clusterMask
1464 << "dedxTPC=" << trkOrig.getdEdx()
1465 << "chi2=" << trkOrig.getChi2()
1466 << "mX=" << trkOrig.getX()
1467 << "mX_ITS=" << mx_ITS
1468 << "nClITS=" << nClITS
1469 << "chi2ITS=" << chi2ITS
1470 << "chi2match_ITSTPC=" << chi2match_ITSTPC
1471 << "PID=" << trkOrig.getPID().getID()
1472 // TPC cov at vertex (without vertex constrained)
1473 << "covTPCAtVertex0=" << covTPCAtVertex0
1474 << "covTPCAtVertex1=" << covTPCAtVertex1
1475 // TPC cov at vertex (with vertex constrained)
1476 << "covTPCConstrVtxP2=" << covTPCConstrVtxP2
1477 << "covTPCConstrVtxP3=" << covTPCConstrVtxP3
1478 << "covTPCConstrVtxP4=" << covTPCConstrVtxP4
1479 // ITS-TPC cov at vertex (with vertex constrained)
1480 << "covITSTPCConstrVtxP2=" << covITSTPCConstrVtxP2
1481 << "covITSTPCConstrVtxP3=" << covITSTPCConstrVtxP3
1482 << "covITSTPCConstrVtxP4=" << covITSTPCConstrVtxP4
1483 // delta Parameter at vertex with TPC track constrained at vertex
1484 << "deltaP2ConstrVtx=" << deltaP2ConstrVtx
1485 << "deltaP3ConstrVtx=" << deltaP3ConstrVtx
1486 << "deltaP4ConstrVtx=" << deltaP4ConstrVtx
1487 //
1488 << "deltaPar0=" << deltaP0
1489 << "deltaPar1=" << deltaP1
1490 << "deltaPar2=" << deltaP2
1491 << "deltaPar3=" << deltaP3
1492 << "deltaPar4=" << deltaP4
1493 << "sigmaY2=" << sigmaY2
1494 << "sigmaZ2=" << sigmaZ2
1495 // meta
1496 << "mult=" << mNTracksWindow[iTrk]
1497 << "time_window_mult=" << mTimeWindowMUS
1498 << "firstTFOrbit=" << mFirstTFOrbit
1499 << "timeMS=" << mTimeMS
1500 << "run=" << mRun
1501 << "mVDrift=" << mVDrift
1502 << "its_flag=" << int(gID)
1503 << "sqrtChi2Match=" << chi2Match
1504 // TOF cluster
1505 << "tpcYDeltaAtTOF=" << tpcYDeltaAtTOF
1506 << "tpcZDeltaAtTOF=" << tpcZDeltaAtTOF
1507 << "mDXatTOF=" << std::get<1>(idxTPCTrackToTOFCluster[iTrk])
1508 << "mDZatTOF=" << std::get<2>(idxTPCTrackToTOFCluster[iTrk])
1509 << "mTOFLength=" << std::get<3>(idxTPCTrackToTOFCluster[iTrk])
1510 << "mTOFSignal=" << std::get<4>(idxTPCTrackToTOFCluster[iTrk])
1511 << "mDeltaTTOFTPC=" << std::get<5>(idxTPCTrackToTOFCluster[iTrk])
1512 << "vertexTime=" << vertexTime
1513 << "trackTime0=" << trackTime0
1514 << "TOFmask=" << std::get<6>(idxTPCTrackToTOFCluster[iTrk])
1515 // TPC delta param
1516 << "deltaTPCParamInOutTgl=" << deltaTPCParamInOutTgl
1517 << "deltaTPCParamInOutQPt=" << deltaTPCParamInOutQPt
1518 // delta parameter between ITS-TPC - TPC at 60 cm
1519 << "deltaP0OuterITS=" << deltaP0OuterITS
1520 << "deltaP1OuterITS=" << deltaP1OuterITS
1521 << "deltaP2OuterITS=" << deltaP2OuterITS
1522 << "deltaP3OuterITS=" << deltaP3OuterITS
1523 << "deltaP4OuterITS=" << deltaP4OuterITS
1524 << "mXOuterMatching=" << mXOuterMatching
1525 // phi of ITS-TPC track at vertex
1526 << "phiITSTPCAtVertex=" << phiITSTPCAtVertex
1527 << "\n";
1528 }
1529 }
1530 }
1531
1532 void sendOutput(ProcessingContext& pc)
1533 {
1534 mBufferDCA.mTSTPC.setStartTime(mTimeMS);
1535 mBufferDCA.mTSITSTPC.setStartTime(mTimeMS);
1536 pc.outputs().snapshot(Output{header::gDataOriginTPC, getDataDescriptionTimeSeries()}, mBufferDCA);
1537 // in case of ROOT output also store the TFinfo in the TTree
1538 if (!mDisableWriter) {
1541 pc.outputs().snapshot(Output{header::gDataOriginTPC, getDataDescriptionTPCTimeSeriesTFId()}, tfinfo);
1542 }
1543 }
1544
1546 void findNearesVertex(const gsl::span<const TrackTPC> tracksTPC, const gsl::span<const o2::dataformats::PrimaryVertex> vertices)
1547 {
1548 // create list of time bins of tracks
1549 const int nVertices = vertices.size();
1550
1551 const int nTracks = tracksTPC.size();
1552 mNearestVtxTPC.clear();
1553 mNearestVtxTPC.resize(nTracks);
1554
1555 // in case no vertices are found
1556 if (!nVertices) {
1557 std::fill(mNearestVtxTPC.begin(), mNearestVtxTPC.end(), -1);
1558 return;
1559 }
1560
1561 // store timestamps of vertices. Assume vertices are already sorted in time!
1562 std::vector<float> times_vtx;
1563 times_vtx.reserve(nVertices);
1564 for (const auto& vtx : vertices) {
1565 times_vtx.emplace_back(vtx.getTimeStamp().getTimeStamp());
1566 }
1567
1568 // loop over tpc tracks and find nearest vertex
1569 auto myThread = [&](int iThread) {
1570 for (int i = iThread; i < nTracks; i += mNThreads) {
1571 const float timeTrack = o2::tpc::ParameterElectronics::Instance().ZbinWidth * tracksTPC[i].getTime0();
1572 const auto lower = std::lower_bound(times_vtx.begin(), times_vtx.end(), timeTrack);
1573 int closestVtx = std::distance(times_vtx.begin(), lower);
1574 // if value is out of bounds use last value
1575 if (closestVtx == nVertices) {
1576 closestVtx -= 1;
1577 } else if (closestVtx > 0) {
1578 // if idx > 0 check preceeding value
1579 double diff1 = std::abs(timeTrack - *lower);
1580 double diff2 = std::abs(timeTrack - *(lower - 1));
1581 if (diff2 < diff1) {
1582 closestVtx -= 1;
1583 }
1584 }
1585 mNearestVtxTPC[i] = closestVtx;
1586 }
1587 };
1588
1589 std::vector<std::thread> threads(mNThreads);
1590 for (int i = 0; i < mNThreads; i++) {
1591 threads[i] = std::thread(myThread, i);
1592 }
1593
1594 // wait for the threads to finish
1595 for (auto& th : threads) {
1596 th.join();
1597 }
1598 }
1599
1601 std::vector<bool> makeClusterBitMask(const TrackTPC& track) const
1602 {
1603 std::vector<bool> tpcClusterMask(Mapper::PADROWS, false);
1604 const int nCl = track.getNClusterReferences();
1605 for (int j = 0; j < nCl; ++j) {
1606 uint8_t sector, padrow;
1607 uint32_t clusterIndexInRow;
1608 track.getClusterReference(mTPCTrackClIdx, j, sector, padrow, clusterIndexInRow);
1609 tpcClusterMask[padrow] = true;
1610 }
1611 return tpcClusterMask;
1612 }
1613
1615 void findNNeighbourTracks(const gsl::span<const TrackTPC> tracksTPC)
1616 {
1617 const float tpcTBinMUS = o2::tpc::ParameterElectronics::Instance().ZbinWidth; // 0.199606f; time bin in MUS
1618 const float windowTimeBins = mTimeWindowMUS / tpcTBinMUS; // number of neighbouring time bins to check
1619
1620 // create list of time bins of tracks
1621 std::vector<float> times;
1622 const int nTracks = tracksTPC.size();
1623 times.reserve(nTracks);
1624 for (const auto& trk : tracksTPC) {
1625 times.emplace_back(trk.getTime0());
1626 }
1627 std::sort(times.begin(), times.end());
1628
1629 mNTracksWindow.clear();
1630 mNTracksWindow.resize(nTracks);
1631
1632 // loop over tpc tracks and count number of neighouring tracks
1633 auto myThread = [&](int iThread) {
1634 for (int i = iThread; i < nTracks; i += mNThreads) {
1635 const float t0 = tracksTPC[i].getTime0();
1636 const auto upperV0 = std::upper_bound(times.begin(), times.end(), t0 + windowTimeBins);
1637 const auto lowerV0 = std::lower_bound(times.begin(), times.end(), t0 - windowTimeBins);
1638 const int nMult = std::distance(times.begin(), upperV0) - std::distance(times.begin(), lowerV0);
1639 mNTracksWindow[i] = nMult;
1640 }
1641 };
1642
1643 std::vector<std::thread> threads(mNThreads);
1644 for (int i = 0; i < mNThreads; i++) {
1645 threads[i] = std::thread(myThread, i);
1646 }
1647
1648 // wait for the threads to finish
1649 for (auto& th : threads) {
1650 th.join();
1651 }
1652 }
1653
1654 std::unordered_map<unsigned int, int> processVertices(const gsl::span<const o2::dataformats::PrimaryVertex> vertices, const gsl::span<const o2::dataformats::VtxTrackIndex> primMatchedTracks, const gsl::span<const o2::dataformats::VtxTrackRef> primMatchedTracksRef, const RecoContainer& recoData)
1655 {
1656 // storing collision vertex to ITS-TPC track index
1657 std::unordered_map<unsigned int, int> indicesITSTPC_vtx; // ITS-TPC track index -> collision vertex ID
1658
1659 std::unordered_map<int, int> nContributors_ITS; // ITS: vertex ID -> n contributors
1660 std::unordered_map<int, int> nContributors_ITSTPC; // ITS-TPC (and ITS-TPC-TRD, ITS-TPC-TOF, ITS-TPC-TRD-TOF): vertex ID -> n contributors
1661
1662 // loop over collisions
1663 if (!vertices.empty()) {
1664 for (const auto& ref : primMatchedTracksRef) {
1665 // loop over ITS and ITS-TPC sources
1666 const std::array<TrkSrc, 5> sources = {TrkSrc::ITSTPC, TrkSrc::ITSTPCTRD, TrkSrc::ITSTPCTOF, TrkSrc::ITSTPCTRDTOF, TrkSrc::ITS};
1667 for (auto source : sources) {
1668 const int vID = ref.getVtxID(); // vertex ID
1669 const int firstEntry = ref.getFirstEntryOfSource(source);
1670 const int nEntries = ref.getEntriesOfSource(source);
1671 // loop over all tracks belonging to the vertex
1672 for (int i = 0; i < nEntries; ++i) {
1673 const auto& matchedTrk = primMatchedTracks[i + firstEntry];
1674 bool pvCont = matchedTrk.isPVContributor();
1675 if (pvCont) {
1676 // store index of ITS-TPC track container and vertex ID
1677 auto refITSTPC = recoData.getSingleDetectorRefs(matchedTrk)[TrkSrc::ITSTPC];
1678 if (refITSTPC.isIndexSet()) {
1679 indicesITSTPC_vtx[refITSTPC] = vID;
1680 ++nContributors_ITSTPC[vID];
1681 } else {
1682 ++nContributors_ITS[vID];
1683 }
1684 }
1685 }
1686 }
1687 }
1688 }
1689
1690 // calculate statistics
1691 std::array<RobustAverage, 4> avgVtxITS;
1692 std::array<RobustAverage, 4> avgVtxITSTPC;
1693 for (int i = 0; i < avgVtxITS.size(); ++i) {
1694 avgVtxITS[i].reserve(vertices.size());
1695 avgVtxITSTPC[i].reserve(vertices.size());
1696 }
1697 for (int ivtx = 0; ivtx < vertices.size(); ++ivtx) {
1698 const auto& vtx = vertices[ivtx];
1699 const float itsFrac = nContributors_ITS[ivtx] / static_cast<float>(vtx.getNContributors());
1700 const float itstpcFrac = (nContributors_ITS[ivtx] + nContributors_ITSTPC[ivtx]) / static_cast<float>(vtx.getNContributors());
1701
1702 const float itsMin = 0.2;
1703 const float itsMax = 0.8;
1704 if ((itsFrac > itsMin) && (itsFrac < itsMax)) {
1705 avgVtxITS[0].addValue(vtx.getX());
1706 avgVtxITS[1].addValue(vtx.getY());
1707 avgVtxITS[2].addValue(vtx.getZ());
1708 avgVtxITS[3].addValue(vtx.getNContributors());
1709 }
1710
1711 const float itstpcMax = 0.95;
1712 if (itstpcFrac < itstpcMax) {
1713 avgVtxITSTPC[0].addValue(vtx.getX());
1714 avgVtxITSTPC[1].addValue(vtx.getY());
1715 avgVtxITSTPC[2].addValue(vtx.getZ());
1716 avgVtxITSTPC[3].addValue(vtx.getNContributors());
1717 }
1718 }
1719
1720 // ITS
1721 mBufferDCA.nPrimVertices_ITS.front() = avgVtxITS[3].getValues().size();
1722 mBufferDCA.nVertexContributors_ITS_Median.front() = avgVtxITS[3].getMedian();
1723 mBufferDCA.nVertexContributors_ITS_RMS.front() = avgVtxITS[3].getStdDev();
1724 mBufferDCA.vertexX_ITS_Median.front() = avgVtxITS[0].getMedian();
1725 mBufferDCA.vertexY_ITS_Median.front() = avgVtxITS[1].getMedian();
1726 mBufferDCA.vertexZ_ITS_Median.front() = avgVtxITS[2].getMedian();
1727 mBufferDCA.vertexX_ITS_RMS.front() = avgVtxITS[0].getStdDev();
1728 mBufferDCA.vertexY_ITS_RMS.front() = avgVtxITS[1].getStdDev();
1729 mBufferDCA.vertexZ_ITS_RMS.front() = avgVtxITS[2].getStdDev();
1730
1731 // ITS-TPC
1732 mBufferDCA.nPrimVertices_ITSTPC.front() = avgVtxITSTPC[3].getValues().size();
1733 mBufferDCA.nVertexContributors_ITSTPC_Median.front() = avgVtxITSTPC[3].getMedian();
1734 mBufferDCA.nVertexContributors_ITSTPC_RMS.front() = avgVtxITSTPC[3].getStdDev();
1735 mBufferDCA.vertexX_ITSTPC_Median.front() = avgVtxITSTPC[0].getMedian();
1736 mBufferDCA.vertexY_ITSTPC_Median.front() = avgVtxITSTPC[1].getMedian();
1737 mBufferDCA.vertexZ_ITSTPC_Median.front() = avgVtxITSTPC[2].getMedian();
1738 mBufferDCA.vertexX_ITSTPC_RMS.front() = avgVtxITSTPC[0].getStdDev();
1739 mBufferDCA.vertexY_ITSTPC_RMS.front() = avgVtxITSTPC[1].getStdDev();
1740 mBufferDCA.vertexZ_ITSTPC_RMS.front() = avgVtxITSTPC[2].getStdDev();
1741
1742 // quantiles and truncated mean
1743 RobustAverage avg(vertices.size(), false);
1744 for (const auto& vtx : vertices) {
1745 if (vtx.getNContributors() > mMinTracksPerVertex) {
1746 // transform by n^0.5 to get more flat distribution
1747 avg.addValue(std::sqrt(vtx.getNContributors()));
1748 }
1749 }
1750
1751 // nPrimVertices_Quantiles
1752 int sizeQ = mBufferDCA.nVertexContributors_Quantiles.size();
1753 const int nBinsQ = 20;
1754 if (sizeQ >= (nBinsQ + 3)) {
1755 for (int iq = 0; iq < nBinsQ; ++iq) {
1756 const float quantile = (iq + 1) / static_cast<float>(nBinsQ);
1757 const float val = avg.getQuantile(quantile, 1);
1758 mBufferDCA.nVertexContributors_Quantiles[iq] = val * val;
1759 }
1760 const float tr0 = avg.getTrunctedMean(0.05, 0.95);
1761 const float tr1 = avg.getTrunctedMean(0.1, 0.9);
1762 const float tr2 = avg.getTrunctedMean(0.2, 0.8);
1763 mBufferDCA.nVertexContributors_Quantiles[sizeQ - 3] = tr0 * tr0;
1764 mBufferDCA.nVertexContributors_Quantiles[sizeQ - 2] = tr1 * tr1;
1765 mBufferDCA.nVertexContributors_Quantiles[sizeQ - 1] = tr2 * tr2;
1766 }
1767 mBufferDCA.nPrimVertices.front() = vertices.size();
1768
1769 return indicesITSTPC_vtx;
1770 }
1771
1773 int getNBins() const { return mBufferDCA.mTSTPC.getNBins(); }
1774
1775 ValsdEdx getdEdxVars(bool useQMax, const TrackTPC& track) const
1776 {
1777 const float dedx = useQMax ? track.getdEdx().dEdxMaxTPC : track.getdEdx().dEdxTotTPC;
1778 const float dedxRatioNorm = mPID.getExpectedSignal(track, o2::track::PID::ID(o2::track::PID::Pion)) / dedx;
1779 float dedxNorm = (dedxRatioNorm > 0) ? std::log(mPID.getExpectedSignal(track, o2::track::PID::ID(o2::track::PID::Pion)) / dedx) : -1;
1780 // restrict to specified range
1781 if (std::abs(dedxNorm) > mMaxdEdxRatio) {
1782 dedxNorm = -1;
1783 }
1784
1785 // get log(dedxRegion / dedx)
1786 float dedxIROC = -1;
1787 float dedxOROC1 = -1;
1788 float dedxOROC2 = -1;
1789 float dedxOROC3 = -1;
1790 if (dedx > 0) {
1791 dedxIROC = useQMax ? track.getdEdx().dEdxMaxIROC : track.getdEdx().dEdxTotIROC;
1792 dedxOROC1 = useQMax ? track.getdEdx().dEdxMaxOROC1 : track.getdEdx().dEdxTotOROC1;
1793 dedxOROC2 = useQMax ? track.getdEdx().dEdxMaxOROC2 : track.getdEdx().dEdxTotOROC2;
1794 dedxOROC3 = useQMax ? track.getdEdx().dEdxMaxOROC3 : track.getdEdx().dEdxTotOROC3;
1795 dedxIROC /= dedx;
1796 dedxOROC1 /= dedx;
1797 dedxOROC2 /= dedx;
1798 dedxOROC3 /= dedx;
1799 dedxIROC = (dedxIROC > 0) ? std::log(dedxIROC) : -1;
1800 dedxOROC1 = (dedxOROC1 > 0) ? std::log(dedxOROC1) : -1;
1801 dedxOROC2 = (dedxOROC2 > 0) ? std::log(dedxOROC2) : -1;
1802 dedxOROC3 = (dedxOROC3 > 0) ? std::log(dedxOROC3) : -1;
1803
1804 // restrict to specified range
1805 if (std::abs(dedxIROC) > mMaxdEdxRegionRatio) {
1806 dedxIROC = -1;
1807 }
1808 if (std::abs(dedxOROC1) > mMaxdEdxRegionRatio) {
1809 dedxOROC1 = -1;
1810 }
1811 if (std::abs(dedxOROC2) > mMaxdEdxRegionRatio) {
1812 dedxOROC2 = -1;
1813 }
1814 if (std::abs(dedxOROC3) > mMaxdEdxRegionRatio) {
1815 dedxOROC3 = -1;
1816 }
1817 }
1818 return ValsdEdx{dedxNorm, dedxIROC, dedxOROC1, dedxOROC2, dedxOROC3};
1819 }
1820};
1821
1822o2::framework::DataProcessorSpec getTPCTimeSeriesSpec(const bool disableWriter, const o2::base::Propagator::MatCorrType matType, const bool enableUnbinnedWriter, GTrackID::mask_t src)
1823{
1824 auto dataRequest = std::make_shared<DataRequest>();
1825 bool useMC = false;
1826 GTrackID::mask_t srcTracks = GTrackID::getSourcesMask("TPC,ITS,ITS-TPC,ITS-TPC-TRD,ITS-TPC-TOF,ITS-TPC-TRD-TOF") & src;
1827 srcTracks.set(GTrackID::TPC); // TPC must be always there
1828 dataRequest->requestTracks(srcTracks, useMC);
1829 dataRequest->requestClusters(GTrackID::getSourcesMask("TPC"), useMC);
1830
1831 bool tpcOnly = srcTracks == GTrackID::getSourcesMask("TPC");
1832 if (!tpcOnly) {
1833 dataRequest->requestFT0RecPoints(useMC);
1834 }
1835 dataRequest->requestPrimaryVertices(useMC);
1836
1837 const bool enableAskMatLUT = matType == o2::base::Propagator::MatCorrType::USEMatCorrLUT;
1838 auto ccdbRequest = std::make_shared<o2::base::GRPGeomRequest>(!disableWriter, // orbitResetTime
1839 false, // GRPECS=true for nHBF per TF
1840 false, // GRPLHCIF
1841 true, // GRPMagField
1842 enableAskMatLUT, // askMatLUT
1844 dataRequest->inputs,
1845 true,
1846 true);
1847
1848 o2::tpc::VDriftHelper::requestCCDBInputs(dataRequest->inputs);
1850 std::vector<OutputSpec> outputs;
1851 outputs.emplace_back(o2::header::gDataOriginTPC, getDataDescriptionTimeSeries(), 0, Lifetime::Sporadic);
1852 if (!disableWriter) {
1853 outputs.emplace_back(o2::header::gDataOriginTPC, getDataDescriptionTPCTimeSeriesTFId(), 0, Lifetime::Sporadic);
1854 }
1855
1856 return DataProcessorSpec{
1857 "tpc-time-series",
1858 dataRequest->inputs,
1859 outputs,
1860 AlgorithmSpec{adaptFromTask<TPCTimeSeries>(ccdbRequest, disableWriter, matType, enableUnbinnedWriter, tpcOnly, dataRequest)},
1861 Options{
1862 {"min-momentum", VariantType::Float, 0.2f, {"Minimum momentum of the tracks"}},
1863 {"min-cluster", VariantType::Int, 80, {"Minimum number of clusters of the tracks"}},
1864 {"max-tgl", VariantType::Float, 1.4f, {"Maximum accepted tgl of the tracks"}},
1865 {"max-qPt", VariantType::Float, 5.f, {"Maximum abs(qPt) bin"}},
1866 {"max-snp", VariantType::Float, 0.85f, {"Maximum sinus(phi) for propagation"}},
1867 {"coarse-step", VariantType::Float, 5.f, {"Coarse step during track propagation"}},
1868 {"fine-step", VariantType::Float, 2.f, {"Fine step during track propagation"}},
1869 {"mX-coarse", VariantType::Float, 40.f, {"Perform coarse propagation up to this mx"}},
1870 {"max-tracks", VariantType::Int, -1, {"Number of maximum tracks to process"}},
1871 {"cut-DCA-median", VariantType::Float, 3.f, {"Cut on the DCA: abs(DCA-medianDCA)<cut-DCA-median"}},
1872 {"cut-DCA-RMS", VariantType::Float, 3.f, {"Sigma cut on the DCA"}},
1873 {"refX-for-sector", VariantType::Float, 108.475f, {"Reference local x position for the sector information (default centre of IROC)"}},
1874 {"refX-for-outer-ITS", VariantType::Float, 60.f, {"Reference local x position for matching at outer ITS"}},
1875 {"tgl-bins", VariantType::Int, 30, {"Number of tgl bins for time series variables"}},
1876 {"phi-bins", VariantType::Int, 54, {"Number of phi bins for time series variables"}},
1877 {"qPt-bins", VariantType::Int, 18, {"Number of qPt bins for time series variables"}},
1878 {"mult-bins", VariantType::Int, 25, {"Number of multiplicity bins for time series variables"}},
1879 {"mult-max", VariantType::Int, 50000, {"MAximum multiplicity bin"}},
1880 {"threads", VariantType::Int, 4, {"Number of parallel threads"}},
1881 {"max-ITS-TPC-DCAr", VariantType::Float, 0.2f, {"Maximum absolut DCAr value for ITS-TPC tracks"}},
1882 {"max-ITS-TPC-DCAz", VariantType::Float, 10.f, {"Maximum absolut DCAz value for ITS-TPC tracks - larger due to vertex spread"}},
1883 {"max-ITS-TPC-DCAr_comb", VariantType::Float, 0.2f, {"Maximum absolut DCAr value for ITS-TPC tracks to vertex for combined DCA"}},
1884 {"max-ITS-TPC-DCAz_comb", VariantType::Float, 0.2f, {"Maximum absolut DCAr value for ITS-TPC tracks to vertex for combined DCA"}},
1885 {"MIP-dedx", VariantType::Float, 50.f, {"MIP dEdx for MIP/dEdx monitoring"}},
1886 {"time-window-mult-mus", VariantType::Float, 50.f, {"Time window in micro s for multiplicity estimate"}},
1887 {"sqrts", VariantType::Float, 13600.f, {"Centre of mass energy used for downsampling"}},
1888 {"min-tracks-per-vertex", VariantType::Int, 6, {"Minimum number of tracks per vertex required"}},
1889 {"max-dedx-ratio", VariantType::Float, 0.3f, {"Maximum absolute log(dedx(pion)/dedx) ratio"}},
1890 {"max-dedx-region-ratio", VariantType::Float, 0.5f, {"Maximum absolute log(dedx(region)/dedx) ratio"}},
1891 {"sample-unbinned-tsallis", VariantType::Bool, false, {"Perform sampling of unbinned data based on Tsallis function"}},
1892 {"sampling-factor", VariantType::Float, 0.001f, {"Sampling factor in case sample-unbinned-tsallis is used"}},
1893 {"disable-min-bias-trigger", VariantType::Bool, false, {"Disable the minimum bias trigger for skimmed data"}},
1894 {"out-file-unbinned", VariantType::String, "time_series_tracks.root", {"name of the output file for the unbinned data"}},
1895 {"max-occupancy-bins", VariantType::Int, 912, {"Maximum number of occupancy bins"}}}};
1896}
1897
1898} // namespace tpc
1899} // end namespace o2
std::vector< unsigned long > times
Wrapper container for different reconstructed object types.
Definition of the TOF cluster.
const auto bins
Definition PID.cxx:49
uint64_t vertex
Definition RawEventData.h:9
int16_t time
Definition RawEventData.h:4
Definition of the FIT RecPoints class.
int32_t i
Helper for geometry and GRP related CCDB requests.
calibrator class for accumulating integrated clusters
Class to store the output of the matching to TOF.
Definition of the parameter class for the detector electronics.
Helper class to extract pressure and temperature.
uint32_t j
Definition RawData.h:0
uint32_t side
Definition RawData.h:0
uint32_t padrow
Definition RawData.h:5
class for performing robust averaging and outlier filtering
Definition of the ITS track.
Result of refitting TPC-ITS matched track.
Helper class to extract VDrift from different sources.
Extention of GlobalTrackID by flags relevant for verter-track association.
Referenc on track indices contributing to the vertex, with possibility chose tracks from specific sou...
auto getOrbitResetTimeMS() const
void checkUpdates(o2::framework::ProcessingContext &pc)
bool finaliseCCDB(o2::framework::ConcreteDataMatcher &matcher, void *obj)
static GRPGeomHelper & instance()
void setRequest(std::shared_ptr< GRPGeomRequest > req)
GPUd() value_type estimateLTFast(o2 static GPUd() float estimateLTIncrement(const o2 PropagatorImpl * Instance(bool uninitialized=false)
Definition Propagator.h:143
static mask_t getSourcesMask(const std::string_view srcList)
Static class with identifiers, bitmasks and names for ALICE detectors.
Definition DetID.h:58
void snapshot(const Output &spec, T const &object)
DataAllocator & outputs()
The data allocator is used to allocate memory for the output data.
Cluster class for TOF.
Definition Cluster.h:37
static constexpr unsigned int PADROWS
total number of pad rows
Definition Mapper.h:528
static void requestCCDBInputs(std::vector< o2::framework::InputSpec > &inputs)
float getPressure(const ULong64_t timestamp) const
get pressure for given time stamp in ms
void extractCCDBInputs(o2::framework::ProcessingContext &pc) const
trigger checking for CCDB objects
float getMeanTemperature(const ULong64_t timestamp) const
get mean temperature over A and C side
bool accountCCDBInputs(const o2::framework::ConcreteDataMatcher &matcher, void *obj)
check for new CCDB objects
void finaliseCCDB(o2::framework::ConcreteDataMatcher &matcher, void *obj) final
void run(ProcessingContext &pc) final
void init(framework::InitContext &ic) final
TPCTimeSeries(std::shared_ptr< o2::base::GRPGeomRequest > req, const bool disableWriter, const o2::base::Propagator::MatCorrType matType, const bool enableUnbinnedWriter, const bool tpcOnly, std::shared_ptr< o2::globaltracking::DataRequest > dr)
\constructor
void endOfStream(EndOfStreamContext &eos) final
static void requestCCDBInputs(std::vector< o2::framework::InputSpec > &inputs, bool laser=true, bool itstpcTgl=true)
void extractCCDBInputs(o2::framework::ProcessingContext &pc, bool laser=true, bool itstpcTgl=true)
const VDriftCorrFact & getVDriftObject() const
bool accountCCDBInputs(const o2::framework::ConcreteDataMatcher &matcher, void *obj)
bool isUpdated() const
pid_constants::ID ID
Definition PID.h:92
static constexpr ID Pion
Definition PID.h:96
GLdouble n
Definition glcorearb.h:1982
GLenum src
Definition glcorearb.h:1767
GLuint buffer
Definition glcorearb.h:655
GLuint GLuint GLfloat weight
Definition glcorearb.h:5477
GLsizei GLsizei GLchar * source
Definition glcorearb.h:798
GLint GLint GLsizei GLint GLenum GLenum type
Definition glcorearb.h:275
GLuint GLfloat * val
Definition glcorearb.h:1582
GLuint GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat t0
Definition glcorearb.h:5034
GLint GLuint mask
Definition glcorearb.h:291
GLsizei GLenum * sources
Definition glcorearb.h:2516
constexpr o2::header::DataOrigin gDataOriginTPC
Definition DataHeader.h:576
uint8_t itsSharedClusterMap uint8_t
constexpr double LHCOrbitMUS
constexpr double LHCBunchSpacingNS
constexpr float TwoPI
Defining PrimaryVertex explicitly as messageable.
Definition TFIDInfo.h:20
std::vector< ConfigParamSpec > Options
o2::tpc::PIDResponse PIDResponse
Definition PIDStudy.cxx:69
int angle2Sector(float phi)
Definition Utils.h:183
float sector2Angle(int sect)
Definition Utils.h:193
uint32_t getFirstTForbit(o2::framework::ProcessingContext &pc)
uint64_t getTimeStamp(o2::framework::ProcessingContext &pc)
uint64_t getRunNumber(o2::framework::ProcessingContext &pc)
constexpr unsigned char SECTORSPERSIDE
Definition Defs.h:40
o2::framework::DataProcessorSpec getTPCTimeSeriesSpec(const bool disableWriter, const o2::base::Propagator::MatCorrType matType, const bool enableUnbinnedWriter, o2::dataformats::GlobalTrackID::mask_t src)
Side
TPC readout sidE.
Definition Defs.h:35
@ A
Definition Defs.h:35
@ C
Definition Defs.h:36
TrackParCovF TrackParCov
Definition Track.h:33
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
static void fillTFIDInfo(o2::framework::ProcessingContext &pc, o2::dataformats::TFIDInfo &ti)
GlobalIDSet getSingleDetectorRefs(GTrackID gidx) const
GTrackID getTPCContributorGID(GTrackID source) const
void collectData(o2::framework::ProcessingContext &pc, const DataRequest &request)
gsl::span< const unsigned int > occupancyMapTPC
externally set TPC clusters occupancy map
static bool downsampleTsallisCharged(float pt, float factorPt, float sqrts, float &weight, float rnd, float mass=0.13957)
Definition Tsallis.cxx:31
std::vector< float > mDCAr_comb_A_RMS
DCAr RMS for ITS-TPC track - A-side.
std::vector< float > mTPCSigmaZ2A_RMS
sigmaZ2 RMS at vertex
std::vector< float > mSqrtITSChi2_Ncl_C_Median
sqrt(ITC chi2 / ncl)
ITSTPC_Matching mITSTPCAll
ITS-TPC matching efficiency for ITS standalone + afterburner.
std::vector< float > mITSTPCDeltaP4_C_RMS
RMS of track param TPC - track param ITS-TPC for param 4 - A-side.
std::vector< float > mTPCSigmaZ2C_RMS
sigmaZ2 RMS at vertex
TimeSeries mTSTPC
TPC standalone DCAs.
std::vector< float > mITSTPCDeltaP4_A_Median
track param TPC - track param ITS-TPC for param 4 - A-side
ITSTPC_Matching mITSTPCStandalone
ITS-TPC matching efficiency for ITS standalone.
std::vector< float > mITSTPCDeltaP4_C_Median
track param TPC - track param ITS-TPC for param 4 - A-side
std::vector< float > mITSTPCDeltaP3_A_RMS
RMS of track param TPC - track param ITS-TPC for param 3 - A-side.
std::vector< float > mITS_A_NCl_Median
its number of clusters
std::vector< float > mDCAr_comb_C_Median
DCAr for ITS-TPC track - C-side.
std::vector< float > mSqrtITSChi2_Ncl_A_Median
sqrt(ITC chi2 / ncl)
std::vector< float > mSqrtITSChi2_Ncl_C_RMS
sqrt(ITC chi2 / ncl)
std::vector< float > mDCAz_comb_A_RMS
DCAz RMS for ITS-TPC track - A-side.
std::vector< float > mITSTPCDeltaP2_C_RMS
RMS of track param TPC - track param ITS-TPC for param 2 - A-side.
std::vector< float > mITSTPCDeltaP3_A_Median
track param TPC - track param ITS-TPC for param 3 - A-side
void resize(const unsigned int nTotal)
resize buffer for accumulated currents
std::vector< float > mITSTPCDeltaP2_A_Median
track param TPC - track param ITS-TPC for param 2 - A-side
std::vector< float > mSqrtITSChi2_Ncl_A_RMS
sqrt(ITC chi2 / ncl)
std::vector< float > mITSTPCDeltaP2_C_Median
track param TPC - track param ITS-TPC for param 2 - A-side
std::vector< unsigned int > mOccupancyMapTPC
cluster occupancy map
std::vector< float > mITSTPCDeltaP4_A_RMS
RMS of track param TPC - track param ITS-TPC for param 4 - A-side.
std::vector< float > mTPCSigmaY2A_Median
sigmaY2 at vertex
std::vector< float > mITSTPCDeltaP3_C_RMS
RMS of track param TPC - track param ITS-TPC for param 3 - A-side.
TimeSeries mTSITSTPC
ITS-TPC standalone DCAs.
std::vector< float > mTPCSigmaY2C_Median
sigmaY2 at vertex
std::vector< float > mTPCSigmaZ2C_Median
sigmaZ2 at vertex
std::vector< float > mITS_C_NCl_RMS
its number of clusters
std::vector< float > mDCAz_comb_A_Median
DCAz for ITS-TPC track - A-side.
std::vector< float > mDCAz_comb_C_RMS
DCAz RMS for ITS-TPC track - C-side.
ITSTPC_Matching mITSTPCAfterburner
ITS-TPC matchin efficiency fir ITS afterburner.
std::vector< float > mTPCSigmaY2C_RMS
sigmaY2 RMS at vertex
std::vector< float > mITSTPCDeltaP3_C_Median
track param TPC - track param ITS-TPC for param 3 - A-side
std::vector< float > mDCAz_comb_C_Median
DCAz for ITS-TPC track - C-side.
std::vector< float > mTPCSigmaZ2A_Median
sigmaZ2 at vertex
std::vector< float > mDCAr_comb_C_RMS
DCAr RMS for ITS-TPC track - C-side.
float mVDrift
drift velocity in cm/us
std::vector< float > mITS_C_NCl_Median
its number of clusters
TimeSeriesdEdx mdEdxQTot
time series for dE/dx qTot monitoring
std::vector< float > mITS_A_NCl_RMS
its number of clusters
std::vector< float > mITSTPCDeltaP2_A_RMS
RMS of track param TPC - track param ITS-TPC for param 2 - A-side.
std::vector< float > mTPCSigmaY2A_RMS
sigmaY2 RMS at vertex
std::vector< float > mDCAr_comb_A_Median
DCAr for ITS-TPC track - A-side.
void setBinning(const int nBinsPhi, const int nBinsTgl, const int qPtBins, const int nBinsMult, float tglMax, float qPtMax, float multMax)
TimeSeriesdEdx mdEdxQMax
time series for dE/dx qMax monitoring
int getIndexInt(int slice=0) const
std::vector< float > mDCAr_C_Median
integrated 1D DCAr for C-side weighted mean in phi/tgl slices
std::vector< float > mDCAz_A_RMS
integrated 1D DCAz for A-side RMS in phi/tgl slices
std::vector< float > mDCAz_C_Median
integrated 1D DCAz for C-side median in phi/tgl slices
int getIndexPhi(const int iPhi, int slice=0) const
std::vector< float > mDCAz_A_Median
integrated 1D DCAz for A-side median in phi/tgl slices
std::vector< float > mDCAr_C_RMS
integrated 1D DCAr for C-side RMS in phi/tgl slices
std::vector< float > mDCAr_A_Median
integrated 1D DCAr for A-side median in phi/tgl slices
int getIndexTgl(const int iTgl, int slice=0) const
std::vector< float > mMIPdEdxRatioQMaxA
ratio of MIP/dEdx - qMax -
std::vector< float > mDCAz_C_RMS
integrated 1D DCAz for C-side RMS in phi/tgl slices
int getIndexqPt(const int iqPt, int slice=0) const
std::vector< float > mDCAr_A_RMS
integrated 1D DCAr for A-side RMS in phi/tgl slices
int getIndexMult(const int iMult, int slice=0) const
std::vector< float > mLogdEdx_A_Median
log(dEdx_exp(pion)/dEdx) - A-side
vec clear()
std::uniform_int_distribution< unsigned long long > distr