CMS 3D CMS Logo

BPhysicsOniaDQM.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * \author S. Bolognesi, Erik - CERN
5  */
6 
8 
12 
16 
22 
24 
25 using namespace std;
26 using namespace edm;
27 using namespace reco;
28 
30  // Muon Collection Label
31  vertex_ = consumes<reco::VertexCollection>(
32  parameters.getParameter<InputTag>("vertex"));
33  theMuonCollectionLabel_ = consumes<reco::MuonCollection>(
34  parameters.getParameter<InputTag>("MuonCollection"));
35  lumiSummaryToken_ = consumes<LumiSummary, edm::InLumi>(
36  parameters.getParameter<InputTag>("lumiSummary"));
37 
38  global_background = nullptr;
39  diMuonMass_global = nullptr;
40  tracker_background = nullptr;
41  diMuonMass_tracker = nullptr;
42  standalone_background = nullptr;
43  diMuonMass_standalone = nullptr;
44 
45  glbSigCut = nullptr;
46  glbSigNoCut = nullptr;
47  glbBkgNoCut = nullptr;
48  staSigCut = nullptr;
49  staSigNoCut = nullptr;
50  staBkgNoCut = nullptr;
51  trkSigCut = nullptr;
52  trkSigNoCut = nullptr;
53  trkBkgNoCut = nullptr;
54 
55  metname = "oniaAnalyzer";
56 }
57 
59 
61  edm::Run const &,
62  edm::EventSetup const &) {
63  iBooker.setCurrentFolder("Physics/BPhysics"); // Use folder with name of PAG
64 
65  global_background = iBooker.book1D(
66  "global_background", "Same-sign global-global dimuon mass", 750, 0, 15);
67  diMuonMass_global =
68  iBooker.book1D("diMuonMass_global",
69  "Opposite-sign global-global dimuon mass", 750, 0, 15);
70  tracker_background = iBooker.book1D(
71  "tracker_background",
72  "Same-sign tracker-tracker (arbitrated) dimuon mass", 750, 0, 15);
73  diMuonMass_tracker = iBooker.book1D(
74  "diMuonMass_tracker",
75  "Opposite-sign tracker-tracker (arbitrated) dimuon mass", 750, 0, 15);
76  standalone_background =
77  iBooker.book1D("standalone_background",
78  "Same-sign standalone-standalone dimuon mass", 500, 0, 15);
79  diMuonMass_standalone = iBooker.book1D(
80  "diMuonMass_standalone",
81  "Opposite-sign standalone-standalone dimuon mass", 500, 0, 15);
82 
83  glbSigCut = iBooker.book1D("glbSigCut", "Opposite-sign glb-glb dimuon mass",
84  650, 0, 130);
85  glbSigNoCut = iBooker.book1D(
86  "glbSigNoCut", "Opposite-sign glb-glb dimuon mass (no cut)", 650, 0, 130);
87  glbBkgNoCut = iBooker.book1D(
88  "glbBkgNoCut", "Same-sign glb-glb dimuon mass (no cut)", 650, 0, 130);
89  staSigCut = iBooker.book1D("staSigCut", "Opposite-sign sta-sta dimuon mass",
90  430, 0, 129);
91  staSigNoCut = iBooker.book1D(
92  "staSigNoCut", "Opposite-sign sta-sta dimuon mass (no cut)", 430, 0, 129);
93  staBkgNoCut = iBooker.book1D(
94  "staBkgNoCut", "Same-sign sta-sta dimuon mass (no cut)", 430, 0, 129);
95  trkSigCut = iBooker.book1D("trkSigCut", "Opposite-sign trk-trk dimuon mass",
96  650, 0, 130);
97  trkSigNoCut = iBooker.book1D(
98  "trkSigNoCut", "Opposite-sign trk-trk dimuon mass (no cut)", 650, 0, 130);
99  trkBkgNoCut = iBooker.book1D(
100  "trkBkgNoCutt", "Same-sign trk-trk dimuon mass (no cut)", 650, 0, 130);
101 }
102 
103 void BPhysicsOniaDQM::analyze(const Event &iEvent, const EventSetup &iSetup) {
104  LogTrace(metname) << "[BPhysicsOniaDQM] Analysis of event # ";
105 
106  // Take the STA muon container
108  iEvent.getByToken(theMuonCollectionLabel_, muons);
109 
111  iEvent.getByToken(vertex_, privtxs);
112  VertexCollection::const_iterator privtx;
113 
114  if (privtxs->begin() != privtxs->end()) {
115  privtx = privtxs->begin();
116  RefVtx = privtx->position();
117  } else {
118  RefVtx.SetXYZ(0., 0., 0.);
119  }
120 
121  if (muons.isValid()) {
122  for (MuonCollection::const_iterator recoMu1 = muons->begin();
123  recoMu1 != muons->end(); ++recoMu1) {
124  // only loop over the remaining muons if recoMu1 is one of the following
125  if (recoMu1->isGlobalMuon() || recoMu1->isTrackerMuon() ||
126  recoMu1->isStandAloneMuon()) {
127  for (MuonCollection::const_iterator recoMu2 = recoMu1 + 1;
128  recoMu2 != muons->end(); ++recoMu2) {
129  // fill the relevant histograms if recoMu2 satisfies one of the
130  // following
131  if (recoMu1->isGlobalMuon() && recoMu2->isGlobalMuon()) {
132  math::XYZVector vec1 = recoMu1->globalTrack()->momentum();
133  math::XYZVector vec2 = recoMu2->globalTrack()->momentum();
134  float massJPsi = computeMass(vec1, vec2);
135 
136  // if opposite charges, fill glbSig, else fill glbBkg
137  if (((*recoMu1).charge() * (*recoMu2).charge()) < 0) {
138  if (diMuonMass_global != nullptr) { // BPhysicsOniaDQM original one
139  diMuonMass_global->Fill(massJPsi);
140  }
141 
142  if (glbSigNoCut != nullptr) {
143  glbSigNoCut->Fill(massJPsi);
144  if (selGlobalMuon(*recoMu1) && selGlobalMuon(*recoMu2)) {
145  if (glbSigCut != nullptr) glbSigCut->Fill(massJPsi);
146  }
147  }
148  } else {
149  if (global_background != nullptr) { // BPhysicsOniaDQM original one
150  global_background->Fill(massJPsi);
151  }
152 
153  if (glbBkgNoCut != nullptr) {
154  glbBkgNoCut->Fill(massJPsi);
155  }
156  }
157  }
158 
159  if (recoMu1->isStandAloneMuon() && recoMu2->isStandAloneMuon() &&
160  fabs(recoMu1->outerTrack()->d0()) < 5 &&
161  fabs(recoMu1->outerTrack()->dz()) < 30 &&
162  fabs(recoMu2->outerTrack()->d0()) < 5 &&
163  fabs(recoMu2->outerTrack()->dz()) < 30) {
164  math::XYZVector vec1 = recoMu1->outerTrack()->momentum();
165  math::XYZVector vec2 = recoMu2->outerTrack()->momentum();
166  float massJPsi = computeMass(vec1, vec2);
167 
168  // if opposite charges, fill staSig, else fill staBkg
169  if (((*recoMu1).charge() * (*recoMu2).charge()) < 0) {
170  if (diMuonMass_standalone != nullptr) {
171  diMuonMass_standalone->Fill(massJPsi);
172  }
173 
174  if (staSigNoCut != nullptr) {
175  staSigNoCut->Fill(massJPsi);
176  }
177  } else {
178  if (standalone_background != nullptr) {
179  standalone_background->Fill(massJPsi);
180  }
181 
182  if (staBkgNoCut != nullptr) {
183  staBkgNoCut->Fill(massJPsi);
184  }
185  }
186  }
187 
188  if (recoMu1->isTrackerMuon() && recoMu2->isTrackerMuon() &&
191  math::XYZVector vec1 = recoMu1->innerTrack()->momentum();
192  math::XYZVector vec2 = recoMu2->innerTrack()->momentum();
193  float massJPsi = computeMass(vec1, vec2);
194 
195  // if opposite charges, fill trkSig, else fill trkBkg
196  if (((*recoMu1).charge() * (*recoMu2).charge()) < 0) {
197  if (diMuonMass_tracker != nullptr) {
198  diMuonMass_tracker->Fill(massJPsi);
199  }
200 
201  if (trkSigNoCut != nullptr) {
202  trkSigNoCut->Fill(massJPsi);
203  if (selTrackerMuon(*recoMu1) && selTrackerMuon(*recoMu2)) {
204  if (trkSigCut != nullptr) trkSigCut->Fill(massJPsi);
205  }
206  }
207  } else {
208  if (tracker_background != nullptr) {
209  tracker_background->Fill(massJPsi);
210  }
211 
212  if (trkBkgNoCut != nullptr) {
213  trkBkgNoCut->Fill(massJPsi);
214  }
215  }
216  }
217  } // end of 2nd MuonCollection
218  } // end of GLB,STA,TRK muon check
219  } // end of 1st MuonCollection
220  } // Is this MuonCollection vaild?
221 }
222 
224  const math::XYZVector &vec2) {
225  // mass of muon
226  float massMu = 0.10566;
227  float eMu1 = -999;
228  if (massMu * massMu + vec1.Mag2() > 0)
229  eMu1 = sqrt(massMu * massMu + vec1.Mag2());
230  float eMu2 = -999;
231  if (massMu * massMu + vec2.Mag2() > 0)
232  eMu2 = sqrt(massMu * massMu + vec2.Mag2());
233 
234  float pJPsi = -999;
235  if ((vec1 + vec2).Mag2() > 0) pJPsi = sqrt((vec1 + vec2).Mag2());
236  float eJPsi = eMu1 + eMu2;
237 
238  float massJPsi = -999;
239  if ((eJPsi * eJPsi - pJPsi * pJPsi) > 0)
240  massJPsi = sqrt(eJPsi * eJPsi - pJPsi * pJPsi);
241 
242  return massJPsi;
243 }
244 
246  return (fabs(recoMu.eta()) < 2.4 &&
247  ((fabs(recoMu.eta()) < 1.3 && recoMu.pt() > 3.3) ||
248  (fabs(recoMu.eta()) > 1.3 && fabs(recoMu.eta()) < 2.2 &&
249  recoMu.p() > 2.9) ||
250  (fabs(recoMu.eta()) > 2.2 && recoMu.pt() > 0.8)));
251 }
252 
254  TrackRef iTrack = recoMu.innerTrack();
255  const reco::HitPattern &p = iTrack->hitPattern();
256 
257  TrackRef gTrack = recoMu.globalTrack();
258  const reco::HitPattern &q = gTrack->hitPattern();
259 
260  return (isMuonInAccept(recoMu) && iTrack->found() > 11 &&
261  gTrack->chi2() / gTrack->ndof() < 20.0 &&
262  q.numberOfValidMuonHits() > 0 &&
263  iTrack->chi2() / iTrack->ndof() < 4.0 &&
264  // recoMu.muonID("TrackerMuonArbitrated") &&
265  // recoMu.muonID("TMLastStationAngTight") &&
266  p.pixelLayersWithMeasurement() > 1 &&
267  fabs(iTrack->dxy(RefVtx)) < 3.0 && fabs(iTrack->dz(RefVtx)) < 15.0);
268 }
269 
271  TrackRef iTrack = recoMu.innerTrack();
272  const reco::HitPattern &p = iTrack->hitPattern();
273 
274  return (isMuonInAccept(recoMu) && iTrack->found() > 11 &&
275  iTrack->chi2() / iTrack->ndof() < 4.0 &&
276  // recoMu.muonID("TrackerMuonArbitrated") &&
277  // recoMu.muonID("TMLastStationAngTight") &&
278  p.pixelLayersWithMeasurement() > 1 &&
279  fabs(iTrack->dxy(RefVtx)) < 3.0 && fabs(iTrack->dz(RefVtx)) < 15.0);
280 }
281 
282 // Local Variables:
283 // show-trailing-whitespace: t
284 // truncate-lines: t
285 // End:
T getParameter(std::string const &) const
bool selGlobalMuon(const reco::Muon &recoMu)
double eta() const final
momentum pseudorapidity
virtual TrackRef innerTrack() const
Definition: Muon.h:48
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:508
const std::string metname
double pt() const final
transverse momentum
int pixelLayersWithMeasurement() const
Definition: HitPattern.cc:499
void analyze(const edm::Event &, const edm::EventSetup &) override
Get the analysis.
int iEvent
Definition: GenABIO.cc:230
bool isMuonInAccept(const reco::Muon &recoMu)
std::vector< double > vec1
Definition: HCALResponse.h:15
T sqrt(T t)
Definition: SSEVec.h:18
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
BPhysicsOniaDQM(const edm::ParameterSet &)
Constructor.
~BPhysicsOniaDQM() override
Destructor.
bool isValid() const
Definition: HandleBase.h:74
#define LogTrace(id)
bool isGoodMuon(const reco::Muon &muon, SelectionType type, reco::Muon::ArbitrationType arbitrationType=reco::Muon::SegmentAndTrackArbitration)
main GoodMuon wrapper call
double p() const final
magnitude of momentum vector
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
fixed size matrix
HLT enums.
bool selTrackerMuon(const reco::Muon &recoMu)
float computeMass(const math::XYZVector &vec1, const math::XYZVector &vec2)
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
int numberOfValidMuonHits() const
Definition: HitPattern.h:833
std::vector< vec1 > vec2
Definition: HCALResponse.h:16
Definition: Run.h:43
virtual TrackRef globalTrack() const
reference to Track reconstructed in both tracked and muon detector
Definition: Muon.h:54