File indexing completed on 2025-08-05 08:18:02
0001 #include "TruthRecoTrackMatching.h"
0002 #include "TrkrClusterIsMatcher.h"
0003 #include "ClusKeyIter.h"
0004
0005 #include <g4detectors/PHG4TpcCylinderGeom.h>
0006 #include <g4detectors/PHG4TpcCylinderGeomContainer.h>
0007
0008 #include <g4main/PHG4TruthInfoContainer.h>
0009
0010 #include <g4tracking/EmbRecoMatch.h>
0011 #include <g4tracking/EmbRecoMatchContainer.h>
0012 #include <g4tracking/EmbRecoMatchContainerv1.h>
0013 #include <g4tracking/EmbRecoMatchv1.h>
0014 #include <g4tracking/TrkrTruthTrack.h>
0015 #include <g4tracking/TrkrTruthTrackContainer.h>
0016
0017 #include <trackbase/TpcDefs.h>
0018 #include <trackbase/TrkrCluster.h>
0019 #include <trackbase/TrkrClusterContainer.h>
0020 #include <trackbase/TrkrDefs.h>
0021
0022 #include <trackbase_historic/SvtxTrack.h>
0023 #include <trackbase_historic/SvtxTrackMap.h>
0024
0025 #include <fun4all/Fun4AllReturnCodes.h>
0026
0027 #include <phool/PHCompositeNode.h>
0028 #include <phool/PHDataNode.h> // for PHDataNode
0029 #include <phool/PHIODataNode.h>
0030 #include <phool/PHNode.h> // for PHNode
0031 #include <phool/PHNodeIterator.h>
0032 #include <phool/PHObject.h> // for PHObject
0033 #include <phool/PHRandomSeed.h>
0034 #include <phool/getClass.h>
0035 #include <phool/phool.h> // for PHWHERE
0036
0037 #include <TFile.h>
0038 #include <TSystem.h>
0039 #include <TTree.h>
0040
0041 #include <boost/format.hpp>
0042
0043 #include <algorithm>
0044
0045
0046
0047
0048 TruthRecoTrackMatching::TruthRecoTrackMatching(
0049 TrkrClusterIsMatcher* _ismatcher
0050 , const unsigned short _nmincluster_match
0051 , const float _nmincluster_ratio
0052 , const float _cutoff_dphi
0053 , const float _same_dphi
0054 , const float _cutoff_deta
0055 , const float _same_deta
0056 , const unsigned short _max_nreco_per_truth
0057 , const unsigned short _max_ntruth_per_reco)
0058 : m_ismatcher { _ismatcher }
0059 , m_nmincluster_match{_nmincluster_match}
0060 , m_nmincluster_ratio{_nmincluster_ratio}
0061
0062 , m_cutoff_dphi{_cutoff_dphi}
0063 , m_same_dphi{_same_dphi}
0064 , m_cutoff_deta{_cutoff_deta}
0065 , m_same_deta{_same_deta}
0066
0067
0068 , m_max_nreco_per_truth{_max_nreco_per_truth}
0069 , m_max_ntruth_per_reco{_max_ntruth_per_reco}
0070 {
0071 m_TCEval.ismatcher = m_ismatcher;
0072 if (Verbosity() > 50)
0073 {
0074 std::cout << " Starting TruthRecoTrackMatching.cc " << std::endl;
0075 }
0076 }
0077
0078 int TruthRecoTrackMatching::InitRun(PHCompositeNode* topNode)
0079 {
0080 if (Verbosity() > 10)
0081 {
0082 topNode->print();
0083 }
0084 auto init_status = m_ismatcher->init(topNode);
0085 if (init_status == Fun4AllReturnCodes::ABORTRUN)
0086 {
0087 return init_status;
0088 }
0089
0090 if (createNodes(topNode) != Fun4AllReturnCodes::EVENT_OK)
0091 {
0092 return Fun4AllReturnCodes::ABORTRUN;
0093 }
0094 m_nmatched_id_reco = &(m_EmbRecoMatchContainer->map_nTruthPerReco());
0095 m_nmatched_id_true = &(m_EmbRecoMatchContainer->map_nRecoPerTruth());
0096 return Fun4AllReturnCodes::EVENT_OK;
0097 }
0098
0099 int TruthRecoTrackMatching::process_event(PHCompositeNode* topnode)
0100 {
0101 if (topnode == nullptr)
0102 {
0103 return Fun4AllReturnCodes::ABORTRUN;
0104 }
0105 if (Verbosity() > 1000)
0106 {
0107 topnode->print();
0108 }
0109
0110 m_nmatched_index_true.clear();
0111
0112
0113
0114
0115
0116
0117
0118
0119
0120
0121
0122
0123
0124
0125
0126 if (Verbosity() > 60)
0127 {
0128 std::cout << "reco tracks size: " << m_SvtxTrackMap->size() << std::endl;
0129 }
0130
0131 recoData.clear();
0132
0133 for (auto& reco : *m_SvtxTrackMap)
0134 {
0135 auto index_reco = reco.first;
0136 auto track = reco.second;
0137 recoData.push_back({track->get_phi(), track->get_eta(), track->get_pt(), index_reco});
0138 }
0139
0140 std::sort(recoData.begin(), recoData.end(), CompRECOtoPhi());
0141
0142
0143
0144
0145
0146
0147
0148
0149
0150
0151
0152
0153
0154
0155 RECOvec wraps{};
0156 auto wrap_from_start = std::upper_bound(recoData.begin(),
0157 recoData.end(), (-M_PI + m_cutoff_dphi), CompRECOtoPhi());
0158 for (auto iter = recoData.begin(); iter != wrap_from_start; ++iter)
0159 {
0160 auto entry = *iter;
0161 std::get<RECOphi>(entry) = std::get<RECOphi>(entry) + 2 * M_PI;
0162 wraps.push_back(entry);
0163 }
0164 for (auto E : wraps)
0165 {
0166 recoData.push_back(E);
0167 }
0168
0169 if (Verbosity() > 70)
0170 {
0171 std::cout << " ****************************************** " << std::endl;
0172 std::cout << " This is the RECO map of tracks to match to " << std::endl;
0173 std::cout << " ****************************************** " << std::endl;
0174 for (auto& E : recoData)
0175 {
0176 std::cout << (boost::format(" id:%2i (phi:eta:pt) (%5.2f:%5.2f:%5.2f)") %(std::get<RECOid>(E))
0177 %(std::get<RECOphi>(E)) %(std::get<RECOeta>(E)) %(std::get<RECOpt>(E))).str()
0178 << std::endl;
0179 }
0180 std::cout << " ****end*listing*********************** " << std::endl;
0181 }
0182
0183
0184
0185
0186
0187
0188
0189
0190
0191
0192
0193
0194
0195
0196 std::vector<std::pair<unsigned short, unsigned short>> outerbox_pairs{};
0197 std::vector<std::pair<unsigned short, unsigned short>> innerbox_pairs{};
0198
0199 if (topnode == nullptr)
0200 {
0201 return Fun4AllReturnCodes::ABORTRUN;
0202 }
0203 if (Verbosity() > 70)
0204 {
0205 std::cout << "Number of truth tracks: " << m_TrkrTruthTrackContainer->getMap().size() << std::endl;
0206 for (auto& _pair : m_TrkrTruthTrackContainer->getMap())
0207 {
0208 auto& track = _pair.second;
0209 std::cout << (boost::format(" id:%2i (phi:eta:pt) (%5.2f:%5.2f:%5.2f nclus: %i)") %track->getTrackid()
0210 %track->getPhi() %track->getPseudoRapidity() %track->getPt()
0211 %((int) track->getClusters().size())).str()
0212 << std::endl;
0213 }
0214 std::cout << " ----end-listing-truth-tracks---------- " << std::endl;
0215 }
0216
0217 for (auto _pair : m_TrkrTruthTrackContainer->getMap())
0218 {
0219 auto id_true = _pair.first;
0220 auto track = _pair.second;
0221
0222 auto match_indices = find_box_matches(track->getPhi(),
0223 track->getPseudoRapidity(), track->getPt());
0224
0225
0226 for (auto& id_reco : match_indices.first)
0227 {
0228 innerbox_pairs.emplace_back(id_true, id_reco);
0229 }
0230 for (auto& id_reco : match_indices.second)
0231 {
0232 outerbox_pairs.emplace_back(id_true, id_reco);
0233 }
0234
0235 if (Verbosity() > 80)
0236 {
0237 std::cout << (boost::format(" T(%i) find box(%5.2f:%5.2f:%5.2f)")
0238 %((int) track->getTrackid()) %track->getPhi()
0239 %track->getPseudoRapidity() %track->getPt()).str();
0240 for (auto& id_reco : match_indices.first)
0241 {
0242 std::cout << "->IB(" << id_reco << ")";
0243 }
0244 for (auto& id_reco : match_indices.second)
0245 {
0246 std::cout << "->OB(" << id_reco << ")";
0247 }
0248 std::cout << std::endl;
0249 }
0250 }
0251
0252 if (Verbosity() > 100)
0253 {
0254 std::cout << "Innerbox pairs" << std::endl;
0255 }
0256 match_tracks_in_box(innerbox_pairs);
0257 if (Verbosity() > 100)
0258 {
0259 std::cout << "Outerbox pairs" << std::endl;
0260 }
0261 match_tracks_in_box(outerbox_pairs);
0262
0263
0264 for (auto _pair : m_TrkrTruthTrackContainer->getMap())
0265 {
0266 auto id_true = _pair.first;
0267 m_EmbRecoMatchContainer->checkfill_idsTruthUnmatched(id_true);
0268 }
0269
0270 if (Verbosity() > 50)
0271 {
0272 std::cout << " --START--print-all-stored-matches--" << std::endl;
0273 std::cout << " --0-- Printing all matches stored (start)" << std::endl;
0274
0275 for (auto match : m_EmbRecoMatchContainer->getMatches())
0276 {
0277 std::cout << (boost::format(" Match id(%2i->%2i) nClusMatch-nClusTrue-nClusReco (%2i:%2i:%2i)")
0278 %match->idTruthTrack() %match->idRecoTrack()
0279 %match->nClustersMatched() %match->nClustersTruth() %match->nClustersReco()).str()
0280 << std::endl;
0281 }
0282 std::cout << " --END--print-all-stored-matches----" << std::endl;
0283 }
0284
0285 if (m_write_diag)
0286 {
0287 fill_tree();
0288 }
0289
0290 return Fun4AllReturnCodes::EVENT_OK;
0291 }
0292
0293 int TruthRecoTrackMatching::End(PHCompositeNode* )
0294 {
0295 TFile* s_current = gDirectory->GetFile();
0296 if (m_write_diag)
0297 {
0298 m_diag_file->cd();
0299
0300
0301
0302
0303
0304
0305
0306
0307
0308
0309
0310
0311 m_diag_tree->Write();
0312 m_diag_file->Save();
0313 m_diag_file->Close();
0314 if (s_current != nullptr)
0315 {
0316 s_current->cd();
0317 }
0318 }
0319
0320 return Fun4AllReturnCodes::EVENT_OK;
0321 }
0322
0323
0324
0325
0326
0327 int TruthRecoTrackMatching::createNodes(PHCompositeNode* topNode)
0328 {
0329 PHNodeIterator iter(topNode);
0330 PHCompositeNode* dstNode = dynamic_cast<PHCompositeNode*>(iter.findFirst("PHCompositeNode", "DST"));
0331 if (!dstNode)
0332 {
0333 std::cout << PHWHERE << "DST Node missing, doing nothing." << std::endl;
0334 exit(1);
0335 }
0336
0337
0338 m_PHG4TruthInfoContainer = findNode::getClass<PHG4TruthInfoContainer>(topNode, "G4TruthInfo");
0339 if (!m_PHG4TruthInfoContainer)
0340 {
0341 std::cout << "Could not locate G4TruthInfo node when running "
0342 << "\"TruthRecoTrackMatching\" module." << std::endl;
0343 return Fun4AllReturnCodes::ABORTEVENT;
0344 }
0345
0346 m_SvtxTrackMap = findNode::getClass<SvtxTrackMap>(topNode, "SvtxTrackMap");
0347 if (!m_SvtxTrackMap)
0348 {
0349 std::cout << "Could not locate SvtxTrackMap node when running "
0350 << "\"TruthRecoTrackMatching\" module." << std::endl;
0351 return Fun4AllReturnCodes::ABORTEVENT;
0352 }
0353
0354 m_TruthClusterContainer = findNode::getClass<TrkrClusterContainer>(topNode,
0355 "TRKR_TRUTHCLUSTERCONTAINER");
0356 if (!m_TruthClusterContainer)
0357 {
0358 std::cout << "Could not locate TRKR_TRUTHCLUSTERCONTAINER node when "
0359 << "running \"TruthRecoTrackMatching\" module." << std::endl;
0360 return Fun4AllReturnCodes::ABORTEVENT;
0361 }
0362
0363 m_RecoClusterContainer = findNode::getClass<TrkrClusterContainer>(topNode, "TRKR_CLUSTER");
0364 if (!m_RecoClusterContainer)
0365 {
0366 std::cout << "Could not locate TRKR_CLUSTER node when running "
0367 << "\"TruthRecoTrackMatching\" module." << std::endl;
0368 return Fun4AllReturnCodes::ABORTEVENT;
0369 }
0370
0371 m_TrkrTruthTrackContainer = findNode::getClass<TrkrTruthTrackContainer>(topNode,
0372 "TRKR_TRUTHTRACKCONTAINER");
0373 if (!m_TrkrTruthTrackContainer)
0374 {
0375 std::cout << "Could not locate TRKR_TRUTHTRACKCONTAINER node when "
0376 << "running \"TruthRecoTrackMatching\" module." << std::endl;
0377 return Fun4AllReturnCodes::ABORTEVENT;
0378 }
0379
0380 m_PHG4TpcCylinderGeomContainer =
0381 findNode::getClass<PHG4TpcCylinderGeomContainer>(topNode, "CYLINDERCELLGEOM_SVTX");
0382 if (!m_PHG4TpcCylinderGeomContainer)
0383 {
0384 std::cout << "Could not locate CYLINDERCELLGEOM_SVTX node when "
0385 << "running \"TruthRecoTrackMatching\" module." << std::endl;
0386
0387 return Fun4AllReturnCodes::ABORTEVENT;
0388 }
0389
0390
0391
0392
0393
0394
0395
0396 m_EmbRecoMatchContainer = findNode::getClass<EmbRecoMatchContainer>(topNode,
0397 "TRKR_EMBRECOMATCHCONTAINER");
0398 if (!m_EmbRecoMatchContainer)
0399 {
0400 PHNodeIterator dstiter(dstNode);
0401
0402 auto DetNode = dynamic_cast<PHCompositeNode*>(dstiter.findFirst("PHCompositeNode", "TRKR"));
0403 if (!DetNode)
0404 {
0405 DetNode = new PHCompositeNode("TRKR");
0406 dstNode->addNode(DetNode);
0407 }
0408
0409 m_EmbRecoMatchContainer = new EmbRecoMatchContainerv1;
0410 auto newNode = new PHIODataNode<PHObject>(m_EmbRecoMatchContainer,
0411 "TRKR_EMBRECOMATCHCONTAINER", "PHObject");
0412 DetNode->addNode(newNode);
0413 }
0414
0415 m_ActsGeometry = findNode::getClass<ActsGeometry>(topNode, "ActsGeometry");
0416
0417 return Fun4AllReturnCodes::EVENT_OK;
0418 }
0419
0420 std::pair<std::vector<unsigned short>, std::vector<unsigned short>>
0421 TruthRecoTrackMatching::find_box_matches(float truth_phi, float truth_eta, float truth_pt)
0422 {
0423
0424
0425
0426
0427
0428
0429 RECO_pair_iter mix_pair{recoData.begin(), recoData.end()};
0430
0431 mix_pair.first = std::lower_bound(mix_pair.first, mix_pair.second,
0432 truth_phi - m_cutoff_dphi, CompRECOtoPhi());
0433
0434 mix_pair.second = std::upper_bound(mix_pair.first, mix_pair.second,
0435 truth_phi + m_cutoff_dphi, CompRECOtoPhi());
0436
0437 if (mix_pair.first == mix_pair.second)
0438 {
0439
0440 return {{}, {}};
0441 }
0442
0443
0444
0445
0446 std::sort(mix_pair.first, mix_pair.second, CompRECOtoEta());
0447
0448 RECO_pair_iter outer_box = mix_pair;
0449 outer_box.first = std::lower_bound(mix_pair.first, mix_pair.second,
0450 truth_eta - m_cutoff_deta, CompRECOtoEta());
0451
0452 outer_box.second = std::lower_bound(outer_box.first, outer_box.second,
0453 truth_eta + m_cutoff_deta, CompRECOtoEta());
0454
0455 if (outer_box.first == outer_box.second)
0456 {
0457
0458 std::sort(mix_pair.first, mix_pair.second, CompRECOtoPhi());
0459 return {{}, {}};
0460 }
0461
0462
0463 RECO_pair_iter inner_box = outer_box;
0464 const float _delta_outer_pt = delta_outer_pt(truth_pt);
0465 const float _delta_inner_pt = delta_inner_pt(truth_pt);
0466 if (_delta_outer_pt > 0)
0467 {
0468 std::sort(outer_box.first, outer_box.second, CompRECOtoPt());
0469 outer_box.first = std::lower_bound(outer_box.first, outer_box.second, truth_pt - _delta_outer_pt, CompRECOtoPt());
0470 outer_box.second = std::upper_bound(outer_box.first, outer_box.second, truth_pt + _delta_outer_pt, CompRECOtoPt());
0471
0472
0473 if (outer_box.first == outer_box.second)
0474 {
0475
0476 std::sort(mix_pair.first, mix_pair.second, CompRECOtoPhi());
0477 return {{}, {}};
0478 }
0479
0480
0481 inner_box = outer_box;
0482 if (_delta_inner_pt > 0)
0483 {
0484
0485 inner_box.first = std::lower_bound(inner_box.first,
0486 inner_box.second, truth_pt - _delta_inner_pt, CompRECOtoPt());
0487 inner_box.second = std::upper_bound(inner_box.first,
0488 inner_box.second, truth_pt + _delta_inner_pt, CompRECOtoPt());
0489 }
0490
0491
0492 std::sort(inner_box.first, inner_box.second, CompRECOtoEta());
0493 }
0494
0495
0496
0497
0498
0499 if (inner_box.first != inner_box.second)
0500 {
0501 inner_box.first = std::lower_bound(inner_box.first, inner_box.second,
0502 truth_eta - m_same_deta, CompRECOtoEta());
0503
0504 inner_box.second = std::lower_bound(inner_box.first, inner_box.second,
0505 truth_eta + m_same_deta, CompRECOtoEta());
0506 }
0507
0508
0509 if (inner_box.first != inner_box.second)
0510 {
0511 std::sort(inner_box.first, inner_box.second, CompRECOtoPhi());
0512
0513 inner_box.first = std::lower_bound(inner_box.first, inner_box.second,
0514 truth_phi - m_same_dphi, CompRECOtoPhi());
0515
0516 inner_box.second = std::lower_bound(inner_box.first, inner_box.second,
0517 truth_phi + m_same_dphi, CompRECOtoPhi());
0518 }
0519
0520
0521
0522
0523 std::vector<unsigned short> inner_box_matches, outer_box_matches;
0524
0525
0526 for (auto iter = inner_box.first; iter != inner_box.second; ++iter)
0527 {
0528 inner_box_matches.push_back(std::get<RECOid>(*iter));
0529 }
0530
0531 for (auto iter = outer_box.first; iter != inner_box.first; ++iter)
0532 {
0533 outer_box_matches.push_back(std::get<RECOid>(*iter));
0534 }
0535 for (auto iter = inner_box.second; iter != outer_box.second; ++iter)
0536 {
0537 outer_box_matches.push_back(std::get<RECOid>(*iter));
0538 }
0539
0540
0541 std::sort(mix_pair.first, mix_pair.second, CompRECOtoPhi());
0542
0543
0544 return std::make_pair(inner_box_matches, outer_box_matches);
0545 }
0546
0547 void TruthRecoTrackMatching::match_tracks_in_box(
0548 std::vector<std::pair<unsigned short, unsigned short>>& box_pairs
0549 )
0550 {
0551 if (box_pairs.size() == 0)
0552 {
0553 return;
0554 }
0555
0556 std::sort(box_pairs.begin(), box_pairs.end());
0557 std::vector<PossibleMatch> poss_matches;
0558
0559 auto ipair = box_pairs.begin();
0560 while (ipair != box_pairs.end())
0561 {
0562 auto id_true = ipair->first;
0563
0564
0565 if (at_nmax_index_true(id_true))
0566 {
0567 while (ipair != box_pairs.end())
0568 {
0569 ++ipair;
0570 if (ipair->first != id_true)
0571 {
0572 break;
0573 }
0574 }
0575 continue;
0576 }
0577
0578
0579
0580 m_TCEval.addClusKeys(m_TrkrTruthTrackContainer->getTruthTrack(id_true));
0581
0582
0583
0584
0585
0586
0587
0588 while (ipair != box_pairs.end())
0589 {
0590
0591 if (ipair->first != id_true)
0592 {
0593 break;
0594 }
0595 if (at_nmax_id_reco(ipair->second))
0596 {
0597 ++ipair;
0598 continue;
0599 }
0600
0601 SvtxTrack* reco_track = m_SvtxTrackMap->get(ipair->second);
0602 m_TCEval.addClusKeys(reco_track);
0603 m_TCEval.find_matches();
0604
0605
0606
0607
0608
0609
0610
0611
0612
0613
0614
0615
0616
0617
0618
0619
0620
0621
0622
0623 unsigned short nclus_match = m_TCEval.phg4_n_matched();
0624 unsigned short nclus_true = m_TCEval.phg4_nclus();
0625 unsigned short nclus_reco = m_TCEval.svtx_nclus();
0626 unsigned short nclus_nomatch = (int) (m_TCEval.svtx_keys.size() - nclus_match);
0627
0628 if (Verbosity() > 100)
0629 {
0630 auto truth_track = m_TrkrTruthTrackContainer->getTruthTrack(id_true);
0631 std::cout << (boost::format(
0632 "possmatch:(phi,eta,pT:id) true(%5.2f,%5.2f,%4.2f:%2i) reco(%5.2f,%5.2f,%4.2f:%2i) "
0633 "nCl(match:true:reco:nomatch)(%2i-%2i-%2i-%2i)")
0634 %truth_track->getPhi() %truth_track->getPseudoRapidity() %truth_track->getPt()
0635 %((int) truth_track->getTrackid())
0636 %reco_track->get_phi() %reco_track->get_eta() %reco_track->get_pt()
0637 %ipair->second
0638 %((int) nclus_match) %((int) nclus_true) %((int) nclus_reco) %((int) nclus_nomatch)).str()
0639 << std::endl;
0640 }
0641 if (nclus_match >= m_nmincluster_match && (static_cast<float>(nclus_match) / nclus_true >= m_nmincluster_ratio))
0642 {
0643 poss_matches.push_back(
0644 {nclus_match, nclus_true, nclus_reco,
0645 ipair->first, ipair->second});
0646 }
0647 ++ipair;
0648 }
0649 }
0650
0651
0652
0653
0654
0655 std::sort(poss_matches.begin(), poss_matches.end(), SortPossibleMatch());
0656
0657 if (Verbosity() > 200)
0658 {
0659 std::cout << " All chosen possible matches (" << poss_matches.size() << ") track pairs (nClMatched-nClTrue-nClReco : idTrue-idReco)" << std::endl;
0660 int i{0};
0661 for (auto match : poss_matches)
0662 {
0663
0664
0665 std::cout << (boost::format(" pair(%2i): %2i-%2i-%2i-<%2i>-%2i ") %(i++) %((int) std::get<PM_nmatch>(match)) %((int) std::get<PM_ntrue>(match)) %((int) std::get<PM_nreco>(match)) %((int) std::get<PM_idtrue>(match)) %((int) std::get<PM_idreco>(match))).str() << std::endl;
0666 }
0667 }
0668
0669
0670 auto iter = poss_matches.begin();
0671 while (iter != poss_matches.end())
0672 {
0673 if (skip_match(*iter))
0674 {
0675 ++iter;
0676 continue;
0677 }
0678 std::vector<std::pair<float, int>> sigma_metric = {{0., 0}};
0679 int n_sigma = 0;
0680 auto iter_same = iter + 1;
0681 while (
0682 iter_same != poss_matches.end() && (*iter_same)[PM_nmatch] == (*iter)[PM_nmatch] && (*iter_same)[PM_ntrue] == (*iter)[PM_ntrue] && (*iter_same)[PM_nreco] == (*iter)[PM_nreco])
0683 {
0684 ++n_sigma;
0685 if (n_sigma == 1)
0686 {
0687 sigma_metric[0].first = sigma_CompMatchClusters(*iter);
0688 }
0689 sigma_metric.emplace_back(sigma_CompMatchClusters(*iter_same), n_sigma);
0690 ++iter_same;
0691 }
0692 std::sort(sigma_metric.begin(), sigma_metric.end());
0693
0694 bool first = true;
0695 for (auto& sigma : sigma_metric)
0696 {
0697 if (first)
0698 {
0699 first = false;
0700 }
0701 else if (skip_match(*(iter + sigma.second)))
0702 {
0703 continue;
0704 }
0705 auto match = *(iter + sigma.second);
0706 auto id_true = match[PM_idtrue];
0707 auto id_reco = match[PM_idreco];
0708 auto save_match = new EmbRecoMatchv1(id_true, id_reco,
0709 match[PM_ntrue], match[PM_nreco], match[PM_nmatch]);
0710 m_EmbRecoMatchContainer->addMatch(save_match);
0711
0712 if (m_nmatched_index_true.find(id_true) == m_nmatched_index_true.end())
0713 {
0714 m_nmatched_index_true[id_true] = 1;
0715 }
0716 else
0717 {
0718 m_nmatched_index_true[id_true] += 1;
0719 }
0720 }
0721 iter += sigma_metric.size();
0722 }
0723 }
0724
0725
0726
0727
0728 inline float TruthRecoTrackMatching::abs_dphi(float phi0, float phi1)
0729 {
0730 float dphi = std::fabs(phi0 - phi1);
0731 while (dphi > M_PI)
0732 {
0733 dphi = std::fabs(dphi - 2 * M_PI);
0734 }
0735 return dphi;
0736 }
0737
0738 float TruthRecoTrackMatching::sigma_CompMatchClusters(PossibleMatch& match)
0739 {
0740 auto id_true = match[PM_idtrue];
0741 auto id_reco = match[PM_idreco];
0742
0743 m_TCEval.addClusKeys(m_TrkrTruthTrackContainer->getTruthTrack(id_true));
0744 m_TCEval.addClusKeys(m_SvtxTrackMap->get(id_reco));
0745
0746 m_TCEval.find_matches();
0747
0748 int n_matched = m_TCEval.phg4_n_matched();
0749
0750 if (n_matched)
0751 {
0752 return std::numeric_limits<float>::max();
0753 }
0754 else
0755 {
0756 return m_TCEval.match_stat / n_matched;
0757 }
0758 }
0759
0760
0761
0762
0763
0764
0765
0766
0767
0768
0769
0770
0771
0772
0773
0774
0775
0776
0777
0778
0779
0780
0781
0782
0783
0784
0785
0786
0787 inline bool TruthRecoTrackMatching::skip_match(PossibleMatch& match)
0788 {
0789 return at_nmax_id_reco(match[PM_idreco]) || at_nmax_index_true(match[PM_idreco]);
0790 }
0791
0792
0793
0794
0795
0796
0797
0798
0799 float TruthRecoTrackMatching::delta_outer_pt(float pt) const
0800 {
0801 return -1. + pt * 0.;
0802 }
0803 float TruthRecoTrackMatching::delta_inner_pt(float pt) const
0804 {
0805 return -1. + pt * 0.;
0806 }
0807
0808 bool TruthRecoTrackMatching::at_nmax_index_true(unsigned short index_true)
0809 {
0810 if (m_nmatched_index_true.find(index_true) == m_nmatched_index_true.end())
0811 {
0812 return false;
0813 }
0814 return m_nmatched_index_true[index_true] >= m_max_nreco_per_truth;
0815 }
0816
0817 bool TruthRecoTrackMatching::at_nmax_id_reco(unsigned short id_reco)
0818 {
0819 if (m_nmatched_id_reco->find(id_reco) == m_nmatched_id_reco->end())
0820 {
0821 return false;
0822 }
0823 return (*m_nmatched_id_reco)[id_reco] >= m_max_ntruth_per_reco;
0824 }
0825
0826 void TruthRecoTrackMatching::set_diagnostic_file(const std::string& file_name)
0827 {
0828 m_write_diag = true;
0829 TFile* s_current = gDirectory->GetFile();
0830 m_diag_file = new TFile(file_name.c_str(), "recreate");
0831
0832 m_diag_file->cd();
0833 m_diag_tree = new TTree("T", "Tree of Reco and True Clusters");
0834
0835 m_diag_tree->Branch("event", &m_event);
0836
0837 m_diag_tree->Branch("trkid_reco_matched", &m_trkid_reco_matched);
0838 m_diag_tree->Branch("itrk_reco_matched", &m_cnt_reco_matched);
0839 m_diag_tree->Branch("i0_reco_matched", &m_i0_reco_matched);
0840 m_diag_tree->Branch("i1_reco_matched", &m_i1_reco_matched);
0841 m_diag_tree->Branch("layer_reco_matched", &m_layer_reco_matched);
0842 m_diag_tree->Branch("x_reco_matched", &m_x_reco_matched);
0843 m_diag_tree->Branch("y_reco_matched", &m_y_reco_matched);
0844 m_diag_tree->Branch("z_reco_matched", &m_z_reco_matched);
0845
0846 m_diag_tree->Branch("trkid_reco_notmatched", &m_trkid_reco_notmatched);
0847 m_diag_tree->Branch("itrk_reco_notmatched", &m_cnt_reco_notmatched);
0848 m_diag_tree->Branch("i0_reco_notmatched", &m_i0_reco_notmatched);
0849 m_diag_tree->Branch("i1_reco_notmatched", &m_i1_reco_notmatched);
0850 m_diag_tree->Branch("layer_reco_notmatched", &m_layer_reco_notmatched);
0851 m_diag_tree->Branch("x_reco_notmatched", &m_x_reco_notmatched);
0852 m_diag_tree->Branch("y_reco_notmatched", &m_y_reco_notmatched);
0853 m_diag_tree->Branch("z_reco_notmatched", &m_z_reco_notmatched);
0854
0855 m_diag_tree->Branch("trkid_true_matched", &m_trkid_true_matched);
0856 m_diag_tree->Branch("itrk_true_matched", &m_cnt_true_matched);
0857 m_diag_tree->Branch("i0_true_matched", &m_i0_true_matched);
0858 m_diag_tree->Branch("i1_true_matched", &m_i1_true_matched);
0859 m_diag_tree->Branch("layer_true_matched", &m_layer_true_matched);
0860 m_diag_tree->Branch("x_true_matched", &m_x_true_matched);
0861 m_diag_tree->Branch("y_true_matched", &m_y_true_matched);
0862 m_diag_tree->Branch("z_true_matched", &m_z_true_matched);
0863
0864 m_diag_tree->Branch("trkid_true_notmatched", &m_trkid_true_notmatched);
0865 m_diag_tree->Branch("itrk_true_notmatched", &m_cnt_true_notmatched);
0866 m_diag_tree->Branch("i0_true_notmatched", &m_i0_true_notmatched);
0867 m_diag_tree->Branch("i1_true_notmatched", &m_i1_true_notmatched);
0868 m_diag_tree->Branch("layer_true_notmatched", &m_layer_true_notmatched);
0869 m_diag_tree->Branch("x_true_notmatched", &m_x_true_notmatched);
0870 m_diag_tree->Branch("y_true_notmatched", &m_y_true_notmatched);
0871 m_diag_tree->Branch("z_true_notmatched", &m_z_true_notmatched);
0872
0873 if (s_current != nullptr)
0874 {
0875 s_current->cd();
0876 }
0877 }
0878
0879 void TruthRecoTrackMatching::clear_branch_vectors()
0880 {
0881 m_trkid_reco_matched.clear();
0882 m_i0_reco_matched.clear();
0883 m_i1_reco_matched.clear();
0884 m_layer_reco_matched.clear();
0885 m_x_reco_matched.clear();
0886 m_y_reco_matched.clear();
0887 m_z_reco_matched.clear();
0888
0889 m_trkid_reco_notmatched.clear();
0890 m_i0_reco_notmatched.clear();
0891 m_i1_reco_notmatched.clear();
0892 m_layer_reco_notmatched.clear();
0893 m_x_reco_notmatched.clear();
0894 m_y_reco_notmatched.clear();
0895 m_z_reco_notmatched.clear();
0896
0897 m_trkid_true_matched.clear();
0898 m_i0_true_matched.clear();
0899 m_i1_true_matched.clear();
0900 m_layer_true_matched.clear();
0901 m_x_true_matched.clear();
0902 m_y_true_matched.clear();
0903 m_z_true_matched.clear();
0904
0905 m_trkid_true_notmatched.clear();
0906 m_i0_true_notmatched.clear();
0907 m_i1_true_notmatched.clear();
0908 m_layer_true_notmatched.clear();
0909 m_x_true_notmatched.clear();
0910 m_y_true_notmatched.clear();
0911 m_z_true_notmatched.clear();
0912 }
0913
0914 void TruthRecoTrackMatching::fill_tree()
0915 {
0916
0917 int cnt = 0;
0918 int itrk = 0;
0919 for (auto& trkid : m_EmbRecoMatchContainer->ids_TruthUnmatched())
0920 {
0921 m_trkid_true_notmatched.push_back(trkid);
0922 m_i0_true_notmatched.push_back(cnt);
0923 auto track = m_TrkrTruthTrackContainer->getTruthTrack(trkid);
0924 for (auto& ckey : track->getClusters())
0925 {
0926 auto cluster = m_TruthClusterContainer->findCluster(ckey);
0927 m_cnt_true_notmatched.push_back(itrk);
0928 Eigen::Vector3d gloc = m_ActsGeometry->getGlobalPosition(ckey, cluster);
0929 m_layer_true_notmatched.push_back(TrkrDefs::getLayer(ckey));
0930 m_x_true_notmatched.push_back(gloc[0]);
0931 m_y_true_notmatched.push_back(gloc[1]);
0932 m_z_true_notmatched.push_back(gloc[2]);
0933 ++cnt;
0934 }
0935 m_i1_true_notmatched.push_back(cnt);
0936 ++itrk;
0937 }
0938
0939
0940 cnt = 0;
0941 itrk = 0;
0942 for (auto& trkid : m_EmbRecoMatchContainer->ids_TruthMatched())
0943 {
0944 m_trkid_true_matched.push_back(trkid);
0945 m_i0_true_matched.push_back(cnt);
0946 auto track = m_TrkrTruthTrackContainer->getTruthTrack(trkid);
0947 for (auto& ckey : track->getClusters())
0948 {
0949 auto cluster = m_TruthClusterContainer->findCluster(ckey);
0950 m_cnt_true_matched.push_back(itrk);
0951 Eigen::Vector3d gloc = m_ActsGeometry->getGlobalPosition(ckey, cluster);
0952 m_layer_true_matched.push_back(TrkrDefs::getLayer(ckey));
0953 m_x_true_matched.push_back(gloc[0]);
0954 m_y_true_matched.push_back(gloc[1]);
0955 m_z_true_matched.push_back(gloc[2]);
0956 ++cnt;
0957 }
0958 m_i1_true_matched.push_back(cnt);
0959 ++itrk;
0960 }
0961
0962
0963 std::set<unsigned int> set_reco_matched;
0964 cnt = 0;
0965 itrk = 0;
0966 for (auto& trkid : m_EmbRecoMatchContainer->ids_RecoMatched())
0967 {
0968 set_reco_matched.insert(trkid);
0969 m_trkid_reco_matched.push_back(trkid);
0970 SvtxTrack* reco_track = m_SvtxTrackMap->get(trkid);
0971 m_i0_reco_matched.push_back(cnt);
0972
0973 for (auto reco_ckey : ClusKeyIter(reco_track))
0974 {
0975 auto cluster = m_RecoClusterContainer->findCluster(reco_ckey);
0976 Eigen::Vector3d gloc = m_ActsGeometry->getGlobalPosition(reco_ckey, cluster);
0977 m_layer_reco_matched.push_back(TrkrDefs::getLayer(reco_ckey));
0978 m_cnt_reco_matched.push_back(itrk);
0979 m_x_reco_matched.push_back(gloc[0]);
0980 m_y_reco_matched.push_back(gloc[1]);
0981 m_z_reco_matched.push_back(gloc[2]);
0982 ++cnt;
0983 }
0984 m_i1_reco_matched.push_back(cnt);
0985 ++itrk;
0986 }
0987
0988
0989 cnt = 0;
0990 itrk = 0;
0991 for (auto& reco : *m_SvtxTrackMap)
0992 {
0993 auto trkid = reco.first;
0994 if (set_reco_matched.count(trkid))
0995 {
0996 continue;
0997 }
0998 m_trkid_reco_notmatched.push_back(trkid);
0999 SvtxTrack* reco_track = m_SvtxTrackMap->get(trkid);
1000 m_i0_reco_notmatched.push_back(cnt);
1001 for (auto reco_ckey : ClusKeyIter(reco_track))
1002 {
1003 auto cluster = m_RecoClusterContainer->findCluster(reco_ckey);
1004 Eigen::Vector3d gloc = m_ActsGeometry->getGlobalPosition(reco_ckey, cluster);
1005 m_layer_reco_notmatched.push_back(TrkrDefs::getLayer(reco_ckey));
1006 m_cnt_reco_notmatched.push_back(itrk);
1007 m_x_reco_notmatched.push_back(gloc[0]);
1008 m_y_reco_notmatched.push_back(gloc[1]);
1009 m_z_reco_notmatched.push_back(gloc[2]);
1010 ++cnt;
1011 }
1012 m_i1_reco_notmatched.push_back(cnt);
1013 ++itrk;
1014 }
1015 ++m_event;
1016 m_diag_tree->Fill();
1017 clear_branch_vectors();
1018 return;
1019 }