Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:26:42

0001 // -*- C++ -*-
0002 //
0003 // Package:    METAlgorithms
0004 // Class:      SignPFSpecificAlgo
0005 //
0006 // Authors: A. Khukhunaishvili (Cornell), L. Gibbons (Cornell)
0007 // First Implementation: November 11, 2011
0008 //
0009 //
0010 
0011 //____________________________________________________________________________||
0012 #include "RecoMET/METAlgorithms/interface/SignPFSpecificAlgo.h"
0013 
0014 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0015 
0016 //____________________________________________________________________________||
0017 metsig::SignPFSpecificAlgo::SignPFSpecificAlgo() : resolutions_(nullptr), algo_() { clusteredParticlePtrs_.clear(); }
0018 
0019 //____________________________________________________________________________||
0020 void metsig::SignPFSpecificAlgo::setResolutions(metsig::SignAlgoResolutions* resolutions) {
0021   resolutions_ = resolutions;
0022 }
0023 
0024 //____________________________________________________________________________||
0025 void metsig::SignPFSpecificAlgo::addPFJets(const edm::View<reco::PFJet>* PFJets) {
0026   std::vector<metsig::SigInputObj> vobj;
0027   for (edm::View<reco::PFJet>::const_iterator jet = PFJets->begin(); jet != PFJets->end(); ++jet) {
0028     vobj.push_back(resolutions_->evalPFJet(&(*jet)));
0029     std::vector<reco::PFCandidatePtr> pfs = jet->getPFConstituents();
0030     for (std::vector<reco::PFCandidatePtr>::const_iterator it = pfs.begin(); it != pfs.end(); ++it) {
0031       reco::CandidatePtr ptr(*it);
0032       clusteredParticlePtrs_.insert(ptr);
0033     }
0034   }
0035   algo_.addObjects(vobj);
0036 }
0037 
0038 //____________________________________________________________________________||
0039 void metsig::SignPFSpecificAlgo::useOriginalPtrs(const edm::ProductID& productID) {
0040   std::set<reco::CandidatePtr>::const_iterator it = clusteredParticlePtrs_.begin();
0041   reco::CandidatePtr ptr(*it);
0042   if (ptr.id() == productID)
0043     return;  //If the first element is from the right product, return
0044 
0045   std::set<reco::CandidatePtr> temp;
0046   for (; it != clusteredParticlePtrs_.end(); ++it) {
0047     reco::CandidatePtr ptr(*it);
0048     while (ptr.id() != productID) {
0049       ptr = ptr->sourceCandidatePtr(0);
0050       if (ptr.isNull())
0051         return;  //if it does not get to the correct product, return
0052     }
0053     temp.insert(ptr);
0054   }
0055   clusteredParticlePtrs_.clear();
0056   clusteredParticlePtrs_ = temp;
0057 }
0058 
0059 //____________________________________________________________________________||
0060 void metsig::SignPFSpecificAlgo::addPFCandidate(reco::PFCandidatePtr pf) {
0061   if (clusteredParticlePtrs_.find(pf) != clusteredParticlePtrs_.end()) {
0062     return;  //pf candidate already added in jet collection
0063   }
0064   std::vector<metsig::SigInputObj> vobj;
0065   vobj.push_back(resolutions_->evalPF(&(*pf)));
0066   algo_.addObjects(vobj);
0067 }
0068 
0069 //____________________________________________________________________________||
0070 reco::METCovMatrix metsig::SignPFSpecificAlgo::mkSignifMatrix(edm::Handle<edm::View<reco::Candidate> >& PFCandidates) {
0071   useOriginalPtrs(PFCandidates.id());
0072   for (edm::View<reco::Candidate>::const_iterator iParticle = (PFCandidates.product())->begin();
0073        iParticle != (PFCandidates.product())->end();
0074        ++iParticle) {
0075     const reco::PFCandidate* pfCandidate = dynamic_cast<const reco::PFCandidate*>(&(*iParticle));
0076     if (!pfCandidate)
0077       continue;
0078     reco::CandidatePtr dau(PFCandidates, iParticle - PFCandidates->begin());
0079     if (dau.isNull())
0080       continue;
0081     if (!dau.isAvailable())
0082       continue;
0083     reco::PFCandidatePtr pf(dau.id(), pfCandidate, dau.key());
0084     addPFCandidate(pf);
0085   }
0086   return getSignifMatrix();
0087 }
0088 
0089 //____________________________________________________________________________||