1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
|
#include "DQM/SiStripCommissioningAnalysis/interface/PedsOnlyAlgorithm.h"
#include "CondFormats/SiStripObjects/interface/PedsOnlyAnalysis.h"
#include "DataFormats/SiStripCommon/interface/SiStripHistoTitle.h"
#include "DataFormats/SiStripCommon/interface/SiStripEnumsAndStrings.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "TProfile.h"
#include "TH1.h"
#include <iostream>
#include <iomanip>
#include <cmath>
using namespace sistrip;
// ----------------------------------------------------------------------------
//
PedsOnlyAlgorithm::PedsOnlyAlgorithm(const edm::ParameterSet& pset, PedsOnlyAnalysis* const anal)
: CommissioningAlgorithm(anal), hPeds_(nullptr, ""), hNoise_(nullptr, "") {}
// ----------------------------------------------------------------------------
//
void PedsOnlyAlgorithm::extract(const std::vector<TH1*>& histos) {
if (!anal()) {
edm::LogWarning(mlCommissioning_) << "[PedsOnlyAlgorithm::" << __func__ << "]"
<< " NULL pointer to Analysis object!";
return;
}
// Check number of histograms
if (histos.size() != 2) {
anal()->addErrorCode(sistrip::numberOfHistos_);
}
// Extract FED key from histo title
if (!histos.empty()) {
anal()->fedKey(extractFedKey(histos.front()));
}
// Extract histograms
std::vector<TH1*>::const_iterator ihis = histos.begin();
for (; ihis != histos.end(); ihis++) {
// Check for NULL pointer
if (!(*ihis)) {
continue;
}
// Check run type
SiStripHistoTitle title((*ihis)->GetName());
if (title.runType() != sistrip::PEDS_ONLY) {
anal()->addErrorCode(sistrip::unexpectedTask_);
continue;
}
// Extract peds and raw noise histos (check for legacy names first!)
if (title.extraInfo().find(sistrip::extrainfo::pedsAndRawNoise_) != std::string::npos) {
hPeds_.first = *ihis;
hPeds_.second = (*ihis)->GetName();
hNoise_.first = *ihis;
hNoise_.second = (*ihis)->GetName();
PedsOnlyAnalysis* a = dynamic_cast<PedsOnlyAnalysis*>(const_cast<CommissioningAnalysis*>(anal()));
if (a) {
a->legacy_ = true;
}
} else if (title.extraInfo().find(sistrip::extrainfo::pedestals_) != std::string::npos) {
hPeds_.first = *ihis;
hPeds_.second = (*ihis)->GetName();
} else if (title.extraInfo().find(sistrip::extrainfo::rawNoise_) != std::string::npos) {
hNoise_.first = *ihis;
hNoise_.second = (*ihis)->GetName();
} else {
anal()->addErrorCode(sistrip::unexpectedExtraInfo_);
}
}
}
// -----------------------------------------------------------------------------
//
void PedsOnlyAlgorithm::analyse() {
if (!anal()) {
edm::LogWarning(mlCommissioning_) << "[PedsOnlyAlgorithm::" << __func__ << "]"
<< " NULL pointer to base Analysis object!";
return;
}
CommissioningAnalysis* tmp = const_cast<CommissioningAnalysis*>(anal());
PedsOnlyAnalysis* anal = dynamic_cast<PedsOnlyAnalysis*>(tmp);
if (!anal) {
edm::LogWarning(mlCommissioning_) << "[PedsOnlyAlgorithm::" << __func__ << "]"
<< " NULL pointer to derived Analysis object!";
return;
}
if (!hPeds_.first) {
anal->addErrorCode(sistrip::nullPtr_);
return;
}
if (!hNoise_.first) {
anal->addErrorCode(sistrip::nullPtr_);
return;
}
TProfile* peds_histo = dynamic_cast<TProfile*>(hPeds_.first);
TProfile* raw_histo = dynamic_cast<TProfile*>(hNoise_.first);
if (!peds_histo) {
anal->addErrorCode(sistrip::nullPtr_);
return;
}
if (!raw_histo) {
anal->addErrorCode(sistrip::nullPtr_);
return;
}
if (peds_histo->GetNbinsX() != 256) {
anal->addErrorCode(sistrip::numberOfBins_);
return;
}
if (raw_histo->GetNbinsX() != 256) {
anal->addErrorCode(sistrip::numberOfBins_);
return;
}
// Iterate through APVs
for (uint16_t iapv = 0; iapv < 2; iapv++) {
// Used to calc mean and rms for peds and noise
float p_sum = 0., p_sum2 = 0., p_max = -1. * sistrip::invalid_, p_min = sistrip::invalid_;
float r_sum = 0., r_sum2 = 0., r_max = -1. * sistrip::invalid_, r_min = sistrip::invalid_;
// Iterate through strips of APV
for (uint16_t istr = 0; istr < 128; istr++) {
uint16_t strip = iapv * 128 + istr;
// Pedestals
if (peds_histo) {
if (peds_histo->GetBinEntries(strip + 1)) {
anal->peds_[iapv][istr] = peds_histo->GetBinContent(strip + 1);
p_sum += anal->peds_[iapv][istr];
p_sum2 += (anal->peds_[iapv][istr] * anal->peds_[iapv][istr]);
if (anal->peds_[iapv][istr] > p_max) {
p_max = anal->peds_[iapv][istr];
}
if (anal->peds_[iapv][istr] < p_min) {
p_min = anal->peds_[iapv][istr];
}
}
}
// Raw noise
if (!anal->legacy_) {
if (raw_histo) {
if (raw_histo->GetBinEntries(strip + 1)) {
anal->raw_[iapv][istr] = raw_histo->GetBinContent(strip + 1);
r_sum += anal->raw_[iapv][istr];
r_sum2 += (anal->raw_[iapv][istr] * anal->raw_[iapv][istr]);
if (anal->raw_[iapv][istr] > r_max) {
r_max = anal->raw_[iapv][istr];
}
if (anal->raw_[iapv][istr] < r_min) {
r_min = anal->raw_[iapv][istr];
}
}
}
} else {
if (peds_histo) {
if (peds_histo->GetBinEntries(strip + 1)) {
anal->raw_[iapv][istr] = raw_histo->GetBinError(strip + 1);
r_sum += anal->raw_[iapv][istr];
r_sum2 += (anal->raw_[iapv][istr] * anal->raw_[iapv][istr]);
if (anal->raw_[iapv][istr] > r_max) {
r_max = anal->raw_[iapv][istr];
}
if (anal->raw_[iapv][istr] < r_min) {
r_min = anal->raw_[iapv][istr];
}
}
}
}
} // strip loop
// Calc mean and rms for peds
if (!anal->peds_[iapv].empty()) {
p_sum /= static_cast<float>(anal->peds_[iapv].size());
p_sum2 /= static_cast<float>(anal->peds_[iapv].size());
anal->pedsMean_[iapv] = p_sum;
anal->pedsSpread_[iapv] = sqrt(fabs(p_sum2 - p_sum * p_sum));
}
// Calc mean and rms for raw noise
if (!anal->raw_[iapv].empty()) {
r_sum /= static_cast<float>(anal->raw_[iapv].size());
r_sum2 /= static_cast<float>(anal->raw_[iapv].size());
anal->rawMean_[iapv] = r_sum;
anal->rawSpread_[iapv] = sqrt(fabs(r_sum2 - r_sum * r_sum));
}
// Set max and min values for peds and raw noise
if (p_max > -1. * sistrip::maximum_) {
anal->pedsMax_[iapv] = p_max;
}
if (p_min < 1. * sistrip::maximum_) {
anal->pedsMin_[iapv] = p_min;
}
if (r_max > -1. * sistrip::maximum_) {
anal->rawMax_[iapv] = r_max;
}
if (r_min < 1. * sistrip::maximum_) {
anal->rawMin_[iapv] = r_min;
}
} // apv loop
}
|