Project
Loading...
Searching...
No Matches
PVertexer.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
18#include "Math/SMatrix.h"
19#include "Math/SVector.h"
20#include "MathUtils/fit.h"
21#include <unordered_map>
23#include <TH2F.h>
24
25using namespace o2::vertexing;
27constexpr float PVertexer::kAlmost0F;
28constexpr double PVertexer::kAlmost0D;
29constexpr float PVertexer::kHugeF;
30
31//___________________________________________________________________
33{
34 mMeanVertex.setSigma({0.5f, 0.5f, 6.0f});
35}
36
37//___________________________________________________________________
38int PVertexer::runVertexing(const gsl::span<o2d::GlobalTrackID> gids, const gsl::span<InteractionCandidate> intCand,
39 std::vector<PVertex>& vertices, std::vector<o2d::VtxTrackIndex>& vertexTrackIDs, std::vector<V2TRef>& v2tRefs,
40 gsl::span<const o2::MCCompLabel> lblTracks, std::vector<o2::MCEventLabel>& lblVtx)
41{
42#ifdef _PV_DEBUG_TREE_
43 doDBGPoolDump(lblTracks);
44#endif
45 mTimeDBScan.Start();
46 dbscan_clusterize();
47 mTimeDBScan.Stop();
48 mTotTrials = 0;
49 mMaxTrialPerCluster = 0;
50 mLongestClusterTimeMS = 0;
51 mLongestClusterMult = 0;
52 mNIniFound = 0;
53 mNKilledBCValid = 0;
54 mNKilledIntCand = 0;
55 mNKilledDebris = 0;
56 mNKilledQuality = 0;
57 mNKilledITSOnly = 0;
58 mPoolDumpProduced = false;
59
60 std::vector<PVertex> verticesLoc;
61 std::vector<uint32_t> trackIDs;
62 std::vector<V2TRef> v2tRefsLoc;
63 std::vector<float> validationTimes;
64 std::vector<o2::MCEventLabel> lblVtxLoc;
65 mTimeVertexing.Start();
66 int cntTZ = 0;
67 for (auto tc : mTimeZClusters) {
68 size_t nnold = verticesLoc.size();
70 inp.idRange = gsl::span<int>(tc.trackIDs);
71 inp.scaleSigma2 = mPVParams->iniScale2;
72 inp.timeEst = tc.timeEst;
73#ifdef _PV_DEBUG_TREE_
74 doDBScanDump(inp, lblTracks);
75#endif
76 findVertices(inp, verticesLoc, trackIDs, v2tRefsLoc);
77 }
78 mTimeVertexing.Stop();
79 // sort in time
80 std::vector<int> vtTimeSortID(verticesLoc.size());
81 std::iota(vtTimeSortID.begin(), vtTimeSortID.end(), 0);
82 std::sort(vtTimeSortID.begin(), vtTimeSortID.end(), [&verticesLoc](int i, int j) {
83 return verticesLoc[i].getTimeStamp().getTimeStamp() < verticesLoc[j].getTimeStamp().getTimeStamp();
84 });
85
86#ifdef _PV_DEBUG_TREE_
87 if (lblTracks.size()) { // at this stage labels are needed just for the debug output
88 createMCLabels(lblTracks, trackIDs, v2tRefsLoc, lblVtxLoc);
89 }
90#endif
91 mNIniFound = verticesLoc.size();
92
93 if (mValidateWithIR && mPVParams->minNContributorsForIRcutIni >= 0) {
94 applInteractionValidation(verticesLoc, vtTimeSortID, intCand, mPVParams->minNContributorsForIRcutIni);
95 }
96
97 if ((mPVParams->minITSOnlyFraction > 0.f || mPVParams->maxITSOnlyFraction < 1.0f) && !mITSOnly) {
98 applITSOnlyFractionCut(verticesLoc, vtTimeSortID, v2tRefsLoc, trackIDs);
99 }
100
101 if (mPVParams->applyDebrisReduction) {
102 mTimeDebris.Start();
103 reduceDebris(verticesLoc, vtTimeSortID, lblVtxLoc);
104 mTimeDebris.Stop();
105 }
106
107 if (mPVParams->applyReattachment) {
108 mTimeReAttach.Start();
109 reAttach(verticesLoc, vtTimeSortID, trackIDs, v2tRefsLoc);
110 mTimeReAttach.Stop();
111 }
112
113 if (lblTracks.size()) { // at this stage labels are needed just for the debug output
114 createMCLabels(lblTracks, trackIDs, v2tRefsLoc, lblVtxLoc);
115 }
116
117#ifdef _PV_DEBUG_TREE_
118 doVtxDump(verticesLoc, trackIDs, v2tRefsLoc, lblTracks);
119#endif
120
121 vertices.clear();
122 v2tRefs.clear();
123 vertexTrackIDs.clear();
124 vertices.reserve(verticesLoc.size());
125 v2tRefs.reserve(v2tRefsLoc.size());
126 vertexTrackIDs.reserve(trackIDs.size());
127
128 int trCopied = 0, count = 0, vtimeID = 0;
129 for (auto& i : vtTimeSortID) {
130 if (i < 0) {
131 continue; // vertex was suppressed
132 }
133 auto& vtx = verticesLoc[i];
134 if (!setCompatibleIR(vtx)) {
135 i = -1;
136 mNKilledBCValid++;
137 }
138 }
139 // do we need to validate by Int. records ?
140 if (mValidateWithIR) {
141 applInteractionValidation(verticesLoc, vtTimeSortID, intCand, mPVParams->minNContributorsForIRcut);
142 }
143 //
144 if (mPVParams->maxTMAD > 0 || mPVParams->maxZMAD > 0) {
145 mTimeMADSel.Start();
146 applyMADSelection(verticesLoc, vtTimeSortID, v2tRefsLoc, trackIDs);
147 mTimeMADSel.Stop();
148 }
149 //
150 for (auto i : vtTimeSortID) {
151 if (i < 0) {
152 continue; // vertex was suppressed
153 }
154 auto& vtx = verticesLoc[i];
155 vertices.push_back(vtx);
156 if (lblVtxLoc.size()) {
157 lblVtx.push_back(lblVtxLoc[i]);
158 }
159 int it = v2tRefsLoc[i].getFirstEntry(), itEnd = it + v2tRefsLoc[i].getEntries(), dest0 = vertexTrackIDs.size();
160 for (; it < itEnd; it++) {
161 const auto& trc = mTracksPool[trackIDs[it]];
162 auto& gid = vertexTrackIDs.emplace_back(trc.gid); // assign global track ID
163 gid.setPVContributor();
164 }
165 v2tRefs.emplace_back(dest0, v2tRefsLoc[i].getEntries());
166 LOG(debug) << "#" << count++ << " " << vertices.back() << " | " << v2tRefs.back().getEntries() << " indices from " << v2tRefs.back().getFirstEntry(); // RS REM
167 }
168
169 return vertices.size();
170}
171
172//______________________________________________
173int PVertexer::findVertices(const VertexingInput& input, std::vector<PVertex>& vertices, std::vector<uint32_t>& trackIDs, std::vector<V2TRef>& v2tRefs)
174{
175 // find vertices using tracks with indices (sorted in time) from idRange from "tracks" pool. The pool may containt arbitrary number of tracks,
176 // only those which are in the idRange and have canUse()==true, will be used.
177 // Results are placed in vertices and v2tRefs vectors
178
179 int nfound = 0, ntr = 0;
180 auto seedHistoTZ = buildHistoTZ(input); // histo for seeding peak finding
181
182#ifdef _PV_DEBUG_TREE_
183 static int dbsCount = -1;
184 dbsCount++;
185 int trialCount = 0; // TODO REM
186 auto hh = seedHistoTZ.createTH2F(o2::utils::Str::concat_string("htz", std::to_string(dbsCount)));
187 hh->SetDirectory(nullptr);
188 mDebugDumpFile->cd();
189 hh->Write();
190#endif
191 long tStart = std::chrono::time_point_cast<std::chrono::milliseconds>(std::chrono::system_clock::now()).time_since_epoch().count(), tCurr = tStart;
192 long mult = input.idRange.size();
193 int nTrials = 0;
194 while (nfound < mPVParams->maxVerticesPerCluster && nTrials < mPVParams->maxTrialsPerCluster) {
195 int peakBin = seedHistoTZ.findPeakBin();
196 if (!seedHistoTZ.isValidBin(peakBin)) {
197 break;
198 }
199 int peakBinT = seedHistoTZ.getXBin(peakBin), peakBinZ = seedHistoTZ.getYBin(peakBin);
200 float tv = seedHistoTZ.getBinXCenter(peakBinT);
201 float zv = seedHistoTZ.getBinYCenter(peakBinZ);
202 LOG(debug) << "Seeding with T=" << tv << " Z=" << zv << " bin " << peakBin << " on trial " << nTrials << " for vertex " << nfound << " mult=" << mult;
203
204 PVertex vtx;
205 mMeanVertex.setMeanXYVertexAtZ(vtx, zv);
206 vtx.setTimeStamp({tv, 0.f});
207 if (findVertex(input, vtx)) {
208 finalizeVertex(input, vtx, vertices, v2tRefs, trackIDs, &seedHistoTZ);
209 nfound++;
210 nTrials = 0;
211 } else { // suppress failed seeding bin and its proximities
212 seedHistoTZ.setBinContent(peakBin, -1);
213 }
214 nTrials++;
215
216#ifdef _PV_DEBUG_TREE_
217 auto hh1 = seedHistoTZ.createTH2F(o2::utils::Str::concat_string("htz", std::to_string(dbsCount), "_", std::to_string(trialCount++)));
218 hh1->SetDirectory(nullptr);
219 mDebugDumpFile->cd();
220 hh1->Write();
221#endif
222 tCurr = std::chrono::time_point_cast<std::chrono::milliseconds>(std::chrono::system_clock::now()).time_since_epoch().count();
223 auto clTime = tCurr - tStart;
224 if (clTime > mPVParams->maxTimeMSPerCluster) {
225 LOGP(warn, "Time per TZ-cluster ({}ms) of {} tracks exceeded limit after {} trials, abandon", clTime, mult, nTrials);
226 if (!mPoolDumpProduced) {
227 dumpPool();
228 }
229 break;
230 }
231 }
232 mTotTrials += nTrials;
233 if (size_t(nTrials) > mMaxTrialPerCluster) {
234 mMaxTrialPerCluster = nTrials;
235 }
236 if (tCurr - tStart > mLongestClusterTimeMS) {
237 mLongestClusterTimeMS = tCurr - tStart;
238 mLongestClusterMult = mult;
239 }
240 return nfound;
241}
242
243//______________________________________________
245{
246 // fit vertex taking provided vertex as a seed
247 // tracks pool may contain arbitrary number of tracks, only those which are in
248 // the idRange (indices of tracks sorted in time) will be used.
249
250 int ntr = input.idRange.size(); // RSREM
251
252 VertexSeed vtxSeed(vtx);
253 vtxSeed.setScale(input.scaleSigma2, mTukey2I);
254 vtxSeed.scaleSigma2Prev = input.scaleSigma2;
255 // vtxSeed.setTimeStamp( timeEstimate(input) );
256 LOG(debug) << "Start time guess: " << vtxSeed.getTimeStamp();
257 vtx.setChi2(1.e30);
258 //
260 bool found = false;
262 vtxSeed.resetForNewIteration();
263 vtxSeed.nIterations++;
264 LOG(debug) << "iter " << vtxSeed.nIterations << " with scale=" << vtxSeed.scaleSigma2 << " prevScale=" << vtxSeed.scaleSigma2Prev
265 << " ntr=" << ntr << " Zv=" << vtxSeed.getZ() << " Tv=" << vtxSeed.getTimeStamp().getTimeStamp();
266 result = fitIteration(input, vtxSeed);
267
268 if (result == FitStatus::OK) {
269 result = evalIterations(vtxSeed, vtx);
270 } else if (result == FitStatus::NotEnoughTracks) {
271 if (vtxSeed.nIterations <= mPVParams->maxIterations && upscaleSigma(vtxSeed)) {
272 LOG(debug) << "Upscaling scale to " << vtxSeed.scaleSigma2;
274 continue; // redo with stronger rescaling
275 } else {
276 break;
277 }
279 break;
280 } else {
281 LOG(fatal) << "Unknown fit status " << int(result);
282 }
283 }
284 LOG(debug) << "Stopped with scale=" << vtxSeed.scaleSigma2 << " prevScale=" << vtxSeed.scaleSigma2Prev << " result = " << int(result) << " ntr=" << ntr;
285
286 if (result != FitStatus::OK) {
287 vtx.setChi2(vtxSeed.maxScaleSigma2Tested);
288 return false;
289 } else {
290 return true;
291 }
292}
293
294//___________________________________________________________________
295PVertexer::FitStatus PVertexer::fitIteration(const VertexingInput& input, VertexSeed& vtxSeed)
296{
297 int nTested = 0;
298 for (int i : input.idRange) {
299 if (mTracksPool[i].canUse()) {
300 accountTrack(mTracksPool[i], vtxSeed);
301 // printf("#%d z:%f t:%f te:%f w:%f wh:%f\n", nTested, mTracksPool[i].z, mTracksPool[i].timeEst.getTimeStamp(),
302 // mTracksPool[i].timeEst.getTimeStampError(), mTracksPool[i].wgh, mTracksPool[i].wghHisto);
303 nTested++;
304 }
305 }
306
307 vtxSeed.maxScaleSigma2Tested = vtxSeed.scaleSigma2;
308 if (vtxSeed.getNContributors() < mPVParams->minTracksPerVtx) {
310 }
311 if (mPVParams->useMeanVertexConstraint) {
312 applyConstraint(vtxSeed);
313 }
314 if (!solveVertex(vtxSeed)) {
315 return FitStatus::Failure;
316 }
317 return FitStatus::OK;
318}
319
320//___________________________________________________________________
321void PVertexer::accountTrack(TrackVF& trc, VertexSeed& vtxSeed) const
322{
323 // deltas defined as track - vertex
324 bool useTime = vtxSeed.getTimeStamp().getTimeStampError() >= 0.f;
325 auto chi2T = trc.evalChi2ToVertex(vtxSeed, useTime && mPVParams->useTimeInChi2);
326 float wghT = (1.f - chi2T * vtxSeed.scaleSig2ITuk2I); // weighted distance to vertex
327 if (wghT < kAlmost0F) {
328 trc.wgh = 0.f;
329 return;
330 }
331 wghT *= wghT;
332 float syyI(trc.sig2YI), szzI(trc.sig2ZI), syzI(trc.sigYZI);
333
334 auto timeErrorFromTB = [&trc]() {
335 // decide if the time error is from the time bracket rather than gaussian error
336 return trc.gid.getSource() == GTrackID::ITS;
337 };
338
339 //
340 vtxSeed.wghSum += wghT;
341 vtxSeed.wghChi2 += wghT * chi2T;
342 //
343 syyI *= wghT;
344 syzI *= wghT;
345 szzI *= wghT;
346 trc.wgh = wghT;
347 //
348 // aux variables
349 double tmpSP = trc.sinAlp * trc.tgP, tmpCP = trc.cosAlp * trc.tgP,
350 tmpSC = trc.sinAlp + tmpCP, tmpCS = -trc.cosAlp + tmpSP,
351 tmpCL = trc.cosAlp * trc.tgL, tmpSL = trc.sinAlp * trc.tgL,
352 tmpYXP = trc.y - trc.tgP * trc.x, tmpZXL = trc.z - trc.tgL * trc.x,
353 tmpCLzz = tmpCL * szzI, tmpSLzz = tmpSL * szzI, tmpSCyz = tmpSC * syzI,
354 tmpCSyz = tmpCS * syzI, tmpCSyy = tmpCS * syyI, tmpSCyy = tmpSC * syyI,
355 tmpSLyz = tmpSL * syzI, tmpCLyz = tmpCL * syzI;
356 //
357 // symmetric matrix equation
358 vtxSeed.cxx += tmpCL * (tmpCLzz + tmpSCyz + tmpSCyz) + tmpSC * tmpSCyy; // dchi^2/dx/dx
359 vtxSeed.cxy += tmpCL * (tmpSLzz + tmpCSyz) + tmpSL * tmpSCyz + tmpSC * tmpCSyy; // dchi^2/dx/dy
360 vtxSeed.cxz += -trc.sinAlp * syzI - tmpCLzz - tmpCP * syzI; // dchi^2/dx/dz
361 vtxSeed.cx0 += -(tmpCLyz + tmpSCyy) * tmpYXP - (tmpCLzz + tmpSCyz) * tmpZXL; // RHS
362 //
363 vtxSeed.cyy += tmpSL * (tmpSLzz + tmpCSyz + tmpCSyz) + tmpCS * tmpCSyy; // dchi^2/dy/dy
364 vtxSeed.cyz += -(tmpCSyz + tmpSLzz); // dchi^2/dy/dz
365 vtxSeed.cy0 += -tmpYXP * (tmpCSyy + tmpSLyz) - tmpZXL * (tmpCSyz + tmpSLzz); // RHS
366 //
367 vtxSeed.czz += szzI; // dchi^2/dz/dz
368 vtxSeed.cz0 += tmpZXL * szzI + tmpYXP * syzI; // RHS
369 //
370 if (useTime) {
371 float trErr2I = wghT / (trc.timeEst.getTimeStampError() * trc.timeEst.getTimeStampError());
372 if (timeErrorFromTB()) {
373 vtxSeed.tMeanAccTB += trc.timeEst.getTimeStamp() * trErr2I;
374 vtxSeed.tMeanAccErrTB += trErr2I;
375 vtxSeed.nContributorsTB++;
376 vtxSeed.wghSumTB += wghT;
377 } else {
378 vtxSeed.tMeanAcc += trc.timeEst.getTimeStamp() * trErr2I;
379 vtxSeed.tMeanAccErr += trErr2I;
380 }
381 }
382 vtxSeed.addContributor();
383}
384
385//___________________________________________________________________
386bool PVertexer::solveVertex(VertexSeed& vtxSeed) const
387{
389 mat(0, 0) = vtxSeed.cxx;
390 mat(0, 1) = vtxSeed.cxy;
391 mat(0, 2) = vtxSeed.cxz;
392 mat(1, 1) = vtxSeed.cyy;
393 mat(1, 2) = vtxSeed.cyz;
394 mat(2, 2) = vtxSeed.czz;
395 if (!mat.InvertFast()) {
396 LOG(error) << "Failed to invert matrix" << mat;
397 return false;
398 }
399 ROOT::Math::SVector<double, 3> rhs(vtxSeed.cx0, vtxSeed.cy0, vtxSeed.cz0);
400 auto sol = mat * rhs;
401 vtxSeed.setXYZ(sol(0), sol(1), sol(2));
402 vtxSeed.setCov(mat(0, 0), mat(1, 0), mat(1, 1), mat(2, 0), mat(2, 1), mat(2, 2));
403 if (vtxSeed.tMeanAccErr + vtxSeed.tMeanAccErrTB > 0.) {
404 // since the time error from the ITS measurements does not improve with statistics, we downscale it with number of such tracks
405 auto t = vtxSeed.tMeanAcc;
406 auto e2i = vtxSeed.tMeanAccErr;
407 if (vtxSeed.wghSumTB > 0.) {
408 t += vtxSeed.tMeanAccTB / vtxSeed.wghSumTB;
409 e2i += vtxSeed.tMeanAccErrTB / vtxSeed.wghSumTB;
410 }
411 auto err2 = 1. / e2i;
412 vtxSeed.setTimeStamp({float(t * err2), float(std::sqrt(err2))});
413 }
414
415 vtxSeed.setChi2((vtxSeed.getNContributors() - vtxSeed.wghSum) / vtxSeed.scaleSig2ITuk2I); // calculate chi^2
416 auto newScale = vtxSeed.wghChi2 / vtxSeed.wghSum;
417 LOG(debug) << "Solve: wghChi2=" << vtxSeed.wghChi2 << " wghSum=" << vtxSeed.wghSum << " -> scale= "
418 << newScale << " old scale " << vtxSeed.scaleSigma2 << " prevScale: " << vtxSeed.scaleSigma2Prev
419 << " n-contributors=" << vtxSeed.getNContributors();
420 vtxSeed.setScale(newScale < mPVParams->minScale2 ? mPVParams->minScale2 : newScale, mTukey2I);
421 return true;
422}
423
424//___________________________________________________________________
425PVertexer::FitStatus PVertexer::evalIterations(VertexSeed& vtxSeed, PVertex& vtx) const
426{
427 // decide if new iteration should be done, prepare next one if needed
428 // if scaleSigma2 reached its lower limit stop
430
431 if (vtxSeed.nIterations > mPVParams->maxIterations) {
433 return result;
434 } else if (vtxSeed.scaleSigma2Prev <= mPVParams->minScale2 + kAlmost0F) {
436 LOG(debug) << "stop on simga :" << vtxSeed.scaleSigma2 << " prev: " << vtxSeed.scaleSigma2Prev;
437 }
438 if (fair::Logger::Logging(fair::Severity::debug)) {
439 auto dchi = (vtx.getChi2() - vtxSeed.getChi2()) / vtxSeed.getChi2();
440 auto dx = vtxSeed.getX() - vtx.getX(), dy = vtxSeed.getY() - vtx.getY(), dz = vtxSeed.getZ() - vtx.getZ();
441 auto dst = std::sqrt(dx * dx + dy * dy + dz * dz);
442
443 LOG(debug) << "dChi:" << vtx.getChi2() << "->" << vtxSeed.getChi2() << " :-> " << dchi;
444 LOG(debug) << "dx: " << dx << " dy: " << dy << " dz: " << dz << " -> " << dst;
445 }
446
447 vtx = reinterpret_cast<const PVertex&>(vtxSeed);
448
450 auto chi2Mean = vtxSeed.getChi2() / vtxSeed.getNContributors();
451 if (chi2Mean > mPVParams->maxChi2Mean) {
453 LOG(debug) << "Rejecting at iteration " << vtxSeed.nIterations << " and ScalePrev " << vtxSeed.scaleSigma2Prev << " with meanChi2 = " << chi2Mean;
454 } else {
455 return result;
456 }
457 }
458
459 if (vtxSeed.scaleSigma2 > vtxSeed.scaleSigma2Prev) {
460 if (++vtxSeed.nScaleIncrease > mPVParams->maxNScaleIncreased) {
462 LOG(debug) << "Rejecting at iteration " << vtxSeed.nIterations << " with NScaleIncreased " << vtxSeed.nScaleIncrease;
463 }
464 } else if (vtxSeed.scaleSigma2 > mPVParams->slowConvergenceFactor * vtxSeed.scaleSigma2Prev) {
465 if (++vtxSeed.nScaleSlowConvergence > mPVParams->maxNScaleSlowConvergence) {
466 if (vtxSeed.scaleSigma2 < mPVParams->acceptableScale2) {
467 vtxSeed.setScale(mPVParams->minScale2, mTukey2I);
468 LOG(debug) << "Forcing scale2 to " << mPVParams->minScale2;
470 } else {
472 LOG(debug) << "Rejecting at iteration " << vtxSeed.nIterations << " with NScaleSlowConvergence " << vtxSeed.nScaleSlowConvergence;
473 }
474 }
475 } else {
476 vtxSeed.nScaleSlowConvergence = 0;
477 }
478
479 return result;
480}
481
482//___________________________________________________________________
483void PVertexer::reAttach(std::vector<PVertex>& vertices, std::vector<int>& timeSort, std::vector<uint32_t>& trackIDs, std::vector<V2TRef>& v2tRefs)
484{
485 float tRange = 0.5 * std::max(mITSROFrameLengthMUS, mDBScanDeltaT) + mPVParams->timeMarginReattach; // consider only vertices in this proximity to tracks
486 std::vector<std::pair<int, TimeEst>> vtvec; // valid vertex times and indices
487 int nvtOrig = vertices.size();
488 vtvec.reserve(nvtOrig);
489 mTimeZClusters.resize(nvtOrig);
490 for (int ivt = 0; ivt < nvtOrig; ivt++) {
491 mTimeZClusters[ivt].trackIDs.clear();
492 mTimeZClusters[ivt].trackIDs.reserve(int(vertices[ivt].getNContributors() * 1.2));
493 }
494 for (auto ivt : timeSort) {
495 if (ivt >= 0) {
496 vtvec.emplace_back(ivt, vertices[ivt].getTimeStamp());
497 }
498 }
499 int ntr = mTracksPool.size(), nvt = vtvec.size();
500 int vtStart = 0;
501 for (int itr = 0; itr < ntr; itr++) {
502 auto& trc = mTracksPool[itr];
504 trc.wgh = kAlmost0F;
505 for (int ivt = vtStart; ivt < nvt; ivt++) {
506 auto dt = vtvec[ivt].second.getTimeStamp() - trc.timeEst.getTimeStamp();
507 if (dt < -tRange) { // all following tracks will have higher time than the vertex vtStart, move it
508 vtStart = ivt + 1;
509 continue;
510 }
511 if (dt > tRange) { // all following vertices will be have higher time than this track, stop checking
512 break;
513 }
514 bool useTime = (trc.gid.getSource() != GTrackID::ITS) && mPVParams->useTimeInChi2; // TODO Check if it is ok to not account time error for ITS tracks
515 float wgh = 1.f - mTukey2I * trc.evalChi2ToVertex(vertices[vtvec[ivt].first], useTime);
516 if (wgh > trc.wgh) {
517 trc.wgh = wgh;
518 trc.vtxID = vtvec[ivt].first;
519 }
520 }
521 if (trc.vtxID != TrackVF::kNoVtx) {
522 mTimeZClusters[trc.vtxID].trackIDs.push_back(itr);
523 trc.vtxID = TrackVF::kNoVtx; // to enable for using in the fit
525 }
526 }
527 // refit vertices with reattached tracks
528 v2tRefs.clear();
529 trackIDs.clear();
530 std::vector<PVertex> verticesUpd;
531 for (int ivt = 0; ivt < nvtOrig; ivt++) {
532 auto& clusZT = mTimeZClusters[ivt];
533 auto& vtx = vertices[ivt];
534 if (clusZT.trackIDs.size() < mPVParams->minTracksPerVtx) {
535 continue;
536 }
537 VertexingInput inp;
538 inp.idRange = gsl::span<int>(clusZT.trackIDs);
539 inp.scaleSigma2 = 1.;
540 inp.timeEst = vtx.getTimeStamp();
541 if (!findVertex(inp, vtx)) {
542 vtx.setNContributors(0);
543 continue;
544 }
545 finalizeVertex(inp, vtx, verticesUpd, v2tRefs, trackIDs);
546 }
547 // reorder in time since the time-stamp of vertices might have been changed
548 vertices.swap(verticesUpd);
549 timeSort.resize(vertices.size());
550 std::iota(timeSort.begin(), timeSort.end(), 0);
551 std::sort(timeSort.begin(), timeSort.end(), [&vertices](int i, int j) {
552 return vertices[i].getTimeStamp().getTimeStamp() < vertices[j].getTimeStamp().getTimeStamp();
553 });
554}
555
556//___________________________________________________________________
557void PVertexer::applyMADSelection(std::vector<PVertex>& vertices, std::vector<int>& timeSort, const std::vector<V2TRef>& v2tRefs, const std::vector<uint32_t>& trackIDs)
558{
559 int nv = timeSort.size(), nkill = 0;
560 std::vector<float> dvec;
561 for (int ivt = 0; ivt < nv; ivt++) {
562 int iv = timeSort[ivt];
563 if (iv < 0) { // already rejected?
564 continue;
565 }
566 auto& pv = vertices[iv];
567 const auto trefs = v2tRefs[iv];
568 int idMin = trefs.getFirstEntry(), idMax = idMin + trefs.getEntries();
569 dvec.clear();
570 dvec.reserve(pv.getNContributors());
571 if (mPVParams->maxTMAD > 0) {
572 int nAdjusted = 0;
573 float tVtx = pv.getTimeStamp().getTimeStamp();
574 for (int i = idMin; i < idMax; i++) {
575 const auto& trc = mTracksPool[trackIDs[i]];
576 if (trc.gid.getSource() == GTrackID::ITS || // ITS should not be accounted in the mad MAD as it has no real errors
577 (trc.gid.getSource() == GTrackID::ITSTPC && trc.isITSTPCAdjusted() && (++nAdjusted) > 1)) { // adjusted tracks are counted only once since they are all the same
578 continue;
579 } else {
580 dvec.push_back(trc.timeEst.getTimeStamp() - tVtx);
581 }
582 }
583 float tmad = o2::math_utils::MAD2Sigma(dvec.size(), dvec.data());
584 dvec.clear();
585 if (tmad > mPVParams->maxTMAD || tmad < mPVParams->minTMAD) {
586 timeSort[ivt] = -1; // disable vertex
587 mNKilledQuality++;
588 LOGP(debug, "Killing vertex {} with TMAD {}, {} of {} killed", iv, tmad, ++nkill, nv);
589 continue;
590 }
591 pv.setTMAD(tmad);
592 }
593 if (mPVParams->maxZMAD > 0) {
594 float zVtx = pv.getZ();
595 for (int i = idMin; i < idMax; i++) {
596 const auto& trc = mTracksPool[trackIDs[i]];
597 dvec.push_back(trc.z - zVtx);
598 }
599 float zmad = o2::math_utils::MAD2Sigma(dvec.size(), dvec.data());
600 if (zmad > mPVParams->maxZMAD || zmad < mPVParams->minZMAD) {
601 timeSort[ivt] = -1; // disable vertex
602 mNKilledQuality++;
603 LOGP(debug, "Killing vertex {} with ZMAD {}, {} of {} killed", iv, zmad, ++nkill, nv);
604 continue;
605 }
606 pv.setZMAD(zmad);
607 }
608 }
609}
610
611//___________________________________________________________________
612void PVertexer::applITSOnlyFractionCut(std::vector<PVertex>& vertices, std::vector<int>& timeSort, const std::vector<V2TRef>& v2tRefs, const std::vector<uint32_t>& trackIDs)
613{
614 int nv = timeSort.size();
615 for (int ivt = 0; ivt < nv; ivt++) {
616 int iv = timeSort[ivt];
617 if (iv < 0) { // already rejected?
618 continue;
619 }
620 const auto& pv = vertices[iv];
621 float tVtx = pv.getTimeStamp().getTimeStamp();
622 const auto trefs = v2tRefs[iv];
623 int idMin = trefs.getFirstEntry(), idMax = idMin + trefs.getEntries();
624 int nITS = 0;
625 for (int i = idMin; i < idMax; i++) {
626 const auto& trc = mTracksPool[trackIDs[i]];
627 if (trc.gid.getSource() == GTrackID::Source::ITS) { // ITS should not be accounted in the mad MAD as it has no real errors
628 nITS++;
629 }
630 }
631 float frac = nITS / float(trefs.getEntries());
632 if (frac > mPVParams->maxITSOnlyFraction || frac < mPVParams->minITSOnlyFraction) {
633 timeSort[ivt] = -1; // disable vertex
634 mNKilledITSOnly++;
635 }
636 }
637}
638
639//___________________________________________________________________
640void PVertexer::applInteractionValidation(std::vector<PVertex>& vertices, std::vector<int>& timeSort, const gsl::span<InteractionCandidate> intCand, int minContrib)
641{
642 // get max error
643 int nv = timeSort.size(), nkill = 0;
644 float maxErr = std::max(mPVParams->minTError, mPVParams->nSigmaTimeCut * mPVParams->maxTError) + mPVParams->timeMarginVertexTime;
645
646 int intCandMin = 0, nInt = intCand.size();
647 for (int ivt = 0; ivt < nv; ivt++) {
648 int iv = timeSort[ivt];
649 if (iv < 0) { // already rejected?
650 continue;
651 }
652 auto& pv = vertices[iv];
653 while (intCandMin < nInt && intCand[intCandMin].time < (pv.getTimeStamp().getTimeStamp() + mPVParams->timeBiasMS - maxErr)) {
654 intCandMin++; // skip all times which have no chance to be matched
655 }
656 int i = intCandMin, nCompatible = 0, best = -1;
657 float bestDf = 1e12;
658 while (i < nInt) {
659 float tv = pv.getTimeStamp().getTimeStamp() + mPVParams->timeBiasMS;
660 float tvE = std::max(mPVParams->minTError, mPVParams->nSigmaTimeCut * std::min(mPVParams->maxTError, pv.getTimeStamp().getTimeStampError())) + mPVParams->timeMarginVertexTime;
661 if (intCand[i].time > (tv + tvE)) {
662 break;
663 }
664 if (intCand[i].time < (tv - tvE)) {
665 i++;
666 continue;
667 }
668 nCompatible++;
669 auto dfa = std::abs(tv - intCand[intCandMin].time);
670 if (dfa < bestDf) {
671 bestDf = dfa;
672 best = i;
673 }
674 i++;
675 }
676 if (best >= 0) {
677 pv.setFlags(PVertex::TimeValidated);
678 if (nCompatible == 1) {
679 pv.setIR(intCand[best]);
680 }
681 } else if (pv.getNContributors() >= minContrib) {
682 timeSort[ivt] = -1; // discard
683 mNKilledIntCand++;
684 }
685 }
686}
687
688//___________________________________________________________________
689void PVertexer::reduceDebris(std::vector<PVertex>& vertices, std::vector<int>& timeSort, const std::vector<o2::MCEventLabel>& lblVtx)
690{
691 // eliminate low multiplicity vertices in the close proximity of high mult ones, assuming that these are their debries
692 // The timeSort vector indicates the time ordering of the vertices
693 int nv = vertices.size();
694 std::vector<int> multSort(nv); // sort time indices in multiplicity
695 std::iota(multSort.begin(), multSort.end(), 0);
696 std::sort(multSort.begin(), multSort.end(), [&timeSort, vertices](int i, int j) {
697 return timeSort[i] < 0 ? false : (timeSort[j] < 0 ? true : vertices[timeSort[i]].getNContributors() > vertices[timeSort[j]].getNContributors());
698 });
699
700 // suppress vertex pointed by j if needed, if time difference between i and j is too large, return true to stop checking
701 // pairs starting with i.
702 auto checkPair = [&vertices, &timeSort, &lblVtx, this](int i, int j) {
703 auto &vtI = vertices[timeSort[i]], &vtJ = vertices[timeSort[j]];
704 auto tDiff = std::abs(vtI.getTimeStamp().getTimeStamp() - vtJ.getTimeStamp().getTimeStamp());
705
706 if (tDiff > this->mMaxTDiffDebrisFiducial) {
707 return true; // don't continue checking other neighbours in time
708 }
709 if (vtI.getNContributors() < vtJ.getNContributors()) {
710 return false; // comparison goes from higher to lower mult vtx
711 }
712 bool rej = false;
713 float zDiff = std::abs(vtI.getZ() - vtJ.getZ());
714 float chi2z = 0.f, chi2t = 0.f, chi2zE = 0.f, chi2tE = 0.f;
715 if (zDiff < this->mMaxZDiffDebrisFiducial && float(vtJ.getNContributors() < float(vtI.getNContributors()) * this->mMaxMultRatDebrisFiducial)) {
716 if (this->mMaxTDiffDebrisExtra <= 0 || // if no extra cut requested, no need to recheck the differences
717 (zDiff < this->mPVParams->maxZDiffDebris && tDiff < this->mMaxTDiffDebris && float(vtJ.getNContributors()) < float(vtI.getNContributors()) * this->mPVParams->maxMultRatDebris)) {
718 chi2z = zDiff * zDiff / (vtI.getSigmaZ2() + vtJ.getSigmaZ2() + this->mPVParams->addZSigma2Debris);
719 chi2t = tDiff * tDiff / (vtI.getTimeStamp().getTimeStampError2() + vtJ.getTimeStamp().getTimeStampError2() + this->mPVParams->addTimeSigma2Debris);
720 rej = (chi2z + chi2t) < this->mPVParams->maxChi2TZDebris;
721 }
722 if (!rej && this->mMaxTDiffDebrisExtra > 0 &&
723 (zDiff < this->mPVParams->maxZDiffDebrisExtra && tDiff < this->mMaxTDiffDebrisExtra && float(vtJ.getNContributors()) < float(vtI.getNContributors()) * this->mPVParams->maxMultRatDebrisExtra)) {
724 chi2zE = zDiff * zDiff / (vtI.getSigmaZ2() + vtJ.getSigmaZ2() + this->mPVParams->addZSigma2DebrisExtra);
725 chi2tE = tDiff * tDiff / (vtI.getTimeStamp().getTimeStampError2() + vtJ.getTimeStamp().getTimeStampError2() + this->mPVParams->addTimeSigma2DebrisExtra);
726 rej = (chi2zE + chi2tE) < this->mPVParams->maxChi2TZDebrisExtra;
727 }
728#ifdef _PV_DEBUG_TREE_
729 o2::MCEventLabel dummyLbl;
730 this->mDebugDumpPVComp.emplace_back(vtI, vtJ, chi2z, chi2t, chi2zE, chi2tE, rej);
731 if (!lblVtx.empty()) {
732 this->mDebugDumpPVCompLbl0.push_back(lblVtx[timeSort[i]]);
733 this->mDebugDumpPVCompLbl1.push_back(lblVtx[timeSort[j]]);
734 }
735#endif
736 }
737 if (rej) {
738 timeSort[j] = -1;
739 mNKilledDebris++;
740 vtJ.setNContributors(0);
741 }
742 return false;
743 };
744
745 for (int im = 0; im < nv; im++) { // loop from highest multiplicity to lowest one
746 int it = multSort[im];
747 if (timeSort[it] < 0) { // if <0, the vertex was already discarded
748 continue;
749 }
750
751 int itL = it; // look for vertices with smaller time
752 while (itL) {
753 if (timeSort[--itL] >= 0) { // if <0, the vertex was already discarded
754 if (checkPair(it, itL)) { // if too far in time, don't compare further
755 break;
756 }
757 }
758 } // itL loop
759 int itH = it; // look for vertices with higher time
760 while (++itH < nv) {
761 if (timeSort[itH] >= 0) { // if <0, the vertex was already discarded
762 if (checkPair(it, itH)) { // if too far in time, don't compare further
763 break;
764 }
765 }
766 } // itH loop
767 }
768#ifdef _PV_DEBUG_TREE_
769 mDebugVtxCompTree->Fill();
770 mDebugDumpPVComp.clear();
771 mDebugDumpPVCompLbl0.clear();
772 mDebugDumpPVCompLbl1.clear();
773#endif
774}
775
776//___________________________________________________________________
777void PVertexer::initMeanVertexConstraint()
778{
779 // set mean vertex constraint and its errors
780 float extraErr2 = mPVParams->meanVertexExtraErrConstraint * mPVParams->meanVertexExtraErrConstraint;
781 float ex2 = mMeanVertex.getSigmaX2() + extraErr2, ey2 = mMeanVertex.getSigmaY2() + extraErr2, exy = mMeanVertex.getSigmaXY();
782 double det = ex2 * ey2 - exy * exy;
783 if (det <= kAlmost0D || ex2 < kAlmost0D || ey2 < kAlmost0D) {
784 throw std::runtime_error(fmt::format("Singular matrix for vertex constraint: sxx={:+.4e} syy={:+.4e} sxy={:+.4e}",
785 ex2, ey2, exy));
786 }
787 mXYConstraintInvErr[0] = ey2 / det;
788 mXYConstraintInvErr[1] = -exy / det;
789 mXYConstraintInvErr[2] = ex2 / det;
790}
791
792//______________________________________________
794{
795 // convert 1/tukey^2 to tukey
796 return sqrtf(1. / mTukey2I);
797}
798
799//___________________________________________________________________
800TimeEst PVertexer::timeEstimate(const VertexingInput& input) const
801{
803 for (int i : input.idRange) {
804 if (mTracksPool[i].canUse()) {
805 const auto& timeT = mTracksPool[i].timeEst;
806 auto trErr2I = 1. / (timeT.getTimeStampError() * timeT.getTimeStampError());
807 test.add(timeT.getTimeStamp(), trErr2I);
808 }
809 }
810
811 const auto [t, te2] = test.getMeanRMS2<float>();
812 return {t, te2};
813}
814
815//___________________________________________________________________
817{
818 mPVParams = &PVertexerParams::Instance();
819 initMeanVertexConstraint();
820 setTukey(mPVParams->tukey);
821 auto* prop = o2::base::Propagator::Instance();
822 setBz(prop->getNominalBz());
823 mDBScanDeltaT = mPVParams->dbscanDeltaT > 0.f ? mPVParams->dbscanDeltaT : mITSROFrameLengthMUS * (-mPVParams->dbscanDeltaT);
824 mDBSMaxZ2InvCorePoint = mPVParams->dbscanMaxSigZCorPoint > 0 ? 1. / (mPVParams->dbscanMaxSigZCorPoint * mPVParams->dbscanMaxSigZCorPoint) : 1e6;
825
826 mMaxTDiffDebris = mPVParams->maxTDiffDebris < 0 ? mITSROFrameLengthMUS * (-mPVParams->maxTDiffDebris) : mPVParams->maxTDiffDebris;
827 mMaxTDiffDebrisExtra = mPVParams->maxTDiffDebrisExtra == 0 ? -1 : (mPVParams->maxTDiffDebrisExtra < 0 ? mITSROFrameLengthMUS * (-mPVParams->maxTDiffDebrisExtra) : mPVParams->maxTDiffDebrisExtra);
828 mMaxTDiffDebrisFiducial = mMaxTDiffDebrisExtra > 0 ? std::max(mMaxTDiffDebrisExtra, mMaxTDiffDebris) : mMaxTDiffDebris;
829 mMaxZDiffDebrisFiducial = mMaxTDiffDebrisExtra > 0 ? std::max(mPVParams->maxZDiffDebrisExtra, mPVParams->maxZDiffDebris) : mPVParams->maxZDiffDebris;
830 mMaxMultRatDebrisFiducial = mMaxTDiffDebrisExtra > 0 ? std::max(mPVParams->maxMultRatDebrisExtra, mPVParams->maxMultRatDebris) : mPVParams->maxMultRatDebris;
831
832#ifdef _PV_DEBUG_TREE_
833 mDebugDumpFile = std::make_unique<TFile>("pvtxDebug.root", "recreate");
834
835 mDebugPoolTree = std::make_unique<TTree>("pvtxTrackPool", "PVertexer tracks pool debug output");
836 mDebugPoolTree->Branch("trc", &mDebugDumpDBSTrc);
837 mDebugPoolTree->Branch("gid", &mDebugDumpDBSGID);
838 mDebugPoolTree->Branch("mc", &mDebugDumpDBSTrcMC);
839
840 mDebugDBScanTree = std::make_unique<TTree>("pvtxDBScan", "PVertexer DBScan debug output");
841 mDebugDBScanTree->Branch("trc", &mDebugDumpDBSTrc);
842 mDebugDBScanTree->Branch("gid", &mDebugDumpDBSGID);
843 mDebugDBScanTree->Branch("mc", &mDebugDumpDBSTrcMC);
844
845 mDebugVtxTree = std::make_unique<TTree>("pvtx", "final PVertexer debug output");
846 mDebugVtxTree->Branch("vtx", &mDebugDumpVtx);
847 mDebugVtxTree->Branch("trc", &mDebugDumpVtxTrc);
848 mDebugVtxTree->Branch("gid", &mDebugDumpVtxGID);
849 mDebugVtxTree->Branch("mc", &mDebugDumpVtxTrcMC);
850
851 mDebugVtxCompTree = std::make_unique<TTree>("pvtxComp", "PVertexer neighbouring vertices debud output");
852 mDebugVtxCompTree->Branch("vtxComp", &mDebugDumpPVComp);
853 mDebugVtxCompTree->Branch("vtxCompLbl0", &mDebugDumpPVCompLbl0);
854 mDebugVtxCompTree->Branch("vtxCompLbl1", &mDebugDumpPVCompLbl1);
855#endif
856}
857
858//___________________________________________________________________
860{
861#ifdef _PV_DEBUG_TREE_
862 if (mDebugDumpFile) {
863 mDebugDumpFile->cd();
864 mDebugPoolTree->Write();
865 mDebugDBScanTree->Write();
866 mDebugVtxCompTree->Write();
867 mDebugVtxTree->Write();
868
869 mDebugPoolTree.reset();
870 mDebugDBScanTree.reset();
871 mDebugVtxCompTree.reset();
872 mDebugVtxTree.reset();
873
874 mDebugDumpFile->Close();
875 mDebugDumpFile.reset();
876 }
877#endif
878}
879
880//___________________________________________________________________
881void PVertexer::finalizeVertex(const VertexingInput& input, const PVertex& vtx,
882 std::vector<PVertex>& vertices, std::vector<V2TRef>& v2tRefs, std::vector<uint32_t>& trackIDs,
883 SeedHistoTZ* histo)
884{
885 int lastID = vertices.size();
886 vertices.emplace_back(vtx);
887 auto& ref = v2tRefs.emplace_back(trackIDs.size(), 0);
888 for (int i : input.idRange) {
889 auto& trc = mTracksPool[i];
890 if (trc.canAssign()) {
891 trackIDs.push_back(i);
892 trc.vtxID = lastID;
893 if (trc.bin != TrackVF::kDummyHBin && histo) {
894 histo->setBinContent(trc.bin, -1.f); // discard used bin
895 }
896 }
897 }
898 ref.setEntries(trackIDs.size() - ref.getFirstEntry());
899}
900
901//___________________________________________________________________
902void PVertexer::createMCLabels(gsl::span<const o2::MCCompLabel> lblTracks,
903 const std::vector<uint32_t>& trackIDs, const std::vector<V2TRef>& v2tRefs,
904 std::vector<o2::MCEventLabel>& lblVtx)
905{
906 lblVtx.clear();
907 if (!lblTracks.size()) {
908 LOG(error) << "Track labels are not provided";
909 return;
910 }
911 std::unordered_map<o2::MCEventLabel, int> labelOccurence;
912
913 auto bestLbl = [](std::unordered_map<o2::MCEventLabel, int> mp, int norm) -> o2::MCEventLabel {
914 o2::MCEventLabel best;
915 int bestCount = 0;
916 for (auto [lbl, cnt] : mp) {
917 if (cnt > bestCount && lbl.isSet()) {
918 bestCount = cnt;
919 best = lbl;
920 }
921 }
922 if (bestCount && norm) {
923 best.setCorrWeight(float(bestCount) / norm);
924 }
925 return best;
926 };
927
928 for (const auto& v2t : v2tRefs) {
929 int tref = v2t.getFirstEntry(), last = tref + v2t.getEntries();
930 labelOccurence.clear();
931 o2::MCEventLabel winner; // unset at the moment
932 for (; tref < last; tref++) {
933 const auto& lbl = lblTracks[mTracksPool[trackIDs[tref]].entry];
934 if (!lbl.isSet()) {
935 break;
936 }
937 // if (!lbl.isFake()) { // RS account all labels, not only correct ones
938 labelOccurence[{lbl.getEventID(), lbl.getSourceID(), 0.}]++;
939 // }
940 }
941 if (labelOccurence.size()) {
942 winner = bestLbl(labelOccurence, v2t.getEntries());
943 }
944 lblVtx.push_back(winner);
945 }
946}
947
948
949//___________________________________________________________________
951{
952 mBunchFilling = bf;
953 // find closest (from above) filled bunch
954 int minBC = bf.getFirstFilledBC(), maxBC = bf.getLastFilledBC();
955 if (minBC < 0) {
956 LOG(error) << "Empty bunch filling is provided, all vertices will be rejected";
957 mClosestBunchAbove[0] = mClosestBunchBelow[0] = -1;
958 return;
959 }
960 int bcAbove = minBC;
961 for (int i = o2::constants::lhc::LHCMaxBunches; i--;) {
962 if (bf.testBC(i)) {
963 bcAbove = i;
964 }
965 mClosestBunchAbove[i] = bcAbove;
966 }
967 int bcBelow = maxBC;
968 for (int i = 0; i < o2::constants::lhc::LHCMaxBunches; i++) {
969 if (bf.testBC(i)) {
970 bcBelow = i;
971 }
972 mClosestBunchBelow[i] = bcBelow;
973 }
974}
975
976//___________________________________________________________________
978{
979 // assign compatible IRs accounting for the bunch filling scheme
980 if (mClosestBunchAbove[0] < 0 && mPVParams->doBCValidation) { // empty or no BF was provided
981 return false;
982 }
983 const auto& vtxT = vtx.getTimeStamp();
984 o2::InteractionRecord irMin(mStartIR), irMax(mStartIR);
985 auto rangeT = std::max(mPVParams->minTError, mPVParams->nSigmaTimeCut * std::min(mPVParams->maxTError, vtxT.getTimeStampError())) + mPVParams->timeMarginVertexTime;
986 float t = vtxT.getTimeStamp() + mPVParams->timeBiasMS;
987 if (t > rangeT) {
988 irMin += o2::InteractionRecord(1.e3 * (t - rangeT));
989 }
990 if (t < -rangeT) {
991 return false; // discard vertex at negative time
992 }
993 irMax += o2::InteractionRecord(1.e3 * (t + rangeT)); // RS TODO: make sure that irMax does not exceed TF edge
994 irMax++; // to account for rounding
995 if (mPVParams->doBCValidation) {
996 // restrict using bunch filling
997 int bc = mClosestBunchAbove[irMin.bc];
998 LOG(debug) << "irMin.bc = " << irMin.bc << " bcAbove = " << bc;
999 if (bc < irMin.bc) {
1000 irMin.orbit++;
1001 }
1002 irMin.bc = bc;
1003 bc = mClosestBunchBelow[irMax.bc];
1004 LOG(debug) << "irMax.bc = " << irMax.bc << " bcBelow = " << bc;
1005 if (bc > irMax.bc) {
1006 if (irMax.orbit == 0) {
1007 return false;
1008 }
1009 irMax.orbit--;
1010 }
1011 irMax.bc = bc;
1012 }
1013 vtx.setIRMin(irMin);
1014 vtx.setIRMax(irMax);
1015 if (irMin > irMax) {
1016 LOG(debug) << "Reject VTX " << vtx.asString() << " trange = " << rangeT;
1017 }
1018 return irMax >= irMin;
1019}
1020
1021//___________________________________________________________________
1022int PVertexer::dbscan_RangeQuery(int id, std::vector<int>& cand, std::vector<int>& status)
1023{
1024 // find neighbours for dbscan cluster core point candidate
1025 // Since we use asymmetric distance definition, is it bit more complex than simple search within chi2 proximity
1026 int nFound = 0;
1027 const auto& tI = mTracksPool[id];
1028 if (tI.sig2ZI < mDBSMaxZ2InvCorePoint) {
1029 return nFound;
1030 }
1031 int ntr = mTracksPool.size();
1032
1033 auto procPnt = [this, &tI, &status, &cand, &nFound, id](int idN) {
1034 const auto& tL = this->mTracksPool[idN];
1035 if (std::abs(tI.timeEst.getTimeStamp() - tL.timeEst.getTimeStamp()) > this->mDBScanDeltaT) {
1036 return -1;
1037 }
1038 auto statN = status[idN], stat = status[id];
1039 if (statN >= 0 && (stat < 0 || (stat >= 0 && statN != stat))) { // do not consider as a neighbour if already added to other cluster
1040 return 0;
1041 }
1042 auto dist2 = tL.getDist2(tI);
1043 if (dist2 < this->mPVParams->dbscanMaxDist2) {
1044 nFound++;
1045 if (statN < 0 && statN > DBS_INCHECK) { // no point in adding for check already assigned point, or which is already in the list (i.e. < INCHECK)
1046 cand.push_back(idN);
1047 status[idN] += DBS_INCHECK; // flag that the track is in the candidates list (i.e. DBS_UDEF-10 = -12 or DPB_NOISE-10 = -11).
1048 }
1049 }
1050 return 1;
1051 };
1052 int idL = id;
1053 while (--idL >= 0) { // index in time decreasing direction
1054 if (procPnt(idL) < 0) {
1055 break;
1056 }
1057 }
1058 int idU = id;
1059 while (++idU < ntr) { // index in time increasing direction
1060 if (procPnt(idU) < 0) {
1061 break;
1062 }
1063 }
1064 return nFound;
1065}
1066
1067//_____________________________________________________
1068void PVertexer::dbscan_clusterize()
1069{
1070 mTimeZClusters.clear();
1071 int ntr = mTracksPool.size();
1072 std::vector<int> status(ntr, DBS_UNDEF);
1073 int clID = -1;
1074
1075 std::vector<int> nbVec;
1076 for (int it = 0; it < ntr; it++) {
1077 if (status[it] != DBS_UNDEF) {
1078 continue;
1079 }
1080 nbVec.clear();
1081 auto nnb0 = dbscan_RangeQuery(it, nbVec, status);
1082 int minNeighbours = mPVParams->minTracksPerVtx - 1;
1083 if (nnb0 < minNeighbours) {
1084 status[it] = DBS_NOISE; // noise
1085 continue;
1086 }
1087 if (nnb0 > minNeighbours) {
1088 minNeighbours = std::max(minNeighbours, int(nnb0 * mPVParams->dbscanAdaptCoef));
1089 }
1090 status[it] = ++clID;
1091 auto& clusVec = mTimeZClusters.emplace_back().trackIDs; // new cluster
1092 clusVec.push_back(it);
1093
1094 for (int j = 0; j < nnb0; j++) {
1095 int jt = nbVec[j];
1096 auto statjt = status[jt];
1097 if (statjt >= 0) {
1098 LOG(error) << "assigned track " << jt << " with status " << statjt << " head is " << it << " clID= " << clID;
1099 continue;
1100 }
1101 status[jt] = clID;
1102 clusVec.push_back(jt);
1103 if (statjt == DBS_NOISE + DBS_INCHECK) { // was border point, no check for being core point is needed
1104 continue;
1105 }
1106 int ncurr = nbVec.size();
1107 if (clusVec.size() > minNeighbours) {
1108 minNeighbours = std::max(minNeighbours, int(clusVec.size() * mPVParams->dbscanAdaptCoef));
1109 }
1110 auto nnb1 = dbscan_RangeQuery(jt, nbVec, status);
1111 if (nnb1 < minNeighbours) {
1112 for (unsigned k = ncurr; k < nbVec.size(); k++) {
1113 if (status[nbVec[k]] < DBS_INCHECK) {
1114 status[nbVec[k]] -= DBS_INCHECK; // remove from checks
1115 }
1116 }
1117 nbVec.resize(ncurr); // not a core point, reset the seeds pool to the state before RangeQuery
1118 } else {
1119 nnb0 = nbVec.size(); // core point, its neighbours need to be checked
1120 }
1121 }
1122 }
1123
1124 for (auto& clus : mTimeZClusters) {
1125 if (clus.trackIDs.size() < mPVParams->minTracksPerVtx) {
1126 clus.trackIDs.clear();
1127 continue;
1128 }
1129 float tMean = 0;
1130 for (const auto tid : clus.trackIDs) {
1131 tMean += mTracksPool[tid].timeEst.getTimeStamp();
1132 }
1133 clus.timeEst.setTimeStamp(tMean / clus.trackIDs.size());
1134 }
1135 mNTZClustersIni = mTimeZClusters.size();
1136}
1137
1138//___________________________________________________________________
1139std::pair<int, int> PVertexer::getBestIR(const PVertex& vtx, const gsl::span<InteractionCandidate> intCand, int& currEntry) const
1140{
1141 // select best matching interaction record
1142 int best = -1, n = intCand.size();
1143 while (currEntry < n && intCand[currEntry] < vtx.getIRMin()) {
1144 currEntry++; // skip all times which have no chance to be matched
1145 }
1146 int i = currEntry, nCompatible = 0;
1147 float bestDf = 1e12;
1148 auto tVtxNS = (vtx.getTimeStamp().getTimeStamp() + mPVParams->timeBiasMS) * 1e3; // time in ns
1149 while (i < n) {
1150 if (intCand[i] > vtx.getIRMax()) {
1151 break;
1152 }
1153 nCompatible++;
1154 auto dfa = std::abs(intCand[i].differenceInBCNS(mStartIR) - tVtxNS);
1155 if (dfa <= bestDf) {
1156 bestDf = dfa;
1157 best = i;
1158 }
1159 i++;
1160 }
1161 return {best, nCompatible};
1162}
1163
1164//___________________________________________________________________
1165SeedHistoTZ PVertexer::buildHistoTZ(const VertexingInput& input)
1166{
1167 // build histo for tracks time / z entries weigthed by their inverse error, to be used for seeding peak finding
1168 // estimat the range of TZ histo
1169
1170 float hZMin = 1e9, hZMax = -1e8, hTMin = 1e9, hTMax = -1e9;
1171 for (int i : input.idRange) {
1172 const auto& trc = mTracksPool[i];
1173 if (trc.canUse()) {
1174 if (trc.z > hZMax) {
1175 hZMax = trc.z;
1176 }
1177 if (trc.z < hZMin) {
1178 hZMin = trc.z;
1179 }
1180 if (trc.timeEst.getTimeStamp() > hTMax) {
1181 hTMax = trc.timeEst.getTimeStamp();
1182 }
1183 if (trc.timeEst.getTimeStamp() < hTMin) {
1184 hTMin = trc.timeEst.getTimeStamp();
1185 }
1186 }
1187 }
1188
1189 float dz = hZMax - hZMin, dt = hTMax - hTMin;
1190 int nbz = 1 + int((dz) / mPVParams->histoBinZSize), nbt = 1 + int((dt) / mPVParams->histoBinTSize);
1191 if (nbz * nbt > 0x7fff) {
1192 float factor = float(nbz * nbt / 0x7fff);
1193 if (nbt > nbz) {
1194 nbt = (0x7fff - 1) / nbz;
1195 } else {
1196 nbz = (0x7fff - 1) / nbt;
1197 }
1198 }
1199 float dzh = 0.5f * (nbz * mPVParams->histoBinZSize - dz), dth = 0.5f * (nbt * mPVParams->histoBinTSize - dt);
1200 SeedHistoTZ seedHistoTZ(nbt, hTMin - dth, hTMax + dth, nbz, hZMin - dzh, hZMax + dzh);
1201
1202 for (int i : input.idRange) {
1203 auto& trc = mTracksPool[i];
1204 if (trc.canUse()) {
1205 trc.bin = seedHistoTZ.fillAndFlagBin(trc.timeEst.getTimeStamp(), trc.z, trc.wghHisto);
1206 }
1207 }
1208
1209 return std::move(seedHistoTZ);
1210}
1211
1212//______________________________________________
1213bool PVertexer::relateTrackToMeanVertex(o2::track::TrackParCov& trc, float vtxErr2)
1214{
1215 o2d::DCA dca;
1216 auto z = trc.getZAt(0., mBz);
1217 if (z < -998.) {
1218 z = mMeanVertex.getZ();
1219 }
1220 mMeanVertex.setMeanXYVertexAtZ(mMeanVertexSeed, z);
1221 if (!o2::base::Propagator::Instance()->propagateToDCA(mMeanVertex, trc, mBz, 2.0f, o2::base::Propagator::MatCorrType::USEMatCorrLUT, &dca, nullptr, 0, mPVParams->dcaTolerance)) {
1222 return false;
1223 }
1224 return dca.getY() * dca.getY() / (dca.getSigmaY2() + vtxErr2) < mPVParams->pullIniCut;
1225}
1226
1227//______________________________________________
1228bool PVertexer::relateTrackToVertex(o2::track::TrackParCov& trc, const o2d::VertexBase& vtxSeed) const
1229{
1230 return o2::base::Propagator::Instance()->propagateToDCA(vtxSeed, trc, mBz, 2.0f, o2::base::Propagator::MatCorrType::USEMatCorrLUT);
1231}
1232
1233//______________________________________________
1234void PVertexer::doDBScanDump(const VertexingInput& input, gsl::span<const o2::MCCompLabel> lblTracks)
1235{
1236 // dump tracks for T-Z clusters identified by the DBScan
1237#ifdef _PV_DEBUG_TREE_
1238 for (int i : input.idRange) {
1239 const auto& trc = mTracksPool[i];
1240 if (trc.canUse()) {
1241 mDebugDumpDBSTrc.emplace_back(trc.z, trc.sig2ZI, trc.timeEst.getTimeStamp(), trc.timeEst.getTimeStampError(), trc.wghHisto);
1242 mDebugDumpDBSGID.push_back(trc.gid);
1243 if (lblTracks.size()) {
1244 mDebugDumpDBSTrcMC.push_back(lblTracks[trc.entry]);
1245 }
1246 }
1247 }
1248 mDebugDBScanTree->Fill();
1249 mDebugDumpDBSTrc.clear();
1250 mDebugDumpDBSGID.clear();
1251 mDebugDumpDBSTrcMC.clear();
1252#endif
1253}
1254
1255//______________________________________________
1256void PVertexer::doDBGPoolDump(gsl::span<const o2::MCCompLabel> lblTracks)
1257{
1258 // dump tracks of the pool
1259#ifdef _PV_DEBUG_TREE_
1260 for (const auto& trc : mTracksPool) {
1261 mDebugDumpDBSTrc.emplace_back(trc.z, trc.sig2ZI, trc.timeEst.getTimeStamp(), trc.timeEst.getTimeStampError(), trc.wghHisto);
1262 mDebugDumpDBSGID.push_back(trc.gid);
1263 if (lblTracks.size()) {
1264 mDebugDumpDBSTrcMC.push_back(lblTracks[trc.entry]);
1265 }
1266 }
1267 mDebugPoolTree->Fill();
1268 mDebugDumpDBSTrc.clear();
1269 mDebugDumpDBSGID.clear();
1270 mDebugDumpDBSTrcMC.clear();
1271#endif
1272}
1273
1274//______________________________________________
1275void PVertexer::doVtxDump(std::vector<PVertex>& vertices, std::vector<uint32_t> trackIDsLoc, std::vector<V2TRef>& v2tRefsLoc,
1276 gsl::span<const o2::MCCompLabel> lblTracks)
1277{
1278 // dump tracks for T-Z clusters identified by the DBScan
1279#ifdef _PV_DEBUG_TREE_
1280 int nv = vertices.size();
1281 for (int iv = 0; iv < nv; iv++) {
1282 mDebugDumpVtx = vertices[iv];
1283 if (mDebugDumpVtx.getNContributors() == 0) { // discarded
1284 continue;
1285 }
1286 int start = v2tRefsLoc[iv].getFirstEntry(), stop = start + v2tRefsLoc[iv].getEntries();
1287 for (int it = start; it < stop; it++) {
1288 const auto& trc = mTracksPool[trackIDsLoc[it]];
1289 mDebugDumpVtxTrc.emplace_back(trc.z, trc.sig2ZI, trc.timeEst.getTimeStamp(), trc.timeEst.getTimeStampError(), trc.wgh);
1290 mDebugDumpVtxGID.push_back(trc.gid);
1291 if (lblTracks.size()) {
1292 mDebugDumpVtxTrcMC.push_back(lblTracks[trc.entry]);
1293 }
1294 }
1295 mDebugVtxTree->Fill();
1296 mDebugDumpVtxTrc.clear();
1297 mDebugDumpVtxGID.clear();
1298 mDebugDumpVtxTrcMC.clear();
1299 }
1300#endif
1301}
1302
1303//______________________________________________
1304PVertex PVertexer::refitVertex(const std::vector<bool> useTrack, const o2d::VertexBase& vtxSeed)
1305{
1306 // Refit the tracks prepared by the successful prepareVertexRefit, possible skipping those tracks wich have useTrack value false
1307 // (useTrack is ignored if empty).
1308 // The vtxSeed is the originally found vertex, must be the same as used for the prepareVertexRefit.
1309 // Refitted PrimaryVertex is returned, negative chi2 means failure of the refit.
1310 // ATTENTION: only the position is refitted, the vertex time and IRMin/IRMax info is dummy.
1311
1312 if (vtxSeed != mVtxRefitOrig) {
1313 throw std::runtime_error("refitVertex must be preceded by successful prepareVertexRefit");
1314 }
1315 VertexingInput inp;
1316 inp.scaleSigma2 = mPVParams->minScale2;
1317 inp.idRange = gsl::span<int>(mRefitTrackIDs);
1318 if (useTrack.size()) {
1319 for (uint32_t i = 0; i < mTracksPool.size(); i++) {
1320 mTracksPool[i].vtxID = useTrack[mTracksPool[i].entry] ? TrackVF::kNoVtx : TrackVF::kDiscarded;
1321 }
1322 }
1323 VertexSeed vtxs;
1324 vtxs.VertexBase::operator=(vtxSeed);
1325 PVertex vtxRes;
1326 vtxs.setScale(inp.scaleSigma2, mTukey2I);
1327 vtxs.setTimeStamp({0.f, -1.}); // time is not refitter
1328 if (fitIteration(inp, vtxs) == FitStatus::OK) {
1329 vtxRes = vtxs;
1330 } else {
1331 vtxRes.setChi2(-1.);
1332 }
1333 return vtxRes;
1334}
1335
1336//______________________________________________
1338{
1339 int cnt = 0;
1340 for (int i : input.idRange) {
1341 if (mTracksPool[i].canUse()) {
1342 LOGP(info, "#{}: {} z:{:.3f}/{:.3f} t:{:.3f}/{:.3f} w:{:.3f} gid:{}", cnt, i, mTracksPool[i].z, std::sqrt(1. / mTracksPool[i].sig2ZI),
1343 mTracksPool[i].timeEst.getTimeStamp(), mTracksPool[i].timeEst.getTimeStampError(), mTracksPool[i].wgh, mTracksPool[i].gid.asString());
1344 }
1345 cnt++;
1346 }
1347}
1348
1349//______________________________________________
1350void PVertexer::dumpPool()
1351{
1352 static int dumpID = 0;
1353 if (mPoolDumpDirectory != "/dev/null") {
1354 TFile dumpFile(fmt::format("{}{}pvtracksPool{}_{}_{}.root", mPoolDumpDirectory, (mPoolDumpDirectory.empty() || mPoolDumpDirectory.back() == '/') ? "" : "/",
1355 dumpID, std::chrono::time_point_cast<std::chrono::milliseconds>(std::chrono::system_clock::now()).time_since_epoch().count(),
1356 o2::utils::Str::getRandomString(5))
1357 .c_str(),
1358 "create");
1359 dumpFile.WriteObjectAny(&mTracksPool, "std::vector<o2::vertexing::TrackVF>", "pool");
1360 LOGP(warn, "Produced tracks pool dump {}", dumpFile.GetName());
1361 }
1362 mPoolDumpProduced = true;
1363}
1364//______________________________________________
1365int PVertexer::processFromExternalPool(const std::vector<TrackVF>& pool, std::vector<PVertex>& vertices, std::vector<o2d::VtxTrackIndex>& vertexTrackIDs, std::vector<V2TRef>& v2tRefs)
1366{
1367 // dummy inputs
1368 std::vector<InteractionCandidate> intCand;
1369 std::vector<o2::MCCompLabel> lblTracks;
1370 std::vector<o2::MCEventLabel> lblVtx;
1371
1372 std::vector<GTrackID> gids;
1373 for (auto tr : pool) {
1374 tr.vtxID = TrackVF::kNoVtx;
1375 tr.wgh = 0.;
1376 mTracksPool.push_back(tr);
1377 gids.push_back(tr.gid);
1378 }
1379 return runVertexing(gids, intCand, vertices, vertexTrackIDs, v2tRefs, lblTracks, lblVtx);
1380}
1381
1382//______________________________________________
1384{
1385 mTrackSrc = s;
1386 // fill vector of sources to consider
1388 int nGloSrc = 0;
1389 for (int is = 0; is < GTrackID::NSources; is++) {
1390 if (mTrackSrc[is]) {
1391 mSrcVec.push_back(is);
1392 if ((GTrackID::getSourceDetectorsMask(is) & ITSTPC) == ITSTPC) {
1393 nGloSrc++;
1394 }
1395 }
1396 }
1397 mITSOnly = nGloSrc == 0;
1398}
uint64_t bc
Definition RawEventData.h:5
int16_t time
Definition RawEventData.h:4
int32_t i
Primary vertex finder.
uint32_t j
Definition RawData.h:0
std::ostringstream debug
int getLastFilledBC(int dir=-1) const
bool testBC(int bcID, int dir=-1) const
int getFirstFilledBC(int dir=-1) const
void setCorrWeight(float corrW)
GPUd() value_type estimateLTFast(o2 static GPUd() float estimateLTIncrement(const o2 PropagatorImpl * Instance(bool uninitialized=false)
Definition Propagator.h:143
void setBinContent(uint32_t bin, T w)
void setSigma(int icoord, float val)
void setMeanXYVertexAtZ(VertexBase &v, float z) const
const InteractionRecord & getIRMin() const
void setIRMin(const InteractionRecord &ir)
const InteractionRecord & getIRMax() const
void setIRMax(const InteractionRecord &ir)
Static class with identifiers, bitmasks and names for ALICE detectors.
Definition DetID.h:58
static constexpr ID ITS
Definition DetID.h:63
static constexpr ID TPC
Definition DetID.h:64
static mask_t getMask(const std::string_view detList)
detector masks from any non-alpha-num delimiter-separated list (empty if NONE is supplied)
Definition DetID.cxx:42
void setTrackSources(GTrackID::mask_t s)
void setBz(float bz)
Definition PVertexer.h:87
void printInpuTracksStatus(const VertexingInput &input) const
bool findVertex(const VertexingInput &input, PVertex &vtx)
PVertex refitVertex(const std::vector< bool > useTrack, const o2d::VertexBase &vtxSeed)
int processFromExternalPool(const std::vector< TrackVF > &pool, std::vector< PVertex > &vertices, std::vector< o2d::VtxTrackIndex > &vertexTrackIDs, std::vector< V2TRef > &v2tRefs)
void setBunchFilling(const o2::BunchFilling &bf)
bool setCompatibleIR(PVertex &vtx)
void setTukey(float t)
Definition PVertexer.h:77
GLdouble n
Definition glcorearb.h:1982
GLint GLsizei count
Definition glcorearb.h:399
GLuint64EXT * result
Definition glcorearb.h:5662
GLenum GLenum dst
Definition glcorearb.h:1767
GLuint start
Definition glcorearb.h:469
GLuint id
Definition glcorearb.h:650
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
Definition A.h:18
constexpr int LHCMaxBunches
D const SVectorGPU< T, D > & rhs
Definition SMatrixGPU.h:191
T MAD2Sigma(int np, T *y)
Definition fit.h:773
return * this
uint64_t getTimeStamp(o2::framework::ProcessingContext &pc)
TrackParCovF TrackParCov
Definition Track.h:33
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
Defining DataPointCompositeObject explicitly as copiable.
std::string to_string(gsl::span< T, Size > span)
Definition common.h:52
FIXME: do not use data model tables.
Common utility functions.
uint32_t orbit
LHC orbit.
uint16_t bc
bunch crossing ID of interaction
static std::string concat_string(Ts const &... ts)
float iniScale2
initial scale to assign
bool applyDebrisReduction
apply algorithm reducing split vertices
float dbscanMaxSigZCorPoint
max sigZ of the track which can be core points in the DBScan
float dbscanMaxDist2
distance^2 cut (eps^2).
float addTimeSigma2DebrisExtra
increment time error^2 by this amount when calculating vertex-to-vertex chi2
int maxNScaleIncreased
max number of scaling-non-decreasing iterations
float pullIniCut
cut on pull (n^2 sigma) on dca to mean vertex
float histoBinZSize
size of the seedTZ histo bin Z
float maxTDiffDebris
when reducing debris, don't consider vertices separated by time > this value in \mus if >0,...
int minNContributorsForIRcutIni
min multiplicity to reject before cleanup if no matching interaction found (if >= 0)
float dbscanDeltaT
abs. time difference cut, should be ~ 0.5 ITS ROF duration if ITS SA tracks used, if < 0 then the val...
float maxChi2TZDebris
don't consider vertices with mutual chi2 exceeding this (for pp should be ~10)
int minTracksPerVtx
min N tracks per vertex
float addZSigma2Debris
increment z error^2 by this amount when calculating vertex-to-vertex chi2
float maxChi2Mean
max mean chi2 of vertex to accept
float maxTError
use min of vertex time error or this for nsigma evaluation
float minScale2
min scaling factor^2
bool doBCValidation
apply validation by interacting BC compatible with the vertex time span
bool applyReattachment
refit vertices reattaching tracks to closest found vertex
float nSigmaTimeCut
eliminate vertex if there is no FT0 or BC signal within this cut
float maxTMAD
max accepted tMAD, not tMAD cleanup if negative
long maxTimeMSPerCluster
max allowed time per TZCluster processing, ms
bool useTimeInChi2
use track-vertex time difference in chi2 calculation
float dcaTolerance
consider tracks within this abs DCA to mean vertex
float maxTDiffDebrisExtra
when reducing debris, don't consider vertices separated by time > this value in \mus if >0,...
float timeMarginVertexTime
additive marginal error to \mus vertex time bracket
float histoBinTSize
size of the seedTZ histo bin T
float dbscanAdaptCoef
adapt dbscan minPts for each cluster as minPts=max(minPts, currentSize*dbscanAdaptCoef).
int minNContributorsForIRcut
do not apply IR cut to vertices below IR tagging efficiency threshold
bool useMeanVertexConstraint
use MeanVertex as extra measured point
float maxMultRatDebris
don't consider vertices with multiplicity ratio above this
float maxITSOnlyFraction
max ITS-only tracks fraction to accept, recommended value for PbPb = 0.85
float meanVertexExtraErrConstraint
extra error to meanvertex sigma used when applying constrant
float tukey
Tukey parameter.
float maxZDiffDebris
don't consider vertices separated by Z > this value in cm
float maxChi2TZDebrisExtra
don't consider vertices with mutual chi2 exceeding this (for pp should be ~10)
float minTError
don't use error smaller than that (~BC/2/minNContributorsForFT0cut)
float timeBiasMS
relative bias in ms to add to TPCITS-based time stamp
float slowConvergenceFactor
consider convergence as slow if ratio new/old scale2 exceeds it
float maxMultRatDebrisExtra
don't consider vertices with multiplicity ratio above this
float maxZMAD
max accepted zMAD, not zMAD cleanup if negative
float timeMarginReattach
safety marging for track time vs vertex time difference during reattachment
float minITSOnlyFraction
min ITS-only tracks fraction to accept
float addTimeSigma2Debris
increment time error^2 by this amount when calculating vertex-to-vertex chi2
int maxIterations
max iterations per vertex fit
int maxNScaleSlowConvergence
max number of weak scaling decrease iterations
float acceptableScale2
if below this factor, try to refit with minScale2
float addZSigma2DebrisExtra
increment z error^2 by this amount when calculating vertex-to-vertex chi2
float maxZDiffDebrisExtra
don't consider vertices separated by Z > this value in cm
float tgL
tangent(lambda)
float tgP
tangent(phi) in tracking frame
float sigYZI
YZ component of inverse cov.matrix.
float evalChi2ToVertex(const PVertex &vtx, bool useTime)
int vtxID
assigned vertex
float wgh
track weight wrt current vertex seed
float sig2YI
YY component of inverse cov.matrix.
float cosAlp
cos of alpha frame
float sig2ZI
ZZ component of inverse cov.matrix.
float sinAlp
sin of alpha frame
void setScale(float scale2, float tukey2I)
weights and scaling params for current vertex
std::chrono::system_clock system_clock
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"