Skip to content

Commit dd1a441

Browse files
committed
add more monitoring histograms to shortened track resolution
1 parent cad258f commit dd1a441

File tree

1 file changed

+105
-11
lines changed

1 file changed

+105
-11
lines changed

DQM/TrackingMonitorSource/plugins/ShortenedTrackResolution.cc

Lines changed: 105 additions & 11 deletions
Original file line numberDiff line numberDiff line change
@@ -39,7 +39,15 @@ class ShortenedTrackResolution : public DQMEDAnalyzer {
3939
const edm::EDGetTokenT<std::vector<reco::Track>> tracksToken_;
4040
const std::vector<edm::EDGetTokenT<std::vector<reco::Track>>> tracksRerecoToken_;
4141

42+
std::vector<MonitorElement *> histsPtRatioAll_;
43+
std::vector<MonitorElement *> histsPtDiffAll_;
44+
std::vector<MonitorElement *> histsEtaDiffAll_;
45+
std::vector<MonitorElement *> histsPhiDiffAll_;
46+
std::vector<MonitorElement *> histsPtRatioVsDeltaRAll_;
47+
std::vector<MonitorElement *> histsDeltaPtOverPtAll_;
4248
std::vector<MonitorElement *> histsPtAll_;
49+
std::vector<MonitorElement *> histsNhitsAll_;
50+
std::vector<MonitorElement *> histsDeltaRAll_;
4351
};
4452

4553
// -----------------------------
@@ -58,25 +66,103 @@ ShortenedTrackResolution::ShortenedTrackResolution(const edm::ParameterSet &ps)
5866
tracksToken_(consumes<std::vector<reco::Track>>(tracksTag_)),
5967
tracksRerecoToken_(edm::vector_transform(
6068
tracksRerecoTag_, [this](edm::InputTag const &tag) { return consumes<std::vector<reco::Track>>(tag); })) {
69+
histsPtRatioAll_.clear();
70+
histsPtDiffAll_.clear();
71+
histsEtaDiffAll_.clear();
72+
histsPhiDiffAll_.clear();
73+
histsPtRatioVsDeltaRAll_.clear();
74+
histsDeltaPtOverPtAll_.clear();
6175
histsPtAll_.clear();
76+
histsNhitsAll_.clear();
77+
histsDeltaRAll_.clear();
78+
79+
const size_t n = hitsRemain_.size();
80+
histsPtRatioAll_.reserve(n);
81+
histsPtDiffAll_.reserve(n);
82+
histsEtaDiffAll_.reserve(n);
83+
histsPhiDiffAll_.reserve(n);
84+
histsPtRatioVsDeltaRAll_.reserve(n);
85+
histsDeltaPtOverPtAll_.reserve(n);
86+
histsPtAll_.reserve(n);
87+
histsNhitsAll_.reserve(n);
88+
histsDeltaRAll_.reserve(n);
6289
}
6390

6491
//__________________________________________________________________________________
6592
void ShortenedTrackResolution::bookHistograms(DQMStore::IBooker &iBook,
6693
edm::Run const &iRun,
6794
edm::EventSetup const &iSetup) {
68-
std::string currentFolder = folderName_ + "/";
95+
auto book1D = [&](const std::string &name, const std::string &title, int bins, double min, double max) {
96+
return iBook.book1D(name.c_str(), title.c_str(), bins, min, max);
97+
};
98+
99+
auto book2D = [&](const std::string &name,
100+
const std::string &title,
101+
int binsX,
102+
double minX,
103+
double maxX,
104+
int binsY,
105+
double minY,
106+
double maxY) {
107+
return iBook.book2D(name.c_str(), title.c_str(), binsX, minX, maxX, binsY, minY, maxY);
108+
};
109+
110+
std::string currentFolder = folderName_ + "/Resolutions";
69111
iBook.setCurrentFolder(currentFolder);
70112

71-
for (int i = 0; i < int(hitsRemain_.size()); ++i) {
72-
histsPtAll_.push_back(iBook.book1D(
73-
fmt::sprintf("trackPtRatio_%s", hitsRemain_[i]).c_str(),
74-
fmt::sprintf("Short Track p_{T} / Full Track p_{T} - %s layers;p_{T}^{short}/p_{T}^{full};n. tracks",
75-
hitsRemain_[i])
76-
.c_str(),
77-
101,
78-
-0.05,
79-
2.05));
113+
for (const auto &label : hitsRemain_) {
114+
std::string name, title;
115+
116+
name = fmt::sprintf("trackPtRatio_%s", label);
117+
title =
118+
fmt::sprintf("Short Track p_{T} / Full Track p_{T} - %s layers;p_{T}^{short}/p_{T}^{full};n. tracks", label);
119+
histsPtRatioAll_.push_back(book1D(name, title, 100, 0.5, 1.5));
120+
121+
name = fmt::sprintf("trackPtDiff_%s", label);
122+
title = fmt::sprintf(
123+
"Short Track p_{T} - Full Track p_{T} - %s layers;p_{T}^{short} - p_{T}^{full} [GeV];n. tracks", label);
124+
histsPtDiffAll_.push_back(book1D(name, title, 100, -10., 10.));
125+
126+
name = fmt::sprintf("trackEtaDiff_%s", label);
127+
title = fmt::sprintf("Short Track #eta - Full Track #eta - %s layers;#eta^{short} - #eta^{full};n. tracks", label);
128+
histsEtaDiffAll_.push_back(book1D(name, title, 100, -0.001, 0.001));
129+
130+
name = fmt::sprintf("trackPhiDiff_%s", label);
131+
title = fmt::sprintf("Short Track #phi - Full Track #phi - %s layers;#phi^{short} - #phi^{full};n. tracks", label);
132+
histsPhiDiffAll_.push_back(book1D(name, title, 100, -0.001, 0.001));
133+
134+
name = fmt::sprintf("trackPtRatioVsDeltaR_%s", label);
135+
title = fmt::sprintf(
136+
"Short Track p_{T} / Full Track p_{T} - %s layers vs "
137+
"#DeltaR;#DeltaR(short,full);p_{T}^{short}/p_{T}^{full} [GeV];n. tracks",
138+
label);
139+
histsPtRatioVsDeltaRAll_.push_back(book2D(name, title, 100, 0., 0.01, 101, -0.05, 2.05));
140+
141+
name = fmt::sprintf("trackDeltaPtOverPt_%s", label);
142+
title = fmt::sprintf(
143+
"Short Track p_{T} - Full Track p_{T} / Full Track p_{T} - %s layers;"
144+
"p_{T}^{short} - p_{T}^{full} / p^{full}_{T};n. tracks",
145+
label);
146+
histsDeltaPtOverPtAll_.push_back(book1D(name, title, 101, -5., 5.));
147+
}
148+
149+
currentFolder = folderName_ + "/TrackProperties";
150+
iBook.setCurrentFolder(currentFolder);
151+
152+
for (const auto &label : hitsRemain_) {
153+
std::string name, title;
154+
155+
name = fmt::sprintf("trackPt_%s", label);
156+
title = fmt::sprintf("Short Track p_{T} - %s layers;p_{T}^{short} [GeV];n. tracks", label);
157+
histsPtAll_.push_back(book1D(name, title, 100, 0., 100.));
158+
159+
name = fmt::sprintf("trackNhits_%s", label);
160+
title = fmt::sprintf("Short Track n. hits - %s layers;n. hits per track;n. tracks", label);
161+
histsNhitsAll_.push_back(book1D(name, title, 20, -0.5, 19.5));
162+
163+
name = fmt::sprintf("trackDeltaR_%s", label);
164+
title = fmt::sprintf("Short Track / Full Track #DeltaR - %s layers;#DeltaR(short,full);n. tracks", label);
165+
histsDeltaRAll_.push_back(book1D(name, title, 100, 0., 0.005));
80166
}
81167
}
82168

@@ -110,7 +196,15 @@ void ShortenedTrackResolution::analyze(edm::Event const &iEvent, edm::EventSetup
110196
if (deltaR < maxDr_) {
111197
if (track_rereco.pt() >= minTracksPt_ && track_rereco.pt() <= maxTracksPt_ &&
112198
std::abs(track_rereco.eta()) >= minTracksEta_ && std::abs(track_rereco.eta()) <= maxTracksEta_) {
113-
histsPtAll_[i]->Fill(1.0 * track_rereco.pt() / track.pt());
199+
histsPtRatioAll_[i]->Fill(1.0 * track_rereco.pt() / track.pt());
200+
histsPtDiffAll_[i]->Fill(track_rereco.pt() - track.pt());
201+
histsDeltaPtOverPtAll_[i]->Fill((track_rereco.pt() - track.pt()) / track.pt());
202+
histsEtaDiffAll_[i]->Fill(track_rereco.eta() - track.eta());
203+
histsPhiDiffAll_[i]->Fill(track_rereco.phi() - track.phi());
204+
histsPtRatioVsDeltaRAll_[i]->Fill(deltaR, track_rereco.pt() / track.pt());
205+
histsPtAll_[i]->Fill(track_rereco.pt());
206+
histsNhitsAll_[i]->Fill(track_rereco.numberOfValidHits());
207+
histsDeltaRAll_[i]->Fill(deltaR);
114208
}
115209
}
116210
}

0 commit comments

Comments
 (0)