Project
Loading...
Searching...
No Matches
SVertexer.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
15
27
28#ifdef WITH_OPENMP
29#include <omp.h>
30#endif
31
33
34using namespace o2::vertexing;
35namespace o2f = o2::framework;
40
41//__________________________________________________________________
43{
44 mRecoCont = &recoData;
45 mNV0s = mNCascades = mN3Bodies = 0;
46 updateTimeDependentParams(); // TODO RS: strictly speaking, one should do this only in case of the CCDB objects update
47 mPVertices = recoData.getPrimaryVertices();
48 buildT2V(recoData); // build track->vertex refs from vertex->track (if other workflow will need this, consider producing a message in the VertexTrackMatcher)
49 int ntrP = mTracksPool[POS].size(), ntrN = mTracksPool[NEG].size();
50 if (mStrTracker) {
51 mStrTracker->loadData(recoData);
52 mStrTracker->prepareITStracks();
53 }
54#ifdef WITH_OPENMP
55 int dynGrp = std::min(4, std::max(1, mNThreads / 2));
56#pragma omp parallel for schedule(dynamic, dynGrp) num_threads(mNThreads)
57#endif
58 for (int itp = 0; itp < ntrP; itp++) {
59 auto& seedP = mTracksPool[POS][itp];
60 const int firstN = mVtxFirstTrack[NEG][seedP.vBracket.getMin()];
61 if (firstN < 0) {
62 LOG(debug) << "No partner is found for pos.track " << itp << " out of " << ntrP;
63 continue;
64 }
65 for (int itn = firstN; itn < ntrN; itn++) { // start from the 1st negative track of lowest-ID vertex of positive
66 auto& seedN = mTracksPool[NEG][itn];
67 if (seedN.vBracket > seedP.vBracket) { // all vertices compatible with seedN are in future wrt that of seedP
68 LOG(debug) << "Brackets do not match";
69 break;
70 }
71 if (mSVParams->maxPVContributors < 2 && seedP.gid.isPVContributor() + seedN.gid.isPVContributor() > mSVParams->maxPVContributors) {
72 continue;
73 }
74#ifdef WITH_OPENMP
75 int iThread = omp_get_thread_num();
76#else
77 int iThread = 0;
78#endif
79 checkV0(seedP, seedN, itp, itn, iThread);
80 }
81 }
82
83 produceOutput(pc);
84}
85
86//__________________________________________________________________
88{
89 // sort V0s and Cascades in vertex id
90 struct vid {
91 int thrID;
92 int entry;
93 int vtxID;
94 };
95 for (int ith = 0; ith < mNThreads; ith++) {
96 mNV0s += mV0sIdxTmp[ith].size();
97 mNCascades += mCascadesIdxTmp[ith].size();
98 mN3Bodies += m3bodyIdxTmp[ith].size();
99 }
100 std::vector<vid> v0SortID, cascSortID, nbodySortID;
101 v0SortID.reserve(mNV0s);
102 cascSortID.reserve(mNCascades);
103 nbodySortID.reserve(mN3Bodies);
104 for (int ith = 0; ith < mNThreads; ith++) {
105 for (int j = 0; j < (int)mV0sIdxTmp[ith].size(); j++) {
106 v0SortID.emplace_back(vid{ith, j, mV0sIdxTmp[ith][j].getVertexID()});
107 }
108 for (int j = 0; j < (int)mCascadesIdxTmp[ith].size(); j++) {
109 cascSortID.emplace_back(vid{ith, j, mCascadesIdxTmp[ith][j].getVertexID()});
110 }
111 for (int j = 0; j < (int)m3bodyIdxTmp[ith].size(); j++) {
112 nbodySortID.emplace_back(vid{ith, j, m3bodyIdxTmp[ith][j].getVertexID()});
113 }
114 }
115 std::sort(v0SortID.begin(), v0SortID.end(), [](const vid& a, const vid& b) { return a.vtxID < b.vtxID; });
116 std::sort(cascSortID.begin(), cascSortID.end(), [](const vid& a, const vid& b) { return a.vtxID < b.vtxID; });
117 std::sort(nbodySortID.begin(), nbodySortID.end(), [](const vid& a, const vid& b) { return a.vtxID < b.vtxID; });
118
119 // dpl output
120 auto& v0sIdx = pc.outputs().make<std::vector<V0Index>>(o2f::Output{"GLO", "V0S_IDX", 0});
121 auto& cascsIdx = pc.outputs().make<std::vector<CascadeIndex>>(o2f::Output{"GLO", "CASCS_IDX", 0});
122 auto& body3Idx = pc.outputs().make<std::vector<Decay3BodyIndex>>(o2f::Output{"GLO", "DECAYS3BODY_IDX", 0});
123 auto& fullv0s = pc.outputs().make<std::vector<V0>>(o2f::Output{"GLO", "V0S", 0});
124 auto& fullcascs = pc.outputs().make<std::vector<Cascade>>(o2f::Output{"GLO", "CASCS", 0});
125 auto& full3body = pc.outputs().make<std::vector<Decay3Body>>(o2f::Output{"GLO", "DECAYS3BODY", 0});
126 auto& v0Refs = pc.outputs().make<std::vector<RRef>>(o2f::Output{"GLO", "PVTX_V0REFS", 0});
127 auto& cascRefs = pc.outputs().make<std::vector<RRef>>(o2f::Output{"GLO", "PVTX_CASCREFS", 0});
128 auto& vtx3bodyRefs = pc.outputs().make<std::vector<RRef>>(o2f::Output{"GLO", "PVTX_3BODYREFS", 0});
129
130 // sorted V0s
131 v0sIdx.reserve(mNV0s);
132 if (mSVParams->createFullV0s) {
133 fullv0s.reserve(mNV0s);
134 }
135 // sorted Cascades
136 cascsIdx.reserve(mNCascades);
137 if (mSVParams->createFullCascades) {
138 fullcascs.reserve(mNCascades);
139 }
140 // sorted 3 body decays
141 body3Idx.reserve(mN3Bodies);
142 if (mSVParams->createFull3Bodies) {
143 full3body.reserve(mN3Bodies);
144 }
145
146 for (const auto& id : v0SortID) {
147 auto& v0idx = mV0sIdxTmp[id.thrID][id.entry];
148 int pos = v0sIdx.size();
149 v0sIdx.push_back(v0idx);
150 v0idx.setVertexID(pos); // this v0 copy will be discarded, use its vertexID to store the new position of final V0
151 if (mSVParams->createFullV0s) {
152 fullv0s.push_back(mV0sTmp[id.thrID][id.entry]);
153 }
154 }
155 // since V0s were reshuffled, we need to correct the cascade -> V0 reference indices
156 for (int ith = 0; ith < mNThreads; ith++) { // merge results of all threads
157 for (size_t ic = 0; ic < mCascadesIdxTmp[ith].size(); ic++) { // before merging fix cascades references on v0
158 auto& cidx = mCascadesIdxTmp[ith][ic];
159 cidx.setV0ID(mV0sIdxTmp[ith][cidx.getV0ID()].getVertexID());
160 }
161 }
162 int cascCnt = 0;
163 for (const auto& id : cascSortID) {
164 cascsIdx.push_back(mCascadesIdxTmp[id.thrID][id.entry]);
165 mCascadesIdxTmp[id.thrID][id.entry].setVertexID(cascCnt++); // memorize new ID
166 if (mSVParams->createFullCascades) {
167 fullcascs.push_back(mCascadesTmp[id.thrID][id.entry]);
168 }
169 }
170 int b3cnt = 0;
171 for (const auto& id : nbodySortID) {
172 body3Idx.push_back(m3bodyIdxTmp[id.thrID][id.entry]);
173 m3bodyIdxTmp[id.thrID][id.entry].setVertexID(b3cnt++); // memorize new ID
174 if (mSVParams->createFull3Bodies) {
175 full3body.push_back(m3bodyTmp[id.thrID][id.entry]);
176 }
177 }
178 if (mStrTracker) {
179 mNStrangeTracks = 0;
180 for (int ith = 0; ith < mNThreads; ith++) {
181 mNStrangeTracks += mStrTracker->getNTracks(ith);
182 }
183
184 std::vector<o2::dataformats::StrangeTrack> strTracksTmp;
185 std::vector<o2::strangeness_tracking::ClusAttachments> strClusTmp;
186 std::vector<o2::MCCompLabel> mcLabTmp;
187 strTracksTmp.reserve(mNStrangeTracks);
188 strClusTmp.reserve(mNStrangeTracks);
189 if (mStrTracker->getMCTruthOn()) {
190 mcLabTmp.reserve(mNStrangeTracks);
191 }
192
193 for (int ith = 0; ith < mNThreads; ith++) { // merge results of all threads
194 auto& strTracks = mStrTracker->getStrangeTrackVec(ith);
195 auto& strClust = mStrTracker->getClusAttachments(ith);
196 auto& stcTrMCLab = mStrTracker->getStrangeTrackLabels(ith);
197 for (int i = 0; i < (int)strTracks.size(); i++) {
198 auto& t = strTracks[i];
199 if (t.mPartType == o2::dataformats::kStrkV0) {
200 t.mDecayRef = mV0sIdxTmp[ith][t.mDecayRef].getVertexID(); // reassign merged V0 ID
201 } else if (t.mPartType == o2::dataformats::kStrkCascade) {
202 t.mDecayRef = mCascadesIdxTmp[ith][t.mDecayRef].getVertexID(); // reassign merged Cascase ID
203 } else if (t.mPartType == o2::dataformats::kStrkThreeBody) {
204 t.mDecayRef = m3bodyIdxTmp[ith][t.mDecayRef].getVertexID(); // reassign merged Cascase ID
205 } else {
206 LOGP(fatal, "Unknown strange track decay reference type {} for index {}", int(t.mPartType), t.mDecayRef);
207 }
208
209 strTracksTmp.push_back(t);
210 strClusTmp.push_back(strClust[i]);
211 if (mStrTracker->getMCTruthOn()) {
212 mcLabTmp.push_back(stcTrMCLab[i]);
213 }
214 }
215 }
216
217 auto& strTracksOut = pc.outputs().make<std::vector<o2::dataformats::StrangeTrack>>(o2f::Output{"GLO", "STRANGETRACKS", 0});
218 auto& strClustOut = pc.outputs().make<std::vector<o2::strangeness_tracking::ClusAttachments>>(o2f::Output{"GLO", "CLUSUPDATES", 0});
220 strTracksOut.resize(mNStrangeTracks);
221 strClustOut.resize(mNStrangeTracks);
222 if (mStrTracker->getMCTruthOn()) {
223 mcLabsOut.resize(mNStrangeTracks);
224 }
225
226 std::vector<int> sortIdx(strTracksTmp.size());
227 std::iota(sortIdx.begin(), sortIdx.end(), 0);
228 // if mNTreads > 1 we need to sort tracks, clus and MCLabs by their mDecayRef
229 if (mNThreads > 1 && mNStrangeTracks > 1) {
230 std::sort(sortIdx.begin(), sortIdx.end(), [&strTracksTmp](int i1, int i2) { return strTracksTmp[i1].mDecayRef < strTracksTmp[i2].mDecayRef; });
231 }
232
233 for (int i = 0; i < (int)sortIdx.size(); i++) {
234 strTracksOut[i] = strTracksTmp[sortIdx[i]];
235 strClustOut[i] = strClusTmp[sortIdx[i]];
236 if (mStrTracker->getMCTruthOn()) {
237 mcLabsOut[i] = mcLabTmp[sortIdx[i]];
238 }
239 }
240
241 if (mStrTracker->getMCTruthOn()) {
242 auto& strTrMCLableOut = pc.outputs().make<std::vector<o2::MCCompLabel>>(o2f::Output{"GLO", "STRANGETRACKS_MC", 0});
243 strTrMCLableOut.swap(mcLabsOut);
244 }
245 }
246
247 for (int ith = 0; ith < mNThreads; ith++) { // clean unneeded s.vertices
248 mV0sTmp[ith].clear();
249 mCascadesTmp[ith].clear();
250 m3bodyTmp[ith].clear();
251 mV0sIdxTmp[ith].clear();
252 mCascadesIdxTmp[ith].clear();
253 m3bodyIdxTmp[ith].clear();
254 }
255
256 extractPVReferences(v0sIdx, v0Refs, cascsIdx, cascRefs, body3Idx, vtx3bodyRefs);
257}
258
259//__________________________________________________________________
261{
262}
263
264//__________________________________________________________________
265void SVertexer::updateTimeDependentParams()
266{
267 // TODO RS: strictly speaking, one should do this only in case of the CCDB objects update
268 static bool updatedOnce = false;
269 if (!updatedOnce) {
270 updatedOnce = true;
271 mSVParams = &SVertexerParams::Instance();
272 if (mSVParams->mExcludeTPCtracks && !mRecoCont->isTrackSourceLoaded(GIndex::TPC)) {
273 LOGP(fatal, "TPC tracks requested but not provided");
274 }
275 // precalculated selection cuts
276 mMinR2ToMeanVertex = mSVParams->minRToMeanVertex * mSVParams->minRToMeanVertex;
277 mMaxR2ToMeanVertexCascV0 = mSVParams->maxRToMeanVertexCascV0 * mSVParams->maxRToMeanVertexCascV0;
278 mMaxDCAXY2ToMeanVertex = mSVParams->maxDCAXYToMeanVertex * mSVParams->maxDCAXYToMeanVertex;
279 mMaxDCAXY2ToMeanVertexV0Casc = mSVParams->maxDCAXYToMeanVertexV0Casc * mSVParams->maxDCAXYToMeanVertexV0Casc;
280 mMaxDCAXY2ToMeanVertex3bodyV0 = mSVParams->maxDCAXYToMeanVertex3bodyV0 * mSVParams->maxDCAXYToMeanVertex3bodyV0;
281 mMinR2DiffV0Casc = mSVParams->minRDiffV0Casc * mSVParams->minRDiffV0Casc;
282 mMinPt2V0 = mSVParams->minPtV0 * mSVParams->minPtV0;
283 mMaxTgl2V0 = mSVParams->maxTglV0 * mSVParams->maxTglV0;
284 mMinPt2Casc = mSVParams->minPtCasc * mSVParams->minPtCasc;
285 mMaxTgl2Casc = mSVParams->maxTglCasc * mSVParams->maxTglCasc;
286 mMinPt23Body = mSVParams->minPt3Body * mSVParams->minPt3Body;
287 mMaxTgl23Body = mSVParams->maxTgl3Body * mSVParams->maxTgl3Body;
288 setupThreads();
289 }
290 auto bz = o2::base::Propagator::Instance()->getNominalBz();
291 mV0Hyps[HypV0::Photon].set(PID::Photon, PID::Electron, PID::Electron, mSVParams->pidCutsPhoton, bz);
292 mV0Hyps[HypV0::K0].set(PID::K0, PID::Pion, PID::Pion, mSVParams->pidCutsK0, bz);
293 mV0Hyps[HypV0::Lambda].set(PID::Lambda, PID::Proton, PID::Pion, mSVParams->pidCutsLambda, bz);
294 mV0Hyps[HypV0::AntiLambda].set(PID::Lambda, PID::Pion, PID::Proton, mSVParams->pidCutsLambda, bz);
299 mCascHyps[HypCascade::XiMinus].set(PID::XiMinus, PID::Lambda, PID::Pion, mSVParams->pidCutsXiMinus, bz, mSVParams->maximalCascadeWidth);
301
310
311 for (auto& ft : mFitterV0) {
312 ft.setBz(bz);
313 }
314 for (auto& ft : mFitterCasc) {
315 ft.setBz(bz);
316 }
317 for (auto& ft : mFitter3body) {
318 ft.setBz(bz);
319 }
320
321 mPIDresponse.setBetheBlochParams(mSVParams->mBBpars);
322}
323
324//______________________________________________
326{
327 mTPCVDrift = v.refVDrift * v.corrFact;
328 mTPCVDriftCorrFact = v.corrFact;
329 mTPCVDriftRef = v.refVDrift;
330 mTPCDriftTimeOffset = v.getTimeOffset();
331 mTPCBin2Z = mTPCVDrift / mMUS2TPCBin;
332}
333//______________________________________________
335{
336 mTPCCorrMapsHelper = maph;
337}
338
339//__________________________________________________________________
340void SVertexer::setupThreads()
341{
342 if (!mV0sTmp.empty()) {
343 return;
344 }
345 mV0sTmp.resize(mNThreads);
346 mCascadesTmp.resize(mNThreads);
347 m3bodyTmp.resize(mNThreads);
348 mV0sIdxTmp.resize(mNThreads);
349 mCascadesIdxTmp.resize(mNThreads);
350 m3bodyIdxTmp.resize(mNThreads);
351 mFitterV0.resize(mNThreads);
352 mBz = o2::base::Propagator::Instance()->getNominalBz();
353 int fitCounter = 0;
354 for (auto& fitter : mFitterV0) {
355 fitter.setFitterID(fitCounter++);
356 fitter.setBz(mBz);
357 fitter.setUseAbsDCA(mSVParams->useAbsDCA);
358 fitter.setPropagateToPCA(false);
359 fitter.setMaxR(mSVParams->maxRIni);
360 fitter.setMinParamChange(mSVParams->minParamChange);
361 fitter.setMinRelChi2Change(mSVParams->minRelChi2Change);
362 fitter.setMaxDZIni(mSVParams->maxDZIni);
363 fitter.setMaxDXYIni(mSVParams->maxDXYIni);
364 fitter.setMaxChi2(mSVParams->maxChi2);
365 fitter.setMatCorrType(o2::base::Propagator::MatCorrType(mSVParams->matCorr));
366 fitter.setUsePropagator(mSVParams->usePropagator);
367 fitter.setRefitWithMatCorr(mSVParams->refitWithMatCorr);
368 fitter.setMaxStep(mSVParams->maxStep);
369 fitter.setMaxSnp(mSVParams->maxSnp);
370 fitter.setMinXSeed(mSVParams->minXSeed);
371 }
372 mFitterCasc.resize(mNThreads);
373 fitCounter = 1000;
374 for (auto& fitter : mFitterCasc) {
375 fitter.setFitterID(fitCounter++);
376 fitter.setBz(mBz);
377 fitter.setUseAbsDCA(mSVParams->useAbsDCA);
378 fitter.setPropagateToPCA(false);
379 fitter.setMaxR(mSVParams->maxRIniCasc);
380 fitter.setMinParamChange(mSVParams->minParamChange);
381 fitter.setMinRelChi2Change(mSVParams->minRelChi2Change);
382 fitter.setMaxDZIni(mSVParams->maxDZIni);
383 fitter.setMaxDXYIni(mSVParams->maxDXYIni);
384 fitter.setMaxChi2(mSVParams->maxChi2);
385 fitter.setMatCorrType(o2::base::Propagator::MatCorrType(mSVParams->matCorr));
386 fitter.setUsePropagator(mSVParams->usePropagator);
387 fitter.setRefitWithMatCorr(mSVParams->refitWithMatCorr);
388 fitter.setMaxStep(mSVParams->maxStep);
389 fitter.setMaxSnp(mSVParams->maxSnp);
390 fitter.setMinXSeed(mSVParams->minXSeed);
391 }
392
393 mFitter3body.resize(mNThreads);
394 fitCounter = 2000;
395 for (auto& fitter : mFitter3body) {
396 fitter.setFitterID(fitCounter++);
397 fitter.setBz(mBz);
398 fitter.setUseAbsDCA(mSVParams->useAbsDCA);
399 fitter.setPropagateToPCA(false);
400 fitter.setMaxR(mSVParams->maxRIni3body);
401 fitter.setMinParamChange(mSVParams->minParamChange);
402 fitter.setMinRelChi2Change(mSVParams->minRelChi2Change);
403 fitter.setMaxDZIni(mSVParams->maxDZIni);
404 fitter.setMaxDXYIni(mSVParams->maxDXYIni);
405 fitter.setMaxChi2(mSVParams->maxChi2);
406 fitter.setMatCorrType(o2::base::Propagator::MatCorrType(mSVParams->matCorr));
407 fitter.setUsePropagator(mSVParams->usePropagator);
408 fitter.setRefitWithMatCorr(mSVParams->refitWithMatCorr);
409 fitter.setMaxStep(mSVParams->maxStep);
410 fitter.setMaxSnp(mSVParams->maxSnp);
411 fitter.setMinXSeed(mSVParams->minXSeed);
412 }
413}
414
415//__________________________________________________________________
416bool SVertexer::acceptTrack(const GIndex gid, const o2::track::TrackParCov& trc) const
417{
418 if (gid.isPVContributor() && mSVParams->maxPVContributors < 1) {
419 return false;
420 }
421
422 // DCA to mean vertex
423 if (mSVParams->minDCAToPV > 0.f) {
424 o2::track::TrackPar trp(trc);
425 std::array<float, 2> dca;
426 auto* prop = o2::base::Propagator::Instance();
427 if (mSVParams->usePropagator) {
428 if (trp.getX() > mSVParams->minRFor3DField && !prop->PropagateToXBxByBz(trp, mSVParams->minRFor3DField, mSVParams->maxSnp, mSVParams->maxStep, o2::base::Propagator::MatCorrType(mSVParams->matCorr))) {
429 return true; // we don't need actually to propagate to the beam-line
430 }
431 if (!prop->propagateToDCA(mMeanVertex.getXYZ(), trp, prop->getNominalBz(), mSVParams->maxStep, o2::base::Propagator::MatCorrType(mSVParams->matCorr), &dca)) {
432 return true;
433 }
434 } else {
435 if (!trp.propagateParamToDCA(mMeanVertex.getXYZ(), prop->getNominalBz(), &dca)) {
436 return true;
437 }
438 }
439 if (std::abs(dca[0]) < mSVParams->minDCAToPV) {
440 return false;
441 }
442 }
443 return true;
444}
445
446//__________________________________________________________________
447void SVertexer::buildT2V(const o2::globaltracking::RecoContainer& recoData) // accessor to various tracks
448{
449 // build track->vertices from vertices->tracks, rejecting vertex contributors if requested
450 auto trackIndex = recoData.getPrimaryVertexMatchedTracks(); // Global ID's for associated tracks
451 auto vtxRefs = recoData.getPrimaryVertexMatchedTrackRefs(); // references from vertex to these track IDs
452 bool isTPCloaded = recoData.isTrackSourceLoaded(GIndex::TPC);
453 bool isITSloaded = recoData.isTrackSourceLoaded(GIndex::ITS);
454 bool isITSTPCloaded = recoData.isTrackSourceLoaded(GIndex::ITSTPC);
455 if (isTPCloaded && !mSVParams->mExcludeTPCtracks) {
456 mTPCTracksArray = recoData.getTPCTracks();
457 mTPCTrackClusIdx = recoData.getTPCTracksClusterRefs();
458 mTPCClusterIdxStruct = &recoData.inputsTPCclusters->clusterIndex;
459 mTPCRefitterShMap = recoData.clusterShMapTPC;
460 mTPCRefitterOccMap = mRecoCont->occupancyMapTPC;
461 mTPCRefitter = std::make_unique<o2::gpu::GPUO2InterfaceRefit>(mTPCClusterIdxStruct, mTPCCorrMapsHelper, o2::base::Propagator::Instance()->getNominalBz(), mTPCTrackClusIdx.data(), 0, mTPCRefitterShMap.data(), mTPCRefitterOccMap.data(), mTPCRefitterOccMap.size(), nullptr, o2::base::Propagator::Instance());
462 }
463
464 std::unordered_map<GIndex, std::pair<int, int>> tmap;
465 std::unordered_map<GIndex, bool> rejmap;
466 int nv = vtxRefs.size() - 1; // The last entry is for unassigned tracks, ignore them
467 for (int i = 0; i < 2; i++) {
468 mTracksPool[i].clear();
469 mVtxFirstTrack[i].clear();
470 mVtxFirstTrack[i].resize(nv, -1);
471 }
472 for (int iv = 0; iv < nv; iv++) {
473 const auto& vtref = vtxRefs[iv];
474 int it = vtref.getFirstEntry(), itLim = it + vtref.getEntries();
475 for (; it < itLim; it++) {
476 auto tvid = trackIndex[it];
477 if (!recoData.isTrackSourceLoaded(tvid.getSource())) {
478 continue;
479 }
480 if (tvid.getSource() == GIndex::TPC) {
481 if (mSVParams->mExcludeTPCtracks) {
482 continue;
483 }
484 // unconstrained TPC tracks require special treatment: there is no point in checking DCA to mean vertex since it is not precise,
485 // but we need to create a clone of TPC track constrained to this particular vertex time.
486 if (processTPCTrack(mTPCTracksArray[tvid], tvid, iv)) {
487 continue;
488 }
489 }
490
491 if (tvid.isAmbiguous()) { // was this track already processed?
492 auto tref = tmap.find(tvid);
493 if (tref != tmap.end()) {
494 mTracksPool[tref->second.second][tref->second.first].vBracket.setMax(iv); // this track was already processed with other vertex, account the latter
495 continue;
496 }
497 // was it already rejected?
498 if (rejmap.find(tvid) != rejmap.end()) {
499 continue;
500 }
501 }
502 const auto& trc = recoData.getTrackParam(tvid);
503
504 bool hasTPC = false;
505 bool heavyIonisingParticle = false;
506 bool compatibleWithProton = mSVParams->mFractiondEdxforCascBaryons > 0.999f; // if 1 or above, accept all regardless of TPC
507 auto tpcGID = recoData.getTPCContributorGID(tvid);
508 if (tpcGID.isIndexSet() && isTPCloaded) {
509 hasTPC = true;
510 auto& tpcTrack = recoData.getTPCTrack(tpcGID);
511 float dEdxTPC = tpcTrack.getdEdx().dEdxTotTPC;
512 if (dEdxTPC > mSVParams->minTPCdEdx && trc.getP() > mSVParams->minMomTPCdEdx) // accept high dEdx tracks (He3, He4)
513 {
514 heavyIonisingParticle = true;
515 }
516 auto protonId = o2::track::PID::Proton;
517 float dEdxExpected = mPIDresponse.getExpectedSignal(tpcTrack, protonId);
518 float fracDevProton = std::abs((dEdxTPC - dEdxExpected) / dEdxExpected);
519 if (fracDevProton < mSVParams->mFractiondEdxforCascBaryons) {
520 compatibleWithProton = true;
521 }
522 }
523
524 // get Nclusters in the ITS if available
525 int8_t nITSclu = -1;
526 bool shortOBITSOnlyTrack = false;
527 auto itsGID = recoData.getITSContributorGID(tvid);
528 if (itsGID.getSource() == GIndex::ITS) {
529 if (isITSloaded) {
530 auto& itsTrack = recoData.getITSTrack(itsGID);
531 nITSclu = itsTrack.getNumberOfClusters();
532 if (itsTrack.hasHitOnLayer(6) && itsTrack.hasHitOnLayer(5) && itsTrack.hasHitOnLayer(4) && itsTrack.hasHitOnLayer(3)) {
533 shortOBITSOnlyTrack = true;
534 }
535 }
536 } else if (itsGID.getSource() == GIndex::ITSAB) {
537 if (isITSTPCloaded) {
538 auto& itsABTracklet = recoData.getITSABRef(itsGID);
539 nITSclu = itsABTracklet.getNClusters();
540 }
541 }
542 if (!acceptTrack(tvid, trc) && !heavyIonisingParticle) {
543 if (tvid.isAmbiguous()) {
544 rejmap[tvid] = true;
545 }
546 continue;
547 }
548 if ((isTPCloaded && !hasTPC) && (isITSloaded && (nITSclu < mSVParams->mITSSAminNclu && (!shortOBITSOnlyTrack || mSVParams->mRejectITSonlyOBtrack)))) {
549 continue; // reject short ITS-only
550 }
551
552 int posneg = trc.getSign() < 0 ? 1 : 0;
553 float r = std::sqrt(trc.getX() * trc.getX() + trc.getY() * trc.getY());
554 mTracksPool[posneg].emplace_back(TrackCand{trc, tvid, {iv, iv}, r, hasTPC, nITSclu, compatibleWithProton});
555 if (tvid.getSource() == GIndex::TPC) { // constrained TPC track?
556 correctTPCTrack(mTracksPool[posneg].back(), mTPCTracksArray[tvid], -1, -1);
557 }
558 if (tvid.isAmbiguous()) { // track attached to >1 vertex, remember that it was already processed
559 tmap[tvid] = {mTracksPool[posneg].size() - 1, posneg};
560 }
561 }
562 }
563 // register 1st track of each charge for each vertex
564 for (int pn = 0; pn < 2; pn++) {
565 auto& vtxFirstT = mVtxFirstTrack[pn];
566 const auto& tracksPool = mTracksPool[pn];
567 for (unsigned i = 0; i < tracksPool.size(); i++) {
568 const auto& t = tracksPool[i];
569 for (int j{t.vBracket.getMin()}; j <= t.vBracket.getMax(); ++j) {
570 if (vtxFirstT[j] == -1) {
571 vtxFirstT[j] = i;
572 }
573 }
574 }
575 }
576
577 LOG(info) << "Collected " << mTracksPool[POS].size() << " positive and " << mTracksPool[NEG].size() << " negative seeds";
578}
579
580//__________________________________________________________________
581bool SVertexer::checkV0(const TrackCand& seedP, const TrackCand& seedN, int iP, int iN, int ithread)
582{
583 auto& fitterV0 = mFitterV0[ithread];
584 // Fast rough cuts on pairs before feeding to DCAFitter, tracks are not in the same Frame or at same X
585 bool isTPConly = (seedP.gid.getSource() == GIndex::TPC || seedN.gid.getSource() == GIndex::TPC);
586 if (mSVParams->mTPCTrackPhotonTune && isTPConly) {
587 // Check if Tgl is close enough
588 if (std::abs(seedP.getTgl() - seedN.getTgl()) > mSVParams->maxV0TglAbsDiff) {
589 LOG(debug) << "RejTgl";
590 return false;
591 }
592 // Check in transverse plane
593 float sna, csa;
594 o2::math_utils::CircleXYf_t trkPosCircle;
595 seedP.getCircleParams(mBz, trkPosCircle, sna, csa);
596 o2::math_utils::CircleXYf_t trkEleCircle;
597 seedN.getCircleParams(mBz, trkEleCircle, sna, csa);
598 // Does the radius of both tracks compare to their absolute circle center distance
599 float c2c = std::hypot(trkPosCircle.xC - trkEleCircle.xC,
600 trkPosCircle.yC - trkEleCircle.yC);
601 float r2r = trkPosCircle.rC + trkEleCircle.rC;
602 float dcr = c2c - r2r;
603 if (std::abs(dcr) > mSVParams->mTPCTrackD2R) {
604 LOG(debug) << "RejD2R " << c2c << " " << r2r << " " << dcr;
605 return false;
606 }
607 // Will the conversion point look somewhat reasonable
608 float r1_r = trkPosCircle.rC / r2r;
609 float r2_r = trkEleCircle.rC / r2r;
610 float dR = std::hypot(r2_r * trkPosCircle.xC + r1_r * trkEleCircle.xC, r2_r * trkPosCircle.yC + r1_r * trkEleCircle.yC);
611 if (dR > mSVParams->mTPCTrackDR) {
612 LOG(debug) << "RejDR" << dR;
613 return false;
614 }
615
616 // Setup looser cuts for the DCAFitter
617 fitterV0.setMaxDZIni(mSVParams->mTPCTrackMaxDZIni);
618 fitterV0.setMaxDXYIni(mSVParams->mTPCTrackMaxDXYIni);
619 fitterV0.setMaxChi2(mSVParams->mTPCTrackMaxChi2);
620 fitterV0.setCollinear(true);
621 }
622
623 // feed DCAFitter
624 int nCand = fitterV0.process(seedP, seedN);
625 if (mSVParams->mTPCTrackPhotonTune && isTPConly) {
626 // Reset immediately to the defaults
627 fitterV0.setMaxDZIni(mSVParams->maxDZIni);
628 fitterV0.setMaxDXYIni(mSVParams->maxDXYIni);
629 fitterV0.setMaxChi2(mSVParams->maxChi2);
630 fitterV0.setCollinear(false);
631 }
632
633 if (nCand == 0) { // discard this pair
634 LOG(debug) << "RejDCAFitter";
635 return false;
636 }
637 const auto& v0XYZ = fitterV0.getPCACandidate();
638 // validate V0 radial position
639 // check closeness to the beam-line
640 float dxv0 = v0XYZ[0] - mMeanVertex.getX(), dyv0 = v0XYZ[1] - mMeanVertex.getY(), r2v0 = dxv0 * dxv0 + dyv0 * dyv0;
641 if (r2v0 < mMinR2ToMeanVertex) {
642 LOG(debug) << "RejMinR2ToMeanVertex";
643 return false;
644 }
645 float rv0 = std::sqrt(r2v0), drv0P = rv0 - seedP.minR, drv0N = rv0 - seedN.minR;
646 if (drv0P > mSVParams->causalityRTolerance || drv0P < -mSVParams->maxV0ToProngsRDiff ||
647 drv0N > mSVParams->causalityRTolerance || drv0N < -mSVParams->maxV0ToProngsRDiff) {
648 LOG(debug) << "RejCausality " << drv0P << " " << drv0N;
649 return false;
650 }
651 const int cand = 0;
652 if (!fitterV0.isPropagateTracksToVertexDone(cand) && !fitterV0.propagateTracksToVertex(cand)) {
653 LOG(debug) << "RejProp failed";
654 return false;
655 }
656 const auto& trPProp = fitterV0.getTrack(0, cand);
657 const auto& trNProp = fitterV0.getTrack(1, cand);
658 std::array<float, 3> pP{}, pN{};
659 trPProp.getPxPyPzGlo(pP);
660 trNProp.getPxPyPzGlo(pN);
661 // estimate DCA of neutral V0 track to beamline: straight line with parametric equation
662 // x = X0 + pV0[0]*t, y = Y0 + pV0[1]*t reaches DCA to beamline (Xv, Yv) at
663 // t = -[ (x0-Xv)*pV0[0] + (y0-Yv)*pV0[1]) ] / ( pT(pV0)^2 )
664 // Similar equation for 3D distance involving pV0[2]
665 std::array<float, 3> pV0 = {pP[0] + pN[0], pP[1] + pN[1], pP[2] + pN[2]};
666 float pt2V0 = pV0[0] * pV0[0] + pV0[1] * pV0[1], prodXYv0 = dxv0 * pV0[0] + dyv0 * pV0[1], tDCAXY = prodXYv0 / pt2V0;
667 if (pt2V0 < mMinPt2V0) { // pt cut
668 LOG(debug) << "RejPt2 " << pt2V0;
669 return false;
670 }
671 if (pV0[2] * pV0[2] / pt2V0 > mMaxTgl2V0) { // tgLambda cut
672 LOG(debug) << "RejTgL " << pV0[2] * pV0[2] / pt2V0;
673 return false;
674 }
675 float p2V0 = pt2V0 + pV0[2] * pV0[2], ptV0 = std::sqrt(pt2V0);
676 // apply mass selections
677 float p2Pos = pP[0] * pP[0] + pP[1] * pP[1] + pP[2] * pP[2], p2Neg = pN[0] * pN[0] + pN[1] * pN[1] + pN[2] * pN[2];
678
679 bool goodHyp = false, photonOnly = mSVParams->mTPCTrackPhotonTune && isTPConly;
680 std::array<bool, NHypV0> hypCheckStatus{};
681 int nPID = photonOnly ? (Photon + 1) : NHypV0;
682 for (int ipid = 0; (ipid < nPID) && mSVParams->checkV0Hypothesis; ipid++) {
683 if (mV0Hyps[ipid].check(p2Pos, p2Neg, p2V0, ptV0)) {
684 goodHyp = hypCheckStatus[ipid] = true;
685 }
686 }
687 // check tight lambda mass only
688 bool goodLamForCascade = false, goodALamForCascade = false;
689 bool usesTPCOnly = (seedP.hasTPC && !seedP.hasITS()) || (seedN.hasTPC && !seedN.hasITS());
690 bool usesShortITSOnly = (!seedP.hasTPC && seedP.nITSclu < mSVParams->mITSSAminNcluCascades) || (!seedN.hasTPC && seedN.nITSclu < mSVParams->mITSSAminNcluCascades);
691 if (ptV0 > mSVParams->minPtV0FromCascade && (!mSVParams->mSkipTPCOnlyCascade || !usesTPCOnly) && !usesShortITSOnly) {
692 if (mV0Hyps[Lambda].checkTight(p2Pos, p2Neg, p2V0, ptV0) && (!mSVParams->mRequireTPCforCascBaryons || seedP.hasTPC) && seedP.compatibleProton) {
693 goodLamForCascade = true;
694 }
695 if (mV0Hyps[AntiLambda].checkTight(p2Pos, p2Neg, p2V0, ptV0) && (!mSVParams->mRequireTPCforCascBaryons || seedN.hasTPC) && seedN.compatibleProton) {
696 goodALamForCascade = true;
697 }
698 }
699
700 // apply mass selections for 3-body decay
701 bool good3bodyV0Hyp = false;
702 for (int ipid = 2; ipid < 4; ipid++) {
703 float massForLambdaHyp = mV0Hyps[ipid].calcMass(p2Pos, p2Neg, p2V0);
704 if (massForLambdaHyp - mV0Hyps[ipid].getMassV0Hyp() < mV0Hyps[ipid].getMargin(ptV0)) {
705 good3bodyV0Hyp = true;
706 break;
707 }
708 }
709
710 // we want to reconstruct the 3 body decay of hypernuclei starting from the V0 of a proton and a pion (e.g. H3L->d + (p + pi-), or He4L->He3 + (p + pi-)))
711 bool checkFor3BodyDecays = mEnable3BodyDecays &&
712 (!mSVParams->checkV0Hypothesis || good3bodyV0Hyp) &&
713 (pt2V0 > 0.5) &&
714 (!mSVParams->mSkipTPCOnly3Body || !isTPConly);
715 bool rejectAfter3BodyCheck = false; // To reject v0s which can be 3-body decay candidates but not cascade or v0
716 bool checkForCascade = mEnableCascades &&
717 (!mSVParams->mSkipTPCOnlyCascade || !usesTPCOnly) &&
718 r2v0 < mMaxR2ToMeanVertexCascV0 &&
719 (!mSVParams->checkV0Hypothesis || (goodLamForCascade || goodALamForCascade) &&
720 (!isTPConly || !hypCheckStatus[HypV0::Photon]));
721 bool rejectIfNotCascade = false;
722
723 if (!goodHyp && mSVParams->checkV0Hypothesis) {
724 LOG(debug) << "RejHypo";
725 if (!checkFor3BodyDecays && !checkForCascade) {
726 return false;
727 } else {
728 rejectAfter3BodyCheck = true;
729 }
730 }
731
732 float dcaX = dxv0 - pV0[0] * tDCAXY, dcaY = dyv0 - pV0[1] * tDCAXY, dca2 = dcaX * dcaX + dcaY * dcaY;
733 float cosPAXY = prodXYv0 / std::sqrt(r2v0 * pt2V0);
734
735 if (checkForCascade) { // use looser cuts for cascade v0 candidates
736 if (dca2 > mMaxDCAXY2ToMeanVertexV0Casc || cosPAXY < mSVParams->minCosPAXYMeanVertexCascV0) {
737 LOG(debug) << "Rej for cascade DCAXY2: " << dca2 << " << cosPAXY: " << cosPAXY;
738 if (!checkFor3BodyDecays) {
739 return false;
740 } else {
741 rejectAfter3BodyCheck = true;
742 }
743 }
744 }
745 if (checkFor3BodyDecays) { // use looser cuts for 3-body decay candidates
746 if (dca2 > mMaxDCAXY2ToMeanVertex3bodyV0 || cosPAXY < mSVParams->minCosPAXYMeanVertex3bodyV0) {
747 LOG(debug) << "Rej for 3 body decays DCAXY2: " << dca2 << " << cosPAXY: " << cosPAXY;
748 checkFor3BodyDecays = false;
749 }
750 }
751
752 if (dca2 > mMaxDCAXY2ToMeanVertex || cosPAXY < mSVParams->minCosPAXYMeanVertex) {
753 if (checkForCascade) {
754 rejectIfNotCascade = true;
755 } else if (checkFor3BodyDecays) {
756 rejectAfter3BodyCheck = true;
757 } else {
758 if (mSVParams->mTPCTrackPhotonTune && isTPConly) {
759 // Check for looser cut for tpc-only photons only
760 if (dca2 > mSVParams->mTPCTrackMaxDCAXY2ToMeanVertex) {
761 return false;
762 }
763 } else {
764 return false;
765 }
766 }
767 }
768
769 auto vlist = seedP.vBracket.getOverlap(seedN.vBracket); // indices of vertices shared by both seeds
770 bool candFound = false;
771 auto bestCosPA = checkForCascade ? mSVParams->minCosPACascV0 : mSVParams->minCosPA;
772 bestCosPA = checkFor3BodyDecays ? std::min(mSVParams->minCosPA3bodyV0, bestCosPA) : bestCosPA;
773 V0 v0new;
774 V0Index v0Idxnew;
775
776 for (int iv = vlist.getMin(); iv <= vlist.getMax(); iv++) {
777 const auto& pv = mPVertices[iv];
778 const auto v0XYZ = fitterV0.getPCACandidatePos(cand);
779 // check cos of pointing angle
780 float dx = v0XYZ[0] - pv.getX(), dy = v0XYZ[1] - pv.getY(), dz = v0XYZ[2] - pv.getZ(), prodXYZv0 = dx * pV0[0] + dy * pV0[1] + dz * pV0[2];
781 float cosPA = prodXYZv0 / std::sqrt((dx * dx + dy * dy + dz * dz) * p2V0);
782 if (cosPA < bestCosPA) {
783 LOG(debug) << "Rej. cosPA: " << cosPA;
784 continue;
785 }
786 if (!candFound) {
787 new (&v0new) V0(v0XYZ, pV0, fitterV0.calcPCACovMatrixFlat(cand), trPProp, trNProp);
788 new (&v0Idxnew) V0Index(-1, seedP.gid, seedN.gid);
789 v0new.setDCA(fitterV0.getChi2AtPCACandidate(cand));
790 candFound = true;
791 }
792 v0new.setCosPA(cosPA);
793 v0Idxnew.setVertexID(iv);
794 bestCosPA = cosPA;
795 }
796 if (!candFound) {
797 return false;
798 }
799 if (bestCosPA < mSVParams->minCosPACascV0) {
800 rejectAfter3BodyCheck = true;
801 }
802 if (bestCosPA < mSVParams->minCosPA && checkForCascade) {
803 rejectIfNotCascade = true;
804 }
805 int nV0Ini = mV0sIdxTmp[ithread].size();
806 // check 3 body decays
807 if (checkFor3BodyDecays) {
808 int n3bodyDecays = 0;
809 n3bodyDecays += check3bodyDecays(v0Idxnew, v0new, rv0, pV0, p2V0, iN, NEG, vlist, ithread);
810 n3bodyDecays += check3bodyDecays(v0Idxnew, v0new, rv0, pV0, p2V0, iP, POS, vlist, ithread);
811 }
812 if (rejectAfter3BodyCheck) {
813 return false;
814 }
815
816 // check cascades
817 int nCascIni = mCascadesIdxTmp[ithread].size(), nV0Used = 0; // number of times this particular v0 (with assigned PV) was used (not counting using its clones with other PV)
818 if (checkForCascade) {
819 if (goodLamForCascade || !mSVParams->checkCascadeHypothesis) {
820 nV0Used += checkCascades(v0Idxnew, v0new, rv0, pV0, p2V0, iN, NEG, vlist, ithread);
821 }
822 if (goodALamForCascade || !mSVParams->checkCascadeHypothesis) {
823 nV0Used += checkCascades(v0Idxnew, v0new, rv0, pV0, p2V0, iP, POS, vlist, ithread);
824 }
825 }
826
827 if (nV0Used) { // need to fix the index of V0 for the cascades using this v0
828 for (unsigned int ic = nCascIni; ic < mCascadesIdxTmp[ithread].size(); ic++) {
829 if (mCascadesIdxTmp[ithread][ic].getV0ID() == -1) {
830 mCascadesIdxTmp[ithread][ic].setV0ID(nV0Ini);
831 }
832 }
833 }
834
835 if (nV0Used || !rejectIfNotCascade) { // need to add this v0
836 mV0sIdxTmp[ithread].push_back(v0Idxnew);
837 if (!rejectIfNotCascade) {
838 mV0sIdxTmp[ithread].back().setStandaloneV0();
839 }
840 if (photonOnly) {
841 mV0sIdxTmp[ithread].back().setPhotonOnly();
842 mV0sIdxTmp[ithread].back().setCollinear();
843 }
844
845 if (mSVParams->createFullV0s) {
846 mV0sTmp[ithread].push_back(v0new);
847 }
848 }
849
850 if (mStrTracker) {
851 for (int iv = nV0Ini; iv < (int)mV0sIdxTmp[ithread].size(); iv++) {
852 mStrTracker->processV0(iv, v0new, v0Idxnew, ithread);
853 }
854 }
855
856 return mV0sIdxTmp[ithread].size() - nV0Ini != 0;
857}
858
859//__________________________________________________________________
860int SVertexer::checkCascades(const V0Index& v0Idx, const V0& v0, float rv0, std::array<float, 3> pV0, float p2V0, int avoidTrackID, int posneg, VBracket v0vlist, int ithread)
861{
862 // check last added V0 for belonging to cascade
863 auto& fitterCasc = mFitterCasc[ithread];
864 auto& tracks = mTracksPool[posneg];
865 int nCascIni = mCascadesIdxTmp[ithread].size(), nv0use = 0;
866
867 // check if a given PV has already been used in a cascade
868 std::unordered_map<int, int> pvMap;
869
870 // start from the 1st bachelor track compatible with earliest vertex in the v0vlist
871 int firstTr = mVtxFirstTrack[posneg][v0vlist.getMin()], nTr = tracks.size();
872 if (firstTr < 0) {
873 firstTr = nTr;
874 }
875 for (int it = firstTr; it < nTr; it++) {
876 if (it == avoidTrackID) {
877 continue; // skip the track used by V0
878 }
879 auto& bach = tracks[it];
880 if (mSVParams->mSkipTPCOnlyCascade && (bach.gid.getSource() == GIndex::TPC)) {
881 continue; // reject TPC-only bachelors
882 }
883 if (!bach.hasTPC && bach.nITSclu < mSVParams->mITSSAminNcluCascades) {
884 continue; // reject short ITS-only
885 }
886
887 if (bach.vBracket.getMin() > v0vlist.getMax()) {
888 LOG(debug) << "Skipping";
889 break; // all other bachelor candidates will be also not compatible with this PV
890 }
891 auto cascVlist = v0vlist.getOverlap(bach.vBracket); // indices of vertices shared by V0 and bachelor
892 if (mSVParams->selectBestV0) {
893 // select only the best V0 candidate among the compatible ones
894 if (v0Idx.getVertexID() < cascVlist.getMin() || v0Idx.getVertexID() > cascVlist.getMax()) {
895 continue;
896 }
897 cascVlist.setMin(v0Idx.getVertexID());
898 cascVlist.setMax(v0Idx.getVertexID());
899 }
900
901 int nCandC = fitterCasc.process(v0, bach);
902 if (nCandC == 0) { // discard this pair
903 continue;
904 }
905 const int candC = 0;
906 const auto& cascXYZ = fitterCasc.getPCACandidatePos(candC);
907
908 // make sure the cascade radius is smaller than that of the mean vertex
909 float dxc = cascXYZ[0] - mMeanVertex.getX(), dyc = cascXYZ[1] - mMeanVertex.getY(), r2casc = dxc * dxc + dyc * dyc;
910 if (rv0 * rv0 - r2casc < mMinR2DiffV0Casc || r2casc < mMinR2ToMeanVertex) {
911 continue;
912 }
913 // do we want to apply mass cut ?
914 //
915 if (!fitterCasc.isPropagateTracksToVertexDone(candC) && !fitterCasc.propagateTracksToVertex(candC)) {
916 continue;
917 }
918
919 auto& trNeut = fitterCasc.getTrack(0, candC);
920 auto& trBach = fitterCasc.getTrack(1, candC);
921 trNeut.setPID(o2::track::PID::Lambda);
922 trBach.setPID(o2::track::PID::Pion);
923 std::array<float, 3> pNeut, pBach;
924 trNeut.getPxPyPzGlo(pNeut);
925 trBach.getPxPyPzGlo(pBach);
926 std::array<float, 3> pCasc = {pNeut[0] + pBach[0], pNeut[1] + pBach[1], pNeut[2] + pBach[2]};
927
928 float pt2Casc = pCasc[0] * pCasc[0] + pCasc[1] * pCasc[1], p2Casc = pt2Casc + pCasc[2] * pCasc[2];
929 if (pt2Casc < mMinPt2Casc) { // pt cut
930 LOG(debug) << "Casc pt too low";
931 continue;
932 }
933 if (pCasc[2] * pCasc[2] / pt2Casc > mMaxTgl2Casc) { // tgLambda cut
934 LOG(debug) << "Casc tgLambda too high";
935 continue;
936 }
937
938 // compute primary vertex and cosPA of the cascade
939 auto bestCosPA = mSVParams->minCosPACasc;
940 auto cascVtxID = -1;
941
942 for (int iv = cascVlist.getMin(); iv <= cascVlist.getMax(); iv++) {
943 const auto& pv = mPVertices[iv];
944 // check cos of pointing angle
945 float dx = cascXYZ[0] - pv.getX(), dy = cascXYZ[1] - pv.getY(), dz = cascXYZ[2] - pv.getZ(), prodXYZcasc = dx * pCasc[0] + dy * pCasc[1] + dz * pCasc[2];
946 float cosPA = prodXYZcasc / std::sqrt((dx * dx + dy * dy + dz * dz) * p2Casc);
947 if (cosPA < bestCosPA) {
948 LOG(debug) << "Rej. cosPA: " << cosPA;
949 continue;
950 }
951 cascVtxID = iv;
952 bestCosPA = cosPA;
953 }
954 if (cascVtxID == -1) {
955 LOG(debug) << "Casc not compatible with any vertex";
956 continue;
957 }
958
959 const auto& cascPv = mPVertices[cascVtxID];
960 float dxCasc = cascXYZ[0] - cascPv.getX(), dyCasc = cascXYZ[1] - cascPv.getY(), dzCasc = cascXYZ[2] - cascPv.getZ();
961 auto prodPPos = pV0[0] * dxCasc + pV0[1] * dyCasc + pV0[2] * dzCasc;
962 if (prodPPos < 0.) { // causality cut
963 LOG(debug) << "Casc not causally compatible";
964 continue;
965 }
966
967 float p2Bach = pBach[0] * pBach[0] + pBach[1] * pBach[1] + pBach[2] * pBach[2];
968 float ptCasc = std::sqrt(pt2Casc);
969 bool goodHyp = false;
970 for (int ipid = 0; ipid < NHypCascade; ipid++) {
971 if (mCascHyps[ipid].check(p2V0, p2Bach, p2Casc, ptCasc)) {
972 goodHyp = true;
973 break;
974 }
975 }
976 if (!goodHyp) {
977 LOG(debug) << "Casc not compatible with any hypothesis";
978 continue;
979 }
980 // note: at the moment the v0 was not added yet. If some cascade will use v0 (with its PV), the v0 will be added after checkCascades
981 // but not necessarily at the and of current v0s vector, since meanwhile checkCascades may add v0 clones (with PV redefined).
982 Cascade casc(cascXYZ, pCasc, fitterCasc.calcPCACovMatrixFlat(candC), trNeut, trBach);
983 o2::track::TrackParCov trc = casc;
985 if (!trc.propagateToDCA(cascPv, fitterCasc.getBz(), &dca, 5.) ||
986 std::abs(dca.getY()) > mSVParams->maxDCAXYCasc || std::abs(dca.getZ()) > mSVParams->maxDCAZCasc) {
987 LOG(debug) << "Casc not compatible with PV";
988 LOG(debug) << "DCA: " << dca.getY() << " " << dca.getZ();
989 continue;
990 }
991 CascadeIndex cascIdx(cascVtxID, -1, bach.gid); // the v0Idx was not yet added, this will be done after the checkCascades
992
993 LOGP(debug, "cascade successfully validated");
994
995 // clone the V0, set new cosPA and VerteXID, add it to the list of V0s
996 if (cascVtxID != v0Idx.getVertexID()) {
997 auto pvIdx = pvMap.find(cascVtxID);
998 if (pvIdx != pvMap.end()) {
999 cascIdx.setV0ID(pvIdx->second); // V0 already exists, add reference to the cascade
1000 } else { // add V0 clone for this cascade (may be used also by other cascades)
1001 const auto& pv = mPVertices[cascVtxID];
1002 cascIdx.setV0ID(mV0sIdxTmp[ithread].size()); // set the new V0 index in the cascade
1003 pvMap[cascVtxID] = mV0sTmp[ithread].size(); // add the new V0 index to the map
1004 mV0sIdxTmp[ithread].emplace_back(cascVtxID, v0Idx.getProngs());
1005 if (mSVParams->createFullV0s) {
1006 mV0sTmp[ithread].push_back(v0);
1007 float dx = v0.getX() - pv.getX(), dy = v0.getY() - pv.getY(), dz = v0.getZ() - pv.getZ(), prodXYZ = dx * pV0[0] + dy * pV0[1] + dz * pV0[2];
1008 mV0sTmp[ithread].back().setCosPA(prodXYZ / std::sqrt((dx * dx + dy * dy + dz * dz) * p2V0));
1009 }
1010 }
1011 } else {
1012 nv0use++; // original v0 was used
1013 }
1014 mCascadesIdxTmp[ithread].push_back(cascIdx);
1015 if (mSVParams->createFullCascades) {
1016 casc.setCosPA(bestCosPA);
1017 casc.setDCA(fitterCasc.getChi2AtPCACandidate(candC));
1018 mCascadesTmp[ithread].push_back(casc);
1019 }
1020 if (mStrTracker) {
1021 mStrTracker->processCascade(mCascadesIdxTmp[ithread].size() - 1, casc, cascIdx, v0, ithread);
1022 }
1023 }
1024
1025 return nv0use;
1026}
1027
1028//__________________________________________________________________
1029int SVertexer::check3bodyDecays(const V0Index& v0Idx, const V0& v0, float rv0, std::array<float, 3> pV0, float p2V0, int avoidTrackID, int posneg, VBracket v0vlist, int ithread)
1030{
1031 // check last added V0 for belonging to cascade
1032 auto& fitter3body = mFitter3body[ithread];
1033 auto& tracks = mTracksPool[posneg];
1034 int n3BodyIni = m3bodyIdxTmp[ithread].size();
1035
1036 // start from the 1st bachelor track compatible with earliest vertex in the v0vlist
1037 int firstTr = mVtxFirstTrack[posneg][v0vlist.getMin()], nTr = tracks.size();
1038 if (firstTr < 0) {
1039 firstTr = nTr;
1040 }
1041
1042 // If the V0 is a pair of proton and pion, we should pair it with all positive particles, and the positive particle in the V0 is a proton.
1043 // Otherwise, we should pair it with all negative particles, and the negative particle in the V0 is a antiproton.
1044
1045 // start from the 1st track compatible with V0's primary vertex
1046 for (int it = firstTr; it < nTr; it++) {
1047 if (it == avoidTrackID) {
1048 continue; // skip the track used by V0
1049 }
1050 auto& bach = tracks[it];
1051 if (mSVParams->mSkipTPCOnly3Body && (bach.gid.getSource() == GIndex::TPC)) {
1052 continue; // reject TPC-only bachelors
1053 }
1054 if (bach.vBracket > v0vlist.getMax()) {
1055 LOG(debug) << "Skipping";
1056 break; // all other bachelor candidates will be also not compatible with this PV
1057 }
1058 auto decay3bodyVlist = v0vlist.getOverlap(bach.vBracket); // indices of vertices shared by V0 and bachelor
1059 if (mSVParams->selectBestV0) {
1060 // select only the best V0 candidate among the compatible ones
1061 if (v0Idx.getVertexID() < decay3bodyVlist.getMin() || v0Idx.getVertexID() > decay3bodyVlist.getMax()) {
1062 continue;
1063 }
1064 decay3bodyVlist.setMin(v0Idx.getVertexID());
1065 decay3bodyVlist.setMax(v0Idx.getVertexID());
1066 }
1067
1068 if (bach.getPt() < 0.6) {
1069 continue;
1070 }
1071
1072 int n3bodyVtx = fitter3body.process(v0.getProng(0), v0.getProng(1), bach);
1073 if (n3bodyVtx == 0) { // discard this pair
1074 continue;
1075 }
1076 int cand3B = 0;
1077 const auto& vertexXYZ = fitter3body.getPCACandidatePos(cand3B);
1078
1079 // make sure the 3 body vertex radius is close to that of the mean vertex
1080 float dxc = vertexXYZ[0] - mMeanVertex.getX(), dyc = vertexXYZ[1] - mMeanVertex.getY(), dzc = vertexXYZ[2] - mMeanVertex.getZ(), r2vertex = dxc * dxc + dyc * dyc;
1081 if (std::abs(rv0 - std::sqrt(r2vertex)) > mSVParams->maxRDiffV03body || r2vertex < mMinR2ToMeanVertex) {
1082 continue;
1083 }
1084 float drvtxBach = std::sqrt(r2vertex) - bach.minR;
1085 if (drvtxBach > mSVParams->causalityRTolerance || drvtxBach < -mSVParams->maxV0ToProngsRDiff) {
1086 LOG(debug) << "RejCausality " << drvtxBach;
1087 }
1088 //
1089 if (!fitter3body.isPropagateTracksToVertexDone() && !fitter3body.propagateTracksToVertex()) {
1090 continue;
1091 }
1092
1093 auto& tr0 = fitter3body.getTrack(0, cand3B);
1094 auto& tr1 = fitter3body.getTrack(1, cand3B);
1095 auto& tr2 = fitter3body.getTrack(2, cand3B);
1096 std::array<float, 3> p0, p1, p2;
1097 tr0.getPxPyPzGlo(p0);
1098 tr1.getPxPyPzGlo(p1);
1099 tr2.getPxPyPzGlo(p2);
1100
1101 bool goodHyp = false;
1102 o2::track::PID pidHyp = o2::track::PID::Electron; // Update if goodHyp is true
1103 auto decay3bodyVtxID = -1;
1104 auto vtxCosPA = -1;
1105
1106 std::array<float, 3> pbach = {0, 0, 0}, p3B = {0, 0, 0}; // Update during the check of invariant mass
1107 for (int ipid = 0; ipid < NHyp3body; ipid++) {
1108 // check mass based on hypothesis of charge of bachelor (pos and neg expected to be proton/pion)
1109 float bachChargeFactor = m3bodyHyps[ipid].getChargeBachProng() / tr2.getAbsCharge();
1110 pbach = {bachChargeFactor * p2[0], bachChargeFactor * p2[1], bachChargeFactor * p2[2]};
1111 p3B = {p0[0] + p1[0] + pbach[0], p0[1] + p1[1] + pbach[1], p0[2] + p1[2] + pbach[2]};
1112 float sqP0 = p0[0] * p0[0] + p0[1] * p0[1] + p0[2] * p0[2], sqP1 = p1[0] * p1[0] + p1[1] * p1[1] + p1[2] * p1[2], sqPBach = pbach[0] * pbach[0] + pbach[1] * pbach[1] + pbach[2] * pbach[2];
1113 float pt2Candidate = p3B[0] * p3B[0] + p3B[1] * p3B[1], p2Candidate = pt2Candidate + p3B[2] * p3B[2];
1114 float ptCandidate = std::sqrt(pt2Candidate);
1115 if (m3bodyHyps[ipid].check(sqP0, sqP1, sqPBach, p2Candidate, ptCandidate)) {
1116 if (pt2Candidate < mMinPt23Body) { // pt cut
1117 continue;
1118 }
1119 if (p3B[2] * p3B[2] > pt2Candidate * mMaxTgl23Body) { // tgLambda cut
1120 continue;
1121 }
1122
1123 // compute primary vertex and cosPA of the 3-body decay
1124 auto bestCosPA = mSVParams->minCosPA3body;
1125 for (int iv = decay3bodyVlist.getMin(); iv <= decay3bodyVlist.getMax(); iv++) {
1126 const auto& pv = mPVertices[iv];
1127 // check cos of pointing angle
1128 float dx = vertexXYZ[0] - pv.getX(), dy = vertexXYZ[1] - pv.getY(), dz = vertexXYZ[2] - pv.getZ(), prodXYZ3body = dx * p3B[0] + dy * p3B[1] + dz * p3B[2];
1129 float cosPA = prodXYZ3body / std::sqrt((dx * dx + dy * dy + dz * dz) * p2Candidate);
1130 if (cosPA < bestCosPA) {
1131 LOG(debug) << "Rej. cosPA: " << cosPA;
1132 continue;
1133 }
1134 decay3bodyVtxID = iv;
1135 bestCosPA = cosPA;
1136 }
1137 if (decay3bodyVtxID == -1) {
1138 LOG(debug) << "3-body decay not compatible with any vertex";
1139 continue;
1140 }
1141
1142 goodHyp = true;
1143 pidHyp = m3bodyHyps[ipid].getPIDHyp();
1144 vtxCosPA = bestCosPA;
1145 break;
1146 }
1147 }
1148 if (!goodHyp) {
1149 continue;
1150 }
1151
1152 const auto& decay3bodyPv = mPVertices[decay3bodyVtxID];
1153 Decay3Body candidate3B(vertexXYZ, p3B, fitter3body.calcPCACovMatrixFlat(cand3B), tr0, tr1, tr2, pidHyp);
1154 o2::track::TrackParCov trc = candidate3B;
1156 if (!trc.propagateToDCA(decay3bodyPv, fitter3body.getBz(), &dca, 5.) ||
1157 std::abs(dca.getY()) > mSVParams->maxDCAXY3Body || std::abs(dca.getZ()) > mSVParams->maxDCAZ3Body) {
1158 continue;
1159 }
1160 if (mSVParams->createFull3Bodies) {
1161 candidate3B.setCosPA(vtxCosPA);
1162 candidate3B.setDCA(fitter3body.getChi2AtPCACandidate());
1163 m3bodyTmp[ithread].push_back(candidate3B);
1164 }
1165 m3bodyIdxTmp[ithread].emplace_back(decay3bodyVtxID, v0Idx.getProngID(0), v0Idx.getProngID(1), bach.gid);
1166
1167 Decay3BodyIndex decay3bodyIdx(decay3bodyVtxID, v0Idx.getProngID(0), v0Idx.getProngID(1), bach.gid);
1168 if (mStrTracker) {
1169 mStrTracker->process3Body(m3bodyIdxTmp[ithread].size() - 1, candidate3B, decay3bodyIdx, ithread);
1170 }
1171 }
1172 return m3bodyIdxTmp[ithread].size() - n3BodyIni;
1173}
1174
1175//__________________________________________________________________
1176template <class TVI, class TCI, class T3I, class TR>
1177void SVertexer::extractPVReferences(const TVI& v0s, TR& vtx2V0Refs, const TCI& cascades, TR& vtx2CascRefs, const T3I& vtx3, TR& vtx2body3Refs)
1178{
1179 // V0s, cascades and 3bodies are already sorted in PV ID
1180 vtx2V0Refs.clear();
1181 vtx2V0Refs.resize(mPVertices.size());
1182 vtx2CascRefs.clear();
1183 vtx2CascRefs.resize(mPVertices.size());
1184 vtx2body3Refs.clear();
1185 vtx2body3Refs.resize(mPVertices.size());
1186 int nv0 = v0s.size(), nCasc = cascades.size(), n3body = vtx3.size();
1187
1188 // relate V0s to primary vertices
1189 int pvID = -1, nForPV = 0;
1190 for (int iv = 0; iv < nv0; iv++) {
1191 if (pvID < v0s[iv].getVertexID()) {
1192 if (pvID > -1) {
1193 vtx2V0Refs[pvID].setEntries(nForPV);
1194 }
1195 pvID = v0s[iv].getVertexID();
1196 vtx2V0Refs[pvID].setFirstEntry(iv);
1197 nForPV = 0;
1198 }
1199 nForPV++;
1200 }
1201 if (pvID != -1) { // finalize
1202 vtx2V0Refs[pvID].setEntries(nForPV);
1203 // fill empty slots
1204 int ent = nv0;
1205 for (int ip = vtx2V0Refs.size(); ip--;) {
1206 if (vtx2V0Refs[ip].getEntries()) {
1207 ent = vtx2V0Refs[ip].getFirstEntry();
1208 } else {
1209 vtx2V0Refs[ip].setFirstEntry(ent);
1210 }
1211 }
1212 }
1213
1214 // relate Cascades to primary vertices
1215 pvID = -1;
1216 nForPV = 0;
1217 for (int iv = 0; iv < nCasc; iv++) {
1218 if (pvID < cascades[iv].getVertexID()) {
1219 if (pvID > -1) {
1220 vtx2CascRefs[pvID].setEntries(nForPV);
1221 }
1222 pvID = cascades[iv].getVertexID();
1223 vtx2CascRefs[pvID].setFirstEntry(iv);
1224 nForPV = 0;
1225 }
1226 nForPV++;
1227 }
1228 if (pvID != -1) { // finalize
1229 vtx2CascRefs[pvID].setEntries(nForPV);
1230 // fill empty slots
1231 int ent = nCasc;
1232 for (int ip = vtx2CascRefs.size(); ip--;) {
1233 if (vtx2CascRefs[ip].getEntries()) {
1234 ent = vtx2CascRefs[ip].getFirstEntry();
1235 } else {
1236 vtx2CascRefs[ip].setFirstEntry(ent);
1237 }
1238 }
1239 }
1240
1241 // relate 3 body decays to primary vertices
1242 pvID = -1;
1243 nForPV = 0;
1244 for (int iv = 0; iv < n3body; iv++) {
1245 const auto& vertex3body = vtx3[iv];
1246 if (pvID < vertex3body.getVertexID()) {
1247 if (pvID > -1) {
1248 vtx2body3Refs[pvID].setEntries(nForPV);
1249 }
1250 pvID = vertex3body.getVertexID();
1251 vtx2body3Refs[pvID].setFirstEntry(iv);
1252 nForPV = 0;
1253 }
1254 nForPV++;
1255 }
1256 if (pvID != -1) { // finalize
1257 vtx2body3Refs[pvID].setEntries(nForPV);
1258 // fill empty slots
1259 int ent = n3body;
1260 for (int ip = vtx2body3Refs.size(); ip--;) {
1261 if (vtx2body3Refs[ip].getEntries()) {
1262 ent = vtx2body3Refs[ip].getFirstEntry();
1263 } else {
1264 vtx2body3Refs[ip].setFirstEntry(ent);
1265 }
1266 }
1267 }
1268}
1269
1270//__________________________________________________________________
1272{
1273#ifdef WITH_OPENMP
1274 mNThreads = n > 0 ? n : 1;
1275#else
1276 mNThreads = 1;
1277#endif
1278}
1279
1280//______________________________________________
1281bool SVertexer::processTPCTrack(const o2::tpc::TrackTPC& trTPC, GIndex gid, int vtxid)
1282{
1283 if (mSVParams->mTPCTrackMaxX > 0. && trTPC.getX() > mSVParams->mTPCTrackMaxX) {
1284 return true;
1285 }
1286 // if TPC trackis unconstrained, try to create in the tracks pool a clone constrained to vtxid vertex time.
1287 if (trTPC.hasBothSidesClusters()) { // this is effectively constrained track
1288 return false; // let it be processed as such
1289 }
1290 const auto& vtx = mPVertices[vtxid];
1291 auto twe = vtx.getTimeStamp();
1292 int posneg = trTPC.getSign() < 0 ? 1 : 0;
1293
1294 bool compatibleWithProton = false;
1295 if (!(mSVParams->mSkipTPCOnlyCascade)) {
1296 // Cascade retrieve dEdx proton frac
1297 const auto protonId = o2::track::PID::Proton;
1298 float dEdxTPC = trTPC.getdEdx().dEdxTotTPC;
1299 float dEdxExpected = mPIDresponse.getExpectedSignal(trTPC, protonId);
1300 float fracDevProton = std::abs((dEdxTPC - dEdxExpected) / dEdxExpected);
1301 if (fracDevProton < mSVParams->mFractiondEdxforCascBaryons) {
1302 compatibleWithProton = true;
1303 }
1304 }
1305
1306 auto& trLoc = mTracksPool[posneg].emplace_back(TrackCand{trTPC, gid, {vtxid, vtxid}, 0., true, -1, compatibleWithProton});
1307 auto err = correctTPCTrack(trLoc, trTPC, twe.getTimeStamp(), twe.getTimeStampError());
1308 if (err < 0) {
1309 mTracksPool[posneg].pop_back(); // discard
1310 return true;
1311 }
1312
1313 if (mSVParams->mTPCTrackPhotonTune) {
1314 // require minimum of tpc clusters
1315 bool dCls = trTPC.getNClusters() < mSVParams->mTPCTrackMinNClusters;
1316 // check track z cuts
1317 bool dDPV = std::abs(trLoc.getX() * trLoc.getTgl() - trLoc.getZ() + vtx.getZ()) > mSVParams->mTPCTrack2Beam;
1318 // check track transveres cuts
1319 float sna{0}, csa{0};
1321 trLoc.getCircleParams(mBz, trkCircle, sna, csa);
1322 float cR = std::hypot(trkCircle.xC, trkCircle.yC);
1323 float drd2 = std::sqrt(cR * cR - trkCircle.rC * trkCircle.rC);
1324 bool dRD2 = drd2 > mSVParams->mTPCTrackXY2Radius;
1325
1326 if (dCls || dDPV || dRD2) {
1327 mTracksPool[posneg].pop_back();
1328 return true;
1329 }
1330 }
1331
1332 return true;
1333}
1334
1335//______________________________________________
1336float SVertexer::correctTPCTrack(SVertexer::TrackCand& trc, const o2::tpc::TrackTPC& tTPC, float tmus, float tmusErr) const
1337{
1338 // Correct the track copy trc of the TPC track for the assumed interaction time
1339 // return extra uncertainty in Z due to the interaction time uncertainty
1340 // TODO: at the moment, apply simple shift, but with Z-dependent calibration we may
1341 // need to do corrections on TPC cluster level and refit
1342 // This is almosto clone of the MatchTPCITS::correctTPCTrack
1343
1344 float tTB, tTBErr;
1345 if (tmusErr < 0) { // use track data
1346 tTB = tTPC.getTime0();
1347 tTBErr = 0.5 * (tTPC.getDeltaTBwd() + tTPC.getDeltaTFwd());
1348 } else {
1349 tTB = tmus * mMUS2TPCBin;
1350 tTBErr = tmusErr * mMUS2TPCBin;
1351 }
1352 float dDrift = (tTB - tTPC.getTime0()) * mTPCBin2Z;
1353 float driftErr = tTBErr * mTPCBin2Z;
1354 if (driftErr < 0.) { // early return will be discarded anyway
1355 return driftErr;
1356 }
1357 // eventually should be refitted, at the moment we simply shift...
1358 trc.setZ(tTPC.getZ() + (tTPC.hasASideClustersOnly() ? dDrift : -dDrift));
1359 trc.setCov(trc.getSigmaZ2() + driftErr * driftErr, o2::track::kSigZ2);
1360 uint8_t sector, row;
1361 auto cl = &tTPC.getCluster(mTPCTrackClusIdx, tTPC.getNClusters() - 1, *mTPCClusterIdxStruct, sector, row);
1362 float x = 0, y = 0, z = 0;
1363 mTPCCorrMapsHelper->Transform(sector, row, cl->getPad(), cl->getTime(), x, y, z, tTB);
1366 }
1367 trc.minR = std::sqrt(x * x + y * y);
1368 LOGP(debug, "set MinR = {} for row {}, x:{}, y:{}, z:{}", trc.minR, row, x, y, z);
1369 return driftErr;
1370}
1371
1372//______________________________________________
1373std::array<size_t, 3> SVertexer::getNFitterCalls() const
1374{
1375 std::array<size_t, 3> calls{};
1376 for (int i = 0; i < mNThreads; i++) {
1377 calls[0] += mFitterV0[i].getCallID();
1378 calls[1] += mFitterCasc[i].getCallID();
1379 calls[2] += mFitter3body[i].getCallID();
1380 }
1381 return calls;
1382}
Helper class to access correction maps.
std::ostringstream debug
int32_t i
Some ALICE geometry constants of common interest.
Global index for barrel track: provides provenance (detectors combination), index in respective array...
constexpr int p2()
constexpr int p1()
constexpr to accelerate the coordinates changing
uint16_t pos
Definition RawData.h:3
uint32_t j
Definition RawData.h:0
Secondary vertex finder.
class to create TPC fast transformation
Reference on ITS/MFT clusters set.
calibration data from laser track calibration
Helper class to obtain TPC clusters / digits / labels from DPL.
GPUd() value_type estimateLTFast(o2 static GPUd() float estimateLTIncrement(const o2 PropagatorImpl * Instance(bool uninitialized=false)
Definition Propagator.h:178
TO BE DONE: extend to generic N body vertex.
Definition Decay3Body.h:26
const std::array< GIndex, N > & getProngs() const
decltype(auto) make(const Output &spec, Args... args)
DataAllocator & outputs()
The data allocator is used to allocate memory for the output data.
void processV0(int iv0, const V0 &v0, const V0Index &v0Idx, int iThread=0)
std::vector< StrangeTrack > & getStrangeTrackVec(int iThread=0)
void processCascade(int icasc, const Cascade &casc, const CascadeIndex &cascIdx, const V0 &cascV0, int iThread=0)
std::vector< ClusAttachments > & getClusAttachments(int iThread=0)
void process3Body(int i3body, const Decay3Body &dec3body, const Decay3BodyIndex &dec3bodyIdx, int iThread=0)
bool loadData(const o2::globaltracking::RecoContainer &recoData)
std::vector< o2::MCCompLabel > & getStrangeTrackLabels(int iThread=0)
static constexpr ID HyperHelium4
Definition PID.h:117
static constexpr ID Electron
Definition PID.h:94
static constexpr ID HyperTriton
Definition PID.h:113
static constexpr ID Lambda
Definition PID.h:112
static constexpr ID Helium3
Definition PID.h:101
static constexpr ID Kaon
Definition PID.h:97
static constexpr ID Pion
Definition PID.h:96
static constexpr ID HyperHelium5
Definition PID.h:118
static constexpr ID Deuteron
Definition PID.h:99
static constexpr ID K0
Definition PID.h:111
static constexpr ID OmegaMinus
Definition PID.h:116
static constexpr ID Photon
Definition PID.h:110
static constexpr ID Proton
Definition PID.h:98
static constexpr ID Triton
Definition PID.h:100
static constexpr ID XiMinus
Definition PID.h:115
static constexpr ID Hyperhydrog4
Definition PID.h:114
static constexpr ID Alpha
Definition PID.h:102
o2::dataformats::V0 V0
Definition SVertexer.h:61
void setTPCCorrMaps(o2::gpu::CorrectionMapsHelper *maph)
std::array< size_t, 3 > getNFitterCalls() const
void process(const o2::globaltracking::RecoContainer &recoTracks, o2::framework::ProcessingContext &pc)
Definition SVertexer.cxx:42
void setTPCVDrift(const o2::tpc::VDriftCorrFact &v)
static constexpr int POS
Definition SVertexer.h:99
static constexpr int NEG
Definition SVertexer.h:99
o2::dataformats::V0Index V0Index
Definition SVertexer.h:62
void produceOutput(o2::framework::ProcessingContext &pc)
Definition SVertexer.cxx:87
GLdouble n
Definition glcorearb.h:1982
GLint GLenum GLint x
Definition glcorearb.h:403
GLuint entry
Definition glcorearb.h:5735
GLsizeiptr size
Definition glcorearb.h:659
const GLdouble * v
Definition glcorearb.h:832
GLboolean GLboolean GLboolean b
Definition glcorearb.h:1233
GLint y
Definition glcorearb.h:270
GLfloat v0
Definition glcorearb.h:811
GLboolean r
Definition glcorearb.h:1233
GLboolean GLboolean GLboolean GLboolean a
Definition glcorearb.h:1233
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
uint8_t itsSharedClusterMap uint8_t
constexpr float XTPCInnerRef
reference radius at which TPC provides the tracks
Defining PrimaryVertex explicitly as messageable.
void check(const std::vector< std::string > &arguments, const std::vector< ConfigParamSpec > &workflowOptions, const std::vector< DeviceSpec > &deviceSpecs, CheckMatrix &matrix)
std::array< int, 24 > p0
std::vector< T, fair::mq::pmr::polymorphic_allocator< T > > vector
TrackParCovF TrackParCov
Definition Track.h:33
struct o2::upgrades_utils::@453 tracks
structure to keep trigger-related info
GTrackID getITSContributorGID(GTrackID source) const
const o2::tpc::TrackTPC & getTPCTrack(GTrackID id) const
const o2::itsmft::TrkClusRef & getITSABRef(GTrackID gid) const
GTrackID getTPCContributorGID(GTrackID source) const
const o2::track::TrackParCov & getTrackParam(GTrackID gidx) const
gsl::span< const unsigned char > clusterShMapTPC
externally set TPC clusters sharing map
const o2::its::TrackITS & getITSTrack(GTrackID gid) const
std::unique_ptr< o2::tpc::internal::getWorkflowTPCInput_ret > inputsTPCclusters
gsl::span< const unsigned int > occupancyMapTPC
externally set TPC clusters occupancy map
int maxPVContributors
max number PV contributors to allow in V0
float pidCutsHhydrog4[SVertexHypothesis::NPIDParams]
float minCosPA
min cos of PA to PV for prompt V0 candidates
float pidCutsHe5L3body[SVertex3Hypothesis::NPIDParams]
float maxDCAXYToMeanVertex
max DCA of V0 from beam line (mean vertex) for prompt V0 candidates
float pidCutsPhoton[SVertexHypothesis::NPIDParams]
bool useAbsDCA
use abs dca minimization
bool usePropagator
use external propagator
float maxDCAXYToMeanVertex3bodyV0
max DCA of V0 from beam line (mean vertex) for 3body V0 candidates
float maxRIni
don't consider as a seed (circles intersection) if its R exceeds this
bool createFullCascades
fill cascades prongs/kinematics
float maxDCAXYToMeanVertexV0Casc
max DCA of V0 from beam line (mean vertex) for cascade V0 candidates
float minParamChange
stop when tracks X-params being minimized change by less that this value
float causalityRTolerance
V0 radius cannot exceed its contributors minR by more than this value.
float maxSnp
max snp when external propagator is used
float minXSeed
minimal X of seed in prong frame (within the radial resolution track should not go to negative X)
int matCorr
material correction to use
float maxV0TglAbsDiff
max absolute difference in Tgl for V0 for photons only
float maxDZIni
don't consider as a seed (circles intersection) if Z distance exceeds this
float maxStep
max step size when external propagator is used
float pidCutsK0[SVertexHypothesis::NPIDParams]
float pidCutsHe4L3body[SVertex3Hypothesis::NPIDParams]
float minPtV0FromCascade
v0 minimum pT for v0 to be used in cascading (lowest pT Run 2 lambda: 0.4)
float mTPCTrackMaxDCAXY2ToMeanVertex
max DCA^2 of V0 from beam line (mean vertex) for prompt V0 candidates, for photon TPC-only track only
float pidCutsLambda[SVertexHypothesis::NPIDParams]
bool createFull3Bodies
fill 3-body decays prongs/kinematics
float minRelChi2Change
stop when chi2 changes by less than this value
float maxChi2
max dca from prongs to vertex
bool selectBestV0
match only the best v0 for each cascade candidate
float minRFor3DField
above this radius use 3D field
bool refitWithMatCorr
refit V0 applying material corrections
float pidCutsH3L3body[SVertex3Hypothesis::NPIDParams]
float minRDiffV0Casc
cascade should be at least this radial distance below V0
float pidCutsOmegaMinus[SVertexHypothesis::NPIDParams]
float minRToMeanVertex
min radial distance of V0 from beam line (mean vertex)
float pidCutsXiMinus[SVertexHypothesis::NPIDParams]
float maxRDiffV03body
Maximum difference between V0 and 3body radii.
float mTPCTrackMaxDXYIni
don't consider as a seed (circles intersection) if XY distance exceeds this, for photon TPC-only trac...
float minDCAToPV
min DCA to PV of single track to accept
float mTPCTrackMaxDZIni
don't consider as a seed (circles intersection) if Z distance exceeds this, for photon TPC-only track...
float pidCutsHTriton[SVertexHypothesis::NPIDParams]
float maxDXYIni
don't consider as a seed (circles intersection) if XY distance exceeds this
float mTPCTrackMaxChi2
max DCA from prongs to vertex for photon TPC-only track only
float pidCutsH4L3body[SVertex3Hypothesis::NPIDParams]
bool createFullV0s
fill V0s prongs/kinematics
float maxTglV0
maximum tgLambda of V0
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"
std::vector< int > row