CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_9/src/ElectroWeakAnalysis/ZMuMu/plugins/ZMCHistogrammer.cc

Go to the documentation of this file.
00001 #include "FWCore/Framework/interface/EDAnalyzer.h"
00002 #include "FWCore/Utilities/interface/InputTag.h"
00003 #include "TH1.h"
00004 
00005 class ZMCHistogrammer : public edm::EDAnalyzer {
00006 public:
00007   ZMCHistogrammer(const edm::ParameterSet& pset); 
00008 private:
00009   virtual void analyze(const edm::Event& event, const edm::EventSetup& setup);
00010   edm::InputTag  z_, gen_, match_;
00011   unsigned int nbinsMass_, nbinsPt_, nbinsAng_, nbinsMassRes_;
00012   double massMax_, ptMax_, angMax_, massResMax_;
00013   TH1F *h_nZ_, *h_mZ_, *h_ptZ_, *h_phiZ_, *h_thetaZ_, *h_etaZ_, *h_rapidityZ_;
00014   TH1F *h_invmMuMu_;
00015   TH1F *h_nZMC_, *h_mZMC_, *h_ptZMC_, *h_phiZMC_, *h_thetaZMC_, *h_etaZMC_, *h_rapidityZMC_;
00016   TH1F *h_invmMuMuMC_ , * h_weight_histo; 
00017   //TH1F *h_mZ2vs3MC_, *h_ptZ2vs3MC_, *h_phiZ2vs3MC_, *h_thetaZ2vs3MC_, *h_etaZ2vs3MC_, *h_rapidityZ2vs3MC_;
00018   TH1F *h_mResZ_, *h_ptResZ_, *h_phiResZ_, *h_thetaResZ_, *h_etaResZ_, *h_rapidityResZ_;
00019   TH1F *h_mResZMuMu_, *h_mRatioZMuMu_;
00020   TH1F *h_mResZMuMuMC_, *h_mRatioZMuMuMC_;
00021   bool isMCatNLO_; 
00022 };
00023 
00024 #include "DataFormats/Candidate/interface/Candidate.h"
00025 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00026 #include "DataFormats/Candidate/interface/CandMatchMap.h"
00027 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00028 #include "DataFormats/HepMCCandidate/interface/GenParticleFwd.h"
00029 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00030 
00031 #include "HepMC/WeightContainer.h"
00032 #include "HepMC/GenEvent.h"
00033 #include "FWCore/Framework/interface/Event.h"
00034 #include "DataFormats/Common/interface/Handle.h"
00035 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00036 #include "FWCore/ServiceRegistry/interface/Service.h"
00037 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00038 #include <cmath>
00039 #include <iostream>
00040 
00041 using namespace std;
00042 using namespace reco;
00043 using namespace edm;
00044 
00045 ZMCHistogrammer::ZMCHistogrammer(const ParameterSet& pset) :
00046   z_(pset.getParameter<InputTag>("z")),
00047   gen_(pset.getParameter<InputTag>("gen")), 
00048   match_(pset.getParameter<InputTag>("match")), 
00049   nbinsMass_(pset.getUntrackedParameter<unsigned int>("nbinsMass")),
00050   nbinsPt_(pset.getUntrackedParameter<unsigned int>("nbinsPt")),
00051   nbinsAng_(pset.getUntrackedParameter<unsigned int>("nbinsAng")),
00052   nbinsMassRes_(pset.getUntrackedParameter<unsigned int>("nbinsMassRes")),
00053   massMax_(pset.getUntrackedParameter<double>("massMax")),
00054   ptMax_(pset.getUntrackedParameter<double>("ptMax")),
00055   angMax_(pset.getUntrackedParameter<double>("angMax")), 
00056   massResMax_(pset.getUntrackedParameter<double>("massResMax")),
00057   isMCatNLO_(pset.getUntrackedParameter<bool>("isMCatNLO"))  { 
00058   cout << ">>> Z Histogrammer constructor" << endl;
00059   Service<TFileService> fs;
00060   TFileDirectory ZHisto = fs->mkdir( "ZRecoHisto" );
00061   TFileDirectory ZMCHisto = fs->mkdir( "ZMCHisto" );
00062   TFileDirectory ZResHisto = fs->mkdir( "ZResolutionHisto" );
00063   //TFileDirectory Z2vs3MCHisto = fs->mkdir( "Z2vs3MCHisto" );
00064   h_nZ_ = ZHisto.make<TH1F>("ZNumber", "number of Z particles", 11, -0.5, 10.5);
00065   h_mZ_ = ZHisto.make<TH1F>("ZMass", "Z mass (GeV/c^{2})", nbinsMass_,  0, massMax_);
00066   h_ptZ_ = ZHisto.make<TH1F>("ZPt", "Z p_{t} (GeV/c)", nbinsPt_, 0, ptMax_);
00067   h_phiZ_ = ZHisto.make<TH1F>("ZPhi", "Z #phi", nbinsAng_,  -angMax_, angMax_);
00068   h_thetaZ_ = ZHisto.make<TH1F>("Ztheta", "Z #theta", nbinsAng_,  0, angMax_);
00069   h_etaZ_ = ZHisto.make<TH1F>("ZEta", "Z #eta", nbinsAng_,  -angMax_, angMax_);
00070   h_rapidityZ_ = ZHisto.make<TH1F>("ZRapidity", "Z rapidity", nbinsAng_,  -angMax_, angMax_);
00071   h_invmMuMu_ = ZHisto.make<TH1F>("MuMuMass", "#mu #mu invariant mass",nbinsMass_,  0, massMax_);
00072   h_weight_histo  = ZHisto.make<TH1F>("weight_histo","weight_histo",20,-10,10); 
00073   h_nZMC_ = ZMCHisto.make<TH1F>("ZMCNumber", "number of Z MC particles", 11, -0.5, 10.5);
00074   h_mZMC_ = ZMCHisto.make<TH1F>("ZMCMass", "Z MC mass (GeV/c^{2})", nbinsMass_,  0, massMax_);
00075   h_ptZMC_ = ZMCHisto.make<TH1F>("ZMCPt", "Z MC p_{t} (GeV/c)", nbinsPt_, 0, ptMax_);
00076   h_phiZMC_ = ZMCHisto.make<TH1F>("ZMCPhi", "Z MC #phi", nbinsAng_,  -angMax_, angMax_);
00077   h_thetaZMC_ = ZMCHisto.make<TH1F>("ZMCTheta", "Z MC #theta", nbinsAng_,  0, angMax_);
00078   h_etaZMC_ = ZMCHisto.make<TH1F>("ZMCEta", "Z MC #eta", nbinsAng_,  -angMax_, angMax_);
00079   h_rapidityZMC_ = ZMCHisto.make<TH1F>("ZMCRapidity", "Z MC rapidity",  nbinsAng_,  -angMax_, angMax_);
00080   h_invmMuMuMC_ = ZMCHisto.make<TH1F>("MuMuMCMass", "#mu #mu MC invariant mass", 
00081                                       nbinsMass_,  0, massMax_);
00082   /*
00083   h_mZ2vs3MC_ = Z2vs3MCHisto.make<TH1F>("Z2vs3MCMass", "Z MC st 2 vs st 3 mass (GeV/c^{2})", 
00084                                         nbinsMassRes_, -massResMax_, massResMax_);
00085   h_ptZ2vs3MC_ = Z2vs3MCHisto.make<TH1F>("Z2vs3MCPt", "Z MC st 2 vs st 3 p_{t} (GeV/c)", 
00086                                          nbinsPt_, -ptMax_, ptMax_);
00087   h_phiZ2vs3MC_ = Z2vs3MCHisto.make<TH1F>("Z2vs3MCPhi", "Z MC st 2 vs st 3 #phi", 
00088                                           nbinsAng_,  -angMax_, angMax_);
00089   h_thetaZ2vs3MC_ = Z2vs3MCHisto.make<TH1F>("Z2vs3MCTheta", "Z MC st 2 vs st 3 #theta", 
00090                                             nbinsAng_,  -angMax_, angMax_);
00091   h_etaZ2vs3MC_ = Z2vs3MCHisto.make<TH1F>("Z2vs3MCEta", "Z MC st 2 vs st 3 #eta", 
00092                                           nbinsAng_,  -angMax_, angMax_);
00093   h_rapidityZ2vs3MC_ = Z2vs3MCHisto.make<TH1F>("Z2vs3MCRapidity", "Z MC st 2 vs st 3 rapidity", 
00094                                        nbinsAng_,  -angMax_, angMax_);
00095   */
00096   h_mResZ_ = ZResHisto.make<TH1F>("ZMassResolution", "Z mass Resolution (GeV/c^{2})", 
00097                                   nbinsMassRes_, -massResMax_, massResMax_);
00098   h_ptResZ_ = ZResHisto.make<TH1F>("ZPtResolution", "Z p_{t} Resolution (GeV/c)", 
00099                                    nbinsPt_, -ptMax_, ptMax_);
00100   h_phiResZ_ = ZResHisto.make<TH1F>("ZPhiResolution", "Z #phi Resolution", 
00101                                     nbinsAng_,  -angMax_, angMax_);
00102   h_thetaResZ_ = ZResHisto.make<TH1F>("ZThetaResolution", "Z #theta Resolution", 
00103                                       nbinsAng_, -angMax_, angMax_);
00104   h_etaResZ_ = ZResHisto.make<TH1F>("ZEtaResolution", "Z #eta Resolution", 
00105                                     nbinsAng_,  -angMax_, angMax_);
00106   h_rapidityResZ_ = ZResHisto.make<TH1F>("ZRapidityResolution", "Z rapidity Resolution", 
00107                                          nbinsAng_,  -angMax_, angMax_);
00108   h_mResZMuMu_ = ZResHisto.make<TH1F>("ZToMuMuRecoMassResolution", 
00109                                       "Z Reco vs matched final state #mu #mu mass Difference (GeV/c^{2})", 
00110                                       nbinsMassRes_, -massResMax_, massResMax_);
00111   h_mRatioZMuMu_ = ZResHisto.make<TH1F>("ZToMuMuRecoMassRatio", 
00112                                         "Z Reco vs matched final state #mu #mu mass Ratio", 
00113                                         4000, 0, 2);
00114   h_mResZMuMuMC_ = ZResHisto.make<TH1F>("ZToMuMuMCMassResolution", 
00115                                         "Z vs final state #mu #mu MC mass Difference (GeV/c^{2})", 
00116                                         nbinsMassRes_/2 + 1, -2*massResMax_/nbinsMassRes_, massResMax_);
00117   h_mRatioZMuMuMC_ = ZResHisto.make<TH1F>("ZToMuMuMCMassRatio", 
00118                                           "Z vs final state #mu #mu MC mass Ratio", 
00119                                           2002, 0.999, 2);
00120 }
00121 
00122 void ZMCHistogrammer::analyze(const edm::Event& event, const edm::EventSetup& setup) { 
00123   cout << ">>> Z Histogrammer analyze" << endl;
00124   Handle<CandidateView> z;
00125   Handle<CandidateView> gen;
00126   Handle<std::vector<GenParticleRef> > match;
00127   event.getByLabel(z_, z);
00128   event.getByLabel(gen_, gen);
00129   event.getByLabel(match_, match);
00130   h_nZ_->Fill(z->size());
00131  
00132 // get HepMC::GenEvent ...
00133    Handle<HepMCProduct> evt_h;
00134    event.getByType(evt_h);
00135    HepMC::GenEvent* evt = new  HepMC::GenEvent(*(evt_h->GetEvent()));
00136 
00137 
00138    // get weight and fill it to histogram
00139   HepMC::WeightContainer weights = evt->weights();   
00140   double weight = weights.front();
00141   if(!weight) weight=1.;
00142   h_weight_histo->Fill(weight); 
00143   
00144   if(isMCatNLO_) {
00145     weight > 0 ?  weight=1. : weight=-1.;
00146   }
00147 
00148   for(unsigned int i = 0; i < z->size(); ++i) {
00149     const Candidate &zCand = (*z)[i];
00150     h_mZ_->Fill(zCand.mass(),weight );
00151     h_ptZ_->Fill(zCand.pt(),weight);
00152     h_phiZ_->Fill(zCand.phi(),weight);
00153     h_thetaZ_->Fill(zCand.theta(),weight);
00154     h_etaZ_->Fill(zCand.eta(),weight);
00155     h_rapidityZ_->Fill(zCand.rapidity(),weight);
00156     CandidateBaseRef zCandRef = z->refAt(i);   
00157    
00158  GenParticleRef zMCMatch = (*match)[i];
00159     if(zMCMatch.isNonnull() && zMCMatch->pdgId()==23) {
00160       h_mResZ_->Fill(zCandRef->mass() - zMCMatch->mass());
00161       h_ptResZ_->Fill(zCandRef->pt() - zMCMatch->pt());
00162       h_phiResZ_->Fill(zCandRef->phi() - zMCMatch->phi());
00163       h_thetaResZ_->Fill(zCandRef->theta() - zMCMatch->theta());
00164       h_etaResZ_->Fill(zCandRef->eta() - zMCMatch->eta());
00165       h_rapidityResZ_->Fill(zCandRef->rapidity() - zMCMatch->rapidity());
00166       const Candidate * dau0 = zMCMatch->daughter(0);
00167       const Candidate * dau1 = zMCMatch->daughter(1);
00168       for(unsigned int i0 = 0; i0 < dau0->numberOfDaughters(); ++i0) {
00169         const Candidate * ddau0 = dau0->daughter(i0);
00170         if(abs(ddau0->pdgId())==13 && ddau0->status()==1) {
00171           dau0 = ddau0; break;
00172         }
00173       }
00174       for(unsigned int i1 = 0; i1 < dau1->numberOfDaughters(); ++i1) {
00175         const Candidate * ddau1 = dau1->daughter(i1);
00176         if(abs(ddau1->pdgId())==13 && ddau1->status()==1) {
00177           dau1 = ddau1; break;
00178         }
00179       }
00180       assert(abs(dau0->pdgId())==13 && dau0->status()==1);
00181       assert(abs(dau1->pdgId())==13 && dau1->status()==1);
00182       double invMass = (dau0->p4()+dau1->p4()).mass();
00183       h_invmMuMu_->Fill(invMass,weight);
00184       h_mResZMuMu_->Fill(zCand.mass() - invMass);
00185       h_mRatioZMuMu_->Fill(zCand.mass()/invMass);
00186     }
00187   }
00188   h_nZMC_->Fill(gen->size());
00189   for(unsigned int i = 0; i < gen->size(); ++i) {
00190     const Candidate &genCand = (*gen)[i];
00191     if((genCand.pdgId() == 23) && (genCand.status() == 2)) //this is an intermediate Z0
00192       cout << ">>> intermediate Z0 found, with " << genCand.numberOfDaughters() 
00193            << " daughters" << endl;
00194     if((genCand.pdgId() == 23)&&(genCand.status() == 3)) { //this is a Z0
00195       cout << ">>> Z0 found, with " << genCand.numberOfDaughters() 
00196            << " daughters" << endl;
00197       h_mZMC_->Fill(genCand.mass(),weight);
00198       h_ptZMC_->Fill(genCand.pt(),weight);
00199       h_phiZMC_->Fill(genCand.phi(),weight);
00200       h_thetaZMC_->Fill(genCand.theta(),weight);
00201       h_etaZMC_->Fill(genCand.eta(),weight);
00202       h_rapidityZMC_->Fill(genCand.rapidity(),weight);
00203       Particle::LorentzVector pZ(0, 0, 0, 0);
00204       int nMu = 0;
00205       for(unsigned int j = 0; j < genCand.numberOfDaughters(); ++j) { 
00206         const Candidate *dauGen = genCand.daughter(j);
00207         /*
00208         if((dauGen->pdgId() == 23) && (dauGen->status() == 2)) { 
00209           h_mZ2vs3MC_->Fill(genCand.mass() - dauGen->mass());
00210           h_ptZ2vs3MC_->Fill(genCand.pt() - dauGen->pt());
00211           h_phiZ2vs3MC_->Fill(genCand.phi() - dauGen->phi());
00212           h_thetaZ2vs3MC_->Fill(genCand.theta() - dauGen->theta());
00213           h_etaZ2vs3MC_->Fill(genCand.eta() - dauGen->eta());
00214           h_rapidityZ2vs3MC_->Fill(genCand.rapidity() - dauGen->rapidity());
00215         }
00216         */
00217         if((abs(dauGen->pdgId()) == 13) && (dauGen->numberOfDaughters() != 0)) {
00218           //we are looking for photons of final state radiation
00219           cout << ">>> The muon " << j 
00220                << " has " << dauGen->numberOfDaughters() << " daughters" <<endl;
00221           for(unsigned int k = 0; k < dauGen->numberOfDaughters(); ++k) {
00222             const Candidate * dauMuGen = dauGen->daughter(k);
00223             cout << ">>> Mu " << j 
00224                  << " daughter MC " << k 
00225                  << " PDG Id " << dauMuGen->pdgId() 
00226                  << ", status " << dauMuGen->status() 
00227                  << ", charge " << dauMuGen->charge() 
00228                  << endl;
00229             if(abs(dauMuGen->pdgId()) == 13 && dauMuGen->status() ==1) {
00230               pZ += dauMuGen->p4();
00231               nMu ++;
00232             }
00233           }
00234         }
00235       }
00236       assert(nMu == 2);
00237       double mZ = pZ.mass();
00238       h_invmMuMuMC_->Fill(mZ,weight);
00239       h_mResZMuMuMC_->Fill(genCand.mass() - mZ);
00240       h_mRatioZMuMuMC_->Fill(genCand.mass()/mZ);
00241     }
00242   }
00243 }
00244 
00245 #include "FWCore/Framework/interface/MakerMacros.h"
00246 
00247 DEFINE_FWK_MODULE(ZMCHistogrammer);
00248