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