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