CMS 3D CMS Logo

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