CMS 3D CMS Logo

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