Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:01:04

0001 #include "CommonTools/CandUtils/interface/zMCLeptonDaughters.h"
0002 #include "DataFormats/Candidate/interface/Candidate.h"
0003 #include "FWCore/Utilities/interface/EDMException.h"
0004 #include <cassert>
0005 using namespace std;
0006 using namespace reco;
0007 
0008 pair<const Candidate*, const Candidate*> zMCLeptonDaughters(const Candidate& z, int leptonPdgId) {
0009   if (z.numberOfDaughters() < 2)
0010     throw cms::Exception("RuntimeError") << "calling helper function reco::zMCLeptonDaughters passing a Z candidate"
0011                                             "with less than 2 daughters ("
0012                                          << z.numberOfDaughters() << ").\n";
0013   const Candidate* dau0 = z.daughter(0);
0014   const Candidate* dau1 = z.daughter(1);
0015   for (size_t i0 = 0; i0 < dau0->numberOfDaughters(); ++i0) {
0016     const Candidate* ddau0 = dau0->daughter(i0);
0017     if (abs(ddau0->pdgId()) == leptonPdgId && ddau0->status() == 1) {
0018       dau0 = ddau0;
0019       break;
0020     }
0021   }
0022   for (size_t i1 = 0; i1 < dau1->numberOfDaughters(); ++i1) {
0023     const Candidate* ddau1 = dau1->daughter(i1);
0024     if (abs(ddau1->pdgId()) == leptonPdgId && ddau1->status() == 1) {
0025       dau1 = ddau1;
0026       break;
0027     }
0028   }
0029   assert(abs(dau0->pdgId()) == leptonPdgId && dau0->status() == 1);
0030   assert(abs(dau1->pdgId()) == leptonPdgId && dau1->status() == 1);
0031   return make_pair(dau0, dau1);
0032 }