CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_1/src/ElectroWeakAnalysis/ZMuMu/plugins/ZMassHistogrammer.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 ZMassHistogrammer : public edm::EDAnalyzer {
00006 public:
00007   ZMassHistogrammer(const edm::ParameterSet& pset); 
00008 private:
00009   virtual void analyze(const edm::Event& event, const edm::EventSetup& setup);
00010   edm::InputTag  z_, gen_;
00011   TH1F *h_mZ_, *h_mZMC_;
00012 };
00013 
00014 #include "DataFormats/Candidate/interface/Candidate.h"
00015 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00016 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00017 #include "DataFormats/HepMCCandidate/interface/GenParticleFwd.h"
00018 #include "FWCore/Framework/interface/Event.h"
00019 #include "DataFormats/Common/interface/Handle.h"
00020 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00021 #include "FWCore/ServiceRegistry/interface/Service.h"
00022 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00023 #include <iostream>
00024 
00025 using namespace std;
00026 using namespace reco;
00027 using namespace edm;
00028 
00029 ZMassHistogrammer::ZMassHistogrammer(const ParameterSet& pset) :
00030   z_(pset.getParameter<InputTag>("z")),
00031   gen_(pset.getParameter<InputTag>("gen")) { 
00032   cout << ">>> Z Mass constructor" << endl;
00033   Service<TFileService> fs;
00034   h_mZ_ = fs->make<TH1F>("ZMass", "Z mass (GeV/c^{2})", 100,  0, 200);
00035   h_mZMC_ = fs->make<TH1F>("ZMCMass", "Z MC mass (GeV/c^{2})", 100,  0, 200);
00036 }
00037 
00038 void ZMassHistogrammer::analyze(const edm::Event& event, const edm::EventSetup& setup) { 
00039   cout << ">>> Z Mass analyze" << endl;
00040   Handle<CandidateCollection> z;
00041   Handle<CandidateCollection> gen;
00042   event.getByLabel(z_, z);
00043   event.getByLabel(gen_, gen);
00044   for(unsigned int i = 0; i < z->size(); ++i) {
00045     const Candidate &zCand = (*z)[i];
00046     h_mZ_->Fill(zCand.mass());
00047   }
00048   for(unsigned int i = 0; i < gen->size(); ++i) {
00049     const Candidate &genCand = (*gen)[i];
00050     if((genCand.pdgId() == 23) && (genCand.status() == 2)) //this is an intermediate Z0
00051       cout << ">>> intermediate Z0 found, with " << genCand.numberOfDaughters() 
00052            << " daughters" << endl;
00053     if((genCand.pdgId() == 23)&&(genCand.status() == 3)) { //this is a Z0
00054       cout << ">>> Z0 found, with " << genCand.numberOfDaughters() 
00055            << " daughters" << endl;
00056       h_mZMC_->Fill(genCand.mass());
00057       if(genCand.numberOfDaughters() == 3) {//Z0 decays in mu+ mu-, the 3rd daughter is the same Z0
00058         const Candidate * dauGen0 = genCand.daughter(0);
00059         const Candidate * dauGen1 = genCand.daughter(1);
00060         const Candidate * dauGen2 = genCand.daughter(2);
00061         cout << ">>> daughter MC 0 PDG Id " << dauGen0->pdgId() 
00062              << ", status " << dauGen0->status() 
00063              << ", charge " << dauGen0->charge() 
00064              << endl;
00065         cout << ">>> daughter MC 1 PDG Id " << dauGen1->pdgId() 
00066              << ", status " << dauGen1->status()
00067              << ", charge " << dauGen1->charge() 
00068              << endl;
00069         cout << ">>> daughter MC 2 PDG Id " << dauGen2->pdgId() 
00070              << ", status " << dauGen2->status()
00071              << ", charge " << dauGen2->charge() << endl;
00072       }
00073     }
00074   }
00075 }
00076 
00077 #include "FWCore/Framework/interface/MakerMacros.h"
00078 
00079 DEFINE_FWK_MODULE(ZMassHistogrammer);
00080