Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-12-16 09:21:44

0001 #include "FillTruthRecoMatchTree.h"
0002 
0003 #include <g4tracking/EmbRecoMatch.h>
0004 #include <g4tracking/EmbRecoMatchContainer.h>
0005 #include <g4tracking/TrkrTruthTrack.h>
0006 #include <g4tracking/TrkrTruthTrackContainer.h>
0007 
0008 #include <trackbase/TrkrCluster.h>
0009 #include <trackbase/TrkrClusterContainer.h>
0010 #include <trackbase/TrkrDefs.h>
0011 
0012 #include <trackbase_historic/SvtxPHG4ParticleMap_v1.h>
0013 #include <trackbase_historic/SvtxTrack.h>
0014 #include <trackbase_historic/SvtxTrackMap.h>
0015 
0016 #include <g4main/PHG4TruthInfoContainer.h>
0017 
0018 #include <fun4all/Fun4AllReturnCodes.h>
0019 #include <fun4all/PHTFileServer.h>
0020 
0021 #include <phool/PHCompositeNode.h>
0022 #include <phool/PHDataNode.h>
0023 #include <phool/PHIODataNode.h>
0024 #include <phool/PHNode.h>
0025 #include <phool/PHNodeIterator.h>
0026 #include <phool/PHObject.h>  // for PHObject
0027 #include <phool/getClass.h>
0028 #include <phool/phool.h>  // for PHWHERE
0029 
0030 #include <TFile.h>
0031 #include <TH2.h>
0032 #include <TTree.h>
0033 
0034 #include <format>
0035 #include <iostream>
0036 
0037 //____________________________________________________________________________..
0038 FillTruthRecoMatchTree::FillTruthRecoMatchTree(
0039     bool _fill_clusters, bool _fill_SvUnMatched, float _cluster_nzwidths, float _cluster_nphiwidths, const std::string& _outfile_name)
0040   : m_cluster_comp{_cluster_nphiwidths, _cluster_nzwidths}
0041   , m_fill_clusters{_fill_clusters}
0042   , m_fill_SvU{_fill_SvUnMatched}
0043   , m_outfile_name{_outfile_name}
0044   , h2_G4_nPixelsPhi(new TH2D("G4_nPixelsPhi", "PHG4 Emb Tracks; cluster pixel width Phi; layer",
0045                               100, -0.5, 99.5, 56, -0.5, 55.5))
0046   , h2_G4_nPixelsZ(new TH2D("G4_nPixelsZ", "PHG4 Emb Tracks; cluster pixel width Z; layer",
0047                             100, -0.5, 99.5, 56, -0.5, 55.5))
0048   , h2_Sv_nPixelsPhi(new TH2D("Sv_nPixelsPhi", "Svtx Reco Tracks; cluster pixel width Phi; layer",
0049                               100, -0.5, 99.5, 56, -0.5, 55.5))
0050   , h2_Sv_nPixelsZ(new TH2D("Sv_nPixelsZ", "Svtx Reco Tracks; cluster pixel width Z; layer",
0051                             100, -0.5, 99.5, 56, -0.5, 55.5))
0052 {
0053   m_cluscntr.set_comparer(&m_cluster_comp);
0054   PHTFileServer::open(m_outfile_name, "RECREATE");
0055 
0056   m_ttree = new TTree("T", "Tracks (and sometimes clusters)");
0057 
0058   m_ttree->Branch("event", &nevent);
0059   m_ttree->Branch("nphg4_part", &nphg4_part);
0060   m_ttree->Branch("centrality", &centrality);
0061   m_ttree->Branch("ntrackmatches", &ntrackmatches);
0062   m_ttree->Branch("nphg4", &nphg4);
0063   m_ttree->Branch("nsvtx", &nsvtx);
0064 
0065   m_ttree->Branch("trackid", &b_trackid);
0066   m_ttree->Branch("is_G4track", &b_is_g4track);
0067   m_ttree->Branch("is_Svtrack", &b_is_Svtrack);
0068   m_ttree->Branch("is_matched", &b_is_matched);
0069 
0070   m_ttree->Branch("trkpt", &b_trkpt);
0071   m_ttree->Branch("trketa", &b_trketa);
0072   m_ttree->Branch("trkphi", &b_trkphi);
0073 
0074   m_ttree->Branch("nclus", &b_nclus);
0075   m_ttree->Branch("nclustpc", &b_nclustpc);
0076   m_ttree->Branch("nclusmvtx", &b_nclusmvtx);
0077   m_ttree->Branch("nclusintt", &b_nclusintt);
0078   m_ttree->Branch("matchrat", &b_matchrat);
0079   m_ttree->Branch("matchrat_intt", &b_matchrat_intt);
0080   m_ttree->Branch("matchrat_mvtx", &b_matchrat_mvtx);
0081   m_ttree->Branch("matchrat_tpc", &b_matchrat_tpc);
0082 
0083   if (m_fill_clusters)
0084   {
0085     m_ttree->Branch("clus_match", &b_clusmatch);
0086     m_ttree->Branch("clus_x", &b_clus_x);
0087     m_ttree->Branch("clus_y", &b_clus_y);
0088     m_ttree->Branch("clus_z", &b_clus_z);
0089     m_ttree->Branch("clus_r", &b_clus_r);
0090     m_ttree->Branch("clus_layer", &b_clus_layer);
0091     m_ttree->Branch("nphibins", &b_clus_nphibins);
0092     m_ttree->Branch("nzbins", &b_clus_ntbins);
0093   }
0094 }
0095 
0096 //____________________________________________________________________________..
0097 int FillTruthRecoMatchTree::Init(PHCompositeNode* topNode)
0098 {
0099   if (Verbosity() > 1)
0100   {
0101     std::cout << " Beginning FillTruthRecoMatchTree " << std::endl;
0102     topNode->print();
0103   }
0104 
0105   return Fun4AllReturnCodes::EVENT_OK;
0106 }
0107 
0108 //____________________________________________________________________________..
0109 int FillTruthRecoMatchTree::InitRun(PHCompositeNode* topNode)
0110 {
0111   auto init_status = m_cluster_comp.init(topNode);
0112   if (init_status == Fun4AllReturnCodes::ABORTRUN)
0113   {
0114     return init_status;
0115   }
0116 
0117   if (createNodes(topNode) != Fun4AllReturnCodes::EVENT_OK)
0118   {
0119     return Fun4AllReturnCodes::ABORTEVENT;
0120   }
0121 
0122   return Fun4AllReturnCodes::EVENT_OK;
0123 }
0124 
0125 int FillTruthRecoMatchTree::createNodes(PHCompositeNode* topNode)
0126 {
0127   PHNodeIterator iter(topNode);
0128 
0129   PHCompositeNode* dstNode = dynamic_cast<PHCompositeNode*>(iter.findFirst("PHCompositeNode", "DST"));
0130 
0131   if (!dstNode)
0132   {
0133     std::cout << PHWHERE << " DST node is missing, quitting" << std::endl;
0134     std::cerr << PHWHERE << " DST node is missing, quitting" << std::endl;
0135     throw std::runtime_error("Failed to find DST node in FillTruthRecoMatchTree::createNodes");
0136   }
0137 
0138   m_EmbRecoMatchContainer = findNode::getClass<EmbRecoMatchContainer>(topNode, "TRKR_EMBRECOMATCHCONTAINER");
0139   if (!m_EmbRecoMatchContainer)
0140   {
0141     std::cout << PHWHERE << " Cannot find node TRKR_EMBRECOMATCHCONTAINER on node tree; quitting " << std::endl;
0142     std::cerr << PHWHERE << " Cannot find node TRKR_EMBRECOMATCHCONTAINER on node tree; quitting " << std::endl;
0143     throw std::runtime_error(" Cannot find node TRKR_EMBRECOMATCHCONTAINER on node tree; quitting");
0144   }
0145 
0146   PHCompositeNode* svtxNode = dynamic_cast<PHCompositeNode*>(iter.findFirst("PHCompositeNode", "SVTX"));
0147   if (!svtxNode)
0148   {
0149     svtxNode = new PHCompositeNode("SVTX");
0150     dstNode->addNode(svtxNode);
0151   }
0152 
0153   m_PHG4TruthInfoContainer = findNode::getClass<PHG4TruthInfoContainer>(topNode, "G4TruthInfo");
0154   if (!m_PHG4TruthInfoContainer)
0155   {
0156     std::cout << "Could not locate G4TruthInfo node when running "
0157               << "\"TruthRecoTrackMatching\" module." << std::endl;
0158     return Fun4AllReturnCodes::ABORTEVENT;
0159   }
0160 
0161   m_SvtxTrackMap = findNode::getClass<SvtxTrackMap>(topNode, "SvtxTrackMap");
0162   if (!m_SvtxTrackMap)
0163   {
0164     std::cout << "Could not locate SvtxTrackMap node when running "
0165               << "\"TruthRecoTrackMatching\" module." << std::endl;
0166     return Fun4AllReturnCodes::ABORTEVENT;
0167   }
0168 
0169   /* m_TruthClusterContainer = findNode::getClass<TrkrClusterContainer>(topNode, */
0170   /*     "TRKR_TRUTHCLUSTERCONTAINER"); */
0171   /* if (!m_TruthClusterContainer) */
0172   /* { */
0173   /*   std::cout << "Could not locate TRKR_TRUTHCLUSTERCONTAINER node when running " */
0174   /*     << "\"TruthRecoTrackMatching\" module." << std::endl; */
0175   /*   return Fun4AllReturnCodes::ABORTEVENT; */
0176   /* } */
0177 
0178   /* m_RecoClusterContainer = findNode::getClass<TrkrClusterContainer>(topNode, "TRKR_CLUSTER"); */
0179   /* if (!m_RecoClusterContainer) */
0180   /* { */
0181   /*   std::cout << "Could not locate TRKR_CLUSTER node when running " */
0182   /*     << "\"TruthRecoTrackMatching\" module." << std::endl; */
0183   /*   return Fun4AllReturnCodes::ABORTEVENT; */
0184   /* } */
0185 
0186   m_TrkrTruthTrackContainer = findNode::getClass<TrkrTruthTrackContainer>(topNode,
0187                                                                           "TRKR_TRUTHTRACKCONTAINER");
0188   if (!m_TrkrTruthTrackContainer)
0189   {
0190     std::cout << "Could not locate TRKR_TRUTHTRACKCONTAINER node when running "
0191               << "\"TruthRecoTrackMatching\" module." << std::endl;
0192     return Fun4AllReturnCodes::ABORTEVENT;
0193   }
0194 
0195   return Fun4AllReturnCodes::EVENT_OK;
0196 }
0197 
0198 int FillTruthRecoMatchTree::process_event(PHCompositeNode* /*topNode*/)
0199 {
0200   if (Verbosity() > 5)
0201   {
0202     std::cout << " FillTruthRecoMatchTree::process_event() " << std::endl;
0203   }
0204 
0205   // fill in the event data
0206   ++nevent;
0207   nphg4 = m_TrkrTruthTrackContainer->getMap().size();
0208   nsvtx = m_SvtxTrackMap->size();
0209   ntrackmatches = m_EmbRecoMatchContainer->getMatches().size();
0210   // get centrality later...
0211 
0212   // fill in pixel widths on truth tracks
0213   for (auto hitsetkey : m_cluscntr.get_PHG4_clusters()->getHitSetKeys())
0214   {
0215     float layer = (float) TrkrDefs::getLayer(hitsetkey);
0216     auto range = m_cluscntr.get_PHG4_clusters()->getClusters(hitsetkey);
0217     for (auto iter = range.first; iter != range.second; ++iter)
0218     {
0219       const auto& cluster = iter->second;
0220       h2_G4_nPixelsPhi->Fill(cluster->getPhiSize(), layer);
0221       h2_G4_nPixelsZ->Fill(cluster->getZSize(), layer);
0222     }
0223   }
0224   // fill in pixel widths on reco tracks
0225   for (auto hitsetkey : m_cluscntr.get_SVTX_clusters()->getHitSetKeys())
0226   {
0227     float layer = (float) TrkrDefs::getLayer(hitsetkey);
0228     auto range = m_cluscntr.get_SVTX_clusters()->getClusters(hitsetkey);
0229     for (auto iter = range.first; iter != range.second; ++iter)
0230     {
0231       const auto& cluster = iter->second;
0232       h2_Sv_nPixelsPhi->Fill(cluster->getPhiSize(), layer);
0233       h2_Sv_nPixelsZ->Fill(cluster->getZSize(), layer);
0234     }
0235   }
0236 
0237   nphg4_part = 0;
0238   const auto range = m_PHG4TruthInfoContainer->GetPrimaryParticleRange();
0239   for (PHG4TruthInfoContainer::ConstIterator iter = range.first;
0240        iter != range.second; ++iter)
0241   {
0242     nphg4_part++;
0243   }
0244 
0245   // unmatched tracks are only entered once
0246   // matches can repeat a given svtx or phg4 track, depending on the
0247   // parameters in teh matching in filltruthrecomatchtree
0248   //
0249   // (1) fill unmatched phg4
0250   // (2) fill unmatched svtx
0251   // (3) fill matched phg4 and svtx
0252   clear_clusvecs(" nothing ");
0253 
0254   if (Verbosity() > 2)
0255   {
0256     std::cout << " getting" << (int) m_EmbRecoMatchContainer->getMatches().size() << std::endl;
0257   }
0258   for (auto& match : m_EmbRecoMatchContainer->getMatches())
0259   {
0260     unsigned int g4_trkid = match->idTruthTrack();
0261     int sv_trkid = match->idRecoTrack();
0262 
0263     auto* g4trk = m_TrkrTruthTrackContainer->getTruthTrack(g4_trkid);
0264     auto* svtrk = m_SvtxTrackMap->get(sv_trkid);
0265 
0266     m_cluscntr.addClusKeys(g4trk);
0267     m_cluscntr.addClusKeys(svtrk);
0268     m_cluscntr.find_matches();
0269 
0270     // <- <- <- <- G4 Matched Tracks
0271     b_is_matched = true;
0272     b_is_g4track = true;
0273     b_is_Svtrack = false;
0274 
0275     b_trackid = g4_trkid;
0276     b_trkpt = g4trk->getPt();
0277     b_trketa = g4trk->getPseudoRapidity();
0278     b_trkphi = g4trk->getPhi();
0279 
0280     auto cnt = m_cluscntr.phg4_cntclus();
0281     auto cnt_match = m_cluscntr.phg4_cnt_matchedclus();
0282 
0283     b_nclus = cnt[4];
0284     b_nclusmvtx = cnt[0];
0285     b_nclusintt = cnt[1];
0286     b_nclustpc = cnt[2];
0287 
0288     b_matchrat = (float) cnt_match[4] / cnt[4];
0289     b_matchrat_mvtx = (float) cnt_match[0] / cnt[0];
0290     b_matchrat_intt = (float) cnt_match[1] / cnt[1];
0291     b_matchrat_tpc = (float) cnt_match[2] / cnt[2];
0292 
0293     if (m_fill_clusters)
0294     {
0295       auto clusters = m_cluscntr.phg4_clusloc_unmatched();
0296       for (auto& loc : clusters)
0297       {
0298         b_clusmatch.push_back(false);
0299         b_clus_layer.push_back(std::get<0>(loc));
0300         auto x = std::get<1>(loc)[0];
0301         auto y = std::get<1>(loc)[1];
0302         b_clus_x.push_back(x);
0303         b_clus_y.push_back(y);
0304         b_clus_z.push_back(std::get<1>(loc)[2]);
0305         b_clus_r.push_back(pow(x * x + y * y, 0.5));
0306         b_clus_nphibins.push_back(std::get<2>(loc));
0307         b_clus_ntbins.push_back(std::get<3>(loc));
0308       }
0309 
0310       clusters = m_cluscntr.clusloc_matched();
0311       for (auto& loc : clusters)
0312       {
0313         b_clusmatch.push_back(true);
0314         b_clus_layer.push_back(std::get<0>(loc));
0315         auto x = std::get<1>(loc)[0];
0316         auto y = std::get<1>(loc)[1];
0317         b_clus_x.push_back(x);
0318         b_clus_y.push_back(y);
0319         b_clus_z.push_back(std::get<1>(loc)[2]);
0320         b_clus_r.push_back(pow(x * x + y * y, 0.5));
0321         b_clus_nphibins.push_back(std::get<2>(loc));
0322         b_clus_ntbins.push_back(std::get<3>(loc));
0323       }
0324     }
0325     m_ttree->Fill();
0326     clear_clusvecs();
0327     /* clear_clusvecs("apple0 g4_matched"); */
0328 
0329     // <- <- <- <- Svtx Matched Tracks
0330     b_is_g4track = false;
0331     b_is_Svtrack = true;
0332     b_trackid = sv_trkid;
0333     b_trkpt = svtrk->get_pt();
0334     b_trketa = svtrk->get_eta();
0335     b_trkphi = svtrk->get_phi();
0336 
0337     cnt = m_cluscntr.svtx_cntclus();
0338     b_nclus = cnt[4];
0339     b_nclusmvtx = cnt[0];
0340     b_nclusintt = cnt[1];
0341     b_nclustpc = cnt[2];
0342 
0343     b_matchrat = (float) cnt_match[4] / cnt[4];
0344     b_matchrat_mvtx = (float) cnt_match[0] / cnt[0];
0345     b_matchrat_intt = (float) cnt_match[1] / cnt[1];
0346     b_matchrat_tpc = (float) cnt_match[2] / cnt[2];
0347 
0348     /* int _ = 0; */
0349     if (m_fill_clusters)
0350     {
0351       auto clusters = m_cluscntr.svtx_clusloc_unmatched();
0352       for (auto& loc : clusters)
0353       {
0354         b_clusmatch.push_back(false);
0355         b_clus_layer.push_back(std::get<0>(loc));
0356         auto x = std::get<1>(loc)[0];
0357         auto y = std::get<1>(loc)[1];
0358         /* if (_==0) std::cout << " apple x: " << x << " y: " << y << std::endl; */
0359         /* _ += 1; */
0360         b_clus_x.push_back(x);
0361         b_clus_y.push_back(y);
0362         b_clus_z.push_back(std::get<1>(loc)[2]);
0363         b_clus_r.push_back(pow(x * x + y * y, 0.5));
0364         b_clus_nphibins.push_back(std::get<2>(loc));
0365         b_clus_ntbins.push_back(std::get<3>(loc));
0366       }
0367 
0368       clusters = m_cluscntr.clusloc_matched();
0369       for (auto& loc : clusters)
0370       {
0371         b_clusmatch.push_back(true);
0372         b_clus_layer.push_back(std::get<0>(loc));
0373         auto x = std::get<1>(loc)[0];
0374         auto y = std::get<1>(loc)[1];
0375         b_clus_x.push_back(x);
0376         b_clus_y.push_back(y);
0377         b_clus_z.push_back(std::get<1>(loc)[2]);
0378         b_clus_r.push_back(pow(x * x + y * y, 0.5));
0379         b_clus_nphibins.push_back(std::get<2>(loc));
0380         b_clus_ntbins.push_back(std::get<3>(loc));
0381       }
0382     }
0383     m_ttree->Fill();
0384     clear_clusvecs();
0385     /* clear_clusvecs("apple1 s4_matched"); */
0386   }
0387 
0388   // <- <- <- <- G4 un-matched Tracks
0389   b_is_matched = false;
0390   b_is_g4track = true;
0391   b_is_Svtrack = false;
0392   for (auto& g4_trkid : m_EmbRecoMatchContainer->ids_TruthUnmatched())
0393   {
0394     auto* g4trk = m_TrkrTruthTrackContainer->getTruthTrack(g4_trkid);
0395     m_cluscntr.addClusKeys(g4trk);
0396 
0397     b_trackid = g4_trkid;
0398     b_trkpt = g4trk->getPt();
0399     b_trketa = g4trk->getPseudoRapidity();
0400     b_trkphi = g4trk->getPhi();
0401 
0402     auto cnt = m_cluscntr.phg4_cntclus();
0403     b_nclus = cnt[4];
0404     b_nclusmvtx = cnt[0];
0405     b_nclusintt = cnt[1];
0406     b_nclustpc = cnt[2];
0407 
0408     if (m_fill_clusters)
0409     {
0410       auto clusters = m_cluscntr.phg4_clusloc_all();
0411       for (auto& loc : clusters)
0412       {
0413         b_clusmatch.push_back(false);
0414         b_clus_layer.push_back(std::get<0>(loc));
0415         auto x = std::get<1>(loc)[0];
0416         auto y = std::get<1>(loc)[1];
0417         b_clus_x.push_back(x);
0418         b_clus_y.push_back(y);
0419         b_clus_z.push_back(std::get<1>(loc)[2]);
0420         b_clus_r.push_back(pow(x * x + y * y, 0.5));
0421         b_clus_nphibins.push_back(std::get<2>(loc));
0422         b_clus_ntbins.push_back(std::get<3>(loc));
0423       }
0424       // this is an unmatched track, so there are no matched clusters
0425     }
0426     m_ttree->Fill();
0427     clear_clusvecs();
0428     /* clear_clusvecs("apple2 g4_unmatched"); */
0429   }
0430 
0431   // <- <- <- <- Svtx un-matched Tracks
0432   b_is_matched = false;
0433   b_is_matched = false;
0434   b_is_g4track = false;
0435   b_is_Svtrack = true;
0436 
0437   // just put in all svtx tracks, period...
0438   if (m_fill_SvU)
0439   {
0440     for (auto sv_trkid : G4Eval::unmatchedSvtxTrkIds(m_EmbRecoMatchContainer, m_SvtxTrackMap))
0441     {
0442       auto* svtrk = m_SvtxTrackMap->get(sv_trkid);
0443       m_cluscntr.addClusKeys(svtrk);
0444       b_trackid = sv_trkid;
0445       b_trkpt = svtrk->get_pt();
0446       b_trketa = svtrk->get_eta();
0447       b_trkphi = svtrk->get_phi();
0448 
0449       auto cnt = m_cluscntr.svtx_cntclus();
0450       b_nclus = cnt[4];
0451       b_nclusmvtx = cnt[0];
0452       b_nclusintt = cnt[1];
0453       b_nclustpc = cnt[2];
0454 
0455       if (m_fill_clusters)
0456       {
0457         auto clusters = m_cluscntr.svtx_clusloc_all();
0458         for (auto& loc : clusters)
0459         {
0460           b_clusmatch.push_back(false);
0461           b_clus_layer.push_back(std::get<0>(loc));
0462           auto x = std::get<1>(loc)[0];
0463           auto y = std::get<1>(loc)[1];
0464           b_clus_x.push_back(x);
0465           b_clus_y.push_back(y);
0466           b_clus_z.push_back(std::get<1>(loc)[2]);
0467           b_clus_r.push_back(pow(x * x + y * y, 0.5));
0468           b_clus_nphibins.push_back(std::get<2>(loc));
0469           b_clus_ntbins.push_back(std::get<3>(loc));
0470         }
0471       }
0472       m_ttree->Fill();
0473       clear_clusvecs();
0474     }
0475   }
0476 
0477   if (Verbosity() > 100)
0478   {
0479     print_mvtx_diagnostics();
0480   }
0481   return Fun4AllReturnCodes::EVENT_OK;
0482 }
0483 
0484 void FillTruthRecoMatchTree::print_mvtx_diagnostics()
0485 {
0486   std::cout << "To do: "
0487             << " (1)  number of truth tracks and total number of mvtx and ratio " << std::endl
0488             << " (2)  ditto for reco tracks " << std::endl;
0489 
0490   double n_PHG4_tracks = m_TrkrTruthTrackContainer->getMap().size();
0491   // count how many mvtx clusters in the phg4 tracks
0492   double n_in_PHG4_tracks{0.};
0493   for (auto& pair : m_TrkrTruthTrackContainer->getMap())
0494   {
0495     m_cluscntr.addClusKeys(pair.second);
0496     n_in_PHG4_tracks += m_cluscntr.phg4_cntclus()[0];
0497   }
0498   // count how mant mvtx clusters in truth container (should be the same)
0499   double n_in_PHG4_clusters{0.};
0500   for (auto hitsetkey : m_cluscntr.get_PHG4_clusters()->getHitSetKeys())
0501   {
0502     if (TrkrDefs::getLayer(hitsetkey) > 2)
0503     {
0504       continue;
0505     }
0506     auto range = m_cluscntr.get_PHG4_clusters()->getClusters(hitsetkey);
0507     for (auto r = range.first; r != range.second; ++r)
0508     {
0509       n_in_PHG4_clusters += 1.;
0510     }
0511   }
0512 
0513   // count how many svtx tracks
0514   double n_SVTX_tracks = m_SvtxTrackMap->size();
0515   // count how many mvtx clusters in svtx tracks
0516   double n_in_SVTX_tracks{0.};
0517   for (auto& entry : *m_SvtxTrackMap)
0518   {
0519     m_cluscntr.addClusKeys(entry.second);
0520     n_in_SVTX_tracks += m_cluscntr.svtx_cntclus()[0];
0521   }
0522   // count how many mvtx are total in the container
0523   double n_in_SVTX_clusters{0.};
0524   for (auto hitsetkey : m_cluscntr.get_SVTX_clusters()->getHitSetKeys())
0525   {
0526     if (TrkrDefs::getLayer(hitsetkey) > 2)
0527     {
0528       continue;
0529     }
0530     auto range = m_cluscntr.get_SVTX_clusters()->getClusters(hitsetkey);
0531     for (auto r = range.first; r != range.second; ++r)
0532     {
0533       n_in_SVTX_clusters += 1.;
0534     }
0535   }
0536 
0537   std::cout << std::format(
0538                    "MVTX"
0539                    "\nPHG4:  Tracks({:.0f})   Clusters In tracks({:.0f})   Total ({:.0f})"
0540                    "\n       ave. per track: {:6.3f}   ratio in all tracks: {:6.2f}",
0541                    n_PHG4_tracks, n_in_PHG4_tracks, n_in_PHG4_clusters, (n_in_PHG4_tracks / n_PHG4_tracks), (n_in_PHG4_tracks / n_in_PHG4_clusters))
0542             << std::endl;
0543   std::cout << std::format(
0544                    "\nSVTX:  Tracks({:.0f})   Clusters In tracks({:.0f})   Total ({:.0f})"
0545                    "\n       ave. per track: {:6.3f}   ratio in all tracks: {:6.2f}",
0546                    n_SVTX_tracks, n_in_SVTX_tracks, n_in_SVTX_clusters, (n_in_SVTX_tracks / n_SVTX_tracks), (n_in_SVTX_tracks / n_in_SVTX_clusters))
0547             << std::endl;
0548 }
0549 
0550 int FillTruthRecoMatchTree::End(PHCompositeNode* /*unused*/)
0551 {
0552   if (Verbosity() > 2)
0553   {
0554     std::cout << PHWHERE << ": ending FillTruthRecoMatchTree" << std::endl;
0555   }
0556   PHTFileServer::cd(m_outfile_name);
0557 
0558   h2_G4_nPixelsPhi->Write();
0559   h2_G4_nPixelsZ->Write();
0560   h2_Sv_nPixelsPhi->Write();
0561   h2_Sv_nPixelsZ->Write();
0562 
0563   m_ttree->Write();
0564   return Fun4AllReturnCodes::EVENT_OK;
0565 }
0566 
0567 void FillTruthRecoMatchTree::clear_clusvecs(const std::string& tag)
0568 {
0569   /* std::cout << " banana |" << tag << "|"<<std::endl; */
0570   if (!tag.empty())
0571   {
0572     std::cout << std::endl
0573               << " pear printing " << tag << " x(" << b_clus_x.size() << ") ";
0574     for (auto x : b_clus_x)
0575     {
0576       std::cout << x << " ";
0577     }
0578     std::cout << std::endl;
0579   }
0580   // Tracks and clustes
0581   b_clusmatch.clear();
0582   b_clus_x.clear();
0583   b_clus_y.clear();
0584   b_clus_z.clear();
0585   b_clus_r.clear();
0586   b_clus_layer.clear();
0587   b_clus_nphibins.clear();
0588   b_clus_ntbins.clear();
0589 }