CMS 3D CMS Logo

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