CMS 3D CMS Logo

DimuonStatistics.cc
Go to the documentation of this file.
5 #include <iostream>
6 
8 public:
10  void analyze(const edm::Event&, const edm::EventSetup&) override;
11  void endJob() override;
12 private:
15  std::vector<unsigned int> matched_, unMatched_;
17 };
18 
26 using namespace std;
27 using namespace reco;
28 using namespace edm;
29 const unsigned int maxEntries = 10;
30 
32  src_(cfg.getParameter<InputTag>("src")),
34  matched_(maxEntries + 1, 0),
35  unMatched_(maxEntries + 1, 0),
36  ptMin_(cfg.getUntrackedParameter<double>("ptMin")),
37  massMin_(cfg.getUntrackedParameter<double>("massMin")),
38  massMax_(cfg.getUntrackedParameter<double>("massMax")),
39  etaMin_(cfg.getUntrackedParameter<double>("etaMin")),
40  etaMax_(cfg.getUntrackedParameter<double>("etaMax")),
41  trkIso_(cfg.getUntrackedParameter<double>("trkIso")){
42 }
43 
45  cout << src_.encode() << endl ;
46  cout << " == Matched == " << endl;
47  for(unsigned int i = 0; i <= maxEntries; ++i) {
48  cout << i << ": " << matched_[i];
49  if (i < maxEntries) cout << ", ";
50  }
51  cout << endl;
52  cout << " == unMatched == " << endl;
53  for(unsigned int i = 0; i <= maxEntries; ++i) {
54  cout << i << ": " << unMatched_[i];
55  if (i < maxEntries) cout << ", ";
56 }
57  cout << endl;
58 }
59 
62  evt.getByToken(srcToken_, src);
63  double trackIso1=-1;
64  double trackIso2=-1;
65  int j=0;
66  unsigned int matched = 0, unMatched = 0;
67  cout << ">> entries in " << src_ << ": " << src->size() << endl;
68  for(CandidateView::const_iterator i = src->begin(); i != src->end(); ++i) {
69  j++;
70  const Candidate * dau1 = i->daughter(0);
71  const Candidate * dau2 = i->daughter(1);
72  if(dau1 == nullptr|| dau2 == nullptr)
74  "one of the two daughter does not exist\n";
75  const Candidate * c1 = dau1->masterClone().get();
76  GenParticleRef mc1;
77  const pat::Muon * mu1 = dynamic_cast<const pat::Muon*>(c1);
78  if(mu1 != nullptr) {
79  mc1 = mu1->genParticleRef();
80  // if (mc1.isNonnull()) cout << "DimuonStatistics> genParticleRef1 " << mc1->pdgId() << endl;
81  trackIso1=mu1->trackIso();
82  } else {
83  const pat::GenericParticle * gp1 = dynamic_cast<const pat::GenericParticle*>(c1);
84  if(gp1 == nullptr)
86  "first of two daughter is neither a pat::Muon not pat::GenericParticle\n";
87  mc1 = gp1->genParticleRef();
88  }
89  const Candidate * c2 = dau2->masterClone().get();
90  GenParticleRef mc2;
91  const pat::Muon * mu2 = dynamic_cast<const pat::Muon*>(c2);
92  if(mu2 != nullptr) {
93  mc2 = mu2->genParticleRef();
94  // if (mc2.isNonnull()) cout << "DimuonStatistics> genParticleRef2 " << mc2->pdgId() << endl;
95  trackIso2=mu2->trackIso();
96  } else {
97  const pat::GenericParticle * gp2 = dynamic_cast<const pat::GenericParticle*>(c2);
98  if(gp2 == nullptr)
100  "first of two daughter is neither a pat::Muon not pat::GenericParticle\n";
101  mc2 = gp2->genParticleRef();
102  }
103  GenParticleRef dimuonMatch;
104  if(mc1.isNonnull() && mc2.isNonnull()) {
105  cout << "DimuonStatistics> mc1: " << mc1->pdgId() << ", mc2: " << mc2->pdgId() << endl;
106  int k=0;
107  do {
108  k++;
109  mc1 = mc1->numberOfMothers() > 0 ? mc1->motherRef() : GenParticleRef();
110  mc2 = mc2->numberOfMothers() > 0 ? mc2->motherRef() : GenParticleRef();
111  // cout << "DimuonStatistics> do loop: " << k << " id1 " << mc1->pdgId() << " id2 " << mc2->pdgId() << endl;
112  } while(mc1 != mc2 && mc1.isNonnull() && mc2.isNonnull());
113  if(mc1.isNonnull() && mc2.isNonnull() && mc1->pdgId()==23) {
114  dimuonMatch = mc1;
115  }
116  }
117  // cout << "DimuonMatcher> dimuonMatch " << dimuonMatch.isNonnull() << endl;
118  if( (fabs(dau1->eta()) > etaMin_ && fabs(dau1->eta()) < etaMax_) && dau1->pt()>ptMin_ &&
119  ( (fabs(dau2->eta()) > etaMin_ ) && (fabs(dau2->eta()) < etaMax_) ) && dau2->pt()>ptMin_ &&
120  trackIso1 < trkIso_ && trackIso2 < trkIso_ && i->mass() > massMin_ && i->mass() < massMax_
121  ) {
122  cout << "dimuon mass " << i->mass() << endl;
123  if(dimuonMatch.isNonnull())
124  ++matched;
125  else
126  ++unMatched;
127 
128 
129  }
130  }
131  cout << "matched: " << matched << ", unmatched: " << unMatched << endl;
132 
133  if(matched > maxEntries) matched = maxEntries;
134  if(unMatched > maxEntries) unMatched = maxEntries;
135  ++matched_[matched];
136  ++unMatched_[unMatched];
137 }
138 
139 
141 
value_type const * get() const
Definition: RefToBase.h:234
edm::EDGetTokenT< reco::CandidateView > srcToken_
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:253
reco::GenParticleRef genParticleRef(size_t idx=0) const
Definition: PATObject.h:225
void endJob() override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
DimuonStatistics(const edm::ParameterSet &cfg)
void analyze(const edm::Event &, const edm::EventSetup &) override
size_type size() const
const unsigned int maxEntries
std::vector< unsigned int > unMatched_
std::string encode() const
Definition: InputTag.cc:166
float trackIso() const
Definition: Muon.h:171
std::vector< unsigned int > matched_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
Analysis-level Generic Particle class (e.g. for hadron or muon not fully reconstructed) ...
edm::Ref< edm::HepMCProduct, HepMC::GenParticle > GenParticleRef
const_iterator begin() const
int k[5][pyjets_maxn]
virtual const CandidateBaseRef & masterClone() const =0
virtual double eta() const =0
momentum pseudorapidity
virtual double pt() const =0
transverse momentum
fixed size matrix
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
edm::InputTag src_
const_iterator end() const
Analysis-level muon class.
Definition: Muon.h:50