17#include <TStopwatch.h>
66 : mDataRequest(dr), mGGCCDBRequest(gr), mTracksSrc(
src), mUseMC(useMC)
83 bool refitPV(
o2::dataformats::PrimaryVertex& pv,
int vid);
85 bool processITSTrack(const
o2::its::
TrackITS& iTrack, const
o2::dataformats::PrimaryVertex& pv,
o2::checkresid::
Track& resTrack);
87 o2::globaltracking::RecoContainer* mRecoData =
nullptr;
89 bool mMeanVertexUpdated = false;
90 float mITSROFrameLengthMUS = 0.
f;
91 o2::dataformats::MeanVertexObject mMeanVtx{};
92 std::vector<o2::BaseCluster<float>> mITSClustersArray;
95 std::shared_ptr<DataRequest> mDataRequest;
96 std::shared_ptr<o2::base::GRPGeomRequest> mGGCCDBRequest;
98 std::unique_ptr<o2::utils::TreeStreamRedirector> mDBGOut;
107 std::string dbgnm = maxLanes == 1 ?
"checkResid.root" : fmt::format(
"checkResid_t{}.root", lane);
108 mDBGOut = std::make_unique<o2::utils::TreeStreamRedirector>(dbgnm.c_str(),
"recreate");
109 mNThreads = ic.
options().
get<
int>(
"nthreads");
112 LOGP(warn,
"No OpenMP");
122 mRecoData = &recoData;
124 mRecoData = &recoData;
125 updateTimeDependentParams(pc);
136 static bool initOnceDone =
false;
143 if (!grp->isDetContinuousReadOut(
DetID::ITS)) {
144 mITSROFrameLengthMUS = alpParams.roFrameLengthTrig / 1.e3;
153 if (mMeanVertexUpdated) {
154 mMeanVertexUpdated =
false;
158 bool updateMaps =
false;
181 LOGP(fatal,
"ITS data is not loaded");
189 auto pattIt = patterns.begin();
190 mITSClustersArray.clear();
191 mITSClustersArray.reserve(clusITS.size());
199 static int TFCount = 0;
200 int nv = vtxRefs.size() - 1;
201 std::vector<std::vector<checkresid::Track>> slots;
202 slots.resize(mNThreads);
203 int nvGood = 0, nvUse = 0, nvRefFail = 0;
204 long pvFitDuration{};
205 for (
int iv = 0; iv < nv; iv++) {
206 const auto& vtref = vtxRefs[iv];
207 auto pve = pvvec[iv];
208 if (pve.getNContributors() <
params.minPVContributors) {
213 LOGP(
debug,
"Refitting PV#{} of {} tracks", iv, pve.getNContributors());
214 auto tStartPVF = std::chrono::time_point_cast<std::chrono::microseconds>(std::chrono::system_clock::now()).time_since_epoch().count();
215 bool res = refitPV(pve, iv);
216 pvFitDuration += std::chrono::time_point_cast<std::chrono::microseconds>(std::chrono::system_clock::now()).time_since_epoch().count() - tStartPVF;
227 int idMin = vtref.getFirstEntryOfSource(is), idMax = idMin + vtref.getEntriesOfSource(is);
233 LOGP(fatal,
"Cut on TPC tracks is requested by they are not loaded");
236#pragma omp parallel for schedule(dynamic) num_threads(mNThreads)
238 for (
int i = idMin;
i < idMax;
i++) {
239 auto vid = trackIndex[
i];
240 bool pvCont = vid.isPVContributor();
241 if (!pvCont &&
params.pvcontribOnly) {
252 auto pt = trc.getPt();
253 if (pt < params.minPt || pt >
params.maxPt) {
256 const auto& itsTrack = mRecoData->
getITSTrack(gidITS);
257 if (itsTrack.getNClusters() <
params.minITSCl) {
261 auto& accum = slots[omp_get_thread_num()];
263 auto& accum = slots[0];
265 auto& resTrack = accum.emplace_back();
267 if (!processITSTrack(itsTrack, pve, resTrack)) {
275 for (
const auto& accum : slots) {
276 for (
const auto& tr : accum) {
277 (*mDBGOut) <<
"res" <<
"tr=" << tr <<
"\n";
280 LOGP(info,
"processed {} PVs out of {} good vertices (out of {} in total), PV refits took {} mus, {} refits failed", nvUse, nvGood, nv, pvFitDuration, nvRefFail);
287 auto trFitInw = iTrack.getParamOut();
288 auto trFitOut = iTrack.getParamIn();
292 float bz = prop->getNominalBz();
293 std::array<const o2::BaseCluster<float>*, 8> clArr{};
295 std::array<o2::track::TrackParCov, 8> extrapOut, extrapInw;
298 return refLin ? tr.rotate(
alpha, *refLin, bz) : tr.rotate(
alpha);
303 if (!rotateTrack(tr,
i == 0 ? pvAlpha : geom->getSensorRefAlpha(clArr[
i]->getSensorID()), refLin) ||
304 !prop->propagateTo(tr, refLin, clArr[
i]->getX(), true)) {
308 if (!tr.update(*clArr[
i])) {
312 extrapDest[
i].invalidate();
318 auto inv2d = [](
float s00,
float s11,
float s01) -> std::array<float, 3> {
319 auto det = s00 * s11 - s01 * s01;
321 LOGP(error,
"Singular det {}, input: {} {} {}", det, s00, s11, s01);
322 return {0.f, 0.f, 0.f};
325 return {s11 * det, s00 * det, -s01 * det};
329 if (!prop->propagateToDCA(pv, trFitOut, bz)) {
330 LOGP(
debug,
"Failed to propagateToDCA, {}", trFitOut.asString());
334 if (
params.addPVAsCluster) {
335 float cosAlp, sinAlp;
336 pvAlpha = trFitOut.getAlpha();
337 o2::math_utils::sincos(trFitOut.getAlpha(), sinAlp, cosAlp);
338 bcPV.
setXYZ(pv.getX() * cosAlp + pv.getY() * sinAlp, -pv.getX() * sinAlp + pv.getY() * cosAlp, pv.getZ());
339 bcPV.
setSigmaY2(0.5 * (pv.getSigmaX2() + pv.getSigmaY2()));
345 int nCl = iTrack.getNClusters();
346 for (
int i = 0;
i < nCl;
i++) {
347 const auto& curClu = mITSClustersArray[itsClRefs[iTrack.
getClusterEntry(
i)]];
349 int llr = geom->getLayer(curClu.getSensorID());
350 if (clArr[1 + llr]) {
351 LOGP(error,
"Cluster at lr {} was already assigned, old sens {}, new sens {}", llr, clArr[1 + llr]->
getSensorID(), curClu.getSensorID());
353 clArr[1 + geom->getLayer(curClu.getSensorID())] = &curClu;
356 o2::track::TrackPar refLinIBOut0, refLinOBInw0, *refLinOBInw =
nullptr, *refLinIBOut =
nullptr;
357 if (
params.useStableRef) {
358 refLinOut = &(refLinOut0 = trFitOut);
359 refLinInw = &(refLinInw0 = trFitInw);
361 trFitOut.resetCovariance();
362 trFitOut.setCov(trFitOut.getQ2Pt() * trFitOut.getQ2Pt() * trFitOut.getCov()[14], 14);
363 trFitInw.resetCovariance();
364 trFitInw.setCov(trFitInw.getQ2Pt() * trFitInw.getQ2Pt() * trFitInw.getCov()[14], 14);
366 for (
int i = 0;
i <= 7;
i++) {
369 if (!(resOut = accountCluster(
i, extrapOut, trFitOut, refLinOut)) || !(resInw = accountCluster(7 -
i, extrapInw, trFitInw, refLinInw))) {
375 if (
i == 3 && resOut == 1 && resInw == 1 &&
params.doIBOB && nCl == 7) {
380 refLinIBOut = &(refLinIBOut0 = refLinOut0);
381 refLinOBInw = &(refLinOBInw0 = refLinInw0);
384 if (!resTrack.
trOBInw.getXatLabR(
params.rCompIBOB, xRref, bz) ||
385 !prop->propagateTo(resTrack.
trOBInw, refLinOBInw, xRref,
true) ||
386 !rotateTrack(resTrack.
trOBInw, resTrack.
trOBInw.getPhiPos(), refLinOBInw) ||
387 !rotateTrack(resTrack.
trIBOut, resTrack.
trOBInw.getAlpha(), refLinIBOut) ||
388 !prop->propagateTo(resTrack.
trIBOut, refLinIBOut, resTrack.
trOBInw.getX(),
true)) {
395 bool innerDone =
false;
397 for (
int i = 0;
i <= 7;
i++) {
400 const auto &tInw = extrapInw[
i], &tOut = extrapOut[
i];
401 auto wInw = inv2d(tInw.getSigmaY2(), tInw.getSigmaZ2(), tInw.getSigmaZY());
402 auto wOut = inv2d(tOut.getSigmaY2(), tOut.getSigmaZ2(), tOut.getSigmaZY());
403 if (wInw[0] == 0.f || wOut[0] == 0.f) {
406 std::array<float, 3> wTot = {wInw[0] + wOut[0], wInw[1] + wOut[1], wInw[2] + wOut[2]};
407 auto cTot = inv2d(wTot[0], wTot[1], wTot[2]);
408 auto ywi = wInw[0] * tInw.getY() + wInw[2] * tInw.getZ() + wOut[0] * tOut.getY() + wOut[2] * tOut.getZ();
409 auto zwi = wInw[2] * tInw.getY() + wInw[1] * tInw.getZ() + wOut[2] * tOut.getY() + wOut[1] * tOut.getZ();
410 auto yw = ywi * cTot[0] + zwi * cTot[2];
411 auto zw = ywi * cTot[2] + zwi * cTot[1];
413 auto phi =
i == 0 ? tInw.getPhi() : tInw.getPhiPos();
414 o2::math_utils::bringTo02Pi(phi);
415 resTrack.
points.emplace_back(clArr[
i]->
getY() - yw, clArr[
i]->getZ() - zw, cTot[0] + clArr[
i]->getSigmaY2(), cTot[1] + clArr[
i]->getSigmaZ2(), phi, clArr[
i]->getZ(), clArr[
i]->
getSensorID(),
i - 1);
417 resTrack.
track = tInw;
421 LOGP(
debug,
"No cluster on lr {}",
i);
431 std::vector<o2::track::TrackParCov>
tracks;
432 std::vector<bool> useTrack;
433 std::vector<GTrackID> gidsITS;
434 int ntr = pv.getNContributors(), ntrIni = ntr;
436 useTrack.reserve(ntr);
437 gidsITS.reserve(ntr);
440 int itr = vtref.getFirstEntry(), itLim = itr + vtref.getEntries();
441 for (; itr < itLim; itr++) {
442 auto vid = trackIndex[itr];
443 if (vid.isPVContributor()) {
449 useTrack.resize(ntr);
451#pragma omp parallel for schedule(dynamic) num_threads(mNThreads)
453 for (
int itr = 0; itr < ntr; itr++) {
454 if (!(useTrack[itr] = refitITStrack(tracks[itr], gidsITS[itr]))) {
459 for (
auto v : useTrack) {
463 LOGP(warn,
"Abandon vertex refit: NcontribNew = {} vs NcontribOld = {}", ntr, ntrIni);
466 LOGP(
debug,
"Original vtx: Nc:{} {}, chi2={}", pv.getNContributors(), pv.
asString(), pv.getChi2());
469 LOGP(
debug,
"Refitted vtx: Nc:{} {}, chi2={}", ntr, pv.
asString(), pv.getChi2());
470 if (pv.getChi2() < 0.f) {
471 LOGP(warn,
"Failed to refit PV {}", pvSave.asString());
483 auto pid = track.getPID();
484 track = trkITS.getParamOut();
485 track.resetCovariance();
486 track.setCov(track.getQ2Pt() * track.getQ2Pt() * track.getCov()[14], 14);
488 auto nCl = trkITS.getNumberOfClusters();
491 float bz = prop->getNominalBz();
494 for (
int iCl = 0; iCl < nCl; iCl++) {
495 const auto& cls = mITSClustersArray[itsClRefs[trkITS.getClusterEntry(iCl)]];
496 auto alpha = geom->getSensorRefAlpha(cls.getSensorID());
497 if (!(
params.useStableRef ? track.rotate(
alpha, refLin, bz) : track.rotate(
alpha)) ||
498 !prop->propagateTo(track,
params.useStableRef ? &refLin : nullptr, cls.
getX(), true)) {
499 LOGP(
debug,
"refitITStrack failed on propagation to cl#{}, alpha={}, x={} | {}", iCl,
alpha, cls.getX(), track.asString());
502 if (!track.update(cls)) {
503 LOGP(
debug,
"refitITStrack failed on update with cl#{}, | {}", iCl, track.asString());
531 mMeanVertexUpdated =
true;
535 LOG(info) <<
"cluster dictionary updated";
543 std::vector<OutputSpec> outputs;
544 auto dataRequest = std::make_shared<DataRequest>();
545 dataRequest->requestTracks(srcTracks, useMC);
546 dataRequest->requestClusters(srcClusters, useMC);
547 dataRequest->requestPrimaryVertices(useMC);
548 auto ggRequest = std::make_shared<o2::base::GRPGeomRequest>(
false,
556 dataRequest->inputs.emplace_back(
"meanvtx",
"GLO",
"MEANVERTEX", 0, Lifetime::Condition,
ccdbParamSpec(
"GLO/Calib/MeanVertex", {}, 1));
558 {
"nthreads", VariantType::Int, 1, {
"number of threads"}},
567 AlgorithmSpec{adaptFromTask<CheckResidSpec>(dataRequest, ggRequest, srcTracks, useMC )},
Definition of the GeometryManager class.
Helper for geometry and GRP related CCDB requests.
Global index for barrel track: provides provenance (detectors combination), index in respective array...
Definition of the GeometryTGeo class.
Utility functions for MC particles.
Definition of the Names Generator class.
Wrapper container for different reconstructed object types.
o2::track::TrackParCov TrackParCov
Result of refitting TPC-ITS matched track.
Reference on ITS/MFT clusters set.
Referenc on track indices contributing to the vertex, with possibility chose tracks from specific sou...
void setSensorID(std::int16_t sid)
void setXYZ(T x, T y, T z)
void checkUpdates(o2::framework::ProcessingContext &pc)
static GRPGeomHelper & instance()
void setRequest(std::shared_ptr< GRPGeomRequest > req)
GPUd() value_type estimateLTFast(o2 static GPUd() float estimateLTIncrement(const o2 PropagatorImpl * Instance(bool uninitialized=false)
void endOfStream(EndOfStreamContext &ec) final
This is invoked whenever we have an EndOfStream event.
CheckResidSpec(std::shared_ptr< DataRequest > dr, std::shared_ptr< o2::base::GRPGeomRequest > gr, GTrackID::mask_t src, bool useMC)
~CheckResidSpec() final=default
void run(ProcessingContext &pc) final
void finaliseCCDB(ConcreteDataMatcher &matcher, void *obj) final
void init(InitContext &ic) final
static const CheckResidConfig & Instance()
static void updateFromString(std::string const &)
Static class with identifiers, bitmasks and names for ALICE detectors.
T get(const char *key) const
ServiceRegistryRef services()
ConfigParamRegistry const & options()
InputRecord & inputs()
The inputs associated with this processing context.
static GeometryTGeo * Instance()
int getClusterEntry(int i) const
void initMeanVertexConstraint()
PVertex refitVertexFull(const std::vector< bool > useTrack, const o2d::VertexBase &vtxSeed)
bool prepareVertexRefit(const TR &tracks, const o2d::VertexBase &vtxSeed)
void setMeanVertex(const o2d::MeanVertexObject *v)
GLfloat GLfloat GLfloat alpha
GLenum const GLfloat * params
o2::framework::DataProcessorSpec getCheckResidSpec(o2::dataformats::GlobalTrackID::mask_t srcTracks, o2::dataformats::GlobalTrackID::mask_t srcClus, bool useMC)
create a processor spec
constexpr double LHCBunchSpacingNS
Defining PrimaryVertex explicitly as messageable.
std::vector< ConfigParamSpec > ccdbParamSpec(std::string const &path, int runDependent, std::vector< CCDBMetadata > metadata={}, int qrate=0)
std::vector< ConfigParamSpec > Options
void convertCompactClusters(gsl::span< const itsmft::CompClusterExt > clusters, gsl::span< const unsigned char >::iterator &pattIt, std::vector< o2::BaseCluster< float > > &output, const itsmft::TopologyDictionary *dict)
convert compact clusters to 3D spacepoints
double * getX(double *xyDxy, int N)
double * getY(double *xyDxy, int N)
struct o2::upgrades_utils::@453 tracks
structure to keep trigger-related info
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
o2::track::TrackPar track
o2::track::TrackParCov trIBOut
o2::track::TrackParCov trOBInw
std::vector< Point > points
auto getITSTracks() const
GTrackID getITSContributorGID(GTrackID source) const
bool isTrackSourceLoaded(int src) const
const o2::tpc::TrackTPC & getTPCTrack(GTrackID id) const
auto getITSTracksClusterRefs() const
auto getPrimaryVertices() const
auto getPrimaryVertexMatchedTracks() const
auto getPrimaryVertexMatchedTrackRefs() const
auto getITSClustersPatterns() const
GTrackID getTPCContributorGID(GTrackID source) const
const o2::track::TrackParCov & getTrackParam(GTrackID gidx) const
void collectData(o2::framework::ProcessingContext &pc, const DataRequest &request)
const o2::its::TrackITS & getITSTrack(GTrackID gid) const
auto getITSClusters() const
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"