CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/ElectroWeakAnalysis/ZMuMu/plugins/DimuonMCMatcher.cc

Go to the documentation of this file.
00001 #include "FWCore/Framework/interface/EDProducer.h"
00002 #include "FWCore/Utilities/interface/InputTag.h"
00003 #include <iostream>
00004 
00005 class DimuonMCMatcher : public edm::EDProducer {
00006 public:
00007   DimuonMCMatcher(const edm::ParameterSet & cfg);
00008   virtual void produce(edm::Event&, const edm::EventSetup&);
00009 private:
00010   edm::InputTag src_;
00011 };
00012 
00013 #include "DataFormats/Candidate/interface/Candidate.h"
00014 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00015 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00016 #include "DataFormats/HepMCCandidate/interface/GenParticleFwd.h"
00017 #include "FWCore/Framework/interface/Event.h"
00018 #include "DataFormats/Common/interface/Handle.h"
00019 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00020 #include "DataFormats/PatCandidates/interface/Muon.h"
00021 #include "DataFormats/PatCandidates/interface/GenericParticle.h"
00022 using namespace std;
00023 using namespace reco;
00024 using namespace edm;
00025 
00026 DimuonMCMatcher::DimuonMCMatcher(const edm::ParameterSet & cfg) : 
00027   src_(cfg.getParameter<InputTag>("src")) {
00028   produces<vector<GenParticleRef> >();
00029 }
00030 
00031 void DimuonMCMatcher::produce(edm::Event& evt, const edm::EventSetup&) {
00032   Handle<CandidateView> src;
00033   evt.getByLabel(src_, src);
00034   auto_ptr<vector<GenParticleRef> > matched(new vector<GenParticleRef>);
00035   matched->reserve(src->size());
00036   int j=0;
00037   for(CandidateView::const_iterator i = src->begin(); i != src->end(); ++i) {
00038     j++;
00039     const Candidate * dau1 = i->daughter(0);
00040     const Candidate * dau2 = i->daughter(1);
00041     if(dau1 == 0|| dau2 == 0) 
00042       throw Exception(errors::InvalidReference) <<
00043         "one of the two daughter does not exist\n";
00044     const Candidate * c1 = dau1->masterClone().get();
00045     GenParticleRef mc1;
00046     const pat::Muon * mu1 = dynamic_cast<const pat::Muon*>(c1);
00047     if(mu1 != 0) {
00048       mc1 = mu1->genParticleRef();
00049       //     if (mc1.isNonnull()) cout << "DimuonMCMatcher> genParticleRef1 " << mc1->pdgId() << endl;
00050     } else {
00051       const pat::GenericParticle * gp1 = dynamic_cast<const pat::GenericParticle*>(c1);
00052       if(gp1 == 0) 
00053         throw Exception(errors::InvalidReference) <<
00054           "first of two daughter is neither a pat::Muon not pat::GenericParticle\n";
00055       mc1 = gp1->genParticleRef();
00056     }
00057     const Candidate * c2 = dau2->masterClone().get();
00058     GenParticleRef mc2;
00059     const pat::Muon * mu2 = dynamic_cast<const pat::Muon*>(c2);
00060     if(mu2 != 0) {
00061       mc2 = mu2->genParticleRef();
00062       //      if (mc2.isNonnull()) cout << "DimuonMCMatcher> genParticleRef2 " << mc2->pdgId() << endl;
00063     } else {
00064       const pat::GenericParticle * gp2 = dynamic_cast<const pat::GenericParticle*>(c2);
00065       if(gp2 == 0) 
00066         throw Exception(errors::InvalidReference) <<
00067           "first of two daughter is neither a pat::Muon not pat::GenericParticle\n";
00068       mc2 = gp2->genParticleRef();
00069     }
00070     GenParticleRef dimuonMatch;
00071     //    cout << "DimuonMatcher>  mc1 " << mc1.isNonnull() << "   mc2 " << mc2.isNonnull() << endl;
00072     if(mc1.isNonnull() && mc2.isNonnull()) {
00073       int k=0;
00074       do {
00075         k++;
00076         mc1 = mc1->numberOfMothers() > 0 ? mc1->motherRef() : GenParticleRef();
00077         mc2 = mc2->numberOfMothers() > 0 ? mc2->motherRef() : GenParticleRef();
00078         //      cout << "DimuonMCMatcher> do loop: " << k << "  id1 " << mc1->pdgId() << "    id2 " << mc2->pdgId() << endl;
00079       } while(mc1 != mc2 && mc1.isNonnull() && mc2.isNonnull());
00080       if(mc1.isNonnull() && mc2.isNonnull() && mc1->pdgId()==23) {
00081         dimuonMatch = mc1;
00082       }
00083     }
00084     //    cout << "DimuonMatcher> dimuonMatch " << dimuonMatch.isNonnull() << endl;
00085     matched->push_back(dimuonMatch);
00086   }
00087 
00088   evt.put(matched);
00089 }
00090 
00091 
00092 #include "FWCore/Framework/interface/MakerMacros.h"
00093 
00094 DEFINE_FWK_MODULE(DimuonMCMatcher);