00001 #include "FWCore/Framework/interface/EDAnalyzer.h"
00002 #include "FWCore/Utilities/interface/InputTag.h"
00003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00004 #include "FWCore/Framework/interface/Event.h"
00005 #include "DataFormats/Common/interface/Handle.h"
00006 #include "DataFormats/Candidate/interface/Candidate.h"
00007 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00008 #include <iostream>
00009
00010 using namespace edm;
00011 using namespace reco;
00012 using namespace std;
00013
00014 const Candidate * mcMuDaughter(const Candidate * c) {
00015 unsigned int n = c->numberOfDaughters();
00016 for(unsigned int i = 0; i < n; ++i) {
00017 const Candidate * d = c->daughter(i);
00018 if(fabs(d->pdgId())==13) return d;
00019 }
00020 return 0;
00021 }
00022
00023 struct ZSelector {
00024 ZSelector(double ptMin, double etaDau0Min, double etaDau0Max,double etaDau1Min, double etaDau1Max, double massMin, double massMax, double massMinZMC, double massMaxZMC) :
00025 ptMin_(ptMin), etaDau0Min_(etaDau0Min),etaDau0Max_(etaDau0Max), etaDau1Min_(etaDau1Min),etaDau1Max_(etaDau1Max),
00026 massMin_(massMin), massMax_(massMax), massMinZMC_(massMinZMC), massMaxZMC_(massMaxZMC) { }
00027 int operator()(const Candidate& c) const {
00028
00029 if (c.numberOfDaughters()<2) return 0;
00030 if (c.numberOfDaughters()>=6) return 0;
00031 const Candidate * d0 = c.daughter(0);
00032 const Candidate * d1 = c.daughter(1);
00033 if(c.numberOfDaughters()>2) {
00034 if (d0->numberOfDaughters()>0) d0 = mcMuDaughter(d0);
00035 if (d1->numberOfDaughters()>0) d1 = mcMuDaughter(d1);
00036 }
00037 int temp_cut= 0;
00039 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_ ) ) {
00040 temp_cut=1;
00041 if(d0->pt() > ptMin_ && d1->pt() > ptMin_) {
00042 temp_cut=2;
00043 double m = (d0->p4() + d1->p4()).mass();
00044 if(m > massMin_ && m < massMax_) temp_cut=3;
00045 if (c.mass()> massMinZMC_ && c.mass() < massMaxZMC_) temp_cut =4;
00046 }
00047 }
00048
00049 return temp_cut;
00050 }
00051 double ptMin_, etaDau0Min_, etaDau0Max_, etaDau1Min_, etaDau1Max_, massMin_, massMax_, massMinZMC_, massMaxZMC_;
00052 };
00053
00054 class MCAcceptanceAnalyzer : public EDAnalyzer {
00055 public:
00056 MCAcceptanceAnalyzer(const ParameterSet& cfg);
00057 private:
00058 void analyze(const Event&, const EventSetup&);
00059 void endJob();
00060 InputTag zToMuMu_, zToMuMuMC_, zToMuMuMatched_;
00061 long nZToMuMu_, selZToMuMu_, nZToMuMuMC_, selZToMuMuMC_, nZToMuMuMCMatched_, selZToMuMuMCMatched_, nZToMuMuMCDen_;
00062 ZSelector select_, select_OnlyMassCut_;
00063 };
00064
00065 MCAcceptanceAnalyzer::MCAcceptanceAnalyzer(const ParameterSet& cfg) :
00066 zToMuMu_(cfg.getParameter<InputTag>("zToMuMu")),
00067 zToMuMuMC_(cfg.getParameter<InputTag>("zToMuMuMC")),
00068 zToMuMuMatched_(cfg.getParameter<InputTag>("zToMuMuMatched")),
00069 nZToMuMu_(0), selZToMuMu_(0),
00070 nZToMuMuMC_(0), selZToMuMuMC_(0),
00071 nZToMuMuMCMatched_(0), selZToMuMuMCMatched_(0), nZToMuMuMCDen_(0),
00072 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") ),
00073 select_OnlyMassCut_(-1, -9999,9999 , -9999, 9999,0, 0, cfg.getParameter<double>("massMinZMC"), cfg.getParameter<double>("massMaxZMC") )
00074 {
00075 }
00076
00077 void MCAcceptanceAnalyzer::analyze(const Event& evt, const EventSetup&) {
00078 Handle<CandidateView> zToMuMu;
00079 evt.getByLabel(zToMuMu_, zToMuMu);
00080 Handle<CandidateView> zToMuMuMC;
00081 evt.getByLabel(zToMuMuMC_, zToMuMuMC);
00082 Handle<GenParticleMatch > zToMuMuMatched;
00083 evt.getByLabel(zToMuMuMatched_, zToMuMuMatched);
00084
00085 long nZToMuMuMC = zToMuMuMC->size();
00086 long nZToMuMuMatched = zToMuMuMatched->size();
00087
00088
00089
00090
00091
00092
00093
00094
00095 nZToMuMuMC_ += nZToMuMuMC;
00096 for(long i = 0; i < nZToMuMuMC; ++i) {
00097 const Candidate & z = (*zToMuMuMC)[i];
00098 if(select_(z)==4) ++selZToMuMuMC_;
00099 if(select_OnlyMassCut_(z)==4) ++nZToMuMuMCDen_;
00100
00101 }
00102
00103
00104 for(long i = 0; i < nZToMuMuMatched; ++i) {
00105
00106 const Candidate & z = (*zToMuMu)[i];
00107 CandidateBaseRef zRef = zToMuMu->refAt(i);
00108 GenParticleRef mcRef = (*zToMuMuMatched)[zRef];
00109
00110
00111 if(mcRef.isNonnull()) {
00112 ++nZToMuMu_;
00113 ++nZToMuMuMCMatched_;
00114
00115 int selectZ = select_(z);
00116 if(selectZ==4) ++selZToMuMu_;
00117
00118
00119
00120 int selectMC = select_(*mcRef);
00121
00122 if(selectMC==4) ++selZToMuMuMCMatched_;
00123
00124 if(selectZ != selectMC) {
00125 cout << ">>> select reco: " << selectZ << ", select mc: " << selectMC << endl;
00126 if ((selectZ * selectMC) ==0 ) break;
00127 if (z.numberOfDaughters()> 1){
00128 const Candidate * d0 = z.daughter(0), * d1 = z.daughter(1);
00129 if (mcRef->numberOfDaughters()>1){
00130 const Candidate * mcd0 = mcMuDaughter(mcRef->daughter(0)),
00131 * mcd1 = mcMuDaughter(mcRef->daughter(1));
00132 double m = z.mass(), mcm = (mcd0->p4()+mcd1->p4()).mass();
00133 cout << ">>> reco pt1, eta1: " << d0->pt() <<", " << d0->eta()
00134 << ", 2: " << d1->pt() << ", " << d1->eta()
00135 << ", mass = " << m << endl;
00136 cout << ">>> mc pt1, eta1: " << mcd0->pt() <<", " << mcd0->eta()
00137 << ", 2: " << mcd1->pt() << ", " << mcd1->eta()
00138 << ", mass = " << mcm << endl;
00139 }
00140 }
00141
00142 }
00143
00144 if ((selectZ==3) && (selectMC==3)) break;
00145 }
00146 }
00147
00148 }
00149
00150 void MCAcceptanceAnalyzer::endJob() {
00151 double effZToMuMu = double(selZToMuMu_)/double(nZToMuMu_);
00152 double errZToMuMu = sqrt(effZToMuMu*(1. - effZToMuMu)/nZToMuMu_);
00153 double effZToMuMuMC = double(selZToMuMuMC_)/double(nZToMuMuMC_);
00154 double errZToMuMuMC = sqrt(effZToMuMuMC*(1. - effZToMuMuMC)/nZToMuMuMC_);
00155 double effZToMuMuMCDen = double(selZToMuMuMC_)/double(nZToMuMuMCDen_);
00156 double errZToMuMuMCDen = sqrt(effZToMuMuMCDen*(1. - effZToMuMuMCDen)/nZToMuMuMCDen_);
00157 double effZToMuMuMCMatched = double(selZToMuMuMCMatched_)/double(nZToMuMuMCMatched_);
00158 double errZToMuMuMCMatched = sqrt(effZToMuMuMCMatched*(1. - effZToMuMuMCMatched)/nZToMuMuMCMatched_);
00159 cout << ">>> " << zToMuMu_ << ": " << selZToMuMu_ << "/" << nZToMuMu_
00160 << " = " << effZToMuMu << " +/- " << errZToMuMu << endl;
00161 cout << ">>> " << zToMuMuMC_ << " - matched: " << selZToMuMuMCMatched_ << "/" << nZToMuMuMCMatched_
00162 << " = " << effZToMuMuMCMatched << " +/- " << errZToMuMuMCMatched << endl;
00163 cout << " if the two numbers above are the same we can neglete resolution effect and quote the acceptance as the number below.... " << endl;
00164 cout << "********* acceptance m>sampleMCMassCut (usually 20 or 40) ******** " << endl;
00165 cout << ">>> " << zToMuMuMC_ << ": " << selZToMuMuMC_ << "/" << nZToMuMuMC_
00166 << " = " << effZToMuMuMC << " +/- " << errZToMuMuMC << endl;
00167 cout << "********* acceptance in the given mass range ******** " << endl;
00168 cout << ">>> " << zToMuMuMC_ << ": " << selZToMuMuMC_ << "/" << nZToMuMuMCDen_
00169 << " = " << effZToMuMuMCDen << " +/- " << errZToMuMuMCDen << endl;
00170
00171 }
00172
00173 #include "FWCore/Framework/interface/MakerMacros.h"
00174
00175 DEFINE_FWK_MODULE(MCAcceptanceAnalyzer);
00176