CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_6/src/ElectroWeakAnalysis/ZMuMu/plugins/gamma_radiative_analysis.cc

Go to the documentation of this file.
00001 /* \class gamma_radiative_analyzer
00002  * 
00003  * author: Pasquale Noli
00004  *
00005  * Gamma Radiative analyzer
00006  * 
00007  *
00008  */
00009 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00010 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
00011 #include "DataFormats/MuonReco/interface/Muon.h"
00012 #include "DataFormats/Common/interface/AssociationVector.h"
00013 #include "DataFormats/Candidate/interface/CandMatchMap.h" 
00014 #include "DataFormats/Candidate/interface/Particle.h"
00015 #include "DataFormats/Candidate/interface/Candidate.h"
00016 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00017 #include "DataFormats/Candidate/interface/OverlapChecker.h"
00018 #include "DataFormats/PatCandidates/interface/Muon.h" 
00019 #include "DataFormats/PatCandidates/interface/Lepton.h" 
00020 #include "DataFormats/PatCandidates/interface/GenericParticle.h"
00021 #include "DataFormats/PatCandidates/interface/PATObject.h"
00022 #include "DataFormats/RecoCandidate/interface/IsoDepositVetos.h"
00023 #include "DataFormats/RecoCandidate/interface/IsoDepositDirection.h"
00024 #include "DataFormats/RecoCandidate/interface/IsoDeposit.h"
00025 #include "DataFormats/RecoCandidate/interface/IsoDepositFwd.h"
00026 #include "DataFormats/PatCandidates/interface/Isolation.h"
00027 #include "DataFormats/Common/interface/Handle.h"
00028 #include "DataFormats/Common/interface/ValueMap.h"
00029 #include "DataFormats/Candidate/interface/CandAssociation.h"
00030 #include "FWCore/Framework/interface/EDAnalyzer.h"
00031 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00032 #include "FWCore/Framework/interface/Event.h"
00033 #include "FWCore/Framework/interface/EventSetup.h"
00034 #include "FWCore/Utilities/interface/InputTag.h"
00035 #include "FWCore/ServiceRegistry/interface/Service.h"
00036 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00037 #include <iostream>
00038 #include <iterator>
00039 #include <cmath>
00040 #include <vector>
00041 #include "TH1.h"
00042 #include "TH2.h"
00043 //#include "TH3.h"
00044 
00045 
00046 using namespace edm;
00047 using namespace std;
00048 using namespace reco;
00049 using namespace isodeposit;
00050 
00051 class gamma_radiative_analyzer : public edm::EDAnalyzer {
00052 public:
00053   gamma_radiative_analyzer(const edm::ParameterSet& pset);
00054 private:
00055   virtual void analyze(const edm::Event& event, const edm::EventSetup& setup);
00056   virtual void endJob();
00057 
00058   edm::InputTag zMuMu_, zMuMuMatchMap_, zMuTk_, zMuTkMatchMap_,zMuSa_, zMuSaMatchMap_;
00059   double dRVeto_, dRTrk_, ptThreshold_;
00060   //histograms 
00061   TH2D *h_gamma_pt_eta_, *h_mu_pt_eta_FSR_, *h_mu_pt_eta_no_FSR_;
00062 
00063   //boolean 
00064   bool  FSR_mu, FSR_tk, FSR_mu0, FSR_mu1;
00065 
00066   //counter
00067   int  zmmcounter , zmscounter, zmtcounter,numOfEvent,numofGamma;
00068  };
00069 
00070 
00071 typedef edm::ValueMap<float> IsolationCollection;
00072 
00073 gamma_radiative_analyzer::gamma_radiative_analyzer(const ParameterSet& pset) : 
00074   zMuMu_(pset.getParameter<InputTag>("zMuMu")), 
00075   zMuMuMatchMap_(pset.getParameter<InputTag>("zMuMuMatchMap")),
00076   zMuTk_(pset.getParameter<InputTag>("zMuTk")), 
00077   zMuTkMatchMap_(pset.getParameter<InputTag>("zMuTkMatchMap")),
00078   zMuSa_(pset.getParameter<InputTag>("zMuSa")), 
00079   zMuSaMatchMap_(pset.getParameter<InputTag>("zMuSaMatchMap")){
00080   zmmcounter=0;
00081   zmscounter=0;
00082   zmtcounter=0;
00083   numOfEvent=0;
00084   numofGamma=0;
00085   Service<TFileService> fs;
00086    
00087   // general histograms
00088   h_gamma_pt_eta_= fs->make<TH2D>("h_gamma_pt_eta","pt vs eta of gamma",100,20,100,100,-2.0,2.0);
00089   h_mu_pt_eta_FSR_= fs->make<TH2D>("h_mu_pt_eta_FSR","pt vs eta of muon with FSR",100,20,100,100,-2.0,2.0 );
00090   h_mu_pt_eta_no_FSR_= fs->make<TH2D>("h_mu_pt_eta_no_FSR","pt vs eta of of muon withot FSR",100,20,100,100,-2.0,2.0);
00091 }
00092 
00093 void gamma_radiative_analyzer::analyze(const Event& event, const EventSetup& setup) {
00094   Handle<CandidateView> zMuMu;                //Collection of Z made by  Mu global + Mu global 
00095   Handle<GenParticleMatch> zMuMuMatchMap;     //Map of Z made by Mu global + Mu global with MC
00096   event.getByLabel(zMuMu_, zMuMu); 
00097   Handle<CandidateView> zMuTk;                //Collection of Z made by  Mu global + Track 
00098   Handle<GenParticleMatch> zMuTkMatchMap;
00099   event.getByLabel(zMuTk_, zMuTk); 
00100   Handle<CandidateView> zMuSa;                //Collection of Z made by  Mu global + Sa 
00101   Handle<GenParticleMatch> zMuSaMatchMap;
00102   event.getByLabel(zMuSa_, zMuSa); 
00103   numOfEvent++;
00104   // ZMuMu
00105   if (zMuMu->size() > 0 ) {
00106     event.getByLabel(zMuMuMatchMap_, zMuMuMatchMap); 
00107      for(unsigned int i = 0; i < zMuMu->size(); ++i) { //loop on candidates
00108      
00109       const Candidate & zMuMuCand = (*zMuMu)[i]; //the candidate
00110       CandidateBaseRef zMuMuCandRef = zMuMu->refAt(i);
00111 
00112        
00113       CandidateBaseRef dau0 = zMuMuCand.daughter(0)->masterClone();
00114       CandidateBaseRef dau1 = zMuMuCand.daughter(1)->masterClone();
00115       const pat::Muon& mu0 = dynamic_cast<const pat::Muon&>(*dau0);//cast in patMuon
00116       const pat::Muon& mu1 = dynamic_cast<const pat::Muon&>(*dau1);
00117             
00118       double zmass= zMuMuCand.mass();
00119       double pt0 = mu0.pt();
00120       double pt1 = mu1.pt();
00121       double eta0 = mu0.eta();
00122       double eta1 = mu1.eta();
00123       if(pt0>20 && pt1 > 20 && abs(eta0)<2 && abs(eta1)<2 && zmass > 20 && zmass < 200){
00124         GenParticleRef zMuMuMatch = (*zMuMuMatchMap)[zMuMuCandRef];
00125         if(zMuMuMatch.isNonnull()) {  // ZMuMu matched
00126           zmmcounter++;
00127           FSR_mu0 = false;
00128           FSR_mu1 = false;
00129           
00130           //MonteCarlo Study
00131           const reco::GenParticle * muMc0 = mu0.genLepton();
00132           const reco::GenParticle * muMc1 = mu1.genLepton();
00133           const Candidate * motherMu0 =  muMc0->mother();
00134           const Candidate * motherMu1 =  muMc1->mother();
00135           int num_dau_muon0 = motherMu0->numberOfDaughters();
00136           int num_dau_muon1 = motherMu1->numberOfDaughters();
00137           if( num_dau_muon0 > 1 ){
00138             for(int j = 0; j <  num_dau_muon0; ++j){
00139               int id =motherMu0 ->daughter(j)->pdgId();
00140               if(id == 22){
00141                 double etaG = motherMu0 ->daughter(j)->eta();
00142                 double ptG = motherMu0 ->daughter(j)->pt();
00143                 h_gamma_pt_eta_->Fill(ptG,etaG);
00144                 h_mu_pt_eta_FSR_->Fill(pt0,eta0);
00145                 FSR_mu0=true;
00146                 numofGamma++;
00147               }                          
00148             }
00149           }//end check of gamma
00150           if(!FSR_mu0)  h_mu_pt_eta_no_FSR_->Fill(pt0,eta0);
00151           if( num_dau_muon1 > 1 ){
00152             for(int j = 0; j <  num_dau_muon1; ++j){
00153               int id = motherMu1->daughter(j)->pdgId(); 
00154               if(id == 22){
00155                 double etaG = motherMu1 ->daughter(j)->eta();
00156                 double ptG = motherMu1 ->daughter(j)->pt();
00157                 h_gamma_pt_eta_->Fill(ptG,etaG);
00158                 h_mu_pt_eta_FSR_->Fill(pt1,eta1);
00159                 FSR_mu1=true;
00160                 numofGamma++;
00161               }                          
00162             }
00163           }//end check of gamma
00164           if(!FSR_mu1)  h_mu_pt_eta_no_FSR_->Fill(pt1,eta1);
00165         }// end MC match
00166       }//end of cuts
00167      }// end loop on ZMuMu cand
00168   }// end if ZMuMu size > 0
00169   
00170   // ZMuSa
00171   if (zMuSa->size() > 0 ) {
00172     event.getByLabel(zMuSaMatchMap_, zMuSaMatchMap); 
00173      for(unsigned int i = 0; i < zMuSa->size(); ++i) { //loop on candidates
00174      
00175       const Candidate & zMuSaCand = (*zMuSa)[i]; //the candidate
00176       CandidateBaseRef zMuSaCandRef = zMuSa->refAt(i);
00177 
00178        
00179       CandidateBaseRef dau0 = zMuSaCand.daughter(0)->masterClone();
00180       CandidateBaseRef dau1 = zMuSaCand.daughter(1)->masterClone();
00181       const pat::Muon& mu0 = dynamic_cast<const pat::Muon&>(*dau0);//cast in patMuon
00182       const pat::Muon& mu1 = dynamic_cast<const pat::Muon&>(*dau1);
00183             
00184       double zmass= zMuSaCand.mass();
00185       double pt0 = mu0.pt();
00186       double pt1 = mu1.pt();
00187       double eta0 = mu0.eta();
00188       double eta1 = mu1.eta();
00189       if(pt0>20 && pt1 > 20 && abs(eta0)<2 && abs(eta1)<2 && zmass > 20 && zmass < 200){
00190         GenParticleRef zMuSaMatch = (*zMuSaMatchMap)[zMuSaCandRef];
00191         if(zMuSaMatch.isNonnull()) {  // ZMuSa matched
00192           FSR_mu0 = false;
00193           FSR_mu1 = false;
00194           zmscounter++;
00195           //MonteCarlo Study
00196           const reco::GenParticle * muMc0 = mu0.genLepton();
00197           const reco::GenParticle * muMc1 = mu1.genLepton();
00198           const Candidate * motherMu0 =  muMc0->mother();
00199           const Candidate * motherMu1 =  muMc1->mother();
00200           int num_dau_muon0 = motherMu0->numberOfDaughters();
00201           int num_dau_muon1 = motherMu1->numberOfDaughters();
00202           if( num_dau_muon0 > 1 ){
00203             for(int j = 0; j <  num_dau_muon0; ++j){
00204               int id =motherMu0 ->daughter(j)->pdgId();
00205               if(id == 22){
00206                 double etaG = motherMu0 ->daughter(j)->eta();
00207                 double ptG = motherMu0 ->daughter(j)->pt();
00208                 h_gamma_pt_eta_->Fill(ptG,etaG);
00209                 h_mu_pt_eta_FSR_->Fill(pt0,eta0);
00210                 numofGamma++;
00211                 FSR_mu0=true;
00212               }                          
00213             }
00214           }//end check of gamma
00215           if(!FSR_mu0)  h_mu_pt_eta_no_FSR_->Fill(pt0,eta0);
00216           if( num_dau_muon1 > 1 ){
00217             for(int j = 0; j <  num_dau_muon1; ++j){
00218               int id = motherMu1->daughter(j)->pdgId(); 
00219               if(id == 22){
00220                 double etaG = motherMu1 ->daughter(j)->eta();
00221                 double ptG = motherMu1 ->daughter(j)->pt();
00222                 h_gamma_pt_eta_->Fill(ptG,etaG);
00223                 h_mu_pt_eta_FSR_->Fill(pt1,eta1);
00224                 numofGamma++;
00225                 FSR_mu1=true;
00226               } 
00227             }
00228           }//end check of gamma
00229           if(!FSR_mu1)  h_mu_pt_eta_no_FSR_->Fill(pt1,eta1);
00230         }// end MC match
00231       }//end of cuts
00232      }// end loop on ZMuSa cand
00233   }// end if ZMuSa size > 0
00234 
00235 
00236 
00237   //ZMuTk  
00238   if (zMuTk->size() > 0 ) {
00239     event.getByLabel(zMuTkMatchMap_, zMuTkMatchMap); 
00240     for(unsigned int i = 0; i < zMuTk->size(); ++i) { //loop on candidates
00241       const Candidate & zMuTkCand = (*zMuTk)[i]; //the candidate
00242       CandidateBaseRef zMuTkCandRef = zMuTk->refAt(i);
00243       
00244       
00245       CandidateBaseRef dau0 = zMuTkCand.daughter(0)->masterClone();
00246       CandidateBaseRef dau1 = zMuTkCand.daughter(1)->masterClone();
00247       const pat::Muon& mu0 = dynamic_cast<const pat::Muon&>(*dau0);//cast in patMuon
00248       const pat::GenericParticle& mu1 = dynamic_cast<const pat::GenericParticle &>(*dau1);
00249       
00250       
00251       double zmass= zMuTkCand.mass();
00252       double pt0 = mu0.pt();
00253       double pt1 = mu1.pt();
00254       double eta0 = mu0.eta();
00255       double eta1 = mu1.eta();
00256       if(pt0>20 && pt1 > 20 && abs(eta0)<2 && abs(eta1)<2 && zmass > 20 && zmass < 200){//kinematical cuts
00257         GenParticleRef zMuTkMatch = (*zMuTkMatchMap)[zMuTkCandRef];
00258         if(zMuTkMatch.isNonnull()) {  // ZMuTk matched
00259           FSR_mu = false;
00260           FSR_tk = false;
00261           zmtcounter++;
00262           //MonteCarlo Study
00263           const reco::GenParticle * muMc0 = mu0.genLepton();
00264           const reco::GenParticle * muMc1 = mu1.genParticle() ;
00265           const Candidate * motherMu0 =  muMc0->mother();
00266           const Candidate * motherMu1 =  muMc1->mother();
00267           int num_dau_muon0 = motherMu0->numberOfDaughters();
00268           int num_dau_muon1 = motherMu1->numberOfDaughters();
00269           if( num_dau_muon0 > 1 ){
00270             for(int j = 0; j <  num_dau_muon0; ++j){
00271               int id = motherMu0->daughter(j)->pdgId(); 
00272               if(id == 22){
00273                 double etaG = motherMu0 ->daughter(j)->eta();
00274                 double ptG = motherMu0 ->daughter(j)->pt();
00275                 h_gamma_pt_eta_->Fill(ptG,etaG);
00276                 h_mu_pt_eta_FSR_->Fill(pt0,eta0);
00277                 numofGamma++;
00278                 FSR_mu0=true;
00279               } 
00280             }
00281           }//end check of gamma
00282           if(!FSR_mu0)  h_mu_pt_eta_no_FSR_->Fill(pt0,eta0);
00283           if( num_dau_muon1 > 1 ){
00284             for(int j = 0; j <  num_dau_muon1; ++j){
00285               int id = motherMu1->daughter(j)->pdgId(); 
00286               if(id == 22){
00287                 double etaG = motherMu1 ->daughter(j)->eta();
00288                 double ptG = motherMu1 ->daughter(j)->pt();
00289                 h_gamma_pt_eta_->Fill(ptG,etaG);
00290                 h_mu_pt_eta_FSR_->Fill(pt1,eta1);
00291                 numofGamma++;
00292                 FSR_mu1=true;
00293               } 
00294             }
00295           }//end check of gamma
00296           if(!FSR_mu1)  h_mu_pt_eta_no_FSR_->Fill(pt1,eta1);
00297         }// end MC match
00298       }//end Kine-cuts
00299     }// end loop on ZMuTk cand
00300   }// end if ZMuTk size > 0
00301 }// end analyze
00302 
00303 
00304 
00305 void gamma_radiative_analyzer::endJob() {
00306   cout <<" ============= Summary =========="<<endl;
00307   cout <<" Numero di eventi = "<< numOfEvent << endl;
00308   cout <<" 1)Numero di ZMuMu matched dopo i tagli cinematici = "<< zmmcounter << endl;
00309   cout <<" 2)Numero di ZMuSa matched dopo i tagli cinematici = "<< zmscounter << endl;
00310   cout <<" 3)Numero di ZMuTk matched dopo i tagli cinematici = "<< zmtcounter << endl;  
00311   cout <<" 4)Number of gamma = "<< numofGamma << endl;
00312   }
00313   
00314 #include "FWCore/Framework/interface/MakerMacros.h"
00315 
00316 DEFINE_FWK_MODULE(gamma_radiative_analyzer);