File indexing completed on 2024-04-06 12:08:30
0001 #include "DQM/SiStripCommissioningAnalysis/interface/PedsOnlyAlgorithm.h"
0002 #include "CondFormats/SiStripObjects/interface/PedsOnlyAnalysis.h"
0003 #include "DataFormats/SiStripCommon/interface/SiStripHistoTitle.h"
0004 #include "DataFormats/SiStripCommon/interface/SiStripEnumsAndStrings.h"
0005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0006 #include "TProfile.h"
0007 #include "TH1.h"
0008 #include <iostream>
0009 #include <iomanip>
0010 #include <cmath>
0011
0012 using namespace sistrip;
0013
0014
0015
0016 PedsOnlyAlgorithm::PedsOnlyAlgorithm(const edm::ParameterSet& pset, PedsOnlyAnalysis* const anal)
0017 : CommissioningAlgorithm(anal), hPeds_(nullptr, ""), hNoise_(nullptr, "") {}
0018
0019
0020
0021 void PedsOnlyAlgorithm::extract(const std::vector<TH1*>& histos) {
0022 if (!anal()) {
0023 edm::LogWarning(mlCommissioning_) << "[PedsOnlyAlgorithm::" << __func__ << "]"
0024 << " NULL pointer to Analysis object!";
0025 return;
0026 }
0027
0028
0029 if (histos.size() != 2) {
0030 anal()->addErrorCode(sistrip::numberOfHistos_);
0031 }
0032
0033
0034 if (!histos.empty()) {
0035 anal()->fedKey(extractFedKey(histos.front()));
0036 }
0037
0038
0039 std::vector<TH1*>::const_iterator ihis = histos.begin();
0040 for (; ihis != histos.end(); ihis++) {
0041
0042 if (!(*ihis)) {
0043 continue;
0044 }
0045
0046
0047 SiStripHistoTitle title((*ihis)->GetName());
0048 if (title.runType() != sistrip::PEDS_ONLY) {
0049 anal()->addErrorCode(sistrip::unexpectedTask_);
0050 continue;
0051 }
0052
0053
0054 if (title.extraInfo().find(sistrip::extrainfo::pedsAndRawNoise_) != std::string::npos) {
0055 hPeds_.first = *ihis;
0056 hPeds_.second = (*ihis)->GetName();
0057 hNoise_.first = *ihis;
0058 hNoise_.second = (*ihis)->GetName();
0059 PedsOnlyAnalysis* a = dynamic_cast<PedsOnlyAnalysis*>(const_cast<CommissioningAnalysis*>(anal()));
0060 if (a) {
0061 a->legacy_ = true;
0062 }
0063 } else if (title.extraInfo().find(sistrip::extrainfo::pedestals_) != std::string::npos) {
0064 hPeds_.first = *ihis;
0065 hPeds_.second = (*ihis)->GetName();
0066 } else if (title.extraInfo().find(sistrip::extrainfo::rawNoise_) != std::string::npos) {
0067 hNoise_.first = *ihis;
0068 hNoise_.second = (*ihis)->GetName();
0069 } else {
0070 anal()->addErrorCode(sistrip::unexpectedExtraInfo_);
0071 }
0072 }
0073 }
0074
0075
0076
0077 void PedsOnlyAlgorithm::analyse() {
0078 if (!anal()) {
0079 edm::LogWarning(mlCommissioning_) << "[PedsOnlyAlgorithm::" << __func__ << "]"
0080 << " NULL pointer to base Analysis object!";
0081 return;
0082 }
0083
0084 CommissioningAnalysis* tmp = const_cast<CommissioningAnalysis*>(anal());
0085 PedsOnlyAnalysis* anal = dynamic_cast<PedsOnlyAnalysis*>(tmp);
0086 if (!anal) {
0087 edm::LogWarning(mlCommissioning_) << "[PedsOnlyAlgorithm::" << __func__ << "]"
0088 << " NULL pointer to derived Analysis object!";
0089 return;
0090 }
0091
0092 if (!hPeds_.first) {
0093 anal->addErrorCode(sistrip::nullPtr_);
0094 return;
0095 }
0096
0097 if (!hNoise_.first) {
0098 anal->addErrorCode(sistrip::nullPtr_);
0099 return;
0100 }
0101
0102 TProfile* peds_histo = dynamic_cast<TProfile*>(hPeds_.first);
0103 TProfile* raw_histo = dynamic_cast<TProfile*>(hNoise_.first);
0104
0105 if (!peds_histo) {
0106 anal->addErrorCode(sistrip::nullPtr_);
0107 return;
0108 }
0109
0110 if (!raw_histo) {
0111 anal->addErrorCode(sistrip::nullPtr_);
0112 return;
0113 }
0114
0115 if (peds_histo->GetNbinsX() != 256) {
0116 anal->addErrorCode(sistrip::numberOfBins_);
0117 return;
0118 }
0119
0120 if (raw_histo->GetNbinsX() != 256) {
0121 anal->addErrorCode(sistrip::numberOfBins_);
0122 return;
0123 }
0124
0125
0126 for (uint16_t iapv = 0; iapv < 2; iapv++) {
0127
0128 float p_sum = 0., p_sum2 = 0., p_max = -1. * sistrip::invalid_, p_min = sistrip::invalid_;
0129 float r_sum = 0., r_sum2 = 0., r_max = -1. * sistrip::invalid_, r_min = sistrip::invalid_;
0130
0131
0132 for (uint16_t istr = 0; istr < 128; istr++) {
0133 uint16_t strip = iapv * 128 + istr;
0134
0135
0136 if (peds_histo) {
0137 if (peds_histo->GetBinEntries(strip + 1)) {
0138 anal->peds_[iapv][istr] = peds_histo->GetBinContent(strip + 1);
0139 p_sum += anal->peds_[iapv][istr];
0140 p_sum2 += (anal->peds_[iapv][istr] * anal->peds_[iapv][istr]);
0141 if (anal->peds_[iapv][istr] > p_max) {
0142 p_max = anal->peds_[iapv][istr];
0143 }
0144 if (anal->peds_[iapv][istr] < p_min) {
0145 p_min = anal->peds_[iapv][istr];
0146 }
0147 }
0148 }
0149
0150
0151 if (!anal->legacy_) {
0152 if (raw_histo) {
0153 if (raw_histo->GetBinEntries(strip + 1)) {
0154 anal->raw_[iapv][istr] = raw_histo->GetBinContent(strip + 1);
0155 r_sum += anal->raw_[iapv][istr];
0156 r_sum2 += (anal->raw_[iapv][istr] * anal->raw_[iapv][istr]);
0157 if (anal->raw_[iapv][istr] > r_max) {
0158 r_max = anal->raw_[iapv][istr];
0159 }
0160 if (anal->raw_[iapv][istr] < r_min) {
0161 r_min = anal->raw_[iapv][istr];
0162 }
0163 }
0164 }
0165 } else {
0166 if (peds_histo) {
0167 if (peds_histo->GetBinEntries(strip + 1)) {
0168 anal->raw_[iapv][istr] = raw_histo->GetBinError(strip + 1);
0169 r_sum += anal->raw_[iapv][istr];
0170 r_sum2 += (anal->raw_[iapv][istr] * anal->raw_[iapv][istr]);
0171 if (anal->raw_[iapv][istr] > r_max) {
0172 r_max = anal->raw_[iapv][istr];
0173 }
0174 if (anal->raw_[iapv][istr] < r_min) {
0175 r_min = anal->raw_[iapv][istr];
0176 }
0177 }
0178 }
0179 }
0180
0181 }
0182
0183
0184 if (!anal->peds_[iapv].empty()) {
0185 p_sum /= static_cast<float>(anal->peds_[iapv].size());
0186 p_sum2 /= static_cast<float>(anal->peds_[iapv].size());
0187 anal->pedsMean_[iapv] = p_sum;
0188 anal->pedsSpread_[iapv] = sqrt(fabs(p_sum2 - p_sum * p_sum));
0189 }
0190
0191
0192 if (!anal->raw_[iapv].empty()) {
0193 r_sum /= static_cast<float>(anal->raw_[iapv].size());
0194 r_sum2 /= static_cast<float>(anal->raw_[iapv].size());
0195 anal->rawMean_[iapv] = r_sum;
0196 anal->rawSpread_[iapv] = sqrt(fabs(r_sum2 - r_sum * r_sum));
0197 }
0198
0199
0200 if (p_max > -1. * sistrip::maximum_) {
0201 anal->pedsMax_[iapv] = p_max;
0202 }
0203 if (p_min < 1. * sistrip::maximum_) {
0204 anal->pedsMin_[iapv] = p_min;
0205 }
0206 if (r_max > -1. * sistrip::maximum_) {
0207 anal->rawMax_[iapv] = r_max;
0208 }
0209 if (r_min < 1. * sistrip::maximum_) {
0210 anal->rawMin_[iapv] = r_min;
0211 }
0212
0213 }
0214 }