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