CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DiMuonHistograms.cc
Go to the documentation of this file.
1 /* This Class Header */
3 
4 /* Collaborating Class Header */
16 
19 
20 #include "TLorentzVector.h"
21 #include "TFile.h"
22 #include <vector>
23 #include <cmath>
24 
25 /* C++ Headers */
26 #include <iostream>
27 #include <fstream>
28 #include <cmath>
29 using namespace std;
30 using namespace edm;
31 
33 
34  // initialise parameters:
35  parameters = pSet;
36 
37  theDbe = edm::Service<DQMStore>().operator->();
38 
39  // declare consumes:
40  theMuonCollectionLabel_ = consumes<reco::MuonCollection> (parameters.getParameter<edm::InputTag>("MuonCollection"));
41  theBeamSpotLabel_ = mayConsume<reco::BeamSpot> (parameters.getParameter<edm::InputTag>("BeamSpotLabel"));
42  theVertexLabel_ = consumes<reco::VertexCollection>(parameters.getParameter<edm::InputTag>("VertexLabel"));
43 }
44 
46 
48  metname = "DiMuonhistograms";
49  LogTrace(metname)<<"[DiMuonHistograms] beginJob()";
50 
51 }
52 
53 void DiMuonHistograms::beginRun(const edm::Run& iRun, const edm::EventSetup& iSetup) {
54  LogTrace(metname)<<"[DiMuonHistograms] beginRun()";
55 
56  theDbe->cd();
57  theDbe->setCurrentFolder("Muons/DiMuonHistograms");
58 
59  etaBin = parameters.getParameter<int>("etaBin");
60  etaBBin = parameters.getParameter<int>("etaBBin");
61  etaEBin = parameters.getParameter<int>("etaEBin");
62 
63  etaBMin = parameters.getParameter<double>("etaBMin");
64  etaBMax = parameters.getParameter<double>("etaBMax");
65  etaECMin = parameters.getParameter<double>("etaECMin");
66  etaECMax = parameters.getParameter<double>("etaECMax");
67 
68  LowMassMin = parameters.getParameter<double>("LowMassMin");
69  LowMassMax = parameters.getParameter<double>("LowMassMax");
70  HighMassMin = parameters.getParameter<double>("HighMassMin");
71  HighMassMax = parameters.getParameter<double>("HighMassMax");
72 
73  int nBin = 0;
74  for (unsigned int iEtaRegion=0; iEtaRegion<3; iEtaRegion++){
75  if (iEtaRegion==0) { EtaName = ""; nBin = etaBin;}
76  if (iEtaRegion==1) { EtaName = "_Barrel"; nBin = etaBBin;}
77  if (iEtaRegion==2) { EtaName = "_EndCap"; nBin = etaEBin;}
78 
79  GlbGlbMuon_LM.push_back(theDbe->book1D("GlbGlbMuon_LM"+EtaName,"InvMass_{GLB,GLB}"+EtaName,nBin, LowMassMin, LowMassMax));
80  TrkTrkMuon_LM.push_back(theDbe->book1D("TrkTrkMuon_LM"+EtaName,"InvMass_{TRK,TRK}"+EtaName,nBin, LowMassMin, LowMassMax));
81  StaTrkMuon_LM.push_back(theDbe->book1D("StaTrkMuon_LM"+EtaName,"InvMass_{STA,TRK}"+EtaName,nBin, LowMassMin, LowMassMax));
82 
83  GlbGlbMuon_HM.push_back(theDbe->book1D("GlbGlbMuon_HM"+EtaName,"InvMass_{GLB,GLB}"+EtaName,nBin, HighMassMin, HighMassMax));
84  TrkTrkMuon_HM.push_back(theDbe->book1D("TrkTrkMuon_HM"+EtaName,"InvMass_{TRK,TRK}"+EtaName,nBin, HighMassMin, HighMassMax));
85  StaTrkMuon_HM.push_back(theDbe->book1D("StaTrkMuon_HM"+EtaName,"InvMass_{STA,TRK}"+EtaName,nBin, HighMassMin, HighMassMax));
86 
87  // arround the Z peak
88  TightTightMuon.push_back(theDbe->book1D("TightTightMuon"+EtaName,"InvMass_{Tight,Tight}"+EtaName,nBin, 55.0, 125.0));
89 
90  // low-mass resonances
91  SoftSoftMuon.push_back(theDbe->book1D("SoftSoftMuon"+EtaName,"InvMass_{Soft,Soft}"+EtaName,nBin, 5.0, 55.0));
92  }
93 }
94 
96 
97  LogTrace(metname)<<"[DiMuonHistograms] Analyze the mu in different eta regions";
99  iEvent.getByToken(theMuonCollectionLabel_, muons);
100 
101  // =================================================================================
102  // Look for the Primary Vertex (and use the BeamSpot instead, if you can't find it):
103  reco::Vertex::Point posVtx;
104  reco::Vertex::Error errVtx;
105  unsigned int theIndexOfThePrimaryVertex = 999.;
106 
108  iEvent.getByToken(theVertexLabel_, vertex);
109  if (vertex.isValid()){
110  for (unsigned int ind=0; ind<vertex->size(); ++ind) {
111  if ( (*vertex)[ind].isValid() && !((*vertex)[ind].isFake()) ) {
112  theIndexOfThePrimaryVertex = ind;
113  break;
114  }
115  }
116  }
117 
118  if (theIndexOfThePrimaryVertex<100) {
119  posVtx = ((*vertex)[theIndexOfThePrimaryVertex]).position();
120  errVtx = ((*vertex)[theIndexOfThePrimaryVertex]).error();
121  }
122  else {
123  LogInfo("RecoMuonValidator") << "reco::PrimaryVertex not found, use BeamSpot position instead\n";
124 
125  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
126  iEvent.getByToken(theBeamSpotLabel_,recoBeamSpotHandle);
127  reco::BeamSpot bs = *recoBeamSpotHandle;
128 
129  posVtx = bs.position();
130  errVtx(0,0) = bs.BeamWidthX();
131  errVtx(1,1) = bs.BeamWidthY();
132  errVtx(2,2) = bs.sigmaZ();
133  }
134 
135  const reco::Vertex vtx(posVtx,errVtx);
136 
137  if(!muons.isValid()) return;
138 
139  // Loop on muon collection
140  TLorentzVector Mu1, Mu2;
141  float charge = 99.;
142  float InvMass = -99.;
143 
144  for (reco::MuonCollection::const_iterator recoMu1 = muons->begin(); recoMu1!=muons->end(); ++recoMu1) {
145  LogTrace(metname)<<"[DiMuonHistograms] loop over 1st muon"<<endl;
146 
147  // Loop on second muons to fill invariant mass plots
148  for (reco::MuonCollection::const_iterator recoMu2 = recoMu1; recoMu2!=muons->end(); ++recoMu2){
149  LogTrace(metname)<<"[DiMuonHistograms] loop over 2nd muon"<<endl;
150  if (recoMu1==recoMu2) continue;
151 
152  // Global-Global Muon
153  if (recoMu1->isGlobalMuon() && recoMu2->isGlobalMuon()) {
154  LogTrace(metname)<<"[DiMuonHistograms] Glb-Glb pair"<<endl;
155  reco::TrackRef recoCombinedGlbTrack1 = recoMu1->combinedMuon();
156  reco::TrackRef recoCombinedGlbTrack2 = recoMu2->combinedMuon();
157  Mu1.SetPxPyPzE(recoCombinedGlbTrack1->px(), recoCombinedGlbTrack1->py(),recoCombinedGlbTrack1->pz(), recoCombinedGlbTrack1->p());
158  Mu2.SetPxPyPzE(recoCombinedGlbTrack2->px(), recoCombinedGlbTrack2->py(),recoCombinedGlbTrack2->pz(), recoCombinedGlbTrack2->p());
159 
160  charge = recoCombinedGlbTrack1->charge()*recoCombinedGlbTrack2->charge();
161  if (charge < 0) {
162  InvMass = (Mu1+Mu2).M();
163  for (unsigned int iEtaRegion=0; iEtaRegion<3; iEtaRegion++){
164  if (iEtaRegion==0) {EtaCutMin= 0.; EtaCutMax=2.4; }
165  if (iEtaRegion==1) {EtaCutMin= etaBMin; EtaCutMax=etaBMax; }
166  if (iEtaRegion==2) {EtaCutMin= etaECMin; EtaCutMax=etaECMax; }
167 
168  if(fabs(recoCombinedGlbTrack1->eta())>EtaCutMin && fabs(recoCombinedGlbTrack1->eta())<EtaCutMax &&
169  fabs(recoCombinedGlbTrack2->eta())>EtaCutMin && fabs(recoCombinedGlbTrack2->eta())<EtaCutMax){
170  if (InvMass < LowMassMax) GlbGlbMuon_LM[iEtaRegion]->Fill(InvMass);
171  if (InvMass > HighMassMin) GlbGlbMuon_HM[iEtaRegion]->Fill(InvMass);
172  }
173  }
174  }
175  // Also Tight-Tight Muon Selection
176 
177  if ( muon::isTightMuon(*recoMu1, vtx) &&
178  muon::isTightMuon(*recoMu2, vtx) ) {
179 
180  LogTrace(metname)<<"[DiMuonHistograms] Tight-Tight pair"<<endl;
181  for (unsigned int iEtaRegion=0; iEtaRegion<3; iEtaRegion++){
182  if (iEtaRegion==0) {EtaCutMin= 0.; EtaCutMax=2.4; }
183  if (iEtaRegion==1) {EtaCutMin= etaBMin; EtaCutMax=etaBMax; }
184  if (iEtaRegion==2) {EtaCutMin= etaECMin; EtaCutMax=etaECMax; }
185 
186  if(fabs(recoCombinedGlbTrack1->eta())>EtaCutMin && fabs(recoCombinedGlbTrack1->eta())<EtaCutMax &&
187  fabs(recoCombinedGlbTrack2->eta())>EtaCutMin && fabs(recoCombinedGlbTrack2->eta())<EtaCutMax){
188  if (InvMass > 55. && InvMass < 125.) TightTightMuon[iEtaRegion]->Fill(InvMass);
189  }
190  }
191  }
192  }
193 
194  // Now check for STA-TRK
195  if (recoMu2->isStandAloneMuon() && recoMu1->isTrackerMuon()) {
196  LogTrace(metname)<<"[DiMuonHistograms] STA-Trk pair"<<endl;
197  reco::TrackRef recoStaTrack = recoMu2->standAloneMuon();
198  reco::TrackRef recoTrack = recoMu1->track();
199  Mu2.SetPxPyPzE(recoStaTrack->px(), recoStaTrack->py(),recoStaTrack->pz(), recoStaTrack->p());
200  Mu1.SetPxPyPzE(recoTrack->px(), recoTrack->py(),recoTrack->pz(), recoTrack->p());
201 
202  charge = recoStaTrack->charge()*recoTrack->charge();
203  if (charge < 0) {
204  InvMass = (Mu1+Mu2).M();
205  for (unsigned int iEtaRegion=0; iEtaRegion<3; iEtaRegion++){
206  if (iEtaRegion==0) {EtaCutMin= 0.; EtaCutMax=2.4; }
207  if (iEtaRegion==1) {EtaCutMin= etaBMin; EtaCutMax=etaBMax; }
208  if (iEtaRegion==2) {EtaCutMin= etaECMin; EtaCutMax=etaECMax; }
209 
210  if(fabs(recoStaTrack->eta())>EtaCutMin && fabs(recoStaTrack->eta())<EtaCutMax &&
211  fabs(recoTrack->eta())>EtaCutMin && fabs(recoTrack->eta())<EtaCutMax){
212  if (InvMass < LowMassMax) StaTrkMuon_LM[iEtaRegion]->Fill(InvMass);
213  if (InvMass > HighMassMin) StaTrkMuon_HM[iEtaRegion]->Fill(InvMass);
214  }
215  }
216  }
217  }
218  if (recoMu1->isStandAloneMuon() && recoMu2->isTrackerMuon()) {
219  LogTrace(metname)<<"[DiMuonHistograms] STA-Trk pair"<<endl;
220  reco::TrackRef recoStaTrack = recoMu1->standAloneMuon();
221  reco::TrackRef recoTrack = recoMu2->track();
222  Mu1.SetPxPyPzE(recoStaTrack->px(), recoStaTrack->py(),recoStaTrack->pz(), recoStaTrack->p());
223  Mu2.SetPxPyPzE(recoTrack->px(), recoTrack->py(),recoTrack->pz(), recoTrack->p());
224 
225  charge = recoStaTrack->charge()*recoTrack->charge();
226  if (charge < 0) {
227  InvMass = (Mu1+Mu2).M();
228  for (unsigned int iEtaRegion=0; iEtaRegion<3; iEtaRegion++){
229  if (iEtaRegion==0) {EtaCutMin= 0.; EtaCutMax=2.4; }
230  if (iEtaRegion==1) {EtaCutMin= etaBMin; EtaCutMax=etaBMax; }
231  if (iEtaRegion==2) {EtaCutMin= etaECMin; EtaCutMax=etaECMax; }
232 
233  if(fabs(recoStaTrack->eta())>EtaCutMin && fabs(recoStaTrack->eta())<EtaCutMax &&
234  fabs(recoTrack->eta())>EtaCutMin && fabs(recoTrack->eta())<EtaCutMax){
235  if (InvMass < LowMassMax) StaTrkMuon_LM[iEtaRegion]->Fill(InvMass);
236  if (InvMass > HighMassMin) StaTrkMuon_HM[iEtaRegion]->Fill(InvMass);
237  }
238  }
239  }
240  }
241 
242  // TRK-TRK dimuon
243  if (recoMu1->isTrackerMuon() && recoMu2->isTrackerMuon()) {
244  LogTrace(metname)<<"[DiMuonHistograms] Trk-Trk dimuon pair"<<endl;
245  reco::TrackRef recoTrack2 = recoMu2->track();
246  reco::TrackRef recoTrack1 = recoMu1->track();
247  Mu2.SetPxPyPzE(recoTrack2->px(), recoTrack2->py(),recoTrack2->pz(), recoTrack2->p());
248  Mu1.SetPxPyPzE(recoTrack1->px(), recoTrack1->py(),recoTrack1->pz(), recoTrack1->p());
249 
250  charge = recoTrack1->charge()*recoTrack2->charge();
251  if (charge < 0) {
252  InvMass = (Mu1+Mu2).M();
253  for (unsigned int iEtaRegion=0; iEtaRegion<3; iEtaRegion++){
254  if (iEtaRegion==0) {EtaCutMin= 0.; EtaCutMax=2.4; }
255  if (iEtaRegion==1) {EtaCutMin= etaBMin; EtaCutMax=etaBMax; }
256  if (iEtaRegion==2) {EtaCutMin= etaECMin; EtaCutMax=etaECMax; }
257 
258  if(fabs(recoTrack1->eta())>EtaCutMin && fabs(recoTrack1->eta())<EtaCutMax &&
259  fabs(recoTrack2->eta())>EtaCutMin && fabs(recoTrack2->eta())<EtaCutMax){
260  if (InvMass < LowMassMax) TrkTrkMuon_LM[iEtaRegion]->Fill(InvMass);
261  if (InvMass > HighMassMin) TrkTrkMuon_HM[iEtaRegion]->Fill(InvMass);
262  }
263  }
264  }
265 
266 
267  LogTrace(metname)<<"[DiMuonHistograms] Soft-Soft pair"<<endl;
268 
269  if (muon::isSoftMuon(*recoMu1, vtx) &&
270  muon::isSoftMuon(*recoMu2, vtx) ) {
271 
272  if (charge < 0) {
273  InvMass = (Mu1+Mu2).M();
274  for (unsigned int iEtaRegion=0; iEtaRegion<3; iEtaRegion++){
275  if (iEtaRegion==0) {EtaCutMin= 0.; EtaCutMax=2.4; }
276  if (iEtaRegion==1) {EtaCutMin= etaBMin; EtaCutMax=etaBMax; }
277  if (iEtaRegion==2) {EtaCutMin= etaECMin; EtaCutMax=etaECMax; }
278 
279  if(fabs(recoTrack1->eta())>EtaCutMin && fabs(recoTrack1->eta())<EtaCutMax &&
280  fabs(recoTrack2->eta())>EtaCutMin && fabs(recoTrack2->eta())<EtaCutMax){
281  SoftSoftMuon[iEtaRegion]->Fill(InvMass);
282  }
283  }
284  }
285  }
286  }
287  } //muon2
288  } //Muon1
289 }
290 
dictionary parameters
Definition: Parameters.py:2
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
const std::string metname
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:43
double charge(const std::vector< uint8_t > &Ampls)
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
bool isSoftMuon(const reco::Muon &, const reco::Vertex &)
int iEvent
Definition: GenABIO.cc:243
void beginRun(const edm::Run &iRun, const edm::EventSetup &iSetup)
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
double BeamWidthX() const
beam width X
Definition: BeamSpot.h:86
DiMuonHistograms(const edm::ParameterSet &pset)
bool isValid() const
Definition: HandleBase.h:76
#define LogTrace(id)
double sigmaZ() const
sigma z
Definition: BeamSpot.h:80
double BeamWidthY() const
beam width Y
Definition: BeamSpot.h:88
virtual ~DiMuonHistograms()
tuple muons
Definition: patZpeak.py:38
const Point & position() const
position
Definition: BeamSpot.h:62
bool isTightMuon(const reco::Muon &, const reco::Vertex &)
void analyze(const edm::Event &, const edm::EventSetup &)
Definition: Run.h:41