File indexing completed on 2025-08-05 08:17:56
0001 #include "g4evaltools.h"
0002
0003 #include <mvtx/CylinderGeom_Mvtx.h>
0004
0005 #include <intt/CylinderGeomIntt.h>
0006
0007 #include <g4detectors/PHG4CylinderGeomContainer.h>
0008 #include <g4detectors/PHG4TpcCylinderGeom.h>
0009 #include <g4detectors/PHG4TpcCylinderGeomContainer.h>
0010
0011 #include <g4tracking/EmbRecoMatchContainer.h>
0012 #include <g4tracking/TrkrTruthTrack.h>
0013
0014 #include <trackbase_historic/SvtxTrack.h>
0015 #include <trackbase_historic/SvtxTrackMap.h>
0016
0017 #include <trackbase/ActsGeometry.h>
0018 #include <trackbase/TrkrCluster.h>
0019 #include <trackbase/TrkrClusterContainer.h>
0020
0021 #include <phool/PHCompositeNode.h>
0022 #include <phool/PHNode.h> // for PHNode
0023 #include <phool/PHNodeIterator.h>
0024 #include <phool/PHObject.h> // for PHObject
0025 #include <phool/getClass.h>
0026 #include <phool/phool.h> // for PHWHERE
0027
0028 #include <TFile.h>
0029 #include <TObjString.h>
0030
0031 #include <iostream>
0032 #include <numeric> // for std::accumulate
0033
0034 namespace G4Eval
0035 {
0036 std::vector<int> unmatchedSvtxTrkIds(EmbRecoMatchContainer* matches, SvtxTrackMap* m_SvtxTrackMap)
0037 {
0038 std::set<unsigned int> ids_matched{};
0039 for (auto trkid : matches->ids_RecoMatched())
0040 {
0041 ids_matched.insert(trkid);
0042 }
0043
0044 std::set<int> ids_unmatched;
0045 for (auto& reco : *m_SvtxTrackMap)
0046 {
0047 auto trkid = reco.first;
0048 if (ids_matched.count(trkid) == 0)
0049 {
0050 ids_unmatched.insert(trkid);
0051 }
0052 }
0053 std::vector<int> ids_vec;
0054 ids_vec.reserve(ids_unmatched.size());
0055 for (auto id : ids_unmatched)
0056 {
0057 ids_vec.push_back((int) id);
0058 }
0059 std::sort(ids_vec.begin(), ids_vec.end());
0060 return ids_vec;
0061 }
0062
0063
0064
0065
0066 void write_StringToTFile(const std::string& msg_name, const std::string& msg)
0067 {
0068
0069
0070 TFile* s_current = gDirectory->GetFile();
0071 if (s_current == nullptr)
0072 {
0073 std::cout << PHWHERE << " Error no TFile open to which to wrote the "
0074 << std::endl
0075 << " TObjString mesaged." << std::endl;
0076 return;
0077 }
0078 TObjString obj(msg.c_str());
0079 s_current->WriteObject(&obj, msg_name.c_str());
0080 return;
0081 }
0082
0083
0084 int trklayer_0123(TrkrDefs::hitsetkey key)
0085 {
0086 auto layer = TrkrDefs::getLayer(key);
0087 if (layer < 3)
0088 {
0089 return 0;
0090 }
0091 if (layer < 7)
0092 {
0093 return 1;
0094 }
0095 if (layer < 55)
0096 {
0097 return 2;
0098 }
0099 return 3;
0100 }
0101
0102
0103 TrkrClusterComparer::TrkrClusterComparer(float _nphi_widths, float _nz_widths)
0104 : m_nphi_widths{_nphi_widths}
0105 , m_nz_widths{_nz_widths} {};
0106 int TrkrClusterComparer::init(PHCompositeNode* topNode,
0107 const std::string& name_phg4_clusters,
0108 const std::string& name_reco_clusters)
0109 {
0110
0111
0112 PHG4CylinderGeomContainer* geom_container_mvtx = findNode::getClass<PHG4CylinderGeomContainer>(topNode, "CYLINDERGEOM_MVTX");
0113 if (!geom_container_mvtx)
0114 {
0115 std::cout << PHWHERE << " Could not locate CYLINDERGEOM_MVTX " << std::endl;
0116 return Fun4AllReturnCodes::ABORTRUN;
0117 }
0118 for (int this_layer = 0; this_layer < 3; ++this_layer)
0119 {
0120 auto layergeom = dynamic_cast<CylinderGeom_Mvtx*>(geom_container_mvtx->GetLayerGeom(this_layer));
0121 const double pitch = layergeom->get_pixel_x();
0122 const double length = layergeom->get_pixel_z();
0123 m_phistep[this_layer] = pitch;
0124 if (this_layer == 0)
0125 {
0126 m_zstep_mvtx = length;
0127 }
0128 }
0129
0130
0131 PHG4CylinderGeomContainer* geom_container_intt = findNode::getClass<PHG4CylinderGeomContainer>(topNode, "CYLINDERGEOM_INTT");
0132 if (!geom_container_intt)
0133 {
0134 std::cout << PHWHERE << " Could not locate CYLINDERGEOM_INTT " << std::endl;
0135 return Fun4AllReturnCodes::ABORTRUN;
0136 }
0137
0138 for (int this_layer = 3; this_layer < 7; ++this_layer)
0139 {
0140 CylinderGeomIntt* geom =
0141 dynamic_cast<CylinderGeomIntt*>(geom_container_intt->GetLayerGeom(this_layer));
0142 float pitch = geom->get_strip_y_spacing();
0143 m_phistep[this_layer] = pitch;
0144 }
0145
0146
0147 auto geom_tpc =
0148 findNode::getClass<PHG4TpcCylinderGeomContainer>(topNode, "CYLINDERCELLGEOM_SVTX");
0149 if (!geom_tpc)
0150 {
0151 std::cout << PHWHERE << " Could not locate CYLINDERCELLGEOM_SVTX node " << std::endl;
0152 return Fun4AllReturnCodes::ABORTRUN;
0153 }
0154 for (int this_layer = 7; this_layer < 55; ++this_layer)
0155 {
0156 PHG4TpcCylinderGeom* layergeom = geom_tpc->GetLayerCellGeom(this_layer);
0157 if (this_layer == 7)
0158 {
0159 m_zstep_tpc = layergeom->get_zstep();
0160 }
0161 m_phistep[this_layer] = layergeom->get_phistep();
0162 }
0163
0164 m_TruthClusters =
0165 findNode::getClass<TrkrClusterContainer>(topNode, name_phg4_clusters.c_str());
0166 if (!m_TruthClusters)
0167 {
0168 std::cout << PHWHERE << " Could not locate " << name_phg4_clusters << " node" << std::endl;
0169 return Fun4AllReturnCodes::ABORTRUN;
0170 }
0171
0172 m_RecoClusters =
0173 findNode::getClass<TrkrClusterContainer>(topNode, name_reco_clusters.c_str());
0174 if (!m_TruthClusters)
0175 {
0176 std::cout << PHWHERE << " Could not locate " << name_reco_clusters << " node" << std::endl;
0177 return Fun4AllReturnCodes::ABORTRUN;
0178 }
0179
0180 m_ActsGeometry = findNode::getClass<ActsGeometry>(topNode, "ActsGeometry");
0181 if (!m_ActsGeometry)
0182 {
0183 std::cout << PHWHERE << " Could not locate ActsGeometry node" << std::endl;
0184 return Fun4AllReturnCodes::ABORTRUN;
0185 }
0186
0187 return Fun4AllReturnCodes::EVENT_OK;
0188 }
0189
0190
0191
0192
0193
0194
0195
0196 std::pair<bool, float> TrkrClusterComparer::operator()(TrkrDefs::cluskey key_T, TrkrDefs::cluskey key_R)
0197 {
0198
0199 layer = TrkrDefs::getLayer(key_T);
0200 if (layer > 55)
0201 {
0202 std::cout << " Error! Trying to compar cluster in layer > 55, "
0203 << "which is not programmed yet!" << std::endl;
0204 return {false, 0.};
0205 }
0206
0207 in_mvtx = (layer < 3);
0208 in_intt = (layer > 2 && layer < 7);
0209 in_tpc = (layer > 6 && layer < 55);
0210
0211 float phi_step = m_phistep[layer];
0212 float z_step = in_mvtx ? m_zstep_mvtx : m_zstep_tpc;
0213
0214 clus_T = m_TruthClusters->findCluster(key_T);
0215 clus_R = m_RecoClusters->findCluster(key_R);
0216
0217 phi_T = clus_T->getPosition(0);
0218 phi_R = clus_R->getPosition(0);
0219 phisize_R = clus_R->getPhiSize() * m_nphi_widths;
0220 phisize_T = clus_T->getPhiSize() * m_nphi_widths;
0221
0222 z_T = clus_T->getPosition(1);
0223 z_R = clus_R->getPosition(1);
0224
0225 if (!in_intt)
0226 {
0227 zsize_R = clus_R->getZSize() * m_nz_widths;
0228 zsize_T = clus_R->getZSize() * m_nz_widths;
0229 }
0230
0231 phi_delta = fabs(phi_T - phi_R);
0232 while (phi_delta > M_PI)
0233 {
0234 phi_delta = fabs(phi_delta - 2 * M_PI);
0235 }
0236 phi_delta /= phi_step;
0237
0238 z_delta = fabs(z_T - z_R) / z_step;
0239
0240
0241
0242 float phi_stat = std::max(phisize_T, phisize_R);
0243 float fit_statistic = (phi_delta / phi_stat);
0244 is_match = (fit_statistic <= 1.);
0245
0246 float z_stat = 0;
0247 if (!in_intt)
0248 {
0249 z_stat = std::max(zsize_T, zsize_R);
0250
0251 is_match = is_match && (z_delta < z_stat);
0252 fit_statistic += z_delta / z_stat;
0253 }
0254 return {is_match, fit_statistic};
0255 }
0256
0257 ClusLoc TrkrClusterComparer::clusloc_PHG4(
0258 std::pair<TrkrDefs::hitsetkey, TrkrDefs::cluskey> input)
0259 {
0260 auto cluster = m_TruthClusters->findCluster(input.second);
0261 Eigen::Vector3d gloc = m_ActsGeometry->getGlobalPosition(input.second, cluster);
0262 return {TrkrDefs::getLayer(input.first), gloc,
0263 (int) cluster->getPhiSize(), (int) cluster->getZSize()};
0264 }
0265
0266 ClusLoc TrkrClusterComparer::clusloc_SVTX(
0267 std::pair<TrkrDefs::hitsetkey, TrkrDefs::cluskey> input)
0268 {
0269 auto cluster = m_RecoClusters->findCluster(input.second);
0270 Eigen::Vector3d gloc = m_ActsGeometry->getGlobalPosition(input.second, cluster);
0271 return {TrkrDefs::getLayer(input.first), gloc,
0272 (int) cluster->getPhiSize(), (int) cluster->getZSize()};
0273 }
0274
0275
0276
0277
0278
0279
0280 ClusKeyIter::ClusKeyIter(SvtxTrack* _track)
0281 : track{_track}
0282 , in_silicon{_track->get_silicon_seed() != nullptr}
0283 , has_tpc{_track->get_tpc_seed() != nullptr}
0284 , no_data{!in_silicon && !has_tpc}
0285 {
0286 }
0287
0288 ClusKeyIter ClusKeyIter::begin()
0289 {
0290 ClusKeyIter iter0{track};
0291 if (iter0.no_data)
0292 {
0293 return iter0;
0294 }
0295 if (iter0.in_silicon)
0296 {
0297 iter0.iter = track->get_silicon_seed()->begin_cluster_keys();
0298 iter0.iter_end_silicon = track->get_silicon_seed()->end_cluster_keys();
0299 }
0300 else if (has_tpc)
0301 {
0302 iter0.iter = track->get_tpc_seed()->begin_cluster_keys();
0303 }
0304 return iter0;
0305 }
0306
0307 ClusKeyIter ClusKeyIter::end()
0308 {
0309 ClusKeyIter iter0{track};
0310 if (iter0.no_data)
0311 {
0312 return iter0;
0313 }
0314 if (has_tpc)
0315 {
0316 iter0.iter = track->get_tpc_seed()->end_cluster_keys();
0317 }
0318 else if (in_silicon)
0319 {
0320 iter0.iter = track->get_silicon_seed()->end_cluster_keys();
0321 }
0322 return iter0;
0323 }
0324
0325 void ClusKeyIter::operator++()
0326 {
0327 if (no_data)
0328 {
0329 return;
0330 }
0331 ++iter;
0332 if (in_silicon && has_tpc && iter == iter_end_silicon)
0333 {
0334 in_silicon = false;
0335 iter = track->get_tpc_seed()->begin_cluster_keys();
0336 }
0337 }
0338
0339 bool ClusKeyIter::operator!=(const ClusKeyIter& rhs)
0340 {
0341 if (no_data)
0342 {
0343 return false;
0344 }
0345 return iter != rhs.iter;
0346 }
0347
0348 TrkrDefs::cluskey ClusKeyIter::operator*()
0349 {
0350 return *iter;
0351 }
0352
0353 TrkrClusterContainer* ClusCntr::get_PHG4_clusters()
0354 {
0355 if (comp == nullptr)
0356 {
0357 return nullptr;
0358 }
0359 else
0360 {
0361 return comp->m_TruthClusters;
0362 }
0363 }
0364
0365 TrkrClusterContainer* ClusCntr::get_SVTX_clusters()
0366 {
0367 if (comp == nullptr)
0368 {
0369 return nullptr;
0370 }
0371 else
0372 {
0373 return comp->m_RecoClusters;
0374 }
0375 }
0376
0377 std::array<int, 5> ClusCntr::cntclus(Vector& keys)
0378 {
0379 std::array<int, 5> cnt{0, 0, 0, 0, 0};
0380 for (auto& it : keys)
0381 {
0382 cnt[trklayer_0123(it.first)] += 1;
0383 }
0384 for (int i = 0; i < 4; ++i)
0385 {
0386 cnt[4] += cnt[i];
0387 }
0388 return cnt;
0389 }
0390
0391 std::array<int, 5> ClusCntr::cnt_matchedclus(Vector& keys, std::vector<bool>& matches)
0392 {
0393 std::array<int, 5> cnt{0, 0, 0, 0, 0};
0394 if (keys.size() != matches.size())
0395 {
0396 std::cout << PHWHERE << " matching and key vector not the same size. "
0397 << std::endl
0398 << " run find_matches() first." << std::endl;
0399 return cnt;
0400 }
0401 for (unsigned int i = 0; i < keys.size(); ++i)
0402 {
0403 if (matches[i])
0404 {
0405 cnt[trklayer_0123(keys[i].first)] += 1;
0406 }
0407 }
0408 for (int i = 0; i < 4; ++i)
0409 {
0410 cnt[4] += cnt[i];
0411 }
0412 return cnt;
0413 }
0414
0415 int ClusCntr::addClusKeys(SvtxTrack* track)
0416 {
0417 svtx_keys.clear();
0418 for (auto ckey : ClusKeyIter(track))
0419 {
0420 svtx_keys.push_back({TrkrDefs::getHitSetKeyFromClusKey(ckey), ckey});
0421 }
0422 std::sort(svtx_keys.begin(), svtx_keys.end());
0423 return svtx_keys.size();
0424 }
0425
0426 int ClusCntr::addClusKeys(TrkrTruthTrack* track)
0427 {
0428 phg4_keys.clear();
0429 for (auto ckey : track->getClusters())
0430 {
0431 phg4_keys.push_back({TrkrDefs::getHitSetKeyFromClusKey(ckey), ckey});
0432 }
0433 std::sort(phg4_keys.begin(), phg4_keys.end());
0434 return phg4_keys.size();
0435 }
0436
0437 void ClusCntr::reset()
0438 {
0439 phg4_keys.clear();
0440 phg4_matches.clear();
0441 svtx_keys.clear();
0442 svtx_matches.clear();
0443 }
0444
0445 std::array<int, 3> ClusCntr::find_matches()
0446 {
0447 if (comp == nullptr)
0448 {
0449 std::cout << PHWHERE
0450 << " Won't compare tracks because of missing TrkrClusterComparer" << std::endl;
0451 return {0, 0, 0};
0452 }
0453
0454
0455
0456
0457
0458
0459 auto& vA = phg4_keys;
0460 auto& vB = svtx_keys;
0461
0462 auto& matchesA = phg4_matches;
0463 auto& matchesB = svtx_matches;
0464
0465 match_stat = 0.;
0466
0467
0468 matchesA = std::vector<bool>(vA.size(), false);
0469 matchesB = std::vector<bool>(vB.size(), false);
0470
0471
0472 auto iA0 = vA.begin();
0473 auto iA1 = vA.end();
0474
0475 auto iB0 = vB.begin();
0476 auto iB1 = vB.end();
0477
0478 auto iA = iA0;
0479 auto iB = iB0;
0480
0481 int n_match{0};
0482
0483 while (iA != iA1 && iB != iB1)
0484 {
0485 if (iA->first == iB->first)
0486 {
0487 auto hitset = iA->first;
0488
0489
0490 auto sAend = iA + 1;
0491 while (sAend != iA1 && sAend->first == hitset)
0492 {
0493 ++sAend;
0494 }
0495
0496 auto sBend = iB + 1;
0497 while (sBend != iB1 && sBend->first == hitset)
0498 {
0499 ++sBend;
0500 }
0501
0502 for (auto A = iA; A != sAend; ++A)
0503 {
0504 for (auto B = iB; B != sBend; ++B)
0505 {
0506 auto comp_val = comp->operator()(A->second, B->second);
0507 if (comp_val.first)
0508 {
0509 matchesA[A - iA0] = true;
0510 matchesB[B - iB0] = true;
0511 match_stat += comp_val.second;
0512 ++n_match;
0513 }
0514 }
0515 }
0516 iA = sAend;
0517 iB = sBend;
0518 }
0519 else if (iA->first < iB->first)
0520 {
0521 ++iA;
0522 }
0523 else
0524 {
0525 ++iB;
0526 }
0527 }
0528 return {n_match, (int) phg4_keys.size(), (int) svtx_keys.size()};
0529 }
0530
0531 std::array<int, 3> ClusCntr::find_matches(TrkrTruthTrack* g4_track, SvtxTrack* sv_track)
0532 {
0533 addClusKeys(sv_track);
0534 addClusKeys(g4_track);
0535 return find_matches();
0536 }
0537
0538 int ClusCntr::phg4_n_matched()
0539 {
0540 return std::accumulate(phg4_matches.begin(), phg4_matches.end(), 0);
0541 }
0542
0543 int ClusCntr::svtx_n_matched()
0544 {
0545 return std::accumulate(svtx_matches.begin(), svtx_matches.end(), 0);
0546 }
0547
0548 std::vector<ClusLoc> ClusCntr::phg4_clusloc_all()
0549 {
0550 std::vector<ClusLoc> vec{};
0551 for (auto& cluspair : phg4_keys)
0552 {
0553 vec.push_back(comp->clusloc_PHG4(cluspair));
0554 }
0555 return vec;
0556 }
0557
0558 std::vector<ClusLoc> ClusCntr::phg4_clusloc_unmatched()
0559 {
0560 std::vector<ClusLoc> vec{};
0561 auto cnt = phg4_keys.size();
0562 for (unsigned int i = 0; i < cnt; ++i)
0563 {
0564 if (!phg4_matches[i])
0565 {
0566 vec.push_back(comp->clusloc_PHG4(phg4_keys[i]));
0567 }
0568 }
0569 return vec;
0570 }
0571
0572 std::vector<ClusLoc> ClusCntr::svtx_clusloc_all()
0573 {
0574 std::vector<ClusLoc> vec{};
0575 for (auto& cluspair : svtx_keys)
0576 {
0577 vec.push_back(comp->clusloc_SVTX(cluspair));
0578 }
0579 return vec;
0580 }
0581
0582 std::vector<ClusLoc> ClusCntr::svtx_clusloc_unmatched()
0583 {
0584 std::vector<ClusLoc> vec{};
0585 auto cnt = svtx_keys.size();
0586 for (unsigned int i = 0; i < cnt; ++i)
0587 {
0588 if (!svtx_matches[i])
0589 {
0590 vec.push_back(comp->clusloc_SVTX(svtx_keys[i]));
0591 }
0592 }
0593 return vec;
0594 }
0595
0596 std::vector<ClusLoc> ClusCntr::clusloc_matched()
0597 {
0598 std::vector<ClusLoc> vec{};
0599 auto cnt = phg4_keys.size();
0600 for (unsigned int i = 0; i < cnt; ++i)
0601 {
0602 if (phg4_matches[i])
0603 {
0604 vec.push_back(comp->clusloc_PHG4(phg4_keys[i]));
0605 }
0606 }
0607 return vec;
0608 }
0609
0610
0611
0612
0613
0614
0615
0616
0617
0618 }