CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MuonKinVsEtaAnalyzer.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * $Date: 2011/07/14 13:27:35 $
5  * $Revision: 1.3 $
6  * \author S. Goy Lopez, CIEMAT
7  */
8 
10 
15 
19 
21 
22 #include <string>
23 using namespace std;
24 using namespace edm;
25 
26 
27 
29 }
30 
31 
33 
34 
36 
37  metname = "muonKinVsEta";
38 
39  LogTrace(metname)<<"[MuonKinVsEtaAnalyzer] Parameters initialization";
40  dbe->setCurrentFolder("Muons/MuonKinVsEtaAnalyzer");
41 
42  etaBin = parameters.getParameter<int>("etaBin");
43  etaMin = parameters.getParameter<double>("etaMin");
44  etaMax = parameters.getParameter<double>("etaMax");
45 
46  phiBin = parameters.getParameter<int>("phiBin");
47  phiMin = parameters.getParameter<double>("phiMin");
48  phiMax = parameters.getParameter<double>("phiMax");
49 
50  pBin = parameters.getParameter<int>("pBin");
51  pMin = parameters.getParameter<double>("pMin");
52  pMax = parameters.getParameter<double>("pMax");
53 
54  ptBin = parameters.getParameter<int>("ptBin");
55  ptMin = parameters.getParameter<double>("ptMin");
56  ptMax = parameters.getParameter<double>("ptMax");
57 
58  etaBMin = parameters.getParameter<double>("etaBMin");
59  etaBMax = parameters.getParameter<double>("etaBMax");
60  etaECMin = parameters.getParameter<double>("etaECMin");
61  etaECMax = parameters.getParameter<double>("etaECMax");
62  etaOvlpMin = parameters.getParameter<double>("etaOvlpMin");
63  etaOvlpMax = parameters.getParameter<double>("etaOvlpMax");
64 
65  std::string EtaName;
66  for(unsigned int iEtaRegion=0;iEtaRegion<4;iEtaRegion++){
67  if (iEtaRegion==0) EtaName = "Barrel_";
68  if (iEtaRegion==1) EtaName = "EndCap_";
69  if (iEtaRegion==2) EtaName = "Overlap_";
70  if (iEtaRegion==3) EtaName = "_";
71 
72  // monitoring of eta parameter
73  etaGlbTrack.push_back(dbe->book1D("GlbMuon_eta_"+EtaName, "#eta_{GLB} "+EtaName, etaBin, etaMin, etaMax));
74  etaTrack.push_back(dbe->book1D("TkMuon_eta_"+EtaName, "#eta_{TK} "+EtaName, etaBin, etaMin, etaMax));
75  etaStaTrack.push_back(dbe->book1D("StaMuon_eta_"+EtaName, "#eta_{STA} "+EtaName, etaBin, etaMin, etaMax));
76  etaTightTrack.push_back(dbe->book1D("TightMuon_eta_"+EtaName, "#eta_{Tight} "+EtaName, etaBin, etaMin, etaMax));
77 
78  // monitoring of phi paramater
79  phiGlbTrack.push_back(dbe->book1D("GlbMuon_phi_"+EtaName, "#phi_{GLB} "+EtaName+ "(rad)", phiBin, phiMin, phiMax));
80  phiTrack.push_back(dbe->book1D("TkMuon_phi_"+EtaName, "#phi_{TK}" +EtaName +"(rad)", phiBin, phiMin, phiMax));
81  phiStaTrack.push_back(dbe->book1D("StaMuon_phi_"+EtaName, "#phi_{STA}"+EtaName+" (rad)", phiBin, phiMin, phiMax));
82  phiTightTrack.push_back(dbe->book1D("TightMuon_phi_"+EtaName, "#phi_{Tight}_"+EtaName, phiBin, phiMin, phiMax));
83 
84  // monitoring of the momentum
85  pGlbTrack.push_back(dbe->book1D("GlbMuon_p_"+EtaName, "p_{GLB} "+EtaName, pBin, pMin, pMax));
86  pTrack.push_back(dbe->book1D("TkMuon_p"+EtaName, "p_{TK} "+EtaName, pBin, pMin, pMax));
87  pStaTrack.push_back(dbe->book1D("StaMuon_p"+EtaName, "p_{STA} "+EtaName, pBin, pMin, pMax));
88  pTightTrack.push_back(dbe->book1D("TightMuon_p_"+EtaName, "p_{Tight} "+EtaName, pBin, pMin, pMax));
89 
90  // monitoring of the transverse momentum
91  ptGlbTrack.push_back(dbe->book1D("GlbMuon_pt_" +EtaName, "pt_{GLB} "+EtaName, ptBin, ptMin, ptMax));
92  ptTrack.push_back(dbe->book1D("TkMuon_pt_"+EtaName, "pt_{TK} "+EtaName, ptBin, ptMin, ptMax));
93  ptStaTrack.push_back(dbe->book1D("StaMuon_pt_"+EtaName, "pt_{STA} "+EtaName, ptBin, ptMin, pMax));
94  ptTightTrack.push_back(dbe->book1D("TightMuon_pt_"+EtaName, "pt_{Tight} "+EtaName, ptBin, ptMin, ptMax));
95  }
96 }
97 
98 
99 
100 
101 void MuonKinVsEtaAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup, const reco::Muon& recoMu) {
102 
104  Handle<reco::BeamSpot> beamSpotHandle;
105  iEvent.getByLabel("offlineBeamSpot", beamSpotHandle);
106  beamSpot = *beamSpotHandle;
107 
108  LogTrace(metname)<<"[MuonKinVsEtaAnalyzer] Analyze the mu in different eta regions";
109 
110  for(unsigned int iEtaRegion=0;iEtaRegion<4;iEtaRegion++){
111  if (iEtaRegion==0) {EtaCutMin= etaBMin; EtaCutMax=etaBMax;}
112  if (iEtaRegion==1) {EtaCutMin= etaECMin; EtaCutMax=etaECMax;}
113  if (iEtaRegion==2) {EtaCutMin= etaOvlpMin; EtaCutMax=etaOvlpMax;}
114  if (iEtaRegion==3) {EtaCutMin= etaBMin; EtaCutMax=etaECMax;}
115 
116  if(recoMu.isGlobalMuon()) {
117  LogTrace(metname)<<"[MuonKinVsEtaAnalyzer] The mu is global - filling the histos";
118  reco::TrackRef recoCombinedGlbTrack = recoMu.combinedMuon();
119  // get the track combinig the information from both the glb fit"
120  if(fabs(recoCombinedGlbTrack->eta())>EtaCutMin && fabs(recoCombinedGlbTrack->eta())<EtaCutMax){
121  etaGlbTrack[iEtaRegion]->Fill(recoCombinedGlbTrack->eta());
122  phiGlbTrack[iEtaRegion]->Fill(recoCombinedGlbTrack->phi());
123  pGlbTrack[iEtaRegion]->Fill(recoCombinedGlbTrack->p());
124  ptGlbTrack[iEtaRegion]->Fill(recoCombinedGlbTrack->pt());
125  }
126  }
127 
128  if(recoMu.isTrackerMuon()) {
129  LogTrace(metname)<<"[MuonKinVsEtaAnalyzer] The mu is tracker - filling the histos";
130  // get the track using only the tracker data
131  reco::TrackRef recoTrack = recoMu.track();
132  if(fabs(recoTrack->eta())>EtaCutMin && fabs(recoTrack->eta())<EtaCutMax){
133  etaTrack[iEtaRegion]->Fill(recoTrack->eta());
134  phiTrack[iEtaRegion]->Fill(recoTrack->phi());
135  pTrack[iEtaRegion]->Fill(recoTrack->p());
136  ptTrack[iEtaRegion]->Fill(recoTrack->pt());
137  }
138  }
139 
140  if(recoMu.isStandAloneMuon()) {
141  LogTrace(metname)<<"[MuonKinVsEtaAnalyzer] The mu is standalone - filling the histos";
142  // get the track using only the mu spectrometer data
143  reco::TrackRef recoStaTrack = recoMu.standAloneMuon();
144  if(fabs(recoStaTrack->eta())>EtaCutMin && fabs(recoStaTrack->eta())<EtaCutMax){
145  etaStaTrack[iEtaRegion]->Fill(recoStaTrack->eta());
146  phiStaTrack[iEtaRegion]->Fill(recoStaTrack->phi());
147  pStaTrack[iEtaRegion]->Fill(recoStaTrack->p());
148  ptStaTrack[iEtaRegion]->Fill(recoStaTrack->pt());
149  }
150  }
151  if (recoMu.isGlobalMuon() && recoMu.isTrackerMuon() && recoMu.combinedMuon()->normalizedChi2()<10.
152  && recoMu.combinedMuon()->hitPattern().numberOfValidMuonHits()>0 && fabs(recoMu.combinedMuon()->dxy(beamSpot.position()))<0.2
153  && recoMu.combinedMuon()->hitPattern().numberOfValidPixelHits()>0 && recoMu.numberOfMatches() > 1) {
154  LogTrace(metname)<<"[MuonKinVsEtaAnalyzer] The mu is Tight - filling the histos";
155  reco::TrackRef recoTightTrack = recoMu.combinedMuon();
156  if(fabs(recoTightTrack->eta())>EtaCutMin && fabs(recoTightTrack->eta())<EtaCutMax){
157  etaTightTrack[iEtaRegion]->Fill(recoTightTrack->eta());
158  phiTightTrack[iEtaRegion]->Fill(recoTightTrack->phi());
159  pTightTrack[iEtaRegion]->Fill(recoTightTrack->p());
160  ptTightTrack[iEtaRegion]->Fill(recoTightTrack->pt());
161  }
162  }
163  }
164 }
T getParameter(std::string const &) const
std::vector< MonitorElement * > etaTrack
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:717
std::vector< MonitorElement * > phiStaTrack
bool isTrackerMuon() const
Definition: Muon.h:212
virtual TrackRef track() const
reference to a Track
Definition: Muon.h:50
bool isGlobalMuon() const
Definition: Muon.h:211
bool isStandAloneMuon() const
Definition: Muon.h:213
void analyze(const edm::Event &, const edm::EventSetup &, const reco::Muon &recoMu)
Get the analysis.
std::vector< MonitorElement * > ptTrack
int iEvent
Definition: GenABIO.cc:243
void beginJob(DQMStore *dbe)
Inizialize parameters for histo binning.
MuonKinVsEtaAnalyzer(const edm::ParameterSet &, MuonServiceProxy *theService)
Constructor.
std::vector< MonitorElement * > pTrack
std::vector< MonitorElement * > ptTightTrack
std::vector< MonitorElement * > etaTightTrack
edm::ParameterSet parameters
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
#define LogTrace(id)
int numberOfMatches(ArbitrationType type=SegmentAndTrackArbitration) const
get number of chambers with matched segments
Definition: Muon.cc:56
std::vector< MonitorElement * > phiTightTrack
virtual TrackRef combinedMuon() const
reference to a stand-alone muon Track
Definition: Muon.h:56
std::vector< MonitorElement * > pTightTrack
std::vector< MonitorElement * > etaStaTrack
std::vector< MonitorElement * > pGlbTrack
std::vector< MonitorElement * > pStaTrack
std::vector< MonitorElement * > etaGlbTrack
const Point & position() const
position
Definition: BeamSpot.h:63
std::vector< MonitorElement * > phiTrack
std::vector< MonitorElement * > ptStaTrack
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:429
std::vector< MonitorElement * > ptGlbTrack
std::vector< MonitorElement * > phiGlbTrack
virtual ~MuonKinVsEtaAnalyzer()
Destructor.
virtual TrackRef standAloneMuon() const
reference to a stand-alone muon Track
Definition: Muon.h:53