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 // precalculated selection cuts
273 mMinR2ToMeanVertex = mSVParams->minRToMeanVertex * mSVParams->minRToMeanVertex;
274 mMaxR2ToMeanVertexCascV0 = mSVParams->maxRToMeanVertexCascV0 * mSVParams->maxRToMeanVertexCascV0;
275 mMaxDCAXY2ToMeanVertex = mSVParams->maxDCAXYToMeanVertex * mSVParams->maxDCAXYToMeanVertex;
276 mMaxDCAXY2ToMeanVertexV0Casc = mSVParams->maxDCAXYToMeanVertexV0Casc * mSVParams->maxDCAXYToMeanVertexV0Casc;
277 mMaxDCAXY2ToMeanVertex3bodyV0 = mSVParams->maxDCAXYToMeanVertex3bodyV0 * mSVParams->maxDCAXYToMeanVertex3bodyV0;
278 mMinR2DiffV0Casc = mSVParams->minRDiffV0Casc * mSVParams->minRDiffV0Casc;
279 mMinPt2V0 = mSVParams->minPtV0 * mSVParams->minPtV0;
280 mMaxTgl2V0 = mSVParams->maxTglV0 * mSVParams->maxTglV0;
281 mMinPt2Casc = mSVParams->minPtCasc * mSVParams->minPtCasc;
282 mMaxTgl2Casc = mSVParams->maxTglCasc * mSVParams->maxTglCasc;
283 mMinPt23Body = mSVParams->minPt3Body * mSVParams->minPt3Body;
284 mMaxTgl23Body = mSVParams->maxTgl3Body * mSVParams->maxTgl3Body;
285 setupThreads();
286 }
287 auto bz = o2::base::Propagator::Instance()->getNominalBz();
288 mV0Hyps[HypV0::Photon].set(PID::Photon, PID::Electron, PID::Electron, mSVParams->pidCutsPhoton, bz);
289 mV0Hyps[HypV0::K0].set(PID::K0, PID::Pion, PID::Pion, mSVParams->pidCutsK0, bz);
290 mV0Hyps[HypV0::Lambda].set(PID::Lambda, PID::Proton, PID::Pion, mSVParams->pidCutsLambda, bz);
291 mV0Hyps[HypV0::AntiLambda].set(PID::Lambda, PID::Pion, PID::Proton, mSVParams->pidCutsLambda, bz);
296 mCascHyps[HypCascade::XiMinus].set(PID::XiMinus, PID::Lambda, PID::Pion, mSVParams->pidCutsXiMinus, bz, mSVParams->maximalCascadeWidth);
298
307
308 for (auto& ft : mFitterV0) {
309 ft.setBz(bz);
310 }
311 for (auto& ft : mFitterCasc) {
312 ft.setBz(bz);
313 }
314 for (auto& ft : mFitter3body) {
315 ft.setBz(bz);
316 }
317
318 mPIDresponse.setBetheBlochParams(mSVParams->mBBpars);
319}
320
321//______________________________________________
323{
324 mTPCVDrift = v.refVDrift * v.corrFact;
325 mTPCVDriftCorrFact = v.corrFact;
326 mTPCVDriftRef = v.refVDrift;
327 mTPCDriftTimeOffset = v.getTimeOffset();
328 mTPCBin2Z = mTPCVDrift / mMUS2TPCBin;
329}
330//______________________________________________
332{
333 mTPCCorrMapsHelper = maph;
334}
335
336//__________________________________________________________________
337void SVertexer::setupThreads()
338{
339 if (!mV0sTmp.empty()) {
340 return;
341 }
342 mV0sTmp.resize(mNThreads);
343 mCascadesTmp.resize(mNThreads);
344 m3bodyTmp.resize(mNThreads);
345 mV0sIdxTmp.resize(mNThreads);
346 mCascadesIdxTmp.resize(mNThreads);
347 m3bodyIdxTmp.resize(mNThreads);
348 mFitterV0.resize(mNThreads);
349 mBz = o2::base::Propagator::Instance()->getNominalBz();
350 int fitCounter = 0;
351 for (auto& fitter : mFitterV0) {
352 fitter.setFitterID(fitCounter++);
353 fitter.setBz(mBz);
354 fitter.setUseAbsDCA(mSVParams->useAbsDCA);
355 fitter.setPropagateToPCA(false);
356 fitter.setMaxR(mSVParams->maxRIni);
357 fitter.setMinParamChange(mSVParams->minParamChange);
358 fitter.setMinRelChi2Change(mSVParams->minRelChi2Change);
359 fitter.setMaxDZIni(mSVParams->maxDZIni);
360 fitter.setMaxDXYIni(mSVParams->maxDXYIni);
361 fitter.setMaxChi2(mSVParams->maxChi2);
362 fitter.setMatCorrType(o2::base::Propagator::MatCorrType(mSVParams->matCorr));
363 fitter.setUsePropagator(mSVParams->usePropagator);
364 fitter.setRefitWithMatCorr(mSVParams->refitWithMatCorr);
365 fitter.setMaxStep(mSVParams->maxStep);
366 fitter.setMaxSnp(mSVParams->maxSnp);
367 fitter.setMinXSeed(mSVParams->minXSeed);
368 }
369 mFitterCasc.resize(mNThreads);
370 fitCounter = 1000;
371 for (auto& fitter : mFitterCasc) {
372 fitter.setFitterID(fitCounter++);
373 fitter.setBz(mBz);
374 fitter.setUseAbsDCA(mSVParams->useAbsDCA);
375 fitter.setPropagateToPCA(false);
376 fitter.setMaxR(mSVParams->maxRIniCasc);
377 fitter.setMinParamChange(mSVParams->minParamChange);
378 fitter.setMinRelChi2Change(mSVParams->minRelChi2Change);
379 fitter.setMaxDZIni(mSVParams->maxDZIni);
380 fitter.setMaxDXYIni(mSVParams->maxDXYIni);
381 fitter.setMaxChi2(mSVParams->maxChi2);
382 fitter.setMatCorrType(o2::base::Propagator::MatCorrType(mSVParams->matCorr));
383 fitter.setUsePropagator(mSVParams->usePropagator);
384 fitter.setRefitWithMatCorr(mSVParams->refitWithMatCorr);
385 fitter.setMaxStep(mSVParams->maxStep);
386 fitter.setMaxSnp(mSVParams->maxSnp);
387 fitter.setMinXSeed(mSVParams->minXSeed);
388 }
389
390 mFitter3body.resize(mNThreads);
391 fitCounter = 2000;
392 for (auto& fitter : mFitter3body) {
393 fitter.setFitterID(fitCounter++);
394 fitter.setBz(mBz);
395 fitter.setUseAbsDCA(mSVParams->useAbsDCA);
396 fitter.setPropagateToPCA(false);
397 fitter.setMaxR(mSVParams->maxRIni3body);
398 fitter.setMinParamChange(mSVParams->minParamChange);
399 fitter.setMinRelChi2Change(mSVParams->minRelChi2Change);
400 fitter.setMaxDZIni(mSVParams->maxDZIni);
401 fitter.setMaxDXYIni(mSVParams->maxDXYIni);
402 fitter.setMaxChi2(mSVParams->maxChi2);
403 fitter.setMatCorrType(o2::base::Propagator::MatCorrType(mSVParams->matCorr));
404 fitter.setUsePropagator(mSVParams->usePropagator);
405 fitter.setRefitWithMatCorr(mSVParams->refitWithMatCorr);
406 fitter.setMaxStep(mSVParams->maxStep);
407 fitter.setMaxSnp(mSVParams->maxSnp);
408 fitter.setMinXSeed(mSVParams->minXSeed);
409 }
410}
411
412//__________________________________________________________________
413bool SVertexer::acceptTrack(const GIndex gid, const o2::track::TrackParCov& trc) const
414{
415 if (gid.isPVContributor() && mSVParams->maxPVContributors < 1) {
416 return false;
417 }
418
419 // DCA to mean vertex
420 if (mSVParams->minDCAToPV > 0.f) {
421 o2::track::TrackPar trp(trc);
422 std::array<float, 2> dca;
423 auto* prop = o2::base::Propagator::Instance();
424 if (mSVParams->usePropagator) {
425 if (trp.getX() > mSVParams->minRFor3DField && !prop->PropagateToXBxByBz(trp, mSVParams->minRFor3DField, mSVParams->maxSnp, mSVParams->maxStep, o2::base::Propagator::MatCorrType(mSVParams->matCorr))) {
426 return true; // we don't need actually to propagate to the beam-line
427 }
428 if (!prop->propagateToDCA(mMeanVertex.getXYZ(), trp, prop->getNominalBz(), mSVParams->maxStep, o2::base::Propagator::MatCorrType(mSVParams->matCorr), &dca)) {
429 return true;
430 }
431 } else {
432 if (!trp.propagateParamToDCA(mMeanVertex.getXYZ(), prop->getNominalBz(), &dca)) {
433 return true;
434 }
435 }
436 if (std::abs(dca[0]) < mSVParams->minDCAToPV) {
437 return false;
438 }
439 }
440 return true;
441}
442
443//__________________________________________________________________
444void SVertexer::buildT2V(const o2::globaltracking::RecoContainer& recoData) // accessor to various tracks
445{
446 // build track->vertices from vertices->tracks, rejecting vertex contributors if requested
447 auto trackIndex = recoData.getPrimaryVertexMatchedTracks(); // Global ID's for associated tracks
448 auto vtxRefs = recoData.getPrimaryVertexMatchedTrackRefs(); // references from vertex to these track IDs
449 bool isTPCloaded = recoData.isTrackSourceLoaded(GIndex::TPC);
450 bool isITSloaded = recoData.isTrackSourceLoaded(GIndex::ITS);
451 bool isITSTPCloaded = recoData.isTrackSourceLoaded(GIndex::ITSTPC);
452 if (isTPCloaded && !mSVParams->mExcludeTPCtracks) {
453 mTPCTracksArray = recoData.getTPCTracks();
454 mTPCTrackClusIdx = recoData.getTPCTracksClusterRefs();
455 mTPCClusterIdxStruct = &recoData.inputsTPCclusters->clusterIndex;
456 mTPCRefitterShMap = recoData.clusterShMapTPC;
457 mTPCRefitterOccMap = mRecoCont->occupancyMapTPC;
458 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());
459 mTPCRefitter->setTrackReferenceX(900); // disable propagation after refit by setting reference to value > 500
460 }
461
462 std::unordered_map<GIndex, std::pair<int, int>> tmap;
463 std::unordered_map<GIndex, bool> rejmap;
464 int nv = vtxRefs.size() - 1; // The last entry is for unassigned tracks, ignore them
465 for (int i = 0; i < 2; i++) {
466 mTracksPool[i].clear();
467 mVtxFirstTrack[i].clear();
468 mVtxFirstTrack[i].resize(nv, -1);
469 }
470 for (int iv = 0; iv < nv; iv++) {
471 const auto& vtref = vtxRefs[iv];
472 int it = vtref.getFirstEntry(), itLim = it + vtref.getEntries();
473 for (; it < itLim; it++) {
474 auto tvid = trackIndex[it];
475 if (!recoData.isTrackSourceLoaded(tvid.getSource())) {
476 continue;
477 }
478 if (tvid.getSource() == GIndex::TPC) {
479 if (mSVParams->mExcludeTPCtracks) {
480 continue;
481 }
482 // unconstrained TPC tracks require special treatment: there is no point in checking DCA to mean vertex since it is not precise,
483 // but we need to create a clone of TPC track constrained to this particular vertex time.
484 if (processTPCTrack(mTPCTracksArray[tvid], tvid, iv)) {
485 continue;
486 }
487 }
488
489 if (tvid.isAmbiguous()) { // was this track already processed?
490 auto tref = tmap.find(tvid);
491 if (tref != tmap.end()) {
492 mTracksPool[tref->second.second][tref->second.first].vBracket.setMax(iv); // this track was already processed with other vertex, account the latter
493 continue;
494 }
495 // was it already rejected?
496 if (rejmap.find(tvid) != rejmap.end()) {
497 continue;
498 }
499 }
500 const auto& trc = recoData.getTrackParam(tvid);
501
502 bool hasTPC = false;
503 bool heavyIonisingParticle = false;
504 bool compatibleWithProton = mSVParams->mFractiondEdxforCascBaryons > 0.999f; // if 1 or above, accept all regardless of TPC
505 auto tpcGID = recoData.getTPCContributorGID(tvid);
506 if (tpcGID.isIndexSet() && isTPCloaded) {
507 hasTPC = true;
508 auto& tpcTrack = recoData.getTPCTrack(tpcGID);
509 float dEdxTPC = tpcTrack.getdEdx().dEdxTotTPC;
510 if (dEdxTPC > mSVParams->minTPCdEdx && trc.getP() > mSVParams->minMomTPCdEdx) // accept high dEdx tracks (He3, He4)
511 {
512 heavyIonisingParticle = true;
513 }
514 auto protonId = o2::track::PID::Proton;
515 float dEdxExpected = mPIDresponse.getExpectedSignal(tpcTrack, protonId);
516 float fracDevProton = std::abs((dEdxTPC - dEdxExpected) / dEdxExpected);
517 if (fracDevProton < mSVParams->mFractiondEdxforCascBaryons) {
518 compatibleWithProton = true;
519 }
520 }
521
522 // get Nclusters in the ITS if available
523 int8_t nITSclu = -1;
524 bool shortOBITSOnlyTrack = false;
525 auto itsGID = recoData.getITSContributorGID(tvid);
526 if (itsGID.getSource() == GIndex::ITS) {
527 if (isITSloaded) {
528 auto& itsTrack = recoData.getITSTrack(itsGID);
529 nITSclu = itsTrack.getNumberOfClusters();
530 if (itsTrack.hasHitOnLayer(6) && itsTrack.hasHitOnLayer(5) && itsTrack.hasHitOnLayer(4) && itsTrack.hasHitOnLayer(3)) {
531 shortOBITSOnlyTrack = true;
532 }
533 }
534 } else if (itsGID.getSource() == GIndex::ITSAB) {
535 if (isITSTPCloaded) {
536 auto& itsABTracklet = recoData.getITSABRef(itsGID);
537 nITSclu = itsABTracklet.getNClusters();
538 }
539 }
540 if (!acceptTrack(tvid, trc) && !heavyIonisingParticle) {
541 if (tvid.isAmbiguous()) {
542 rejmap[tvid] = true;
543 }
544 continue;
545 }
546
547 if (!hasTPC && nITSclu < mSVParams->mITSSAminNclu && (!shortOBITSOnlyTrack || mSVParams->mRejectITSonlyOBtrack)) {
548 continue; // reject short ITS-only
549 }
550
551 int posneg = trc.getSign() < 0 ? 1 : 0;
552 float r = std::sqrt(trc.getX() * trc.getX() + trc.getY() * trc.getY());
553 mTracksPool[posneg].emplace_back(TrackCand{trc, tvid, {iv, iv}, r, hasTPC, nITSclu, compatibleWithProton});
554 if (tvid.getSource() == GIndex::TPC) { // constrained TPC track?
555 correctTPCTrack(mTracksPool[posneg].back(), mTPCTracksArray[tvid], -1, -1);
556 }
557 if (tvid.isAmbiguous()) { // track attached to >1 vertex, remember that it was already processed
558 tmap[tvid] = {mTracksPool[posneg].size() - 1, posneg};
559 }
560 }
561 }
562 // register 1st track of each charge for each vertex
563 for (int pn = 0; pn < 2; pn++) {
564 auto& vtxFirstT = mVtxFirstTrack[pn];
565 const auto& tracksPool = mTracksPool[pn];
566 for (unsigned i = 0; i < tracksPool.size(); i++) {
567 const auto& t = tracksPool[i];
568 for (int j{t.vBracket.getMin()}; j <= t.vBracket.getMax(); ++j) {
569 if (vtxFirstT[j] == -1) {
570 vtxFirstT[j] = i;
571 }
572 }
573 }
574 }
575
576 LOG(info) << "Collected " << mTracksPool[POS].size() << " positive and " << mTracksPool[NEG].size() << " negative seeds";
577}
578
579//__________________________________________________________________
580bool SVertexer::checkV0(const TrackCand& seedP, const TrackCand& seedN, int iP, int iN, int ithread)
581{
582 auto& fitterV0 = mFitterV0[ithread];
583 // Fast rough cuts on pairs before feeding to DCAFitter, tracks are not in the same Frame or at same X
584 bool isTPConly = (seedP.gid.getSource() == GIndex::TPC || seedN.gid.getSource() == GIndex::TPC);
585 if (mSVParams->mTPCTrackPhotonTune && isTPConly) {
586 // Check if Tgl is close enough
587 if (std::abs(seedP.getTgl() - seedN.getTgl()) > mSVParams->maxV0TglAbsDiff) {
588 LOG(debug) << "RejTgl";
589 return false;
590 }
591 // Check in transverse plane
592 float sna, csa;
593 o2::math_utils::CircleXYf_t trkPosCircle;
594 seedP.getCircleParams(mBz, trkPosCircle, sna, csa);
595 o2::math_utils::CircleXYf_t trkEleCircle;
596 seedN.getCircleParams(mBz, trkEleCircle, sna, csa);
597 // Does the radius of both tracks compare to their absolute circle center distance
598 float c2c = std::hypot(trkPosCircle.xC - trkEleCircle.xC,
599 trkPosCircle.yC - trkEleCircle.yC);
600 float r2r = trkPosCircle.rC + trkEleCircle.rC;
601 float dcr = c2c - r2r;
602 if (std::abs(dcr) > mSVParams->mTPCTrackD2R) {
603 LOG(debug) << "RejD2R " << c2c << " " << r2r << " " << dcr;
604 return false;
605 }
606 // Will the conversion point look somewhat reasonable
607 float r1_r = trkPosCircle.rC / r2r;
608 float r2_r = trkEleCircle.rC / r2r;
609 float dR = std::hypot(r2_r * trkPosCircle.xC + r1_r * trkEleCircle.xC, r2_r * trkPosCircle.yC + r1_r * trkEleCircle.yC);
610 if (dR > mSVParams->mTPCTrackDR) {
611 LOG(debug) << "RejDR" << dR;
612 return false;
613 }
614
615 // Setup looser cuts for the DCAFitter
616 fitterV0.setMaxDZIni(mSVParams->mTPCTrackMaxDZIni);
617 fitterV0.setMaxDXYIni(mSVParams->mTPCTrackMaxDXYIni);
618 fitterV0.setMaxChi2(mSVParams->mTPCTrackMaxChi2);
619 fitterV0.setCollinear(true);
620 }
621
622 // feed DCAFitter
623 int nCand = fitterV0.process(seedP, seedN);
624 if (mSVParams->mTPCTrackPhotonTune && isTPConly) {
625 // Reset immediately to the defaults
626 fitterV0.setMaxDZIni(mSVParams->maxDZIni);
627 fitterV0.setMaxDXYIni(mSVParams->maxDXYIni);
628 fitterV0.setMaxChi2(mSVParams->maxChi2);
629 fitterV0.setCollinear(false);
630 }
631
632 if (nCand == 0) { // discard this pair
633 LOG(debug) << "RejDCAFitter";
634 return false;
635 }
636 const auto& v0XYZ = fitterV0.getPCACandidate();
637 // validate V0 radial position
638 // check closeness to the beam-line
639 float dxv0 = v0XYZ[0] - mMeanVertex.getX(), dyv0 = v0XYZ[1] - mMeanVertex.getY(), r2v0 = dxv0 * dxv0 + dyv0 * dyv0;
640 if (r2v0 < mMinR2ToMeanVertex) {
641 LOG(debug) << "RejMinR2ToMeanVertex";
642 return false;
643 }
644 float rv0 = std::sqrt(r2v0), drv0P = rv0 - seedP.minR, drv0N = rv0 - seedN.minR;
645 if (drv0P > mSVParams->causalityRTolerance || drv0P < -mSVParams->maxV0ToProngsRDiff ||
646 drv0N > mSVParams->causalityRTolerance || drv0N < -mSVParams->maxV0ToProngsRDiff) {
647 LOG(debug) << "RejCausality " << drv0P << " " << drv0N;
648 return false;
649 }
650 const int cand = 0;
651 if (!fitterV0.isPropagateTracksToVertexDone(cand) && !fitterV0.propagateTracksToVertex(cand)) {
652 LOG(debug) << "RejProp failed";
653 return false;
654 }
655 const auto& trPProp = fitterV0.getTrack(0, cand);
656 const auto& trNProp = fitterV0.getTrack(1, cand);
657 std::array<float, 3> pP{}, pN{};
658 trPProp.getPxPyPzGlo(pP);
659 trNProp.getPxPyPzGlo(pN);
660 // estimate DCA of neutral V0 track to beamline: straight line with parametric equation
661 // x = X0 + pV0[0]*t, y = Y0 + pV0[1]*t reaches DCA to beamline (Xv, Yv) at
662 // t = -[ (x0-Xv)*pV0[0] + (y0-Yv)*pV0[1]) ] / ( pT(pV0)^2 )
663 // Similar equation for 3D distance involving pV0[2]
664 std::array<float, 3> pV0 = {pP[0] + pN[0], pP[1] + pN[1], pP[2] + pN[2]};
665 float pt2V0 = pV0[0] * pV0[0] + pV0[1] * pV0[1], prodXYv0 = dxv0 * pV0[0] + dyv0 * pV0[1], tDCAXY = prodXYv0 / pt2V0;
666 if (pt2V0 < mMinPt2V0) { // pt cut
667 LOG(debug) << "RejPt2 " << pt2V0;
668 return false;
669 }
670 if (pV0[2] * pV0[2] / pt2V0 > mMaxTgl2V0) { // tgLambda cut
671 LOG(debug) << "RejTgL " << pV0[2] * pV0[2] / pt2V0;
672 return false;
673 }
674 float p2V0 = pt2V0 + pV0[2] * pV0[2], ptV0 = std::sqrt(pt2V0);
675 // apply mass selections
676 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];
677
678 bool goodHyp = false, photonOnly = mSVParams->mTPCTrackPhotonTune && isTPConly;
679 std::array<bool, NHypV0> hypCheckStatus{};
680 int nPID = photonOnly ? (Photon + 1) : NHypV0;
681 for (int ipid = 0; (ipid < nPID) && mSVParams->checkV0Hypothesis; ipid++) {
682 if (mV0Hyps[ipid].check(p2Pos, p2Neg, p2V0, ptV0)) {
683 goodHyp = hypCheckStatus[ipid] = true;
684 }
685 }
686 // check tight lambda mass only
687 bool goodLamForCascade = false, goodALamForCascade = false;
688 bool usesTPCOnly = (seedP.hasTPC && !seedP.hasITS()) || (seedN.hasTPC && !seedN.hasITS());
689 bool usesShortITSOnly = (!seedP.hasTPC && seedP.nITSclu < mSVParams->mITSSAminNcluCascades) || (!seedN.hasTPC && seedN.nITSclu < mSVParams->mITSSAminNcluCascades);
690 if (ptV0 > mSVParams->minPtV0FromCascade && (!mSVParams->mSkipTPCOnlyCascade || !usesTPCOnly) && !usesShortITSOnly) {
691 if (mV0Hyps[Lambda].checkTight(p2Pos, p2Neg, p2V0, ptV0) && (!mSVParams->mRequireTPCforCascBaryons || seedP.hasTPC) && seedP.compatibleProton) {
692 goodLamForCascade = true;
693 }
694 if (mV0Hyps[AntiLambda].checkTight(p2Pos, p2Neg, p2V0, ptV0) && (!mSVParams->mRequireTPCforCascBaryons || seedN.hasTPC) && seedN.compatibleProton) {
695 goodALamForCascade = true;
696 }
697 }
698
699 // apply mass selections for 3-body decay
700 bool good3bodyV0Hyp = false;
701 for (int ipid = 2; ipid < 4; ipid++) {
702 float massForLambdaHyp = mV0Hyps[ipid].calcMass(p2Pos, p2Neg, p2V0);
703 if (massForLambdaHyp - mV0Hyps[ipid].getMassV0Hyp() < mV0Hyps[ipid].getMargin(ptV0)) {
704 good3bodyV0Hyp = true;
705 break;
706 }
707 }
708
709 // 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-)))
710 bool checkFor3BodyDecays = mEnable3BodyDecays &&
711 (!mSVParams->checkV0Hypothesis || good3bodyV0Hyp) &&
712 (pt2V0 > 0.5) &&
713 (!mSVParams->mSkipTPCOnly3Body || !isTPConly);
714 bool rejectAfter3BodyCheck = false; // To reject v0s which can be 3-body decay candidates but not cascade or v0
715 bool checkForCascade = mEnableCascades &&
716 (!mSVParams->mSkipTPCOnlyCascade || !usesTPCOnly) &&
717 r2v0 < mMaxR2ToMeanVertexCascV0 &&
718 (!mSVParams->checkV0Hypothesis || (goodLamForCascade || goodALamForCascade) &&
719 (!isTPConly || !hypCheckStatus[HypV0::Photon]));
720 bool rejectIfNotCascade = false;
721
722 if (!goodHyp && mSVParams->checkV0Hypothesis) {
723 LOG(debug) << "RejHypo";
724 if (!checkFor3BodyDecays && !checkForCascade) {
725 return false;
726 } else {
727 rejectAfter3BodyCheck = true;
728 }
729 }
730
731 float dcaX = dxv0 - pV0[0] * tDCAXY, dcaY = dyv0 - pV0[1] * tDCAXY, dca2 = dcaX * dcaX + dcaY * dcaY;
732 float cosPAXY = prodXYv0 / std::sqrt(r2v0 * pt2V0);
733
734 if (checkForCascade) { // use looser cuts for cascade v0 candidates
735 if (dca2 > mMaxDCAXY2ToMeanVertexV0Casc || cosPAXY < mSVParams->minCosPAXYMeanVertexCascV0) {
736 LOG(debug) << "Rej for cascade DCAXY2: " << dca2 << " << cosPAXY: " << cosPAXY;
737 if (!checkFor3BodyDecays) {
738 return false;
739 } else {
740 rejectAfter3BodyCheck = true;
741 }
742 }
743 }
744 if (checkFor3BodyDecays) { // use looser cuts for 3-body decay candidates
745 if (dca2 > mMaxDCAXY2ToMeanVertex3bodyV0 || cosPAXY < mSVParams->minCosPAXYMeanVertex3bodyV0) {
746 LOG(debug) << "Rej for 3 body decays DCAXY2: " << dca2 << " << cosPAXY: " << cosPAXY;
747 checkFor3BodyDecays = false;
748 }
749 }
750
751 if (dca2 > mMaxDCAXY2ToMeanVertex || cosPAXY < mSVParams->minCosPAXYMeanVertex) {
752 if (checkForCascade) {
753 rejectIfNotCascade = true;
754 } else if (checkFor3BodyDecays) {
755 rejectAfter3BodyCheck = true;
756 } else {
757 if (mSVParams->mTPCTrackPhotonTune && isTPConly) {
758 // Check for looser cut for tpc-only photons only
759 if (dca2 > mSVParams->mTPCTrackMaxDCAXY2ToMeanVertex) {
760 return false;
761 }
762 } else {
763 return false;
764 }
765 }
766 }
767
768 auto vlist = seedP.vBracket.getOverlap(seedN.vBracket); // indices of vertices shared by both seeds
769 bool candFound = false;
770 auto bestCosPA = checkForCascade ? mSVParams->minCosPACascV0 : mSVParams->minCosPA;
771 bestCosPA = checkFor3BodyDecays ? std::min(mSVParams->minCosPA3bodyV0, bestCosPA) : bestCosPA;
772 V0 v0new;
773 V0Index v0Idxnew;
774
775 for (int iv = vlist.getMin(); iv <= vlist.getMax(); iv++) {
776 const auto& pv = mPVertices[iv];
777 const auto v0XYZ = fitterV0.getPCACandidatePos(cand);
778 // check cos of pointing angle
779 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];
780 float cosPA = prodXYZv0 / std::sqrt((dx * dx + dy * dy + dz * dz) * p2V0);
781 if (cosPA < bestCosPA) {
782 LOG(debug) << "Rej. cosPA: " << cosPA;
783 continue;
784 }
785 if (!candFound) {
786 new (&v0new) V0(v0XYZ, pV0, fitterV0.calcPCACovMatrixFlat(cand), trPProp, trNProp);
787 new (&v0Idxnew) V0Index(-1, seedP.gid, seedN.gid);
788 v0new.setDCA(fitterV0.getChi2AtPCACandidate(cand));
789 candFound = true;
790 }
791 v0new.setCosPA(cosPA);
792 v0Idxnew.setVertexID(iv);
793 bestCosPA = cosPA;
794 }
795 if (!candFound) {
796 return false;
797 }
798 if (bestCosPA < mSVParams->minCosPACascV0) {
799 rejectAfter3BodyCheck = true;
800 }
801 if (bestCosPA < mSVParams->minCosPA && checkForCascade) {
802 rejectIfNotCascade = true;
803 }
804 int nV0Ini = mV0sIdxTmp[ithread].size();
805 // check 3 body decays
806 if (checkFor3BodyDecays) {
807 int n3bodyDecays = 0;
808 n3bodyDecays += check3bodyDecays(v0Idxnew, v0new, rv0, pV0, p2V0, iN, NEG, vlist, ithread);
809 n3bodyDecays += check3bodyDecays(v0Idxnew, v0new, rv0, pV0, p2V0, iP, POS, vlist, ithread);
810 }
811 if (rejectAfter3BodyCheck) {
812 return false;
813 }
814
815 // check cascades
816 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)
817 if (checkForCascade) {
818 if (goodLamForCascade || !mSVParams->checkCascadeHypothesis) {
819 nV0Used += checkCascades(v0Idxnew, v0new, rv0, pV0, p2V0, iN, NEG, vlist, ithread);
820 }
821 if (goodALamForCascade || !mSVParams->checkCascadeHypothesis) {
822 nV0Used += checkCascades(v0Idxnew, v0new, rv0, pV0, p2V0, iP, POS, vlist, ithread);
823 }
824 }
825
826 if (nV0Used) { // need to fix the index of V0 for the cascades using this v0
827 for (unsigned int ic = nCascIni; ic < mCascadesIdxTmp[ithread].size(); ic++) {
828 if (mCascadesIdxTmp[ithread][ic].getV0ID() == -1) {
829 mCascadesIdxTmp[ithread][ic].setV0ID(nV0Ini);
830 }
831 }
832 }
833
834 if (nV0Used || !rejectIfNotCascade) { // need to add this v0
835 mV0sIdxTmp[ithread].push_back(v0Idxnew);
836 if (!rejectIfNotCascade) {
837 mV0sIdxTmp[ithread].back().setStandaloneV0();
838 }
839 if (photonOnly) {
840 mV0sIdxTmp[ithread].back().setPhotonOnly();
841 mV0sIdxTmp[ithread].back().setCollinear();
842 }
843
844 if (mSVParams->createFullV0s) {
845 mV0sTmp[ithread].push_back(v0new);
846 }
847 }
848
849 if (mStrTracker) {
850 for (int iv = nV0Ini; iv < (int)mV0sIdxTmp[ithread].size(); iv++) {
851 mStrTracker->processV0(iv, v0new, v0Idxnew, ithread);
852 }
853 }
854
855 return mV0sIdxTmp[ithread].size() - nV0Ini != 0;
856}
857
858//__________________________________________________________________
859int 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)
860{
861 // check last added V0 for belonging to cascade
862 auto& fitterCasc = mFitterCasc[ithread];
863 auto& tracks = mTracksPool[posneg];
864 int nCascIni = mCascadesIdxTmp[ithread].size(), nv0use = 0;
865
866 // check if a given PV has already been used in a cascade
867 std::unordered_map<int, int> pvMap;
868
869 // start from the 1st bachelor track compatible with earliest vertex in the v0vlist
870 int firstTr = mVtxFirstTrack[posneg][v0vlist.getMin()], nTr = tracks.size();
871 if (firstTr < 0) {
872 firstTr = nTr;
873 }
874 for (int it = firstTr; it < nTr; it++) {
875 if (it == avoidTrackID) {
876 continue; // skip the track used by V0
877 }
878 auto& bach = tracks[it];
879 if (mSVParams->mSkipTPCOnlyCascade && (bach.gid.getSource() == GIndex::TPC)) {
880 continue; // reject TPC-only bachelors
881 }
882 if (!bach.hasTPC && bach.nITSclu < mSVParams->mITSSAminNcluCascades) {
883 continue; // reject short ITS-only
884 }
885
886 if (bach.vBracket.getMin() > v0vlist.getMax()) {
887 LOG(debug) << "Skipping";
888 break; // all other bachelor candidates will be also not compatible with this PV
889 }
890 auto cascVlist = v0vlist.getOverlap(bach.vBracket); // indices of vertices shared by V0 and bachelor
891 if (mSVParams->selectBestV0) {
892 // select only the best V0 candidate among the compatible ones
893 if (v0Idx.getVertexID() < cascVlist.getMin() || v0Idx.getVertexID() > cascVlist.getMax()) {
894 continue;
895 }
896 cascVlist.setMin(v0Idx.getVertexID());
897 cascVlist.setMax(v0Idx.getVertexID());
898 }
899
900 int nCandC = fitterCasc.process(v0, bach);
901 if (nCandC == 0) { // discard this pair
902 continue;
903 }
904 const int candC = 0;
905 const auto& cascXYZ = fitterCasc.getPCACandidatePos(candC);
906
907 // make sure the cascade radius is smaller than that of the mean vertex
908 float dxc = cascXYZ[0] - mMeanVertex.getX(), dyc = cascXYZ[1] - mMeanVertex.getY(), r2casc = dxc * dxc + dyc * dyc;
909 if (rv0 * rv0 - r2casc < mMinR2DiffV0Casc || r2casc < mMinR2ToMeanVertex) {
910 continue;
911 }
912 // do we want to apply mass cut ?
913 //
914 if (!fitterCasc.isPropagateTracksToVertexDone(candC) && !fitterCasc.propagateTracksToVertex(candC)) {
915 continue;
916 }
917
918 auto& trNeut = fitterCasc.getTrack(0, candC);
919 auto& trBach = fitterCasc.getTrack(1, candC);
920 trNeut.setPID(o2::track::PID::Lambda);
921 trBach.setPID(o2::track::PID::Pion);
922 std::array<float, 3> pNeut, pBach;
923 trNeut.getPxPyPzGlo(pNeut);
924 trBach.getPxPyPzGlo(pBach);
925 std::array<float, 3> pCasc = {pNeut[0] + pBach[0], pNeut[1] + pBach[1], pNeut[2] + pBach[2]};
926
927 float pt2Casc = pCasc[0] * pCasc[0] + pCasc[1] * pCasc[1], p2Casc = pt2Casc + pCasc[2] * pCasc[2];
928 if (pt2Casc < mMinPt2Casc) { // pt cut
929 LOG(debug) << "Casc pt too low";
930 continue;
931 }
932 if (pCasc[2] * pCasc[2] / pt2Casc > mMaxTgl2Casc) { // tgLambda cut
933 LOG(debug) << "Casc tgLambda too high";
934 continue;
935 }
936
937 // compute primary vertex and cosPA of the cascade
938 auto bestCosPA = mSVParams->minCosPACasc;
939 auto cascVtxID = -1;
940
941 for (int iv = cascVlist.getMin(); iv <= cascVlist.getMax(); iv++) {
942 const auto& pv = mPVertices[iv];
943 // check cos of pointing angle
944 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];
945 float cosPA = prodXYZcasc / std::sqrt((dx * dx + dy * dy + dz * dz) * p2Casc);
946 if (cosPA < bestCosPA) {
947 LOG(debug) << "Rej. cosPA: " << cosPA;
948 continue;
949 }
950 cascVtxID = iv;
951 bestCosPA = cosPA;
952 }
953 if (cascVtxID == -1) {
954 LOG(debug) << "Casc not compatible with any vertex";
955 continue;
956 }
957
958 const auto& cascPv = mPVertices[cascVtxID];
959 float dxCasc = cascXYZ[0] - cascPv.getX(), dyCasc = cascXYZ[1] - cascPv.getY(), dzCasc = cascXYZ[2] - cascPv.getZ();
960 auto prodPPos = pV0[0] * dxCasc + pV0[1] * dyCasc + pV0[2] * dzCasc;
961 if (prodPPos < 0.) { // causality cut
962 LOG(debug) << "Casc not causally compatible";
963 continue;
964 }
965
966 float p2Bach = pBach[0] * pBach[0] + pBach[1] * pBach[1] + pBach[2] * pBach[2];
967 float ptCasc = std::sqrt(pt2Casc);
968 bool goodHyp = false;
969 for (int ipid = 0; ipid < NHypCascade; ipid++) {
970 if (mCascHyps[ipid].check(p2V0, p2Bach, p2Casc, ptCasc)) {
971 goodHyp = true;
972 break;
973 }
974 }
975 if (!goodHyp) {
976 LOG(debug) << "Casc not compatible with any hypothesis";
977 continue;
978 }
979 // 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
980 // but not necessarily at the and of current v0s vector, since meanwhile checkCascades may add v0 clones (with PV redefined).
981 Cascade casc(cascXYZ, pCasc, fitterCasc.calcPCACovMatrixFlat(candC), trNeut, trBach);
982 o2::track::TrackParCov trc = casc;
984 if (!trc.propagateToDCA(cascPv, fitterCasc.getBz(), &dca, 5.) ||
985 std::abs(dca.getY()) > mSVParams->maxDCAXYCasc || std::abs(dca.getZ()) > mSVParams->maxDCAZCasc) {
986 LOG(debug) << "Casc not compatible with PV";
987 LOG(debug) << "DCA: " << dca.getY() << " " << dca.getZ();
988 continue;
989 }
990 CascadeIndex cascIdx(cascVtxID, -1, bach.gid); // the v0Idx was not yet added, this will be done after the checkCascades
991
992 LOGP(debug, "cascade successfully validated");
993
994 // clone the V0, set new cosPA and VerteXID, add it to the list of V0s
995 if (cascVtxID != v0Idx.getVertexID()) {
996 auto pvIdx = pvMap.find(cascVtxID);
997 if (pvIdx != pvMap.end()) {
998 cascIdx.setV0ID(pvIdx->second); // V0 already exists, add reference to the cascade
999 } else { // add V0 clone for this cascade (may be used also by other cascades)
1000 const auto& pv = mPVertices[cascVtxID];
1001 cascIdx.setV0ID(mV0sIdxTmp[ithread].size()); // set the new V0 index in the cascade
1002 pvMap[cascVtxID] = mV0sTmp[ithread].size(); // add the new V0 index to the map
1003 mV0sIdxTmp[ithread].emplace_back(cascVtxID, v0Idx.getProngs());
1004 if (mSVParams->createFullV0s) {
1005 mV0sTmp[ithread].push_back(v0);
1006 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];
1007 mV0sTmp[ithread].back().setCosPA(prodXYZ / std::sqrt((dx * dx + dy * dy + dz * dz) * p2V0));
1008 }
1009 }
1010 } else {
1011 nv0use++; // original v0 was used
1012 }
1013 mCascadesIdxTmp[ithread].push_back(cascIdx);
1014 if (mSVParams->createFullCascades) {
1015 casc.setCosPA(bestCosPA);
1016 casc.setDCA(fitterCasc.getChi2AtPCACandidate(candC));
1017 mCascadesTmp[ithread].push_back(casc);
1018 }
1019 if (mStrTracker) {
1020 mStrTracker->processCascade(mCascadesIdxTmp[ithread].size() - 1, casc, cascIdx, v0, ithread);
1021 }
1022 }
1023
1024 return nv0use;
1025}
1026
1027//__________________________________________________________________
1028int 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)
1029{
1030 // check last added V0 for belonging to cascade
1031 auto& fitter3body = mFitter3body[ithread];
1032 auto& tracks = mTracksPool[posneg];
1033 int n3BodyIni = m3bodyIdxTmp[ithread].size();
1034
1035 // start from the 1st bachelor track compatible with earliest vertex in the v0vlist
1036 int firstTr = mVtxFirstTrack[posneg][v0vlist.getMin()], nTr = tracks.size();
1037 if (firstTr < 0) {
1038 firstTr = nTr;
1039 }
1040
1041 // 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.
1042 // Otherwise, we should pair it with all negative particles, and the negative particle in the V0 is a antiproton.
1043
1044 // start from the 1st track compatible with V0's primary vertex
1045 for (int it = firstTr; it < nTr; it++) {
1046 if (it == avoidTrackID) {
1047 continue; // skip the track used by V0
1048 }
1049 auto& bach = tracks[it];
1050 if (mSVParams->mSkipTPCOnly3Body && (bach.gid.getSource() == GIndex::TPC)) {
1051 continue; // reject TPC-only bachelors
1052 }
1053 if (bach.vBracket > v0vlist.getMax()) {
1054 LOG(debug) << "Skipping";
1055 break; // all other bachelor candidates will be also not compatible with this PV
1056 }
1057 auto decay3bodyVlist = v0vlist.getOverlap(bach.vBracket); // indices of vertices shared by V0 and bachelor
1058 if (mSVParams->selectBestV0) {
1059 // select only the best V0 candidate among the compatible ones
1060 if (v0Idx.getVertexID() < decay3bodyVlist.getMin() || v0Idx.getVertexID() > decay3bodyVlist.getMax()) {
1061 continue;
1062 }
1063 decay3bodyVlist.setMin(v0Idx.getVertexID());
1064 decay3bodyVlist.setMax(v0Idx.getVertexID());
1065 }
1066
1067 if (bach.getPt() < 0.6) {
1068 continue;
1069 }
1070
1071 int n3bodyVtx = fitter3body.process(v0.getProng(0), v0.getProng(1), bach);
1072 if (n3bodyVtx == 0) { // discard this pair
1073 continue;
1074 }
1075 int cand3B = 0;
1076 const auto& vertexXYZ = fitter3body.getPCACandidatePos(cand3B);
1077
1078 // make sure the 3 body vertex radius is close to that of the mean vertex
1079 float dxc = vertexXYZ[0] - mMeanVertex.getX(), dyc = vertexXYZ[1] - mMeanVertex.getY(), dzc = vertexXYZ[2] - mMeanVertex.getZ(), r2vertex = dxc * dxc + dyc * dyc;
1080 if (std::abs(rv0 - std::sqrt(r2vertex)) > mSVParams->maxRDiffV03body || r2vertex < mMinR2ToMeanVertex) {
1081 continue;
1082 }
1083 float drvtxBach = std::sqrt(r2vertex) - bach.minR;
1084 if (drvtxBach > mSVParams->causalityRTolerance || drvtxBach < -mSVParams->maxV0ToProngsRDiff) {
1085 LOG(debug) << "RejCausality " << drvtxBach;
1086 }
1087 //
1088 if (!fitter3body.isPropagateTracksToVertexDone() && !fitter3body.propagateTracksToVertex()) {
1089 continue;
1090 }
1091
1092 auto& tr0 = fitter3body.getTrack(0, cand3B);
1093 auto& tr1 = fitter3body.getTrack(1, cand3B);
1094 auto& tr2 = fitter3body.getTrack(2, cand3B);
1095 std::array<float, 3> p0, p1, p2;
1096 tr0.getPxPyPzGlo(p0);
1097 tr1.getPxPyPzGlo(p1);
1098 tr2.getPxPyPzGlo(p2);
1099
1100 bool goodHyp = false;
1101 o2::track::PID pidHyp = o2::track::PID::Electron; // Update if goodHyp is true
1102 auto decay3bodyVtxID = -1;
1103 auto vtxCosPA = -1;
1104
1105 std::array<float, 3> pbach = {0, 0, 0}, p3B = {0, 0, 0}; // Update during the check of invariant mass
1106 for (int ipid = 0; ipid < NHyp3body; ipid++) {
1107 // check mass based on hypothesis of charge of bachelor (pos and neg expected to be proton/pion)
1108 float bachChargeFactor = m3bodyHyps[ipid].getChargeBachProng() / tr2.getAbsCharge();
1109 pbach = {bachChargeFactor * p2[0], bachChargeFactor * p2[1], bachChargeFactor * p2[2]};
1110 p3B = {p0[0] + p1[0] + pbach[0], p0[1] + p1[1] + pbach[1], p0[2] + p1[2] + pbach[2]};
1111 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];
1112 float pt2Candidate = p3B[0] * p3B[0] + p3B[1] * p3B[1], p2Candidate = pt2Candidate + p3B[2] * p3B[2];
1113 float ptCandidate = std::sqrt(pt2Candidate);
1114 if (m3bodyHyps[ipid].check(sqP0, sqP1, sqPBach, p2Candidate, ptCandidate)) {
1115 if (pt2Candidate < mMinPt23Body) { // pt cut
1116 continue;
1117 }
1118 if (p3B[2] * p3B[2] > pt2Candidate * mMaxTgl23Body) { // tgLambda cut
1119 continue;
1120 }
1121
1122 // compute primary vertex and cosPA of the 3-body decay
1123 auto bestCosPA = mSVParams->minCosPA3body;
1124 for (int iv = decay3bodyVlist.getMin(); iv <= decay3bodyVlist.getMax(); iv++) {
1125 const auto& pv = mPVertices[iv];
1126 // check cos of pointing angle
1127 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];
1128 float cosPA = prodXYZ3body / std::sqrt((dx * dx + dy * dy + dz * dz) * p2Candidate);
1129 if (cosPA < bestCosPA) {
1130 LOG(debug) << "Rej. cosPA: " << cosPA;
1131 continue;
1132 }
1133 decay3bodyVtxID = iv;
1134 bestCosPA = cosPA;
1135 }
1136 if (decay3bodyVtxID == -1) {
1137 LOG(debug) << "3-body decay not compatible with any vertex";
1138 continue;
1139 }
1140
1141 goodHyp = true;
1142 pidHyp = m3bodyHyps[ipid].getPIDHyp();
1143 vtxCosPA = bestCosPA;
1144 break;
1145 }
1146 }
1147 if (!goodHyp) {
1148 continue;
1149 }
1150
1151 const auto& decay3bodyPv = mPVertices[decay3bodyVtxID];
1152 Decay3Body candidate3B(vertexXYZ, p3B, fitter3body.calcPCACovMatrixFlat(cand3B), tr0, tr1, tr2, pidHyp);
1153 o2::track::TrackParCov trc = candidate3B;
1155 if (!trc.propagateToDCA(decay3bodyPv, fitter3body.getBz(), &dca, 5.) ||
1156 std::abs(dca.getY()) > mSVParams->maxDCAXY3Body || std::abs(dca.getZ()) > mSVParams->maxDCAZ3Body) {
1157 continue;
1158 }
1159 if (mSVParams->createFull3Bodies) {
1160 candidate3B.setCosPA(vtxCosPA);
1161 candidate3B.setDCA(fitter3body.getChi2AtPCACandidate());
1162 m3bodyTmp[ithread].push_back(candidate3B);
1163 }
1164 m3bodyIdxTmp[ithread].emplace_back(decay3bodyVtxID, v0Idx.getProngID(0), v0Idx.getProngID(1), bach.gid);
1165
1166 Decay3BodyIndex decay3bodyIdx(decay3bodyVtxID, v0Idx.getProngID(0), v0Idx.getProngID(1), bach.gid);
1167 if (mStrTracker) {
1168 mStrTracker->process3Body(m3bodyIdxTmp[ithread].size() - 1, candidate3B, decay3bodyIdx, ithread);
1169 }
1170 }
1171 return m3bodyIdxTmp[ithread].size() - n3BodyIni;
1172}
1173
1174//__________________________________________________________________
1175template <class TVI, class TCI, class T3I, class TR>
1176void SVertexer::extractPVReferences(const TVI& v0s, TR& vtx2V0Refs, const TCI& cascades, TR& vtx2CascRefs, const T3I& vtx3, TR& vtx2body3Refs)
1177{
1178 // V0s, cascades and 3bodies are already sorted in PV ID
1179 vtx2V0Refs.clear();
1180 vtx2V0Refs.resize(mPVertices.size());
1181 vtx2CascRefs.clear();
1182 vtx2CascRefs.resize(mPVertices.size());
1183 vtx2body3Refs.clear();
1184 vtx2body3Refs.resize(mPVertices.size());
1185 int nv0 = v0s.size(), nCasc = cascades.size(), n3body = vtx3.size();
1186
1187 // relate V0s to primary vertices
1188 int pvID = -1, nForPV = 0;
1189 for (int iv = 0; iv < nv0; iv++) {
1190 if (pvID < v0s[iv].getVertexID()) {
1191 if (pvID > -1) {
1192 vtx2V0Refs[pvID].setEntries(nForPV);
1193 }
1194 pvID = v0s[iv].getVertexID();
1195 vtx2V0Refs[pvID].setFirstEntry(iv);
1196 nForPV = 0;
1197 }
1198 nForPV++;
1199 }
1200 if (pvID != -1) { // finalize
1201 vtx2V0Refs[pvID].setEntries(nForPV);
1202 // fill empty slots
1203 int ent = nv0;
1204 for (int ip = vtx2V0Refs.size(); ip--;) {
1205 if (vtx2V0Refs[ip].getEntries()) {
1206 ent = vtx2V0Refs[ip].getFirstEntry();
1207 } else {
1208 vtx2V0Refs[ip].setFirstEntry(ent);
1209 }
1210 }
1211 }
1212
1213 // relate Cascades to primary vertices
1214 pvID = -1;
1215 nForPV = 0;
1216 for (int iv = 0; iv < nCasc; iv++) {
1217 if (pvID < cascades[iv].getVertexID()) {
1218 if (pvID > -1) {
1219 vtx2CascRefs[pvID].setEntries(nForPV);
1220 }
1221 pvID = cascades[iv].getVertexID();
1222 vtx2CascRefs[pvID].setFirstEntry(iv);
1223 nForPV = 0;
1224 }
1225 nForPV++;
1226 }
1227 if (pvID != -1) { // finalize
1228 vtx2CascRefs[pvID].setEntries(nForPV);
1229 // fill empty slots
1230 int ent = nCasc;
1231 for (int ip = vtx2CascRefs.size(); ip--;) {
1232 if (vtx2CascRefs[ip].getEntries()) {
1233 ent = vtx2CascRefs[ip].getFirstEntry();
1234 } else {
1235 vtx2CascRefs[ip].setFirstEntry(ent);
1236 }
1237 }
1238 }
1239
1240 // relate 3 body decays to primary vertices
1241 pvID = -1;
1242 nForPV = 0;
1243 for (int iv = 0; iv < n3body; iv++) {
1244 const auto& vertex3body = vtx3[iv];
1245 if (pvID < vertex3body.getVertexID()) {
1246 if (pvID > -1) {
1247 vtx2body3Refs[pvID].setEntries(nForPV);
1248 }
1249 pvID = vertex3body.getVertexID();
1250 vtx2body3Refs[pvID].setFirstEntry(iv);
1251 nForPV = 0;
1252 }
1253 nForPV++;
1254 }
1255 if (pvID != -1) { // finalize
1256 vtx2body3Refs[pvID].setEntries(nForPV);
1257 // fill empty slots
1258 int ent = n3body;
1259 for (int ip = vtx2body3Refs.size(); ip--;) {
1260 if (vtx2body3Refs[ip].getEntries()) {
1261 ent = vtx2body3Refs[ip].getFirstEntry();
1262 } else {
1263 vtx2body3Refs[ip].setFirstEntry(ent);
1264 }
1265 }
1266 }
1267}
1268
1269//__________________________________________________________________
1271{
1272#ifdef WITH_OPENMP
1273 mNThreads = n > 0 ? n : 1;
1274#else
1275 mNThreads = 1;
1276#endif
1277}
1278
1279//______________________________________________
1280bool SVertexer::processTPCTrack(const o2::tpc::TrackTPC& trTPC, GIndex gid, int vtxid)
1281{
1282 if (mSVParams->mTPCTrackMaxX > 0. && trTPC.getX() > mSVParams->mTPCTrackMaxX) {
1283 return true;
1284 }
1285 // if TPC trackis unconstrained, try to create in the tracks pool a clone constrained to vtxid vertex time.
1286 if (trTPC.hasBothSidesClusters()) { // this is effectively constrained track
1287 return false; // let it be processed as such
1288 }
1289 const auto& vtx = mPVertices[vtxid];
1290 auto twe = vtx.getTimeStamp();
1291 int posneg = trTPC.getSign() < 0 ? 1 : 0;
1292
1293 bool compatibleWithProton = false;
1294 if (!(mSVParams->mSkipTPCOnlyCascade)) {
1295 // Cascade retrieve dEdx proton frac
1296 const auto protonId = o2::track::PID::Proton;
1297 float dEdxTPC = trTPC.getdEdx().dEdxTotTPC;
1298 float dEdxExpected = mPIDresponse.getExpectedSignal(trTPC, protonId);
1299 float fracDevProton = std::abs((dEdxTPC - dEdxExpected) / dEdxExpected);
1300 if (fracDevProton < mSVParams->mFractiondEdxforCascBaryons) {
1301 compatibleWithProton = true;
1302 }
1303 }
1304
1305 auto& trLoc = mTracksPool[posneg].emplace_back(TrackCand{trTPC, gid, {vtxid, vtxid}, 0., true, -1, compatibleWithProton});
1306 auto err = correctTPCTrack(trLoc, trTPC, twe.getTimeStamp(), twe.getTimeStampError());
1307 if (err < 0) {
1308 mTracksPool[posneg].pop_back(); // discard
1309 return true;
1310 }
1311
1312 if (mSVParams->mTPCTrackPhotonTune) {
1313 // require minimum of tpc clusters
1314 bool dCls = trTPC.getNClusters() < mSVParams->mTPCTrackMinNClusters;
1315 // check track z cuts
1316 bool dDPV = std::abs(trLoc.getX() * trLoc.getTgl() - trLoc.getZ() + vtx.getZ()) > mSVParams->mTPCTrack2Beam;
1317 // check track transveres cuts
1318 float sna{0}, csa{0};
1320 trLoc.getCircleParams(mBz, trkCircle, sna, csa);
1321 float cR = std::hypot(trkCircle.xC, trkCircle.yC);
1322 float drd2 = std::sqrt(cR * cR - trkCircle.rC * trkCircle.rC);
1323 bool dRD2 = drd2 > mSVParams->mTPCTrackXY2Radius;
1324
1325 if (dCls || dDPV || dRD2) {
1326 mTracksPool[posneg].pop_back();
1327 return true;
1328 }
1329 }
1330
1331 return true;
1332}
1333
1334//______________________________________________
1335float SVertexer::correctTPCTrack(SVertexer::TrackCand& trc, const o2::tpc::TrackTPC& tTPC, float tmus, float tmusErr) const
1336{
1337 // Correct the track copy trc of the TPC track for the assumed interaction time
1338 // return extra uncertainty in Z due to the interaction time uncertainty
1339 // TODO: at the moment, apply simple shift, but with Z-dependent calibration we may
1340 // need to do corrections on TPC cluster level and refit
1341 // This is almosto clone of the MatchTPCITS::correctTPCTrack
1342
1343 float tTB, tTBErr;
1344 if (tmusErr < 0) { // use track data
1345 tTB = tTPC.getTime0();
1346 tTBErr = 0.5 * (tTPC.getDeltaTBwd() + tTPC.getDeltaTFwd());
1347 } else {
1348 tTB = tmus * mMUS2TPCBin;
1349 tTBErr = tmusErr * mMUS2TPCBin;
1350 }
1351 float dDrift = (tTB - tTPC.getTime0()) * mTPCBin2Z;
1352 float driftErr = tTBErr * mTPCBin2Z;
1353 if (driftErr < 0.) { // early return will be discarded anyway
1354 return driftErr;
1355 }
1356 // eventually should be refitted, at the moment we simply shift...
1357 trc.setZ(tTPC.getZ() + (tTPC.hasASideClustersOnly() ? dDrift : -dDrift));
1358 trc.setCov(trc.getSigmaZ2() + driftErr * driftErr, o2::track::kSigZ2);
1359 uint8_t sector, row;
1360 auto cl = &tTPC.getCluster(mTPCTrackClusIdx, tTPC.getNClusters() - 1, *mTPCClusterIdxStruct, sector, row);
1361 float x = 0, y = 0, z = 0;
1362 mTPCCorrMapsHelper->Transform(sector, row, cl->getPad(), cl->getTime(), x, y, z, tTB);
1365 }
1366 trc.minR = std::sqrt(x * x + y * y);
1367 LOGP(debug, "set MinR = {} for row {}, x:{}, y:{}, z:{}", trc.minR, row, x, y, z);
1368 return driftErr;
1369}
1370
1371//______________________________________________
1372std::array<size_t, 3> SVertexer::getNFitterCalls() const
1373{
1374 std::array<size_t, 3> calls{};
1375 for (int i = 0; i < mNThreads; i++) {
1376 calls[0] += mFitterV0[i].getCallID();
1377 calls[1] += mFitterCasc[i].getCallID();
1378 calls[2] += mFitter3body[i].getCallID();
1379 }
1380 return calls;
1381}
Helper class to access correction maps.
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
std::ostringstream debug
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:143
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.
Definition TFIDInfo.h:20
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, o2::pmr::polymorphic_allocator< T > > vector
TrackParCovF TrackParCov
Definition Track.h:33
struct o2::upgrades_utils::@454 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