Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-06 08:17:39

0001 #include "HFTrackEfficiency.h"
0002 
0003 #include <iomanip>
0004 
0005 //____________________________________________________________________________..
0006 HFTrackEfficiency::HFTrackEfficiency(const std::string &name)
0007   : SubsysReco(name)
0008   , m_triggerOnDecay(false)
0009   , m_input_track_map_node_name("SvtxTrackMap")
0010   , m_output_track_map_node_name("HFSelected")
0011   , m_write_track_map(false)
0012   , m_outfile_name("outputHFTrackEff.root")
0013   , m_outfile(nullptr)
0014   , m_tree(nullptr)
0015   , m_write_nTuple(false)
0016   , m_truthRecoMatchPercent(5.)
0017   , m_nDaughters(2)
0018 {
0019 }
0020 
0021 //____________________________________________________________________________..
0022 HFTrackEfficiency::~HFTrackEfficiency() = default;
0023 
0024 //____________________________________________________________________________..
0025 int HFTrackEfficiency::Init(PHCompositeNode *topNode)
0026 {
0027   m_truthRecoMatchPercent /= 100.;
0028 
0029   if (m_write_track_map)
0030   {
0031     PHNodeIterator nodeIter(topNode);
0032 
0033     PHCompositeNode *dstNode = dynamic_cast<PHCompositeNode *>(nodeIter.findFirst("PHCompositeNode", "DST"));
0034     if (!dstNode)
0035     {
0036       dstNode = new PHCompositeNode("DST");
0037       topNode->addNode(dstNode);
0038       std::cout << "DST node added" << std::endl;
0039     }
0040 
0041     outputNodeName = m_output_track_map_node_name + "_SvtxTrackMap";
0042     m_output_trackMap = new SvtxTrackMap_v2();
0043     PHIODataNode<PHObject> *outputTrackNode = new PHIODataNode<PHObject>(m_output_trackMap, outputNodeName.c_str(), "PHObject");
0044     dstNode->addNode(outputTrackNode);
0045     std::cout << outputNodeName << " node added" << std::endl;
0046   }
0047 
0048   return Fun4AllReturnCodes::EVENT_OK;
0049 }
0050 
0051 //____________________________________________________________________________..
0052 int HFTrackEfficiency::process_event(PHCompositeNode *topNode)
0053 {
0054   PHNodeIterator nodeIter(topNode);
0055   PHCompositeNode *dstNode = static_cast<PHCompositeNode *>(nodeIter.findFirst("PHCompositeNode", "DST"));
0056   assert(dstNode);
0057 
0058   std::string df_node_name = m_df_module_name + "_DecayMap";
0059   m_decayMap = findNode::getClass<DecayFinderContainer_v1>(topNode, df_node_name.c_str());
0060   if (!m_decayMap)
0061   {
0062     std::cout << __FILE__ << ": Missing node " << df_node_name << std::endl;
0063     return Fun4AllReturnCodes::ABORTRUN;
0064   }
0065 
0066   m_input_trackMap = findNode::getClass<SvtxTrackMap>(topNode, m_input_track_map_node_name.c_str());
0067   if (!m_input_trackMap)
0068   {
0069     std::cout << __FILE__ << ": Missing node " << m_input_track_map_node_name << std::endl;
0070     return Fun4AllReturnCodes::ABORTRUN;
0071   }
0072 
0073   m_geneventmap = findNode::getClass<PHHepMCGenEventMap>(topNode, "PHHepMCGenEventMap");
0074   if (!m_geneventmap)
0075   {
0076     std::cout << __FILE__ << ": Missing node PHHepMCGenEventMap" << std::endl;
0077   }
0078 
0079   m_truthInfo = findNode::getClass<PHG4TruthInfoContainer>(topNode, "G4TruthInfo");
0080   if (!m_truthInfo)
0081   {
0082     if (Verbosity() >= VERBOSITY_MORE)
0083     {
0084       std::cout << __FILE__ << ": Missing node G4TruthInfo" << std::endl;
0085     }
0086   }
0087 
0088   m_dst_truth_reco_map = findNode::getClass<PHG4ParticleSvtxMap_v1>(topNode, "PHG4ParticleSvtxMap");
0089   if (m_dst_truth_reco_map)
0090   {
0091     if (Verbosity() >= VERBOSITY_MORE)
0092     {
0093       std::cout << __FILE__ << ": PHG4ParticleSvtxMap found, truth matching will be more accurate" << std::endl;
0094     }
0095   }
0096   else
0097   {
0098     if (Verbosity() >= VERBOSITY_MORE)
0099     {
0100       std::cout << __FILE__ << ": PHG4ParticleSvtxMap not found, reverting to true matching by momentum relations. Truth matching will be less accurate" << std::endl;
0101     }
0102   }
0103 
0104   if (m_decay_descriptor.empty() && !m_decayMap->empty())
0105   {
0106     getDecayDescriptor();
0107     getNDaughters();
0108   }
0109 
0110   bool reconstructableDecay = false;
0111   for (auto &iter : *m_decayMap)
0112   {
0113     Decay decay = iter.second;
0114     ++m_counter_allDecays;
0115 
0116     m_all_tracks_reconstructed = findTracks(topNode, decay);
0117     if (m_all_tracks_reconstructed)
0118     {
0119       ++m_counter_acceptedDecays;
0120       reconstructableDecay = true;
0121     }
0122 
0123     if (m_write_nTuple)
0124     {
0125       if (m_counter_allDecays == 1)
0126       {
0127         initializeBranches();
0128       }
0129       m_tree->Fill();
0130       resetBranches();
0131     }
0132   }
0133 
0134   if (m_triggerOnDecay && !reconstructableDecay)
0135   {
0136     if (Verbosity() >= VERBOSITY_MORE)
0137     {
0138       std::cout << "No decays were reconstructable in this event, skipping" << std::endl;
0139     }
0140     return Fun4AllReturnCodes::ABORTEVENT;
0141   }
0142   else
0143   {
0144     return Fun4AllReturnCodes::EVENT_OK;
0145   }
0146 }
0147 
0148 int HFTrackEfficiency::End(PHCompositeNode * /*topNode*/)
0149 {
0150   PrintEff();
0151 
0152   if (m_write_nTuple && m_counter_allDecays != 0)
0153   {
0154     m_outfile->Write();
0155     m_outfile->Close();
0156     delete m_outfile;
0157   }
0158 
0159   return Fun4AllReturnCodes::EVENT_OK;
0160 }
0161 
0162 bool HFTrackEfficiency::findTracks(PHCompositeNode *topNode, Decay decay)
0163 {
0164   int trackableParticles[] = {11, 13, 211, 321, 2212};
0165   bool recoTrackFound = false;
0166   std::set<SvtxTrack *> selectedTracks;
0167 
0168   CLHEP::HepLorentzVector motherRecoLV;
0169   CLHEP::HepLorentzVector daughterSumTrueLV;
0170   CLHEP::HepLorentzVector *motherTrueLV = new CLHEP::HepLorentzVector();
0171   CLHEP::HepLorentzVector *daughterTrueLV = new CLHEP::HepLorentzVector();
0172 
0173   HepMC::GenEvent *theEvent = nullptr;
0174   if (m_geneventmap)
0175   {
0176     m_genevt = m_geneventmap->get(decay[0].first.first);
0177     assert(m_genevt);
0178 
0179     theEvent = m_genevt->getEvent();
0180     HepMC::GenParticle *mother = theEvent->barcode_to_particle(decay[0].first.second);
0181     assert(mother);
0182     if (Verbosity() >= VERBOSITY_MORE)
0183     {
0184       mother->print();
0185     }
0186 
0187     m_true_mother_pT = mother->momentum().perp();
0188     m_true_mother_p = std::sqrt(std::pow(mother->momentum().px(), 2) + std::pow(mother->momentum().py(), 2) + std::pow(mother->momentum().pz(), 2));  // Must have an old HepMC build, no mag function
0189     m_true_mother_eta = mother->momentum().eta();
0190 
0191     HepMC::GenVertex *thisVtx = mother->production_vertex();
0192     m_primary_vtx_x = thisVtx->point3d().x();
0193     m_primary_vtx_y = thisVtx->point3d().y();
0194     m_primary_vtx_z = thisVtx->point3d().z();
0195   }
0196 
0197   for (unsigned int i = 1; i < decay.size(); ++i)
0198   {
0199     m_dst_track = nullptr;
0200     int truth_ID = -1;
0201     if (std::find(std::begin(trackableParticles), std::end(trackableParticles),
0202                   std::abs(decay[i].second)) != std::end(trackableParticles))
0203     {
0204       if (theEvent && decay[i].first.second > -1)
0205       {
0206         HepMC::GenParticle *daughterHepMC = theEvent->barcode_to_particle(decay[i].first.second);
0207         if (Verbosity() >= VERBOSITY_MORE)
0208         {
0209           daughterHepMC->print();
0210         }
0211 
0212         daughterTrueLV->setVectM(CLHEP::Hep3Vector(daughterHepMC->momentum().px(), daughterHepMC->momentum().py(), daughterHepMC->momentum().pz()), getParticleMass(decay[i].second));
0213         daughterSumTrueLV += *daughterTrueLV;
0214 
0215         m_true_track_PID[i - 1] = daughterHepMC->pdg_id();
0216 
0217         // Now get the decay vertex position
0218         HepMC::GenVertex *thisVtx = daughterHepMC->production_vertex();
0219         m_secondary_vtx_x = thisVtx->point3d().x();
0220         m_secondary_vtx_y = thisVtx->point3d().y();
0221         m_secondary_vtx_z = thisVtx->point3d().z();
0222 
0223         // We need the G4 ID, not the HepMC ID to use the truth/reco map
0224         if (m_dst_truth_reco_map)
0225         {
0226           PHG4TruthInfoContainer::ConstRange range = m_truthInfo->GetParticleRange();
0227 
0228           for (PHG4TruthInfoContainer::ConstIterator iter = range.first; iter != range.second; ++iter)
0229           {
0230             PHG4Particle *daughterG4 = iter->second;
0231 
0232             if (std::abs(daughterG4->get_px() - daughterTrueLV->x()) <= 5e-3 &&
0233                 std::abs(daughterG4->get_py() - daughterTrueLV->y()) <= 5e-3 &&
0234                 std::abs(daughterG4->get_pz() - daughterTrueLV->z()) <= 5e-3 && daughterG4->get_pid() == decay[i].second)
0235             {
0236               truth_ID = daughterG4->get_track_id();
0237               break;
0238             }
0239           }
0240         }
0241       }
0242       else
0243       {
0244         PHG4TruthInfoContainer::ConstRange range = m_truthInfo->GetParticleRange();
0245 
0246         for (PHG4TruthInfoContainer::ConstIterator iter = range.first; iter != range.second; ++iter)
0247         {
0248           PHG4Particle *daughterG4 = iter->second;
0249 
0250           PHG4Particle *motherG4 = nullptr;
0251           if (daughterG4->get_parent_id() != 0)
0252           {
0253             motherG4 = m_truthInfo->GetParticle(daughterG4->get_parent_id());
0254           }
0255           else
0256           {
0257             continue;
0258           }
0259 
0260           if (motherG4->get_pid() == decay[0].second && motherG4->get_barcode() == decay[0].first.second && daughterG4->get_pid() == decay[i].second && daughterG4->get_barcode() == decay[i].first.second)
0261           {
0262             if (Verbosity() >= VERBOSITY_MORE)
0263             {
0264               daughterG4->identify();
0265             }
0266 
0267             CLHEP::Hep3Vector *mother3Vector = new CLHEP::Hep3Vector(motherG4->get_px(), motherG4->get_py(), motherG4->get_pz());
0268             motherTrueLV->setVectM((*mother3Vector), getParticleMass(decay[0].second));
0269             m_true_mother_pT = motherTrueLV->perp();
0270             m_true_mother_p = mother3Vector->mag();
0271             m_true_mother_eta = motherTrueLV->pseudoRapidity();
0272 
0273             PHG4VtxPoint *thisVtx = m_truthInfo->GetVtx(motherG4->get_vtx_id());
0274             m_primary_vtx_x = thisVtx->get_x();
0275             m_primary_vtx_y = thisVtx->get_y();
0276             m_primary_vtx_z = thisVtx->get_z();
0277 
0278             daughterTrueLV->setVectM(CLHEP::Hep3Vector(daughterG4->get_px(), daughterG4->get_py(), daughterG4->get_pz()), getParticleMass(decay[i].second));
0279             daughterSumTrueLV += *daughterTrueLV;
0280 
0281             // Now get the decay vertex position
0282             thisVtx = m_truthInfo->GetVtx(daughterG4->get_vtx_id());
0283             m_secondary_vtx_x = thisVtx->get_x();
0284             m_secondary_vtx_y = thisVtx->get_y();
0285             m_secondary_vtx_z = thisVtx->get_z();
0286 
0287             m_true_track_PID[i - 1] = daughterG4->get_pid();
0288             truth_ID = daughterG4->get_track_id();
0289 
0290             delete mother3Vector;
0291           }
0292         }
0293       }
0294 
0295       m_true_track_pT[i - 1] = (float) daughterTrueLV->perp();
0296       m_true_track_eta[i - 1] = (float) daughterTrueLV->pseudoRapidity();
0297       m_min_true_track_pT = std::min(m_true_track_pT[i - 1], m_min_true_track_pT);
0298       m_max_true_track_pT = std::max(m_true_track_pT[i - 1], m_max_true_track_pT);
0299 
0300       if (m_dst_truth_reco_map && truth_ID >= 0)
0301       {
0302         std::map<float, std::set<unsigned int>> reco_set = m_dst_truth_reco_map->get(truth_ID);
0303         if (reco_set.size() == 0)
0304         {
0305           continue;
0306         }
0307         const auto &best_weight = reco_set.rbegin();
0308         if (best_weight->second.size() == 0)
0309         {
0310           continue;
0311         }
0312         unsigned int best_reco_id = *best_weight->second.rbegin();
0313         m_dst_track = m_input_trackMap->get(best_reco_id);
0314         if (m_dst_track)
0315         {
0316           m_used_truth_reco_map[i - 1] = true;
0317           recoTrackFound = true;
0318         }
0319       }
0320       else
0321       {
0322         for (auto &iter : *m_input_trackMap)
0323         {
0324           m_dst_track = iter.second;
0325           float delta_px = (m_dst_track->get_px() - daughterTrueLV->px()) / daughterTrueLV->px();
0326           float delta_py = (m_dst_track->get_py() - daughterTrueLV->py()) / daughterTrueLV->py();
0327           float delta_pz = (m_dst_track->get_pz() - daughterTrueLV->pz()) / daughterTrueLV->pz();
0328 
0329           if (std::abs(delta_px) <= m_truthRecoMatchPercent && std::abs(delta_py) <= m_truthRecoMatchPercent && std::abs(delta_pz) <= m_truthRecoMatchPercent)
0330           {
0331             recoTrackFound = true;
0332             break;
0333           }
0334         }
0335       }
0336 
0337       if (recoTrackFound)
0338       {
0339         selectedTracks.insert(m_dst_track);
0340         if (Verbosity() >= VERBOSITY_MORE)
0341         {
0342           m_dst_track->identify();
0343         }
0344         m_reco_track_exists[i - 1] = true;
0345         m_reco_track_pT[i - 1] = m_dst_track->get_pt();
0346         m_reco_track_eta[i - 1] = m_dst_track->get_eta();
0347         m_reco_track_chi2nDoF[i - 1] = m_dst_track->get_chisq() / m_dst_track->get_ndf();
0348         if (m_dst_track->get_silicon_seed())
0349         {
0350           m_reco_track_silicon_seeds[i - 1] = static_cast<int>(m_dst_track->get_silicon_seed()->size_cluster_keys());
0351         }
0352         else
0353         {
0354           m_reco_track_silicon_seeds[i - 1] = 0;
0355         }
0356         m_reco_track_tpc_seeds[i - 1] = static_cast<int>(m_dst_track->get_tpc_seed()->size_cluster_keys());
0357         m_min_reco_track_pT = std::min(m_reco_track_pT[i - 1], m_min_reco_track_pT);
0358         m_max_reco_track_pT = std::max(m_reco_track_pT[i - 1], m_max_reco_track_pT);
0359 
0360         CLHEP::HepLorentzVector *daughterRecoLV = new CLHEP::HepLorentzVector();
0361         daughterRecoLV->setVectM(CLHEP::Hep3Vector(m_dst_track->get_px(), m_dst_track->get_py(), m_dst_track->get_pz()), getParticleMass(m_true_track_PID[i - 1]));
0362 
0363         motherRecoLV += *daughterRecoLV;
0364         delete daughterRecoLV;
0365       }
0366 
0367       recoTrackFound = false;
0368     }
0369   }
0370 
0371   m_true_mother_mass = daughterSumTrueLV.m();
0372   bool foundDecay = true;
0373   if (selectedTracks.size() == m_nDaughters)
0374   {
0375     m_reco_mother_mass = motherRecoLV.m();
0376     if (m_write_track_map)
0377     {
0378       m_output_trackMap = findNode::getClass<SvtxTrackMap>(topNode, outputNodeName.c_str());
0379       for (auto &track : selectedTracks)
0380       {
0381         m_output_trackMap->insertWithKey(track, track->get_id());
0382       }
0383     }
0384   }
0385   else
0386   {
0387     foundDecay = false;
0388   }
0389 
0390   selectedTracks.clear();
0391 
0392   delete motherTrueLV;
0393   delete daughterTrueLV;
0394 
0395   return foundDecay;
0396 }
0397 
0398 void HFTrackEfficiency::initializeBranches()
0399 {
0400   m_outfile = new TFile(m_outfile_name.c_str(), "RECREATE");
0401   m_tree = new TTree("HFTrackEfficiency", "HFTrackEfficiency");
0402   m_tree->OptimizeBaskets();
0403   m_tree->SetAutoSave(-5e6);  // Save the output file every 5MB
0404 
0405   m_tree->Branch("all_tracks_reconstructed", &m_all_tracks_reconstructed, "all_tracks_reconstructed/O");
0406   m_tree->Branch("true_mother_mass", &m_true_mother_mass, "true_mother_mass/F");
0407   m_tree->Branch("reco_mother_mass", &m_reco_mother_mass, "reco_mother_mass/F");
0408   m_tree->Branch("true_mother_pT", &m_true_mother_pT, "true_mother_pT/F");
0409   m_tree->Branch("true_mother_p", &m_true_mother_p, "true_mother_p/F");
0410   m_tree->Branch("true_mother_eta", &m_true_mother_eta, "true_mother_eta/F");
0411   m_tree->Branch("min_true_track_pT", &m_min_true_track_pT, "min_true_track_pT/F");
0412   m_tree->Branch("min_reco_track_pT", &m_min_reco_track_pT, "min_reco_track_pT/F");
0413   m_tree->Branch("max_true_track_pT", &m_max_true_track_pT, "max_true_track_pT/F");
0414   m_tree->Branch("max_reco_track_pT", &m_max_reco_track_pT, "max_reco_track_pT/F");
0415 
0416   for (unsigned int iTrack = 0; iTrack < m_nDaughters; ++iTrack)
0417   {
0418     std::string daughter_number = "track_" + std::to_string(iTrack + 1);
0419     m_tree->Branch("reco_" + TString(daughter_number) + "_exists", &m_reco_track_exists[iTrack], "reco_" + TString(daughter_number) + "_exists/O");
0420     m_tree->Branch("reco_" + TString(daughter_number) + "_used_truth_reco_map", &m_used_truth_reco_map[iTrack], "reco_" + TString(daughter_number) + "_used_truth_reco_map/O");
0421     m_tree->Branch("true_" + TString(daughter_number) + "_pT", &m_true_track_pT[iTrack], "true_" + TString(daughter_number) + "_pT/F");
0422     m_tree->Branch("reco_" + TString(daughter_number) + "_pT", &m_reco_track_pT[iTrack], "reco_" + TString(daughter_number) + "_pT/F");
0423     m_tree->Branch("true_" + TString(daughter_number) + "_eta", &m_true_track_eta[iTrack], "true_" + TString(daughter_number) + "_eta/F");
0424     m_tree->Branch("reco_" + TString(daughter_number) + "_eta", &m_reco_track_eta[iTrack], "reco_" + TString(daughter_number) + "_eta/F");
0425     m_tree->Branch("true_" + TString(daughter_number) + "_PID", &m_true_track_PID[iTrack], "true_" + TString(daughter_number) + "_PID/F");
0426     m_tree->Branch("reco_" + TString(daughter_number) + "_chi2nDoF", &m_reco_track_chi2nDoF[iTrack], "reco_" + TString(daughter_number) + "_chi2nDoF/F");
0427     m_tree->Branch("reco_" + TString(daughter_number) + "_silicon_seeds", &m_reco_track_silicon_seeds[iTrack], "reco_" + TString(daughter_number) + "_silicon_seeds/I");
0428     m_tree->Branch("reco_" + TString(daughter_number) + "_tpc_seeds", &m_reco_track_tpc_seeds[iTrack], "reco_" + TString(daughter_number) + "_tpc_seeds/I");
0429   }
0430 
0431   m_tree->Branch("true_primary_vertex_x", &m_primary_vtx_x, "true_primary_vertex_x/F");
0432   m_tree->Branch("true_primary_vertex_y", &m_primary_vtx_y, "true_primary_vertex_y/F");
0433   m_tree->Branch("true_primary_vertex_z", &m_primary_vtx_z, "true_primary_vertex_z/F");
0434   m_tree->Branch("true_secondary_vertex_x", &m_secondary_vtx_x, "true_secondary_vertex_x/F");
0435   m_tree->Branch("true_secondary_vertex_y", &m_secondary_vtx_y, "true_secondary_vertex_y/F");
0436   m_tree->Branch("true_secondary_vertex_z", &m_secondary_vtx_z, "true_secondary_vertex_z/F");
0437 }
0438 
0439 void HFTrackEfficiency::resetBranches()
0440 {
0441   m_all_tracks_reconstructed = false;
0442   m_true_mother_mass = std::numeric_limits<float>::quiet_NaN();
0443   m_reco_mother_mass = std::numeric_limits<float>::quiet_NaN();
0444   m_true_mother_pT = std::numeric_limits<float>::quiet_NaN();
0445   m_true_mother_p = std::numeric_limits<float>::quiet_NaN();
0446   m_true_mother_eta = std::numeric_limits<float>::quiet_NaN();
0447   m_min_true_track_pT = std::numeric_limits<float>::max();
0448   m_min_reco_track_pT = std::numeric_limits<float>::max();
0449   m_max_true_track_pT = -1 * std::numeric_limits<float>::max();
0450   m_max_reco_track_pT = -1 * std::numeric_limits<float>::max();
0451   for (unsigned int iTrack = 0; iTrack < m_nDaughters; ++iTrack)
0452   {
0453     m_reco_track_exists[iTrack] = false;
0454     m_used_truth_reco_map[iTrack] = false;
0455     m_true_track_pT[iTrack] = std::numeric_limits<float>::quiet_NaN();
0456     m_reco_track_pT[iTrack] = std::numeric_limits<float>::quiet_NaN();
0457     m_true_track_eta[iTrack] = std::numeric_limits<float>::quiet_NaN();
0458     m_reco_track_eta[iTrack] = std::numeric_limits<float>::quiet_NaN();
0459     m_true_track_PID[iTrack] = std::numeric_limits<float>::quiet_NaN();
0460     m_reco_track_chi2nDoF[iTrack] = std::numeric_limits<float>::quiet_NaN();
0461     m_reco_track_silicon_seeds[iTrack] = 0;
0462     m_reco_track_tpc_seeds[iTrack] = 0;
0463   }
0464 
0465   m_primary_vtx_x = std::numeric_limits<float>::quiet_NaN();
0466   m_primary_vtx_y = std::numeric_limits<float>::quiet_NaN();
0467   m_primary_vtx_z = std::numeric_limits<float>::quiet_NaN();
0468   m_secondary_vtx_x = std::numeric_limits<float>::quiet_NaN();
0469   m_secondary_vtx_y = std::numeric_limits<float>::quiet_NaN();
0470   m_secondary_vtx_z = std::numeric_limits<float>::quiet_NaN();
0471 }
0472 
0473 void HFTrackEfficiency::getDecayDescriptor()
0474 {
0475   Decay decay = m_decayMap->begin()->second;
0476   m_decay_descriptor = getParticleName(decay[0].second);
0477   m_decay_descriptor += " ->";
0478   for (unsigned int i = 1; i < decay.size(); ++i)
0479   {
0480     m_decay_descriptor += " ";
0481     m_decay_descriptor += getParticleName(decay[i].second);
0482   }
0483 }
0484 
0485 void HFTrackEfficiency::getNDaughters()
0486 {
0487   int trackableParticles[] = {11, 13, 211, 321, 2212};
0488   m_nDaughters = 0;
0489   Decay decay = m_decayMap->begin()->second;
0490   for (unsigned int i = 1; i < decay.size(); ++i)
0491   {
0492     if (std::find(std::begin(trackableParticles), std::end(trackableParticles),
0493                   std::abs(decay[i].second)) != std::end(trackableParticles))
0494     {
0495       ++m_nDaughters;
0496     }
0497   }
0498 }
0499 
0500 std::string HFTrackEfficiency::getParticleName(const int PDGID)
0501 {
0502   return TDatabasePDG::Instance()->GetParticle(PDGID)->GetName();
0503 }
0504 
0505 float HFTrackEfficiency::getParticleMass(const int PDGID)
0506 {
0507   return TDatabasePDG::Instance()->GetParticle(PDGID)->Mass();
0508 }
0509 
0510 //____________________________________________________________________________..
0511 void HFTrackEfficiency::PrintEff()
0512 {
0513   std::cout << "\n--------------- Heavy Flavor Tracking Efficiency ---------------" << std::endl;
0514   std::cout << "Tracking Efficiency for " << m_decay_descriptor << std::endl;
0515   std::cout << "Number of decays fully in acceptance: " << m_counter_allDecays << std::endl;
0516   std::cout << "Number of decays with all tracks reconstructed: " << m_counter_acceptedDecays << std::endl;
0517   std::cout << "Efficiency: " << std::setw(4) << (float) m_counter_acceptedDecays * 100 / float(m_counter_allDecays) << "%" << std::endl;
0518   std::cout << "----------------------------------------------------------------\n"
0519             << std::endl;
0520 }