CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MCAcceptanceAnalyzer.cc
Go to the documentation of this file.
8 #include <iostream>
9 
10 using namespace edm;
11 using namespace reco;
12 using namespace std;
13 
14 const Candidate * mcMuDaughter(const Candidate * c) {
15  unsigned int n = c->numberOfDaughters();
16  for(unsigned int i = 0; i < n; ++i) {
17  const Candidate * d = c->daughter(i);
18  if(fabs(d->pdgId())==13) return d;
19  }
20  return 0;
21 }
22 
23 struct ZSelector { // modify this selector in order to return an integer (0: no eta cut, 1: eta cut only, 2 eta && pt cut, 3: eta, pt and Mass cut, 4: mass cut on the denominator Z MC)
24  ZSelector(double ptMin, double etaDau0Min, double etaDau0Max,double etaDau1Min, double etaDau1Max, double massMin, double massMax, double massMinZMC, double massMaxZMC) :
25  ptMin_(ptMin), etaDau0Min_(etaDau0Min),etaDau0Max_(etaDau0Max), etaDau1Min_(etaDau1Min),etaDau1Max_(etaDau1Max),
26  massMin_(massMin), massMax_(massMax), massMinZMC_(massMinZMC), massMaxZMC_(massMaxZMC) { }
27  int operator()(const Candidate& c) const {
28  // std::cout << "c.numberOfDaughters(): " << c.numberOfDaughters()<< std::endl;
29  if (c.numberOfDaughters()<2) return 0;
30  if (c.numberOfDaughters()>=6) return 0;
31  const Candidate * d0 = c.daughter(0);
32  const Candidate * d1 = c.daughter(1);
33  if(c.numberOfDaughters()>2) {
34  if (d0->numberOfDaughters()>0) d0 = mcMuDaughter(d0);
35  if (d1->numberOfDaughters()>0) d1 = mcMuDaughter(d1);
36  }
37  int temp_cut= 0;
39  if( ( fabs(d0->eta()) > etaDau0Min_ && fabs(d1->eta())>etaDau1Min_ && fabs(d0->eta()) < etaDau0Max_ && fabs(d1->eta()) <etaDau1Max_ ) || ( fabs(d0->eta()) > etaDau1Min_ && fabs(d1->eta())>etaDau0Min_ && fabs(d0->eta()) < etaDau1Max_ && fabs(d1->eta()) <etaDau0Max_ ) ) {
40  temp_cut=1;
41  if(d0->pt() > ptMin_ && d1->pt() > ptMin_) {
42  temp_cut=2;
43  double m = (d0->p4() + d1->p4()).mass();
44  if(m > massMin_ && m < massMax_) temp_cut=3;
45  if (c.mass()> massMinZMC_ && c.mass() < massMaxZMC_) temp_cut =4;
46  }
47  }
48 
49  return temp_cut;
50  }
51  double ptMin_, etaDau0Min_, etaDau0Max_, etaDau1Min_, etaDau1Max_, massMin_, massMax_, massMinZMC_, massMaxZMC_;
52 };
53 
55 public:
57 private:
58  void analyze(const Event&, const EventSetup&) override;
59  void endJob() override;
60  InputTag zToMuMu_, zToMuMuMC_, zToMuMuMatched_;
64  long nZToMuMu_, selZToMuMu_, nZToMuMuMC_, selZToMuMuMC_, nZToMuMuMCMatched_, selZToMuMuMCMatched_, nZToMuMuMCDen_;
66 };
67 
69  zToMuMu_(cfg.getParameter<InputTag>("zToMuMu")),
70  zToMuMuMC_(cfg.getParameter<InputTag>("zToMuMuMC")),
71  zToMuMuMatched_(cfg.getParameter<InputTag>("zToMuMuMatched")),
72  zToMuMuToken_(consumes<CandidateView>(zToMuMu_)),
73  zToMuMuMCToken_(consumes<CandidateView>(zToMuMuMC_)),
74  zToMuMuMatchedToken_(consumes<GenParticleMatch> (zToMuMuMatched_)),
75  nZToMuMu_(0), selZToMuMu_(0),
76  nZToMuMuMC_(0), selZToMuMuMC_(0),
77  nZToMuMuMCMatched_(0), selZToMuMuMCMatched_(0), nZToMuMuMCDen_(0),
78  select_(cfg.getParameter<double>("ptMin"), cfg.getParameter<double>("etaDau0Min"), cfg.getParameter<double>("etaDau0Max"), cfg.getParameter<double>("etaDau1Min"), cfg.getParameter<double>("etaDau1Max"),cfg.getParameter<double>("massMin"), cfg.getParameter<double>("massMax"), cfg.getParameter<double>("massMinZMC"), cfg.getParameter<double>("massMaxZMC") ),
79  select_OnlyMassCut_(-1, -9999,9999 , -9999, 9999,0, 0, cfg.getParameter<double>("massMinZMC"), cfg.getParameter<double>("massMaxZMC") )
80 {
81 }
82 
84  Handle<CandidateView> zToMuMu;
85  evt.getByToken(zToMuMuToken_, zToMuMu);
86  Handle<CandidateView> zToMuMuMC;
87  evt.getByToken(zToMuMuMCToken_, zToMuMuMC);
88  Handle<GenParticleMatch > zToMuMuMatched;
89  evt.getByToken(zToMuMuMatchedToken_, zToMuMuMatched);
90  // long nZToMuMu = zToMuMu->size();
91  long nZToMuMuMC = zToMuMuMC->size();
92  long nZToMuMuMatched = zToMuMuMatched->size();
93 
94 
95  // cout << ">>> " << zToMuMu_ << " has " << nZToMuMu << " entries" << endl;
96  //cout << ">>> " << zToMuMuMC_ << " has " << nZToMuMuMC << " entries" << endl;
97  //cout << ">>> " << zToMuMuMatched_ << " has " << nZToMuMuMatched << " entries" << endl;
98 
99 
100 
101  nZToMuMuMC_ += nZToMuMuMC;
102  for(long i = 0; i < nZToMuMuMC; ++i) {
103  const Candidate & z = (*zToMuMuMC)[i];
104  if(select_(z)==4) ++selZToMuMuMC_;
106 
107  }
108 
109 
110  for(long i = 0; i < nZToMuMuMatched; ++i) {
111 
112  const Candidate & z = (*zToMuMu)[i];
113  CandidateBaseRef zRef = zToMuMu->refAt(i);
114  GenParticleRef mcRef = (*zToMuMuMatched)[zRef];
115 
116 
117  if(mcRef.isNonnull()) { // z candidate matched to Z MC
118  ++nZToMuMu_;
120 
121  int selectZ = select_(z);
122  if(selectZ==4) ++selZToMuMu_;
123 
124 
125 
126  int selectMC = select_(*mcRef);
127 
128  if(selectMC==4) ++selZToMuMuMCMatched_;
129 
130  if(selectZ != selectMC) {
131  cout << ">>> select reco: " << selectZ << ", select mc: " << selectMC << endl;
132  if ((selectZ * selectMC) ==0 ) break;
133  if (z.numberOfDaughters()> 1){
134  const Candidate * d0 = z.daughter(0), * d1 = z.daughter(1);
135  if (mcRef->numberOfDaughters()>1){
136  const Candidate * mcd0 = mcMuDaughter(mcRef->daughter(0)),
137  * mcd1 = mcMuDaughter(mcRef->daughter(1));
138  double m = z.mass(), mcm = (mcd0->p4()+mcd1->p4()).mass();
139  cout << ">>> reco pt1, eta1: " << d0->pt() <<", " << d0->eta()
140  << ", 2: " << d1->pt() << ", " << d1->eta()
141  << ", mass = " << m << endl;
142  cout << ">>> mc pt1, eta1: " << mcd0->pt() <<", " << mcd0->eta()
143  << ", 2: " << mcd1->pt() << ", " << mcd1->eta()
144  << ", mass = " << mcm << endl;
145  }
146  }
147 
148  }
149  // to avoid double counting
150  if ((selectZ==3) && (selectMC==3)) break;
151  }
152  }
153 
154 }
155 
157  double effZToMuMu = double(selZToMuMu_)/double(nZToMuMu_);
158  double errZToMuMu = sqrt(effZToMuMu*(1. - effZToMuMu)/nZToMuMu_);
159  double effZToMuMuMC = double(selZToMuMuMC_)/double(nZToMuMuMC_);
160  double errZToMuMuMC = sqrt(effZToMuMuMC*(1. - effZToMuMuMC)/nZToMuMuMC_);
161  double effZToMuMuMCDen = double(selZToMuMuMC_)/double(nZToMuMuMCDen_);
162  double errZToMuMuMCDen = sqrt(effZToMuMuMCDen*(1. - effZToMuMuMCDen)/nZToMuMuMCDen_);
163  double effZToMuMuMCMatched = double(selZToMuMuMCMatched_)/double(nZToMuMuMCMatched_);
164  double errZToMuMuMCMatched = sqrt(effZToMuMuMCMatched*(1. - effZToMuMuMCMatched)/nZToMuMuMCMatched_);
165  cout << ">>> " << zToMuMu_ << ": " << selZToMuMu_ << "/" << nZToMuMu_
166  << " = " << effZToMuMu << " +/- " << errZToMuMu << endl;
167  cout << ">>> " << zToMuMuMC_ << " - matched: " << selZToMuMuMCMatched_ << "/" << nZToMuMuMCMatched_
168  << " = " << effZToMuMuMCMatched << " +/- " << errZToMuMuMCMatched << endl;
169  cout << " if the two numbers above are the same we can neglete resolution effect and quote the acceptance as the number below.... " << endl;
170  cout << "********* acceptance m>sampleMCMassCut (usually 20 or 40) ******** " << endl;
171  cout << ">>> " << zToMuMuMC_ << ": " << selZToMuMuMC_ << "/" << nZToMuMuMC_
172  << " = " << effZToMuMuMC << " +/- " << errZToMuMuMC << endl;
173  cout << "********* acceptance in the given mass range ******** " << endl;
174  cout << ">>> " << zToMuMuMC_ << ": " << selZToMuMuMC_ << "/" << nZToMuMuMCDen_
175  << " = " << effZToMuMuMCDen << " +/- " << errZToMuMuMCDen << endl;
176 
177 }
178 
180 
182 
int i
Definition: DBlmapReader.cc:9
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
tuple cfg
Definition: looper.py:259
EDGetTokenT< GenParticleMatch > zToMuMuMatchedToken_
EDGetTokenT< CandidateView > zToMuMuToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
virtual double pt() const =0
transverse momentum
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual double mass() const =0
mass
float float float z
virtual size_type numberOfDaughters() const =0
number of daughters
T sqrt(T t)
Definition: SSEVec.h:48
const Candidate * mcMuDaughter(const Candidate *c)
virtual int pdgId() const =0
PDG identifier.
void analyze(const Event &, const EventSetup &) override
int operator()(const Candidate &c) const
EDGetTokenT< CandidateView > zToMuMuMCToken_
tuple cout
Definition: gather_cfg.py:121
MCAcceptanceAnalyzer(const ParameterSet &cfg)
virtual double eta() const =0
momentum pseudorapidity
ZSelector(double ptMin, double etaDau0Min, double etaDau0Max, double etaDau1Min, double etaDau1Max, double massMin, double massMax, double massMinZMC, double massMaxZMC)
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector