Project
Loading...
Searching...
No Matches
MatchGlobalFwd.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
13#include <queue>
14
15using namespace o2::globaltracking;
16
17//_________________________________________________________
19{
20
21 LOG(info) << "Initializing Global Forward Matcher";
22
23 auto& matchingParam = GlobalFwdMatchingParam::Instance();
24
25 setMFTRadLength(matchingParam.MFTRadLength);
26 LOG(info) << "MFT Radiation Length = " << mMFTDiskThicknessInX0 * 5.;
27
28 setAlignResiduals(matchingParam.alignResidual);
29 LOG(info) << "MFT Align residuals = " << mAlignResidual;
30
31 mMatchingPlaneZ = matchingParam.matchPlaneZ;
32 LOG(info) << "MFTMCH matchingPlaneZ = " << mMatchingPlaneZ;
33
34 auto& matchingFcnStr = matchingParam.matchFcn;
35 LOG(info) << "Match function string = " << matchingFcnStr;
36
37 if (matchingParam.isMatchUpstream()) {
38 LOG(info) << " ==> Setting Upstream matching.";
39 mMatchingType = MATCHINGUPSTREAM;
40 } else if (matchingParam.matchingExternalFunction()) {
41 loadExternalMatchingFunction();
42 mMatchingType = MATCHINGFUNC;
43 } else {
44 if (mMatchingFunctionMap.find(matchingFcnStr) != mMatchingFunctionMap.end()) {
45 mMatchFunc = mMatchingFunctionMap[matchingFcnStr];
46 mMatchingType = MATCHINGFUNC;
47 LOG(info) << " Found built-in matching function " << matchingFcnStr;
48 } else {
49 throw std::invalid_argument("Invalid matching function! Aborting...");
50 }
51 }
52
53 auto& cutFcnStr = matchingParam.cutFcn;
54 LOG(info) << "MFTMCH pair candidate cut function string = " << cutFcnStr;
55
56 if (matchingParam.cutExternalFunction()) {
57 loadExternalCutFunction();
58 } else if (mCutFunctionMap.find(cutFcnStr) != mCutFunctionMap.end()) {
59 mCutFunc = mCutFunctionMap[cutFcnStr];
60 LOG(info) << " Found built-in cut function " << cutFcnStr;
61 } else {
62 throw std::invalid_argument("Invalid cut function! Aborting...");
63 }
64
65 mUseMIDMCHMatch = matchingParam.useMIDMatch;
66 LOG(info) << "UseMIDMCH Matching = " << (mUseMIDMCHMatch ? "true" : "false");
67
68 mUseTrackTime = matchingParam.useTrackTime;
69 LOG(info) << "Use track time = " << (mUseTrackTime ? "true" : "false");
70
71 mSaveMode = matchingParam.saveMode;
72 LOG(info) << "Save mode MFTMCH candidates = " << mSaveMode;
73
74 mNCandidates = matchingParam.nCandidates;
75}
76
77//_________________________________________________________
79{
80
81 auto& matchingParam = GlobalFwdMatchingParam::Instance();
82
83 mRecoCont = &inp;
84 mStartIR = inp.startIR;
85
86 clear();
87
88 if (!prepareMFTData() || !prepareMCHData() || !processMCHMIDMatches()) {
89 return;
90 }
91
92 if (matchingParam.MCMatching) { // MC Label matching
93 mMCTruthON ? doMCMatching() : throw std::runtime_error("Label matching requries MC Labels!");
94 } else {
95 switch (mMatchingType) {
96 case MATCHINGFUNC:
97 switch (mSaveMode) {
98 case kBestMatch:
99 doMatching<kBestMatch>();
100 break;
101 case kSaveAll:
102 doMatching<kSaveAll>();
103 break;
105 doMatching<kSaveTrainingData>();
106 break;
107 case kSaveNCandidates:
108 doMatching<kSaveNCandidates>();
109 break;
110 default:
111 LOG(fatal) << "Invalid MFTMCH save mode";
112 }
113 break;
114 case MATCHINGUPSTREAM:
115 loadMatches();
116 break;
117 default:
118 LOG(fatal) << "Invalid MFTMCH matching mode";
119 }
120 }
121
122 fitTracks();
123 finalize();
124}
125
126//_________________________________________________________
128{
129 LOG(info) << " Finalizing GlobalForwardMatch. Pushing " << mMatchedTracks.size() << " matched tracks";
130}
131
132//_________________________________________________________
134{
135 mMCHROFTimes.clear();
136 mMCHWork.clear();
137 mMFTROFTimes.clear();
138 mMFTWork.clear();
139 mMFTClusters.clear();
140 mMatchedTracks.clear();
141 mMatchLabels.clear();
142 mMFTTrackROFContMapping.clear();
143 mMatchingInfo.clear();
144 mCandidates.clear();
145}
146
147//_________________________________________________________
148bool MatchGlobalFwd::prepareMCHData()
149{
150 const auto& inp = *mRecoCont;
151
152 // Load MCH tracks
153 mMCHTracks = inp.getMCHTracks();
154 mMCHTrackROFRec = inp.getMCHTracksROFRecords();
155 if (mMCTruthON) {
156 mMCHTrkLabels = inp.getMCHTracksMCLabels();
157 }
158 int nROFs = mMCHTrackROFRec.size();
159 LOG(info) << "Loaded " << mMCHTracks.size() << " MCH Tracks in " << nROFs << " ROFs";
160 if (mMCHTracks.empty()) {
161 return false;
162 }
163 mMCHWork.reserve(mMCHTracks.size());
164 mMCHID2Work.clear();
165 mMCHID2Work.resize(mMCHTracks.size(), -1);
166 static int BCDiffErrCount = 0;
167 constexpr int MAXBCDiffErrCount = 2;
168
169 for (int irof = 0; irof < nROFs; irof++) {
170 const auto& rofRec = mMCHTrackROFRec[irof];
171
172 int nBC = rofRec.getBCData().differenceInBC(mStartIR);
173 if (nBC < 0) {
174 if (BCDiffErrCount++ < MAXBCDiffErrCount) {
175 LOGP(alarm, "wrong bunches diff. {} for current IR {} wrt 1st TF orbit {} in MCH data", nBC, rofRec.getBCData().asString(), mStartIR.asString());
176 }
177 }
179 float tMax = (nBC + rofRec.getBCWidth()) * o2::constants::lhc::LHCBunchSpacingMUS;
180 auto mchTime = rofRec.getTimeMUS(mStartIR).first;
181
182 mMCHROFTimes.emplace_back(tMin, tMax); // MCH ROF min/max time
183 LOG(debug) << "MCH ROF # " << irof << " " << rofRec.getBCData() << " [tMin;tMax] = [" << tMin << ";" << tMax << "]";
184 int trlim = rofRec.getFirstIdx() + rofRec.getNEntries();
185 for (int it = rofRec.getFirstIdx(); it < trlim; it++) {
186 auto& trcOrig = mMCHTracks[it];
187 int nWorkTracks = mMCHWork.size();
188 mMCHID2Work[it] = nWorkTracks;
189 // working copy MCH track propagated to matching plane and converted to the forward track format
190 o2::mch::TrackParam tempParam(trcOrig.getZ(), trcOrig.getParameters(), trcOrig.getCovariances());
191 if (!o2::mch::TrackExtrap::extrapToVertexWithoutBranson(tempParam, mMatchingPlaneZ)) {
192 LOG(warning) << "MCH track propagation to matching plane failed!";
193 continue;
194 }
195 auto convertedTrack = MCHtoFwd(tempParam);
196 auto& thisMCHTrack = mMCHWork.emplace_back(TrackLocMCH{convertedTrack, {tMin, tMax}});
197 thisMCHTrack.setMCHTrackID(it);
198 thisMCHTrack.setTimeMUS(mchTime);
199 }
200 }
201 return true;
202}
203
204//_________________________________________________________
205bool MatchGlobalFwd::processMCHMIDMatches()
206{
207 if (mUseMIDMCHMatch) {
208 const auto& inp = *mRecoCont;
209
210 // Load MCHMID matches
211 mMCHMIDMatches = inp.getMCHMIDMatches();
212
213 LOG(info) << "Loaded " << mMCHMIDMatches.size() << " MCHMID matches";
214
215 for (const auto& MIDMatch : mMCHMIDMatches) {
216 const auto& MCHId = MIDMatch.getMCHRef().getIndex();
217 const auto& MIDId = MIDMatch.getMIDRef().getIndex();
218 auto& thisMuonTrack = mMCHWork[mMCHID2Work[MCHId]];
219 LOG(debug) << " MCHId: " << MCHId << " --> mMCHID2Work[MCHId]:" << mMCHID2Work[MCHId];
220 const auto& IR = MIDMatch.getIR();
221 int nBC = IR.differenceInBC(mStartIR);
222 float tMin = (nBC - 1) * o2::constants::lhc::LHCBunchSpacingMUS;
223 float tMax = (nBC + 2) * o2::constants::lhc::LHCBunchSpacingMUS;
224 thisMuonTrack.setMIDTrackID(MIDId);
225 thisMuonTrack.setTimeMUS(MIDMatch.getTimeMUS(mStartIR).first);
226 thisMuonTrack.tBracket.set(tMin, tMax);
227 thisMuonTrack.setMIDMatchingChi2(MIDMatch.getMatchChi2OverNDF());
228 }
229 }
230 return true;
231}
232
233//_________________________________________________________
234bool MatchGlobalFwd::prepareMFTData()
235{
236 const auto& inp = *mRecoCont;
237
238 // MFT clusters
239 mMFTClusterROFRec = inp.getMFTClustersROFRecords();
240 mMFTTrackClusIdx = inp.getMFTTracksClusterRefs();
241 const auto clusMFT = inp.getMFTClusters();
242 if (mMFTClusterROFRec.empty() || clusMFT.empty()) {
243 return false;
244 }
245 const auto patterns = inp.getMFTClustersPatterns();
246 auto pattIt = patterns.begin();
247 mMFTClusters.reserve(clusMFT.size());
248 o2::mft::ioutils::convertCompactClusters(clusMFT, pattIt, mMFTClusters, mMFTDict);
249
250 // Load MFT tracks
251 mMFTTracks = inp.getMFTTracks();
252 mMFTTrackROFRec = inp.getMFTTracksROFRecords();
253 if (mMCTruthON) {
254 mMFTTrkLabels = inp.getMFTTracksMCLabels();
255 }
256 int nROFs = mMFTTrackROFRec.size();
257
258 LOG(info) << "Loaded " << mMFTTracks.size() << " MFT Tracks in " << nROFs << " ROFs";
259 if (mMFTTracks.empty()) {
260 return false;
261 }
262 mMFTWork.reserve(mMFTTracks.size());
263 static int BCDiffErrCount = 0;
264 constexpr int MAXBCDiffErrCount = 2;
265
266 for (int irof = 0; irof < nROFs; irof++) {
267 const auto& rofRec = mMFTTrackROFRec[irof];
268 int nBC = rofRec.getBCData().differenceInBC(mStartIR);
269 if (nBC < 0) {
270 if (BCDiffErrCount++ < MAXBCDiffErrCount) {
271 LOGP(alarm, "TF dropped: wrong bunches diff. {} for current IR {} wrt 1st TF orbit {} in MFT data", nBC, rofRec.getBCData().asString(), mStartIR.asString());
272 }
273 return false;
274 }
275 float tMin = (nBC + mMFTROFrameBiasInBC) * o2::constants::lhc::LHCBunchSpacingMUS;
276 float tMax = (nBC + mMFTROFrameLengthInBC + mMFTROFrameBiasInBC) * o2::constants::lhc::LHCBunchSpacingMUS;
277 if (!mMFTTriggered) {
278 auto irofCont = (nBC + mMFTROFrameBiasInBC) / mMFTROFrameLengthInBC;
279 if (mMFTTrackROFContMapping.size() <= irofCont) { // there might be gaps in the non-empty rofs, this will map continuous ROFs index to non empty ones
280 mMFTTrackROFContMapping.resize((1 + irofCont / 128) * 128, 0);
281 }
282 mMFTTrackROFContMapping[irofCont] = irof;
283 }
284 mMFTROFTimes.emplace_back(tMin, tMax); // MFT ROF min/max time
285 LOG(debug) << "MFT ROF # " << irof << " " << rofRec.getBCData() << " [tMin;tMax] = [" << tMin << ";" << tMax << "]";
286
287 int trlim = rofRec.getFirstEntry() + rofRec.getNEntries();
288 for (int it = rofRec.getFirstEntry(); it < trlim; it++) {
289 const auto& trcOrig = mMFTTracks[it];
290
291 int nWorkTracks = mMFTWork.size();
292 // working copy of outer track param
293 auto& trc = mMFTWork.emplace_back(TrackLocMFT{trcOrig, {tMin, tMax}, irof});
294 trc.setParameters(trcOrig.getOutParam().getParameters());
295 trc.setZ(trcOrig.getOutParam().getZ());
296 trc.setCovariances(trcOrig.getOutParam().getCovariances());
297 trc.setTrackChi2(trcOrig.getOutParam().getTrackChi2());
298 // Extrapolate MFT track parameters and covariances matrix to "mMatchingPlaneZ"
299 // Parameters: helix track model; Error propagation: Quadratic
300 // If "mBz" is zero: linear track model
301 trc.propagateToZ(mMatchingPlaneZ, mBz);
302 }
303 }
304
305 return true;
306}
307
308//_________________________________________________________
309void MatchGlobalFwd::loadMatches()
310{
311
312 const auto& inp = *mRecoCont;
313 int nFakes = 0, nTrue = 0;
314
315 // Load MFT-MCH matching info
316 mMatchingInfoUpstream = inp.getMFTMCHMatches();
317
318 LOG(info) << "Loaded " << mMatchingInfoUpstream.size() << " MFTMCH Matches";
319
320 for (const auto& match : mMatchingInfoUpstream) {
321 auto MFTId = match.getMFTTrackID();
322 auto MCHId = match.getMCHTrackID();
323 LOG(debug) << " ==> MFTId = " << MFTId << " MCHId = " << MCHId << std::endl;
324
325 auto& thisMCHTrack = mMCHWork[mMCHID2Work[MCHId]];
326 thisMCHTrack.setMatchInfo(match);
327 mMatchedTracks.emplace_back(thisMCHTrack);
328 if (mMCTruthON) {
329 mMatchLabels.push_back(computeLabel(MCHId, MFTId));
330 mMatchLabels.back().isFake() ? nFakes++ : nTrue++;
331 }
332 }
333
334 LOG(info) << " Done matching from upstream " << mMFTWork.size() << " MFT tracks with " << mMCHWork.size() << " MCH Tracks.";
335 if (mMCTruthON) {
336 LOG(info) << " nFakes = " << nFakes << " nTrue = " << nTrue;
337 }
338}
339
340//_________________________________________________________
341template <Int_t saveAllMode>
342void MatchGlobalFwd::doMatching()
343{
344 // Range of compatible MCH ROFS for the first MFT track
345 int nMCHROFs = mMCHROFTimes.size();
346
347 LOG(info) << "Running MCH-MFT Track Matching.";
348 // ROFrame of first MFT track
349 auto firstMFTTrackIdInROF = 0;
350 auto MFTROFId = mMFTWork.front().roFrame;
351 LOG(debug) << "(*) nMCHROFs: " << nMCHROFs << ", mMFTTracks.size(): " << mMFTTracks.size() << " MFTROFId: " << MFTROFId << ", mMFTTrackROFRec.size(): " << mMFTTrackROFRec.size();
352
353 while ((firstMFTTrackIdInROF < mMFTTracks.size()) && (MFTROFId < mMFTTrackROFRec.size())) {
354 auto MFTROFId = mMFTWork[firstMFTTrackIdInROF].roFrame;
355 const auto& thisMFTBracket = mMFTROFTimes[MFTROFId];
356 auto nMFTTracksInROF = mMFTTrackROFRec[MFTROFId].getNEntries();
357 firstMFTTrackIdInROF = mMFTTrackROFRec[MFTROFId].getFirstEntry();
358 LOG(debug) << "MFT ROF = " << MFTROFId << "; interval: [" << thisMFTBracket.getMin() << "," << thisMFTBracket.getMax() << "]";
359 LOG(debug) << "ROF " << MFTROFId << " : firstMFTTrackIdInROF " << firstMFTTrackIdInROF << " ; nMFTTracksInROF = " << nMFTTracksInROF;
360 firstMFTTrackIdInROF += nMFTTracksInROF;
361
362 int mchROFMatchFirst = -1;
363 int mchROFMatchLast = -1;
364 int mchROF = 0;
365 // loop over MCH ROFs that are not newer than the current MFT ROF
366 while (mchROF < nMCHROFs && !(thisMFTBracket < mMCHROFTimes[mchROF])) {
367 // only consider non-empty MCH ROFs that overlap with the MFT one
368 if (mMCHTrackROFRec[mchROF].getNEntries() > 0 && !(thisMFTBracket.isOutside(mMCHROFTimes[mchROF]))) {
369 // set the index of the first MCH ROF if not yet initialized
370 if (mchROFMatchFirst < 0) {
371 mchROFMatchFirst = mchROF;
372 }
373 // update the index of the last MCH ROF
374 mchROFMatchLast = mchROF;
375 }
376 mchROF++;
377 }
378 // skip if the index of the first MCH ROF is not set
379 if (mchROFMatchFirst < 0) {
380 continue;
381 }
382 LOG(debug) << "FIRST MCH ROF " << mchROFMatchFirst << "; interval: ["
383 << mMCHROFTimes[mchROFMatchFirst].getMin() << ","
384 << mMCHROFTimes[mchROFMatchFirst].getMax() << "] size: " << mMCHTrackROFRec[mchROFMatchFirst].getNEntries();
385 LOG(debug) << "LAST MCH ROF " << mchROFMatchLast << "; interval: ["
386 << mMCHROFTimes[mchROFMatchLast].getMin() << ","
387 << mMCHROFTimes[mchROFMatchLast].getMax() << "] size: " << mMCHTrackROFRec[mchROFMatchLast].getNEntries();
388
389 ROFMatch<saveAllMode>(MFTROFId, mchROFMatchFirst, mchROFMatchLast);
390 }
391
392 if constexpr (saveAllMode == SaveMode::kBestMatch) { // Otherwise output container is filled by ROFMatch()
393 int nFakes = 0, nTrue = 0;
394 for (auto& thisMCHTrack : mMCHWork) {
395 auto bestMFTMatchID = thisMCHTrack.getMFTTrackID();
396 if (bestMFTMatchID >= 0) { // If there is a match, add to output container
397 if (mMCTruthON) {
398 mMatchLabels.push_back(computeLabel(thisMCHTrack.getMCHTrackID(), bestMFTMatchID));
399 mMatchLabels.back().isFake() ? nFakes++ : nTrue++;
400 }
401
402 thisMCHTrack.setMFTTrackID(bestMFTMatchID);
403 LOG(debug) << " thisMCHTrack.getMFTTrackID() = " << thisMCHTrack.getMFTTrackID()
404 << "; thisMCHTrack.getMFTMCHMatchingChi2() = " << thisMCHTrack.getMFTMCHMatchingChi2();
405
406 mMatchedTracks.emplace_back(thisMCHTrack);
407 mMatchingInfo.emplace_back(thisMCHTrack);
408 }
409 }
410 if (mMCTruthON) {
411 LOG(info) << " MFT-MCH Matching: nFakes = " << nFakes << " nTrue = " << nTrue;
412 }
413 } else if constexpr (saveAllMode == SaveMode::kSaveNCandidates) {
414 int nFakes = 0, nTrue = 0;
415 auto& matchAllChi2 = mMatchingFunctionMap["matchALL"];
416 for (auto MCHId = 0; MCHId < mMCHWork.size(); MCHId++) {
417 auto& thisMCHTrack = mMCHWork[MCHId];
418 for (auto& pairCandidate : mCandidates[MCHId]) {
419 thisMCHTrack.setMFTTrackID(pairCandidate.second);
420 auto& thisMFTTrack = mMFTWork[pairCandidate.second];
421 auto chi2 = matchAllChi2(thisMCHTrack, thisMFTTrack); // Matching chi2 is stored independently
422 thisMCHTrack.setMFTMCHMatchingScore(pairCandidate.first);
423 thisMCHTrack.setMFTMCHMatchingChi2(chi2);
424 mMatchedTracks.emplace_back(thisMCHTrack);
425 mMatchingInfo.emplace_back(thisMCHTrack);
426 if (mMCTruthON) {
427 mMatchLabels.push_back(computeLabel(MCHId, pairCandidate.second));
428 mMatchLabels.back().isFake() ? nFakes++ : nTrue++;
429 }
430 }
431 }
432 }
433}
434
435//_________________________________________________________
436template <Int_t saveAllMode>
437void MatchGlobalFwd::ROFMatch(int MFTROFId, int firstMCHROFId, int lastMCHROFId)
438{
440 const auto& thisMFTROF = mMFTTrackROFRec[MFTROFId];
441 const auto& thisMFTBracket = mMFTROFTimes[MFTROFId];
442 const auto& firstMCHROF = mMCHTrackROFRec[firstMCHROFId];
443 const auto& lastMCHROF = mMCHTrackROFRec[lastMCHROFId];
444 int nFakes = 0, nTrue = 0;
445
446 auto compare = [](const std::pair<int, int>& a, const std::pair<int, int>& b) {
447 return a.first < b.first;
448 };
449
450 auto firstMFTTrackID = thisMFTROF.getFirstEntry();
451 auto lastMFTTrackID = firstMFTTrackID + thisMFTROF.getNEntries() - 1;
452
453 auto firstMCHTrackID = firstMCHROF.getFirstIdx();
454 auto lastMCHTrackID = lastMCHROF.getLastIdx();
455
456 auto nMFTTracks = thisMFTROF.getNEntries();
457 auto nMCHTracks = lastMCHTrackID - firstMCHTrackID + 1;
458
459 auto& matchAllChi2 = mMatchingFunctionMap["matchALL"];
460
461 LOG(debug) << "Matching MFT ROF " << MFTROFId << " with MCH ROFs [" << firstMCHROFId << "->" << lastMCHROFId << "]";
462 LOG(debug) << " firstMFTTrackID = " << firstMFTTrackID << " ; lastMFTTrackID = " << lastMFTTrackID;
463 LOG(debug) << " firstMCHTrackID = " << firstMCHTrackID << " ; lastMCHTrackID = " << lastMCHTrackID;
464 LOG(debug) << " thisMFTROF: " << thisMFTROF.getBCData();
465 LOG(debug) << " firstMCHROF: " << firstMCHROF;
466 LOG(debug) << " lastMCHROF: " << lastMCHROF;
467
468 // loop over all MCH tracks
469 for (auto MCHId = firstMCHTrackID; MCHId <= lastMCHTrackID; MCHId++) {
470 auto& thisMCHTrack = mMCHWork[MCHId];
471
472 // If enabled, use the muon track time to check if the track is correlated with the MFT ROF
473 if (mUseTrackTime && (thisMFTBracket.isOutside(thisMCHTrack.tBracket))) {
474 continue;
475 }
476
477 o2::MCCompLabel matchLabel;
478 for (auto MFTId = firstMFTTrackID; MFTId <= lastMFTTrackID; MFTId++) {
479 auto& thisMFTTrack = mMFTWork[MFTId];
480 if (mMCTruthON) {
481 matchLabel = computeLabel(MCHId, MFTId);
482 }
483 if (mCutFunc(thisMCHTrack, thisMFTTrack)) {
484 thisMCHTrack.countMFTCandidate();
485 if (mMCTruthON) {
486 if (matchLabel.isCorrect()) {
487 thisMCHTrack.setCloseMatch();
488 }
489 }
490 auto score = mMatchFunc(thisMCHTrack, thisMFTTrack);
491 if (score < thisMCHTrack.getMFTMCHMatchingScore()) {
492 thisMCHTrack.setMFTTrackID(MFTId);
493 auto chi2 = matchAllChi2(thisMCHTrack, thisMFTTrack); // Matching chi2 is stored independently
494 thisMCHTrack.setMFTMCHMatchingScore(score);
495 thisMCHTrack.setMFTMCHMatchingChi2(chi2);
496 }
497 if constexpr (saveAllMode == SaveMode::kSaveAll) { // In saveAllmode save all pairs to output container
498 thisMCHTrack.setMFTTrackID(MFTId);
499 mMatchedTracks.emplace_back(thisMCHTrack);
500 mMatchingInfo.emplace_back(thisMCHTrack);
501 if (mMCTruthON) {
502 mMatchLabels.push_back(matchLabel);
503 mMatchLabels.back().isFake() ? nFakes++ : nTrue++;
504 }
505 }
506
507 if constexpr (saveAllMode == SaveMode::kSaveNCandidates) { // Save best N matching candidates
508 auto score = mMatchFunc(thisMCHTrack, thisMFTTrack);
509 std::pair<float, int> scoreID = {score, MFTId};
510 mCandidates[MCHId].push_back(scoreID);
511 std::sort(mCandidates[MCHId].begin(), mCandidates[MCHId].end(), compare);
512 if (mCandidates[MCHId].size() > mNCandidates) {
513 mCandidates[MCHId].pop_back();
514 }
515 }
516
517 if constexpr (saveAllMode == SaveMode::kSaveTrainingData) { // In save training data mode store track parameters at matching plane
518 thisMCHTrack.setMFTTrackID(MFTId);
519 mMatchingInfo.emplace_back(thisMCHTrack);
520 mMCHMatchPlaneParams.emplace_back(thisMCHTrack);
521 mMFTMatchPlaneParams.emplace_back(static_cast<o2::mft::TrackMFT>(thisMFTTrack));
522 if (mMCTruthON) {
523 mMatchLabels.push_back(matchLabel);
524 mMatchLabels.back().isFake() ? nFakes++ : nTrue++;
525 }
526 }
527 }
528 }
529 auto bestMFTMatchID = thisMCHTrack.getMFTTrackID();
530 LOG(debug) << " Matching MCHId = " << MCHId << " ==> bestMFTMatchID = " << thisMCHTrack.getMFTTrackID() << " ; thisMCHTrack.getMFTMCHMatchingChi2() = " << thisMCHTrack.getMFTMCHMatchingChi2();
531 LOG(debug) << " MCH COV<X,X> = " << thisMCHTrack.getSigma2X() << " ; COV<Y,Y> = " << thisMCHTrack.getSigma2Y() << " ; pt = " << thisMCHTrack.getPt();
532
533 } // /loop over MCH tracks seeds
534
535 LOG(debug) << "Finished matching MFT ROF " << MFTROFId << ": " << nMFTTracks << " MFT tracks and " << nMCHTracks << " MCH Tracks.";
536 if (mMCTruthON) {
537 LOG(debug) << " nFakes = " << nFakes << " nTrue = " << nTrue;
538 }
539}
540
541//_________________________________________________________
542o2::MCCompLabel MatchGlobalFwd::computeLabel(const int MCHId, const int MFTId)
543{
544 const auto& mchlabel = mMCHTrkLabels[MCHId];
545 const auto& mftlabel = mMFTTrkLabels[MFTId];
546 o2::MCCompLabel matchLabel = mchlabel;
547 matchLabel.setFakeFlag(mftlabel.compare(mchlabel) != 1);
548
549 LOG(debug) << " Computing MFTMCH matching label: MFTTruth = " << mftlabel << " ; MCHTruth = " << mchlabel << " ; Computed label = " << matchLabel;
550
551 return matchLabel;
552}
553
554//_________________________________________________________
555void MatchGlobalFwd::doMCMatching()
556{
557 int nFakes = 0, nTrue = 0;
558
559 // loop over all MCH tracks
560 for (auto MCHId = 0; MCHId < mMCHWork.size(); MCHId++) {
561 auto& thisMCHTrack = mMCHWork[MCHId];
562 const o2::MCCompLabel& thisMCHLabel = mMCHTrkLabels[mMCHID2Work[MCHId]];
563
564 LOG(debug) << " MCH Track # " << MCHId << " Label: " << thisMCHLabel;
565 if (!((thisMCHLabel).isSet())) {
566 continue;
567 }
568 for (auto MFTId = 0; MFTId < mMFTWork.size(); MFTId++) {
569 auto& thisMFTTrack = mMFTWork[MFTId];
570 o2::MCCompLabel matchLabel = computeLabel(MCHId, MFTId);
571
572 if (matchLabel.isCorrect()) {
573 nTrue++;
574 thisMCHTrack.setCloseMatch();
575 auto chi2 = mMatchFunc(thisMCHTrack, thisMFTTrack);
576 thisMCHTrack.setMFTTrackID(MFTId);
577 thisMCHTrack.setMFTMCHMatchingChi2(chi2);
578 mMatchedTracks.emplace_back(thisMCHTrack);
579 mMatchingInfo.emplace_back(thisMCHTrack);
580 mMatchLabels.push_back(matchLabel);
581 auto bestMFTMatchID = thisMCHTrack.getMFTTrackID();
582 LOG(debug) << " Matching MCHId = " << MCHId << " ==> bestMFTMatchID = " << thisMCHTrack.getMFTTrackID() << " ; thisMCHTrack.getMFTMCHMatchingChi2() = " << thisMCHTrack.getMFTMCHMatchingChi2();
583 LOG(debug) << " MCH COV<X,X> = " << thisMCHTrack.getSigma2X() << " ; COV<Y,Y> = " << thisMCHTrack.getSigma2Y() << " ; pt = " << thisMCHTrack.getPt();
584 LOG(debug) << " Label: " << matchLabel;
585 break;
586 }
587 }
588
589 } // /loop over MCH tracks seeds
590
591 auto nMFTTracks = mMFTWork.size();
592 auto nMCHTracks = mMCHWork.size();
593
594 LOG(info) << " Done MC matching of " << nMFTTracks << " MFT tracks with " << nMCHTracks << " MCH Tracks. nFakes = " << nFakes << " nTrue = " << nTrue;
595}
596
597//_________________________________________________________________________________________________
598void MatchGlobalFwd::fitTracks()
599{
600 LOG(info) << "Fitting global muon tracks...";
601
602 auto GTrackID = 0;
603
604 for (auto& track : mMatchedTracks) {
605 LOG(debug) << " ==> Fitting Global Track # " << GTrackID << " with MFT track # " << track.getMFTTrackID() << ":";
606 fitGlobalMuonTrack(track);
607 GTrackID++;
608 }
609
610 LOG(info) << "Finished fitting global muon tracks.";
611}
612
613//_________________________________________________________________________________________________
614void MatchGlobalFwd::fitGlobalMuonTrack(o2::dataformats::GlobalFwdTrack& gTrack)
615{
616 const auto& MFTMatchId = gTrack.getMFTTrackID();
617 const auto& mftTrack = mMFTTracks[MFTMatchId];
618 const auto& mftTrackOut = mMFTWork[MFTMatchId];
619 auto ncls = mftTrack.getNumberOfPoints();
620 auto offset = mftTrack.getExternalClusterIndexOffset();
621
622 LOG(debug) << "***************************** Start Fitting new track *****************************";
623 LOG(debug) << "N Clusters = " << ncls << " Best MFT Track Match ID " << gTrack.getMFTTrackID() << " MCHTrack: X = " << gTrack.getX() << " Y = " << gTrack.getY() << " Z = " << gTrack.getZ() << " Tgl = " << gTrack.getTanl() << " Phi = " << gTrack.getPhi() << " pz = " << gTrack.getPz() << " qpt = " << 1.0 / gTrack.getInvQPt();
624
625 LOG(debug) << "MFTTrack: X = " << mftTrackOut.getX()
626 << " Y = " << mftTrackOut.getY() << " Z = " << mftTrackOut.getZ()
627 << " Tgl = " << mftTrackOut.getTanl()
628 << " Phi = " << mftTrackOut.getPhi() << " pz = " << mftTrackOut.getPz()
629 << " qpt = " << 1.0 / mftTrackOut.getInvQPt();
630 LOG(debug) << " initTrack GlobalTrack: q/pt = " << gTrack.getInvQPt() << std::endl;
631
632 auto lastLayer = mMFTMapping.ChipID2Layer[mMFTClusters[offset + ncls - 1].getSensorID()];
633 LOG(debug) << "Starting by MFTCluster offset " << offset + ncls - 1 << " at lastLayer " << lastLayer;
634
635 for (int icls = ncls - 1; icls > -1; --icls) {
636 auto clsEntry = mMFTTrackClusIdx[offset + icls];
637 auto& thiscluster = mMFTClusters[clsEntry];
638 LOG(debug) << "Computing MFTCluster clsEntry " << clsEntry << " at Z = " << thiscluster.getZ();
639
640 computeCluster(gTrack, thiscluster, lastLayer);
641 }
642}
643
644//_________________________________________________________________________________________________
645bool MatchGlobalFwd::computeCluster(o2::dataformats::GlobalFwdTrack& track, const MFTCluster& cluster, int& startingLayerID)
646{
651
652 const auto& clx = cluster.getX();
653 const auto& cly = cluster.getY();
654 const auto& clz = cluster.getZ();
655 const auto& sigmaX2 = cluster.getSigmaY2() * mAlignResidual * mAlignResidual;
656 ; // ALPIDE local Y coordinate => MFT global X coordinate (ALPIDE rows)
657 const auto& sigmaY2 = cluster.getSigmaZ2() * mAlignResidual * mAlignResidual;
658 ; // ALPIDE local Z coordinate => MFT global Y coordinate (ALPIDE columns)
659
660 const auto& newLayerID = mMFTMapping.ChipID2Layer[cluster.getSensorID()];
661 LOG(debug) << "computeCluster: X = " << clx << " Y = " << cly << " Z = " << clz << " nCluster = " << newLayerID;
662
663 if (!propagateToNextClusterWithMCS(track, clz, startingLayerID, newLayerID)) {
664 return false;
665 }
666
667 LOG(debug) << " AfterExtrap: X = " << track.getX() << " Y = " << track.getY() << " Z = " << track.getZ() << " Tgl = " << track.getTanl() << " Phi = " << track.getPhi() << " pz = " << track.getPz() << " q/pt = " << track.getInvQPt();
668 LOG(debug) << "Track covariances after extrap:" << std::endl
669 << track.getCovariances() << std::endl;
670
671 // recompute parameters
672 const std::array<float, 2>& pos = {clx, cly};
673 const std::array<float, 2>& cov = {sigmaX2, sigmaY2};
674
675 if (track.update(pos, cov)) {
676 LOG(debug) << " New Cluster: X = " << clx << " Y = " << cly << " Z = " << clz;
677 LOG(debug) << " AfterKalman: X = " << track.getX() << " Y = " << track.getY() << " Z = " << track.getZ() << " Tgl = " << track.getTanl() << " Phi = " << track.getPhi() << " pz = " << track.getPz() << " q/pt = " << track.getInvQPt();
678
679 LOG(debug) << "Track covariances after Kalman update: \n"
680 << track.getCovariances() << std::endl;
681
682 return true;
683 }
684 return false;
685}
686
687//_________________________________________________________
689{
690 mMFTROFrameLengthMUS = fums;
691 mMFTROFrameLengthMUSInv = 1. / mMFTROFrameLengthMUS;
692 mMFTROFrameLengthInBC = std::max(1, int(mMFTROFrameLengthMUS / (o2::constants::lhc::LHCBunchSpacingNS * 1e-3)));
693}
694
695//_________________________________________________________
697{
698 mMFTROFrameLengthInBC = nbc;
699 mMFTROFrameLengthMUS = nbc * o2::constants::lhc::LHCBunchSpacingNS * 1e-3;
700 mMFTROFrameLengthMUSInv = 1. / mMFTROFrameLengthMUS;
701}
702
703//_________________________________________________________
705{
706 mMFTROFrameBiasInBC = nbc;
707 mMFTROFrameBiasMUS = nbc * o2::constants::lhc::LHCBunchSpacingNS * 1e-3;
708 mMFTROFrameBiasMUSInv = 1. / mMFTROFrameBiasMUS;
709}
710
711//_________________________________________________________
713{
714 mBunchFilling = bf;
715 // find closest (from above) filled bunch
716 int minBC = bf.getFirstFilledBC(), maxBC = bf.getLastFilledBC();
717 if (minBC < 0) {
718 LOG(error) << "Empty bunch filling is provided to MatchGlobalFwd, checks using it should be ignored";
719 return;
720 }
721 int bcAbove = minBC;
722 for (int i = o2::constants::lhc::LHCMaxBunches; i--;) {
723 if (bf.testBC(i)) {
724 bcAbove = i;
725 }
726 mClosestBunchAbove[i] = bcAbove;
727 }
728 int bcBelow = maxBC;
729 for (int i = 0; i < o2::constants::lhc::LHCMaxBunches; i++) {
730 if (bf.testBC(i)) {
731 bcBelow = i;
732 }
733 mClosestBunchBelow[i] = bcBelow;
734 }
735}
736
737//_________________________________________________________________________________________________
739{
740 // Convert a MCH Track parameters and covariances matrix to the
741 // Forward track format. Must be called after propagation though the absorber
742
743 o2::dataformats::GlobalFwdTrack convertedTrack;
744
745 // Parameter conversion
746 double alpha1, alpha3, alpha4, x2, x3, x4;
747
748 alpha1 = mchParam.getNonBendingSlope();
749 alpha3 = mchParam.getBendingSlope();
750 alpha4 = mchParam.getInverseBendingMomentum();
751
752 x2 = TMath::ATan2(-alpha3, -alpha1);
753 x3 = -1. / TMath::Sqrt(alpha3 * alpha3 + alpha1 * alpha1);
754 x4 = alpha4 * -x3 * TMath::Sqrt(1 + alpha3 * alpha3);
755
756 auto K = alpha1 * alpha1 + alpha3 * alpha3;
757 auto K32 = K * TMath::Sqrt(K);
758 auto L = TMath::Sqrt(alpha3 * alpha3 + 1);
759
760 // Covariances matrix conversion
761 SMatrix55Std jacobian;
762 SMatrix55Sym covariances;
763
764 if (0) {
765
766 std::cout << " MCHtoGlobal - MCH Covariances:\n";
767 std::cout << " mchParam.getCovariances()(0, 0) = "
768 << mchParam.getCovariances()(0, 0)
769 << " ; mchParam.getCovariances()(2, 2) = "
770 << mchParam.getCovariances()(2, 2) << std::endl;
771 }
772 covariances(0, 0) = mchParam.getCovariances()(0, 0);
773 covariances(0, 1) = mchParam.getCovariances()(0, 1);
774 covariances(0, 2) = mchParam.getCovariances()(0, 2);
775 covariances(0, 3) = mchParam.getCovariances()(0, 3);
776 covariances(0, 4) = mchParam.getCovariances()(0, 4);
777
778 covariances(1, 1) = mchParam.getCovariances()(1, 1);
779 covariances(1, 2) = mchParam.getCovariances()(1, 2);
780 covariances(1, 3) = mchParam.getCovariances()(1, 3);
781 covariances(1, 4) = mchParam.getCovariances()(1, 4);
782
783 covariances(2, 2) = mchParam.getCovariances()(2, 2);
784 covariances(2, 3) = mchParam.getCovariances()(2, 3);
785 covariances(2, 4) = mchParam.getCovariances()(2, 4);
786
787 covariances(3, 3) = mchParam.getCovariances()(3, 3);
788 covariances(3, 4) = mchParam.getCovariances()(3, 4);
789
790 covariances(4, 4) = mchParam.getCovariances()(4, 4);
791
792 jacobian(0, 0) = 1;
793
794 jacobian(1, 2) = 1;
795
796 jacobian(2, 1) = -alpha3 / K;
797 jacobian(2, 3) = alpha1 / K;
798
799 jacobian(3, 1) = alpha1 / K32;
800 jacobian(3, 3) = alpha3 / K32;
801
802 jacobian(4, 1) = -alpha1 * alpha4 * L / K32;
803 jacobian(4, 3) = alpha3 * alpha4 * (1 / (TMath::Sqrt(K) * L) - L / K32);
804 jacobian(4, 4) = L / TMath::Sqrt(K);
805
806 // jacobian*covariances*jacobian^T
807 covariances = ROOT::Math::Similarity(jacobian, covariances);
808
809 // Set output
810 convertedTrack.setX(mchParam.getNonBendingCoor());
811 convertedTrack.setY(mchParam.getBendingCoor());
812 convertedTrack.setZ(mchParam.getZ());
813 convertedTrack.setPhi(x2);
814 convertedTrack.setTanl(x3);
815 convertedTrack.setInvQPt(x4);
816 convertedTrack.setCharge(mchParam.getCharge());
817 convertedTrack.setCovariances(covariances);
818
819 return convertedTrack;
820}
821
822//_________________________________________________________________________________________________
824{
825 // Convert Forward Track parameters and covariances matrix to the
826 // MCH track format.
827
828 // Parameter conversion
829 double alpha1, alpha3, alpha4, x2, x3, x4;
830
831 x2 = fwdtrack.getPhi();
832 x3 = fwdtrack.getTanl();
833 x4 = fwdtrack.getInvQPt();
834
835 auto sinx2 = TMath::Sin(x2);
836 auto cosx2 = TMath::Cos(x2);
837
838 alpha1 = cosx2 / x3;
839 alpha3 = sinx2 / x3;
840 alpha4 = x4 / TMath::Sqrt(x3 * x3 + sinx2 * sinx2);
841
842 auto K = TMath::Sqrt(x3 * x3 + sinx2 * sinx2);
843 auto K3 = K * K * K;
844
845 // Covariances matrix conversion
846 SMatrix55Std jacobian;
847 SMatrix55Sym covariances;
848
849 covariances(0, 0) = fwdtrack.getCovariances()(0, 0);
850 covariances(0, 1) = fwdtrack.getCovariances()(0, 1);
851 covariances(0, 2) = fwdtrack.getCovariances()(0, 2);
852 covariances(0, 3) = fwdtrack.getCovariances()(0, 3);
853 covariances(0, 4) = fwdtrack.getCovariances()(0, 4);
854
855 covariances(1, 1) = fwdtrack.getCovariances()(1, 1);
856 covariances(1, 2) = fwdtrack.getCovariances()(1, 2);
857 covariances(1, 3) = fwdtrack.getCovariances()(1, 3);
858 covariances(1, 4) = fwdtrack.getCovariances()(1, 4);
859
860 covariances(2, 2) = fwdtrack.getCovariances()(2, 2);
861 covariances(2, 3) = fwdtrack.getCovariances()(2, 3);
862 covariances(2, 4) = fwdtrack.getCovariances()(2, 4);
863
864 covariances(3, 3) = fwdtrack.getCovariances()(3, 3);
865 covariances(3, 4) = fwdtrack.getCovariances()(3, 4);
866
867 covariances(4, 4) = fwdtrack.getCovariances()(4, 4);
868
869 jacobian(0, 0) = 1;
870
871 jacobian(1, 2) = -sinx2 / x3;
872 jacobian(1, 3) = -cosx2 / (x3 * x3);
873
874 jacobian(2, 1) = 1;
875
876 jacobian(3, 2) = cosx2 / x3;
877 jacobian(3, 3) = -sinx2 / (x3 * x3);
878
879 jacobian(4, 2) = -x4 * sinx2 * cosx2 / K3;
880 jacobian(4, 3) = -x3 * x4 / K3;
881 jacobian(4, 4) = 1 / K;
882 // jacobian*covariances*jacobian^T
883 covariances = ROOT::Math::Similarity(jacobian, covariances);
884
885 double cov[] = {covariances(0, 0), covariances(1, 0), covariances(1, 1), covariances(2, 0), covariances(2, 1), covariances(2, 2), covariances(3, 0), covariances(3, 1), covariances(3, 2), covariances(3, 3), covariances(4, 0), covariances(4, 1), covariances(4, 2), covariances(4, 3), covariances(4, 4)};
886 double param[] = {fwdtrack.getX(), alpha1, fwdtrack.getY(), alpha3, alpha4};
887
888 o2::mch::TrackParam convertedTrack(fwdtrack.getZ(), param, cov);
889 return o2::mch::TrackParam(convertedTrack);
890}
891
892//_________________________________________________________________________________________________
894{
895 mClosestBunchAbove[0] = mClosestBunchAbove[0] = -1;
896
897 // Define built-in matching functions
898 //________________________________________________________________________________
899 mMatchingFunctionMap["matchALL"] = [](const GlobalFwdTrack& mchTrack, const TrackParCovFwd& mftTrack) -> double {
900 // Match two tracks evaluating all parameters: X,Y, phi, tanl & q/pt
901
902 SMatrix55Sym I = ROOT::Math::SMatrixIdentity(), H_k, V_k;
903 SVector5 m_k(mftTrack.getX(), mftTrack.getY(), mftTrack.getPhi(),
904 mftTrack.getTanl(), mftTrack.getInvQPt()),
905 r_k_kminus1;
906 SVector5 GlobalMuonTrackParameters = mchTrack.getParameters();
907 SMatrix55Sym GlobalMuonTrackCovariances = mchTrack.getCovariances();
908 V_k(0, 0) = mftTrack.getCovariances()(0, 0);
909 V_k(1, 1) = mftTrack.getCovariances()(1, 1);
910 V_k(2, 2) = mftTrack.getCovariances()(2, 2);
911 V_k(3, 3) = mftTrack.getCovariances()(3, 3);
912 V_k(4, 4) = mftTrack.getCovariances()(4, 4);
913 H_k(0, 0) = 1.0;
914 H_k(1, 1) = 1.0;
915 H_k(2, 2) = 1.0;
916 H_k(3, 3) = 1.0;
917 H_k(4, 4) = 1.0;
918
919 // Covariance of residuals
920 SMatrix55Std invResCov = (V_k + ROOT::Math::Similarity(H_k, GlobalMuonTrackCovariances));
921 invResCov.Invert();
922
923 // Kalman Gain Matrix
924 SMatrix55Std K_k = GlobalMuonTrackCovariances * ROOT::Math::Transpose(H_k) * invResCov;
925
926 // Update Parameters
927 r_k_kminus1 = m_k - H_k * GlobalMuonTrackParameters; // Residuals of prediction
928
929 auto matchChi2Track = ROOT::Math::Similarity(r_k_kminus1, invResCov);
930
931 return matchChi2Track;
932 };
933
934 //________________________________________________________________________________
935 mMatchingFunctionMap["matchsXYPhiTanl"] = [](const GlobalFwdTrack& mchTrack, const TrackParCovFwd& mftTrack) -> double {
936
937 // Match two tracks evaluating positions & angles
938
939 SMatrix55Sym I = ROOT::Math::SMatrixIdentity();
940 SMatrix45 H_k;
941 SMatrix44 V_k;
942 SVector4 m_k(mftTrack.getX(), mftTrack.getY(), mftTrack.getPhi(),
943 mftTrack.getTanl()),
944 r_k_kminus1;
945 SVector5 GlobalMuonTrackParameters = mchTrack.getParameters();
946 SMatrix55Sym GlobalMuonTrackCovariances = mchTrack.getCovariances();
947 V_k(0, 0) = mftTrack.getCovariances()(0, 0);
948 V_k(1, 1) = mftTrack.getCovariances()(1, 1);
949 V_k(2, 2) = mftTrack.getCovariances()(2, 2);
950 V_k(3, 3) = mftTrack.getCovariances()(3, 3);
951 H_k(0, 0) = 1.0;
952 H_k(1, 1) = 1.0;
953 H_k(2, 2) = 1.0;
954 H_k(3, 3) = 1.0;
955
956 // Covariance of residuals
957 SMatrix44 invResCov = (V_k + ROOT::Math::Similarity(H_k, GlobalMuonTrackCovariances));
958 invResCov.Invert();
959
960 // Kalman Gain Matrix
961 SMatrix54 K_k = GlobalMuonTrackCovariances * ROOT::Math::Transpose(H_k) * invResCov;
962
963 // Residuals of prediction
964 r_k_kminus1 = m_k - H_k * GlobalMuonTrackParameters;
965
966 auto matchChi2Track = ROOT::Math::Similarity(r_k_kminus1, invResCov);
967
968 return matchChi2Track; };
969
970 //________________________________________________________________________________
971 mMatchingFunctionMap["matchXY"] = [](const GlobalFwdTrack& mchTrack, const TrackParCovFwd& mftTrack) -> double {
972
973 // Calculate Matching Chi2 - X and Y positions
974
975 SMatrix55Sym I = ROOT::Math::SMatrixIdentity();
976 SMatrix25 H_k;
977 SMatrix22 V_k;
978 SVector2 m_k(mftTrack.getX(), mftTrack.getY()), r_k_kminus1;
979 SVector5 GlobalMuonTrackParameters = mchTrack.getParameters();
980 SMatrix55Sym GlobalMuonTrackCovariances = mchTrack.getCovariances();
981 V_k(0, 0) = mftTrack.getCovariances()(0, 0);
982 V_k(1, 1) = mftTrack.getCovariances()(1, 1);
983 H_k(0, 0) = 1.0;
984 H_k(1, 1) = 1.0;
985
986 // Covariance of residuals
987 SMatrix22 invResCov = (V_k + ROOT::Math::Similarity(H_k, GlobalMuonTrackCovariances));
988 invResCov.Invert();
989
990 // Kalman Gain Matrix
991 SMatrix52 K_k = GlobalMuonTrackCovariances * ROOT::Math::Transpose(H_k) * invResCov;
992
993 // Residuals of prediction
994 r_k_kminus1 = m_k - H_k * GlobalMuonTrackParameters;
995 auto matchChi2Track = ROOT::Math::Similarity(r_k_kminus1, invResCov);
996
997 return matchChi2Track; };
998
999 //________________________________________________________________________________
1000 mMatchingFunctionMap["matchNeedsName"] = [this](const GlobalFwdTrack& mchTrack, const TrackParCovFwd& mftTrack) -> double {
1001
1002 //Hiroshima's Matching function needs a physics-based name
1003
1004 //Matching constants
1005 Double_t LAbs = 415.; //Absorber Length[cm]
1006 Double_t mumass = 0.106; //mass of muon [GeV/c^2]
1007 Double_t l; //the length that extrapolated MCHtrack passes through absorber
1008
1009 if (mMatchingPlaneZ >= -90.0) {
1010 l = LAbs;
1011 } else {
1012 l = 505.0 + mMatchingPlaneZ;
1013 }
1014
1015 //defference between MFTtrack and MCHtrack
1016
1017 auto dx = mftTrack.getX() - mchTrack.getX();
1018 auto dy = mftTrack.getY() - mchTrack.getY();
1019 auto dthetax = TMath::ATan(mftTrack.getPx() / TMath::Abs(mftTrack.getPz())) - TMath::ATan(mchTrack.getPx() / TMath::Abs(mchTrack.getPz()));
1020 auto dthetay = TMath::ATan(mftTrack.getPy() / TMath::Abs(mftTrack.getPz())) - TMath::ATan(mchTrack.getPy() / TMath::Abs(mchTrack.getPz()));
1021
1022 //Multiple Scattering(=MS)
1023
1024 auto pMCH = mchTrack.getP();
1025 auto lorentzbeta = pMCH / TMath::Sqrt(mumass * mumass + pMCH * pMCH);
1026 auto zMS = copysign(1.0, mchTrack.getCharge());
1027 auto thetaMS = 13.6 / (1000.0 * pMCH * lorentzbeta * 1.0) * zMS * TMath::Sqrt(60.0 * l / LAbs) * (1.0 + 0.038 * TMath::Log(60.0 * l / LAbs));
1028 auto xMS = thetaMS * l / TMath::Sqrt(3.0);
1029
1030 //normalize by theoritical Multiple Coulomb Scattering width to be momentum-independent
1031 //make the dx and dtheta dimensionless
1032
1033 auto dxnorm = dx / xMS;
1034 auto dynorm = dy / xMS;
1035 auto dthetaxnorm = dthetax / thetaMS;
1036 auto dthetaynorm = dthetay / thetaMS;
1037
1038 //rotate distribution
1039
1040 auto dxrot = dxnorm * TMath::Cos(TMath::Pi() / 4.0) - dthetaxnorm * TMath::Sin(TMath::Pi() / 4.0);
1041 auto dthetaxrot = dxnorm * TMath::Sin(TMath::Pi() / 4.0) + dthetaxnorm * TMath::Cos(TMath::Pi() / 4.0);
1042 auto dyrot = dynorm * TMath::Cos(TMath::Pi() / 4.0) - dthetaynorm * TMath::Sin(TMath::Pi() / 4.0);
1043 auto dthetayrot = dynorm * TMath::Sin(TMath::Pi() / 4.0) + dthetaynorm * TMath::Cos(TMath::Pi() / 4.0);
1044
1045 //convert ellipse to circle
1046
1047 auto k = 0.7; //need to optimize!!
1048 auto dxcircle = dxrot;
1049 auto dycircle = dyrot;
1050 auto dthetaxcircle = dthetaxrot / k;
1051 auto dthetaycircle = dthetayrot / k;
1052
1053 //score
1054
1055 auto scoreX = TMath::Sqrt(dxcircle * dxcircle + dthetaxcircle * dthetaxcircle);
1056 auto scoreY = TMath::Sqrt(dycircle * dycircle + dthetaycircle * dthetaycircle);
1057 auto score = TMath::Sqrt(scoreX * scoreX + scoreY * scoreY);
1058
1059 return score; };
1060
1061 // Define built-in candidate cut functions
1062
1063 //________________________________________________________________________________
1064 mCutFunctionMap["cutDisabled"] = [](const GlobalFwdTrack& mchTrack, const TrackParCovFwd& mftTrack) -> bool {
1065 return true;
1066 };
1067
1068 //________________________________________________________________________________
1069 mCutFunctionMap["cut3Sigma"] = [](const GlobalFwdTrack& mchTrack, const TrackParCovFwd& mftTrack) -> bool {
1070 auto dx = mchTrack.getX() - mftTrack.getX();
1071 auto dy = mchTrack.getY() - mftTrack.getY();
1072 auto dPhi = mchTrack.getPhi() - mftTrack.getPhi();
1073 auto dTanl = TMath::Abs(mchTrack.getTanl() - mftTrack.getTanl());
1074 auto dInvQPt = TMath::Abs(mchTrack.getInvQPt() - mftTrack.getInvQPt());
1075 auto distanceSq = dx * dx + dy * dy;
1076 auto cutDistanceSq = 9 * (mchTrack.getSigma2X() + mchTrack.getSigma2Y());
1077 auto cutPhiSq = 9 * (mchTrack.getSigma2Phi() + mftTrack.getSigma2Phi());
1078 auto cutTanlSq = 9 * (mchTrack.getSigma2Tanl() + mftTrack.getSigma2Tanl());
1079 auto cutInvQPtSq = 9 * (mchTrack.getSigma2InvQPt() + mftTrack.getSigma2InvQPt());
1080 return (distanceSq < cutDistanceSq) and (dPhi * dPhi < cutPhiSq) and (dTanl * dTanl < cutTanlSq) and (dInvQPt * dInvQPt < cutInvQPtSq);
1081 };
1082
1083 //________________________________________________________________________________
1084 mCutFunctionMap["cut3SigmaXYAngles"] = [](const GlobalFwdTrack& mchTrack, const TrackParCovFwd& mftTrack) -> bool {
1085 auto dx = mchTrack.getX() - mftTrack.getX();
1086 auto dy = mchTrack.getY() - mftTrack.getY();
1087 auto dPhi = mchTrack.getPhi() - mftTrack.getPhi();
1088 auto dTanl = TMath::Abs(mchTrack.getTanl() - mftTrack.getTanl());
1089 auto distanceSq = dx * dx + dy * dy;
1090 auto cutDistanceSq = 9 * (mchTrack.getSigma2X() + mchTrack.getSigma2Y());
1091 auto cutPhiSq = 9 * (mchTrack.getSigma2Phi() + mftTrack.getSigma2Phi());
1092 auto cutTanlSq = 9 * (mchTrack.getSigma2Tanl() + mftTrack.getSigma2Tanl());
1093 return (distanceSq < cutDistanceSq) and (dPhi * dPhi < cutPhiSq) and (dTanl * dTanl < cutTanlSq);
1094 };
1095}
int32_t i
Class to perform MFT MCH (and MID) matching.
uint16_t pos
Definition RawData.h:3
std::ostringstream debug
T getSigmaZ2() const
Definition BaseCluster.h:66
T getX() const
Definition BaseCluster.h:62
T getSigmaY2() const
Definition BaseCluster.h:65
T getY() const
Definition BaseCluster.h:63
std::int16_t getSensorID() const
Definition BaseCluster.h:81
T getZ() const
Definition BaseCluster.h:64
int getLastFilledBC(int dir=-1) const
bool testBC(int bcID, int dir=-1) const
int getFirstFilledBC(int dir=-1) const
void setFakeFlag(bool v=true)
bool isCorrect() const
Definition MCCompLabel.h:80
const auto & getMFTTrackID() const
void setBunchFilling(const o2::BunchFilling &bf)
set Bunch filling and init helpers for validation by BCs
void setMFTROFrameLengthMUS(float fums)
set MFT ROFrame duration in BC (continuous mode only)
void setMFTROFrameLengthInBC(int nbc)
set MFT ROFrame bias in BC (continuous mode only) or time shift applied already as MFTAlpideParam....
void run(const o2::globaltracking::RecoContainer &inp)
@ MATCHINGFUNC
MFT-MCH matching modes.
@ MATCHINGUPSTREAM
MFT-MCH track matching loaded from input file.
o2::mch::TrackParam FwdtoMCH(const o2::dataformats::GlobalFwdTrack &fwdtrack)
Converts FwdTrack parameters to MCH coordinate system.
o2::dataformats::GlobalFwdTrack MCHtoFwd(const o2::mch::TrackParam &mchTrack)
Converts mchTrack parameters to Forward coordinate system.
static constexpr std::array< int, NChips > ChipID2Layer
static bool extrapToVertexWithoutBranson(TrackParam &trackParam, double zVtx)
Definition TrackExtrap.h:73
track parameters for internal use
Definition TrackParam.h:34
Double_t getInverseBendingMomentum() const
return inverse bending momentum (GeV/c ** -1) times the charge (assumed forward motion)
Definition TrackParam.h:67
Double_t getNonBendingCoor() const
return non bending coordinate (cm)
Definition TrackParam.h:51
Double_t getZ() const
return Z coordinate (cm)
Definition TrackParam.h:47
Double_t getNonBendingSlope() const
return non bending slope (cm ** -1)
Definition TrackParam.h:55
const TMatrixD & getCovariances() const
Double_t getBendingCoor() const
return bending coordinate (cm)
Definition TrackParam.h:59
Double_t getBendingSlope() const
return bending slope (cm ** -1)
Definition TrackParam.h:63
Double_t getCharge() const
return the charge (assumed forward motion)
Definition TrackParam.h:71
void setCovariances(const SMatrix55Sym &covariances)
Definition TrackFwd.h:149
Double_t getSigma2Y() const
Definition TrackFwd.h:153
bool update(const std::array< float, 2 > &p, const std::array< float, 2 > &cov)
Definition TrackFwd.cxx:300
const SMatrix55Sym & getCovariances() const
Definition TrackFwd.h:148
Double_t getSigma2InvQPt() const
Definition TrackFwd.h:157
Double_t getSigma2Phi() const
Definition TrackFwd.h:155
Double_t getSigma2Tanl() const
Definition TrackFwd.h:156
Double_t getSigma2X() const
Definition TrackFwd.h:152
Double_t getPx() const
Definition TrackFwd.h:80
void setCharge(Double_t charge)
set the charge (assumed forward motion)
Definition TrackFwd.h:98
void setTanl(Double_t tanl)
Definition TrackFwd.h:71
void setInvQPt(Double_t invqpt)
Definition TrackFwd.h:76
Double_t getTanl() const
Definition TrackFwd.h:72
Double_t getPhi() const
Definition TrackFwd.h:56
Double_t getZ() const
return Z coordinate (cm)
Definition TrackFwd.h:46
void setPhi(Double_t phi)
Definition TrackFwd.h:55
void setX(Double_t x)
Definition TrackFwd.h:50
Double_t getY() const
Definition TrackFwd.h:52
const SMatrix5 & getParameters() const
return track parameters
Definition TrackFwd.h:106
Double_t getP() const
Definition TrackFwd.h:83
Double_t getCharge() const
return the charge (assumed forward motion)
Definition TrackFwd.h:96
Double_t getX() const
Definition TrackFwd.h:49
void setZ(Double_t z)
set Z coordinate (cm)
Definition TrackFwd.h:48
Double_t getInvQPt() const
Definition TrackFwd.h:77
void setY(Double_t y)
Definition TrackFwd.h:53
Double_t getPz() const
Definition TrackFwd.h:82
Double_t getPy() const
Definition TrackFwd.h:81
bool match(const std::vector< std::string > &queries, const char *pattern)
Definition dcs-ccdb.cxx:229
GLsizeiptr size
Definition glcorearb.h:659
GLuint GLuint end
Definition glcorearb.h:469
GLboolean GLboolean GLboolean b
Definition glcorearb.h:1233
GLintptr offset
Definition glcorearb.h:660
GLenum GLfloat param
Definition glcorearb.h:271
GLboolean GLboolean GLboolean GLboolean a
Definition glcorearb.h:1233
constexpr double LHCBunchSpacingMUS
constexpr int LHCMaxBunches
constexpr double LHCBunchSpacingNS
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 into std::vector<o2::BaseCluster<float>>
Definition IOUtils.cxx:111
unsigned long int nFakes
Enum< T >::Iterator begin(Enum< T >)
Definition Defs.h:173
std::string asString() const
int64_t differenceInBC(const InteractionRecord &other) const
void compare(std::string_view s1, std::string_view s2)
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"