31 vertex_ = consumes<reco::VertexCollection>(
33 theMuonCollectionLabel_ = consumes<reco::MuonCollection>(
35 lumiSummaryToken_ = consumes<LumiSummary, edm::InLumi>(
38 global_background =
NULL;
39 diMuonMass_global =
NULL;
40 tracker_background =
NULL;
41 diMuonMass_tracker =
NULL;
42 standalone_background =
NULL;
43 diMuonMass_standalone =
NULL;
65 global_background = iBooker.
book1D(
66 "global_background",
"Same-sign global-global dimuon mass", 750, 0, 15);
68 iBooker.
book1D(
"diMuonMass_global",
69 "Opposite-sign global-global dimuon mass", 750, 0, 15);
70 tracker_background = iBooker.
book1D(
72 "Same-sign tracker-tracker (arbitrated) dimuon mass", 750, 0, 15);
73 diMuonMass_tracker = iBooker.
book1D(
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);
83 glbSigCut = iBooker.
book1D(
"glbSigCut",
"Opposite-sign glb-glb dimuon mass",
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",
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",
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);
108 iEvent.
getByToken(theMuonCollectionLabel_, muons);
112 VertexCollection::const_iterator privtx;
114 if (privtxs->begin() != privtxs->end()) {
115 privtx = privtxs->begin();
116 RefVtx = privtx->position();
118 RefVtx.SetXYZ(0., 0., 0.);
122 for (MuonCollection::const_iterator recoMu1 = muons->begin();
123 recoMu1 != muons->end(); ++recoMu1) {
125 if (recoMu1->isGlobalMuon() || recoMu1->isTrackerMuon() ||
126 recoMu1->isStandAloneMuon()) {
127 for (MuonCollection::const_iterator recoMu2 = recoMu1 + 1;
128 recoMu2 != muons->end(); ++recoMu2) {
131 if (recoMu1->isGlobalMuon() && recoMu2->isGlobalMuon()) {
134 float massJPsi = computeMass(vec1, vec2);
137 if (((*recoMu1).charge() * (*recoMu2).charge()) < 0) {
138 if (diMuonMass_global !=
NULL) {
139 diMuonMass_global->Fill(massJPsi);
142 if (glbSigNoCut !=
NULL) {
143 glbSigNoCut->Fill(massJPsi);
144 if (selGlobalMuon(*recoMu1) && selGlobalMuon(*recoMu2)) {
145 if (glbSigCut !=
NULL) glbSigCut->Fill(massJPsi);
149 if (global_background !=
NULL) {
150 global_background->Fill(massJPsi);
153 if (glbBkgNoCut !=
NULL) {
154 glbBkgNoCut->Fill(massJPsi);
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) {
166 float massJPsi = computeMass(vec1, vec2);
169 if (((*recoMu1).charge() * (*recoMu2).charge()) < 0) {
170 if (diMuonMass_standalone !=
NULL) {
171 diMuonMass_standalone->Fill(massJPsi);
174 if (staSigNoCut !=
NULL) {
175 staSigNoCut->Fill(massJPsi);
178 if (standalone_background !=
NULL) {
179 standalone_background->Fill(massJPsi);
182 if (staBkgNoCut !=
NULL) {
183 staBkgNoCut->Fill(massJPsi);
188 if (recoMu1->isTrackerMuon() && recoMu2->isTrackerMuon() &&
193 float massJPsi = computeMass(vec1, vec2);
196 if (((*recoMu1).charge() * (*recoMu2).charge()) < 0) {
197 if (diMuonMass_tracker !=
NULL) {
198 diMuonMass_tracker->Fill(massJPsi);
201 if (trkSigNoCut !=
NULL) {
202 trkSigNoCut->Fill(massJPsi);
203 if (selTrackerMuon(*recoMu1) && selTrackerMuon(*recoMu2)) {
204 if (trkSigCut !=
NULL) trkSigCut->Fill(massJPsi);
208 if (tracker_background !=
NULL) {
209 tracker_background->Fill(massJPsi);
212 if (trkBkgNoCut !=
NULL) {
213 trkBkgNoCut->Fill(massJPsi);
226 float massMu = 0.10566;
228 if (massMu * massMu + vec1.Mag2() > 0)
229 eMu1 =
sqrt(massMu * massMu + vec1.Mag2());
231 if (massMu * massMu + vec2.Mag2() > 0)
232 eMu2 =
sqrt(massMu * massMu + vec2.Mag2());
235 if ((vec1 + vec2).Mag2() > 0) pJPsi =
sqrt((vec1 + vec2).Mag2());
236 float eJPsi = eMu1 + eMu2;
238 float massJPsi = -999;
239 if ((eJPsi * eJPsi - pJPsi * pJPsi) > 0)
240 massJPsi =
sqrt(eJPsi * eJPsi - pJPsi * pJPsi);
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 &&
250 (fabs(recoMu.
eta()) > 2.2 && recoMu.
pt() > 0.8)));
260 return (isMuonInAccept(recoMu) && iTrack->found() > 11 &&
261 gTrack->chi2() / gTrack->ndof() < 20.0 &&
263 iTrack->chi2() / iTrack->ndof() < 4.0 &&
267 fabs(iTrack->dxy(RefVtx)) < 3.0 && fabs(iTrack->dz(RefVtx)) < 15.0);
274 return (isMuonInAccept(recoMu) && iTrack->found() > 11 &&
275 iTrack->chi2() / iTrack->ndof() < 4.0 &&
279 fabs(iTrack->dxy(RefVtx)) < 3.0 && fabs(iTrack->dz(RefVtx)) < 15.0);
T getParameter(std::string const &) const
virtual double pt() const final
transverse momentum
bool selGlobalMuon(const reco::Muon &recoMu)
virtual double eta() const final
momentum pseudorapidity
virtual TrackRef innerTrack() const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
const std::string metname
int pixelLayersWithMeasurement() const
void analyze(const edm::Event &, const edm::EventSetup &) override
Get the analysis.
bool isMuonInAccept(const reco::Muon &recoMu)
std::vector< double > vec1
MonitorElement * book1D(Args &&...args)
BPhysicsOniaDQM(const edm::ParameterSet &)
Constructor.
virtual double p() const final
magnitude of momentum vector
bool isGoodMuon(const reco::Muon &muon, SelectionType type, reco::Muon::ArbitrationType arbitrationType=reco::Muon::SegmentAndTrackArbitration)
main GoodMuon wrapper call
void setCurrentFolder(const std::string &fullpath)
XYZVectorD XYZVector
spatial vector with cartesian internal representation
bool selTrackerMuon(const reco::Muon &recoMu)
float computeMass(const math::XYZVector &vec1, const math::XYZVector &vec2)
virtual ~BPhysicsOniaDQM()
Destructor.
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
int numberOfValidMuonHits() const
virtual TrackRef globalTrack() const
reference to Track reconstructed in both tracked and muon detector