CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DimuonMCMatcher.cc
Go to the documentation of this file.
5 #include <iostream>
6 
8 public:
10  virtual void produce(edm::Event&, const edm::EventSetup&) override;
11 private:
13 };
14 
22 using namespace std;
23 using namespace reco;
24 using namespace edm;
25 
27  srcToken_(consumes<CandidateView>(cfg.getParameter<InputTag>("src"))) {
28  produces<vector<GenParticleRef> >();
29 }
30 
33  evt.getByToken(srcToken_, src);
34  auto_ptr<vector<GenParticleRef> > matched(new vector<GenParticleRef>);
35  matched->reserve(src->size());
36  int j=0;
37  for(CandidateView::const_iterator i = src->begin(); i != src->end(); ++i) {
38  j++;
39  const Candidate * dau1 = i->daughter(0);
40  const Candidate * dau2 = i->daughter(1);
41  if(dau1 == 0|| dau2 == 0)
43  "one of the two daughter does not exist\n";
44  const Candidate * c1 = dau1->masterClone().get();
45  GenParticleRef mc1;
46  const pat::Muon * mu1 = dynamic_cast<const pat::Muon*>(c1);
47  if(mu1 != 0) {
48  mc1 = mu1->genParticleRef();
49  // if (mc1.isNonnull()) cout << "DimuonMCMatcher> genParticleRef1 " << mc1->pdgId() << endl;
50  } else {
51  const pat::GenericParticle * gp1 = dynamic_cast<const pat::GenericParticle*>(c1);
52  if(gp1 == 0)
54  "first of two daughter is neither a pat::Muon not pat::GenericParticle\n";
55  mc1 = gp1->genParticleRef();
56  }
57  const Candidate * c2 = dau2->masterClone().get();
58  GenParticleRef mc2;
59  const pat::Muon * mu2 = dynamic_cast<const pat::Muon*>(c2);
60  if(mu2 != 0) {
61  mc2 = mu2->genParticleRef();
62  // if (mc2.isNonnull()) cout << "DimuonMCMatcher> genParticleRef2 " << mc2->pdgId() << endl;
63  } else {
64  const pat::GenericParticle * gp2 = dynamic_cast<const pat::GenericParticle*>(c2);
65  if(gp2 == 0)
67  "first of two daughter is neither a pat::Muon not pat::GenericParticle\n";
68  mc2 = gp2->genParticleRef();
69  }
70  GenParticleRef dimuonMatch;
71  // cout << "DimuonMatcher> mc1 " << mc1.isNonnull() << " mc2 " << mc2.isNonnull() << endl;
72  if(mc1.isNonnull() && mc2.isNonnull()) {
73  int k=0;
74  do {
75  k++;
76  mc1 = mc1->numberOfMothers() > 0 ? mc1->motherRef() : GenParticleRef();
77  mc2 = mc2->numberOfMothers() > 0 ? mc2->motherRef() : GenParticleRef();
78  // cout << "DimuonMCMatcher> do loop: " << k << " id1 " << mc1->pdgId() << " id2 " << mc2->pdgId() << endl;
79  } while(mc1 != mc2 && mc1.isNonnull() && mc2.isNonnull());
80  if(mc1.isNonnull() && mc2.isNonnull() && mc1->pdgId()==23) {
81  dimuonMatch = mc1;
82  }
83  }
84  // cout << "DimuonMatcher> dimuonMatch " << dimuonMatch.isNonnull() << endl;
85  matched->push_back(dimuonMatch);
86  }
87 
88  evt.put(matched);
89 }
90 
91 
93 
value_type const * get() const
Definition: RefToBase.h:213
int i
Definition: DBlmapReader.cc:9
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
reco::GenParticleRef genParticleRef(size_t idx=0) const
Definition: PATObject.h:217
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::Ref< GenParticleCollection > GenParticleRef
persistent reference to a GenParticle
DimuonMCMatcher(const edm::ParameterSet &cfg)
edm::EDGetTokenT< reco::CandidateView > srcToken_
Analysis-level Generic Particle class (e.g. for hadron or muon not fully reconstructed) ...
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
virtual void produce(edm::Event &, const edm::EventSetup &) override
int j
Definition: DBlmapReader.cc:9
int k[5][pyjets_maxn]
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
Analysis-level muon class.
Definition: Muon.h:50
virtual const CandidateBaseRef & masterClone() const =0