CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/DQMOffline/Muon/src/MuonKinVsEtaAnalyzer.cc

Go to the documentation of this file.
00001 /*
00002  *  See header file for a description of this class.
00003  *
00004  *  $Date: 2011/07/14 13:27:35 $
00005  *  $Revision: 1.3 $
00006  *  \author S. Goy Lopez, CIEMAT 
00007  */
00008 
00009 #include "DQMOffline/Muon/src/MuonKinVsEtaAnalyzer.h"
00010 
00011 #include "DataFormats/Common/interface/Handle.h"
00012 #include "DataFormats/MuonReco/interface/Muon.h"
00013 #include "DataFormats/MuonReco/interface/MuonFwd.h" 
00014 #include "DataFormats/MuonReco/interface/MuonEnergy.h"
00015 
00016 #include "DataFormats/TrackReco/interface/Track.h"
00017 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00018 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00019 
00020 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00021 
00022 #include <string>
00023 using namespace std;
00024 using namespace edm;
00025 
00026 
00027 
00028 MuonKinVsEtaAnalyzer::MuonKinVsEtaAnalyzer(const edm::ParameterSet& pSet, MuonServiceProxy *theService):MuonAnalyzerBase(theService){  parameters = pSet;
00029 }
00030 
00031 
00032 MuonKinVsEtaAnalyzer::~MuonKinVsEtaAnalyzer() { }
00033 
00034 
00035 void MuonKinVsEtaAnalyzer::beginJob(DQMStore * dbe) {
00036 
00037   metname = "muonKinVsEta";
00038 
00039   LogTrace(metname)<<"[MuonKinVsEtaAnalyzer] Parameters initialization";
00040   dbe->setCurrentFolder("Muons/MuonKinVsEtaAnalyzer");
00041 
00042   etaBin = parameters.getParameter<int>("etaBin");
00043   etaMin = parameters.getParameter<double>("etaMin");
00044   etaMax = parameters.getParameter<double>("etaMax");
00045 
00046   phiBin = parameters.getParameter<int>("phiBin");
00047   phiMin = parameters.getParameter<double>("phiMin");
00048   phiMax = parameters.getParameter<double>("phiMax");
00049 
00050   pBin = parameters.getParameter<int>("pBin");
00051   pMin = parameters.getParameter<double>("pMin");
00052   pMax = parameters.getParameter<double>("pMax");
00053 
00054   ptBin = parameters.getParameter<int>("ptBin");
00055   ptMin = parameters.getParameter<double>("ptMin");
00056   ptMax = parameters.getParameter<double>("ptMax");
00057 
00058   etaBMin = parameters.getParameter<double>("etaBMin");
00059   etaBMax = parameters.getParameter<double>("etaBMax");
00060   etaECMin = parameters.getParameter<double>("etaECMin");
00061   etaECMax = parameters.getParameter<double>("etaECMax");
00062   etaOvlpMin = parameters.getParameter<double>("etaOvlpMin");
00063   etaOvlpMax = parameters.getParameter<double>("etaOvlpMax");
00064 
00065   std::string EtaName;
00066   for(unsigned int iEtaRegion=0;iEtaRegion<4;iEtaRegion++){
00067     if (iEtaRegion==0) EtaName = "Barrel_";   
00068     if (iEtaRegion==1) EtaName = "EndCap_";   
00069     if (iEtaRegion==2) EtaName = "Overlap_";
00070     if (iEtaRegion==3) EtaName = "_";
00071 
00072     // monitoring of eta parameter
00073     etaGlbTrack.push_back(dbe->book1D("GlbMuon_eta_"+EtaName, "#eta_{GLB} "+EtaName, etaBin, etaMin, etaMax));
00074     etaTrack.push_back(dbe->book1D("TkMuon_eta_"+EtaName, "#eta_{TK} "+EtaName, etaBin, etaMin, etaMax));
00075     etaStaTrack.push_back(dbe->book1D("StaMuon_eta_"+EtaName, "#eta_{STA} "+EtaName, etaBin, etaMin, etaMax));
00076     etaTightTrack.push_back(dbe->book1D("TightMuon_eta_"+EtaName, "#eta_{Tight} "+EtaName, etaBin, etaMin, etaMax));
00077     
00078     // monitoring of phi paramater
00079     phiGlbTrack.push_back(dbe->book1D("GlbMuon_phi_"+EtaName, "#phi_{GLB} "+EtaName+ "(rad)", phiBin, phiMin, phiMax));
00080     phiTrack.push_back(dbe->book1D("TkMuon_phi_"+EtaName, "#phi_{TK}" +EtaName +"(rad)", phiBin, phiMin, phiMax));
00081     phiStaTrack.push_back(dbe->book1D("StaMuon_phi_"+EtaName, "#phi_{STA}"+EtaName+" (rad)", phiBin, phiMin, phiMax));
00082     phiTightTrack.push_back(dbe->book1D("TightMuon_phi_"+EtaName, "#phi_{Tight}_"+EtaName, phiBin, phiMin, phiMax));
00083 
00084     // monitoring of the momentum
00085     pGlbTrack.push_back(dbe->book1D("GlbMuon_p_"+EtaName, "p_{GLB} "+EtaName, pBin, pMin, pMax));
00086     pTrack.push_back(dbe->book1D("TkMuon_p"+EtaName, "p_{TK} "+EtaName, pBin, pMin, pMax));
00087     pStaTrack.push_back(dbe->book1D("StaMuon_p"+EtaName, "p_{STA} "+EtaName, pBin, pMin, pMax));
00088     pTightTrack.push_back(dbe->book1D("TightMuon_p_"+EtaName, "p_{Tight} "+EtaName, pBin, pMin, pMax));
00089 
00090     // monitoring of the transverse momentum
00091     ptGlbTrack.push_back(dbe->book1D("GlbMuon_pt_" +EtaName, "pt_{GLB} "+EtaName, ptBin, ptMin, ptMax));
00092     ptTrack.push_back(dbe->book1D("TkMuon_pt_"+EtaName, "pt_{TK} "+EtaName, ptBin, ptMin, ptMax));
00093     ptStaTrack.push_back(dbe->book1D("StaMuon_pt_"+EtaName, "pt_{STA} "+EtaName, ptBin, ptMin, pMax));
00094     ptTightTrack.push_back(dbe->book1D("TightMuon_pt_"+EtaName, "pt_{Tight} "+EtaName, ptBin, ptMin, ptMax));
00095   }
00096 }
00097 
00098 
00099 
00100 
00101 void MuonKinVsEtaAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup, const reco::Muon& recoMu) {
00102   
00103   reco::BeamSpot beamSpot;
00104   Handle<reco::BeamSpot> beamSpotHandle;
00105   iEvent.getByLabel("offlineBeamSpot", beamSpotHandle);
00106   beamSpot = *beamSpotHandle;
00107   
00108   LogTrace(metname)<<"[MuonKinVsEtaAnalyzer] Analyze the mu in different eta regions";
00109 
00110   for(unsigned int iEtaRegion=0;iEtaRegion<4;iEtaRegion++){
00111     if (iEtaRegion==0) {EtaCutMin= etaBMin; EtaCutMax=etaBMax;}
00112     if (iEtaRegion==1) {EtaCutMin= etaECMin; EtaCutMax=etaECMax;}
00113     if (iEtaRegion==2) {EtaCutMin= etaOvlpMin; EtaCutMax=etaOvlpMax;} 
00114     if (iEtaRegion==3) {EtaCutMin= etaBMin; EtaCutMax=etaECMax;}
00115     
00116     if(recoMu.isGlobalMuon()) {
00117       LogTrace(metname)<<"[MuonKinVsEtaAnalyzer] The mu is global - filling the histos";
00118       reco::TrackRef recoCombinedGlbTrack = recoMu.combinedMuon();
00119       // get the track combinig the information from both the glb fit"
00120       if(fabs(recoCombinedGlbTrack->eta())>EtaCutMin && fabs(recoCombinedGlbTrack->eta())<EtaCutMax){
00121         etaGlbTrack[iEtaRegion]->Fill(recoCombinedGlbTrack->eta());
00122         phiGlbTrack[iEtaRegion]->Fill(recoCombinedGlbTrack->phi());
00123         pGlbTrack[iEtaRegion]->Fill(recoCombinedGlbTrack->p());
00124         ptGlbTrack[iEtaRegion]->Fill(recoCombinedGlbTrack->pt());
00125       }
00126     }
00127 
00128     if(recoMu.isTrackerMuon()) {
00129       LogTrace(metname)<<"[MuonKinVsEtaAnalyzer] The mu is tracker - filling the histos";
00130       // get the track using only the tracker data
00131       reco::TrackRef recoTrack = recoMu.track();
00132       if(fabs(recoTrack->eta())>EtaCutMin && fabs(recoTrack->eta())<EtaCutMax){
00133         etaTrack[iEtaRegion]->Fill(recoTrack->eta());
00134         phiTrack[iEtaRegion]->Fill(recoTrack->phi());
00135         pTrack[iEtaRegion]->Fill(recoTrack->p());
00136         ptTrack[iEtaRegion]->Fill(recoTrack->pt());
00137       }
00138     }
00139 
00140     if(recoMu.isStandAloneMuon()) {
00141       LogTrace(metname)<<"[MuonKinVsEtaAnalyzer] The mu is standalone - filling the histos";
00142       // get the track using only the mu spectrometer data
00143       reco::TrackRef recoStaTrack = recoMu.standAloneMuon();
00144       if(fabs(recoStaTrack->eta())>EtaCutMin && fabs(recoStaTrack->eta())<EtaCutMax){
00145         etaStaTrack[iEtaRegion]->Fill(recoStaTrack->eta());
00146         phiStaTrack[iEtaRegion]->Fill(recoStaTrack->phi());
00147         pStaTrack[iEtaRegion]->Fill(recoStaTrack->p());
00148         ptStaTrack[iEtaRegion]->Fill(recoStaTrack->pt());
00149       }
00150     }
00151     if (recoMu.isGlobalMuon() && recoMu.isTrackerMuon() && recoMu.combinedMuon()->normalizedChi2()<10. 
00152         && recoMu.combinedMuon()->hitPattern().numberOfValidMuonHits()>0 && fabs(recoMu.combinedMuon()->dxy(beamSpot.position()))<0.2 
00153         && recoMu.combinedMuon()->hitPattern().numberOfValidPixelHits()>0 && recoMu.numberOfMatches() > 1) {
00154       LogTrace(metname)<<"[MuonKinVsEtaAnalyzer] The mu is Tight - filling the histos";
00155       reco::TrackRef recoTightTrack = recoMu.combinedMuon();
00156       if(fabs(recoTightTrack->eta())>EtaCutMin && fabs(recoTightTrack->eta())<EtaCutMax){
00157         etaTightTrack[iEtaRegion]->Fill(recoTightTrack->eta());
00158         phiTightTrack[iEtaRegion]->Fill(recoTightTrack->phi());
00159         pTightTrack[iEtaRegion]->Fill(recoTightTrack->p());
00160         ptTightTrack[iEtaRegion]->Fill(recoTightTrack->pt());
00161       }
00162     }
00163   }
00164 }