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
|
#include "DQM/SiStripCommissioningClients/interface/FastFedCablingHistograms.h"
#include "CondFormats/SiStripObjects/interface/FastFedCablingAnalysis.h"
#include "DQM/SiStripCommissioningAnalysis/interface/FastFedCablingAlgorithm.h"
#include "DQM/SiStripCommissioningSummary/interface/FastFedCablingSummaryFactory.h"
#include "DQM/SiStripCommon/interface/ExtractTObject.h"
#include "DataFormats/SiStripCommon/interface/SiStripConstants.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "TProfile.h"
#include <iostream>
#include <memory>
#include <sstream>
#include <iomanip>
using namespace std;
using namespace sistrip;
// -----------------------------------------------------------------------------
/** */
FastFedCablingHistograms::FastFedCablingHistograms(const edm::ParameterSet& pset, DQMStore* bei)
: CommissioningHistograms(
pset.getParameter<edm::ParameterSet>("FastFedCablingParameters"), bei, sistrip::FAST_CABLING) {
factory_ = std::make_unique<FastFedCablingSummaryFactory>();
LogTrace(mlDqmClient_) << "[FastFedCablingHistograms::" << __func__ << "]"
<< " Constructing object...";
}
// -----------------------------------------------------------------------------
/** */
FastFedCablingHistograms::~FastFedCablingHistograms() {
LogTrace(mlDqmClient_) << "[FastFedCablingHistograms::" << __func__ << "]"
<< " Destructing object...";
}
// -----------------------------------------------------------------------------
/** */
void FastFedCablingHistograms::histoAnalysis(bool debug) {
LogTrace(mlDqmClient_) << "[FastFedCablingHistograms::" << __func__ << "]";
// Some initialisation
uint16_t valid = 0;
HistosMap::const_iterator iter;
Analyses::iterator ianal;
std::map<std::string, uint16_t> errors;
// Clear map holding analysis objects
for (ianal = data().begin(); ianal != data().end(); ianal++) {
if (ianal->second) {
delete ianal->second;
}
}
data().clear();
// Iterate through map containing histograms
for (iter = histos().begin(); iter != histos().end(); iter++) {
// Check vector of histos is not empty
if (iter->second.empty()) {
edm::LogWarning(mlDqmClient_) << "[FastFedCablingHistograms::" << __func__ << "]"
<< " Zero histograms found!";
continue;
}
// Retrieve pointers to histos
std::vector<TH1*> profs;
Histos::const_iterator ihis = iter->second.begin();
for (; ihis != iter->second.end(); ihis++) {
TProfile* prof = ExtractTObject<TProfile>().extract((*ihis)->me_);
if (prof) {
profs.push_back(prof);
}
}
// Perform histo analysis
FastFedCablingAnalysis* anal = new FastFedCablingAnalysis(iter->first);
FastFedCablingAlgorithm algo(this->pset(), anal);
FedToFecMap::const_iterator ifed = mapping().find(iter->first);
if (ifed != mapping().end()) {
anal->fecKey(ifed->second);
}
algo.analysis(profs);
data()[iter->first] = anal;
if (anal->isValid()) {
valid++;
}
if (!anal->getErrorCodes().empty()) {
errors[anal->getErrorCodes()[0]]++;
}
}
if (!histos().empty()) {
edm::LogVerbatim(mlDqmClient_) << "[FastFedCablingHistograms::" << __func__ << "]"
<< " Analyzed histograms for " << histos().size() << " FED channels, of which "
<< valid << " (" << 100 * valid / histos().size() << "%) are valid.";
if (!errors.empty()) {
uint16_t count = 0;
std::stringstream ss;
ss << std::endl;
std::map<std::string, uint16_t>::const_iterator ii;
for (ii = errors.begin(); ii != errors.end(); ++ii) {
ss << " " << ii->first << ": " << ii->second << std::endl;
count += ii->second;
}
edm::LogWarning(mlDqmClient_) << "[FastFedCablingHistograms::" << __func__ << "]"
<< " Found " << count << " error strings: " << ss.str();
}
} else {
edm::LogWarning(mlDqmClient_) << "[FastFedCablingHistograms::" << __func__ << "]"
<< " No histograms to analyze!";
}
}
// -----------------------------------------------------------------------------
/** */
void FastFedCablingHistograms::printAnalyses() {
Analyses::iterator ianal = data().begin();
Analyses::iterator janal = data().end();
for (; ianal != janal; ++ianal) {
FastFedCablingAnalysis* anal = dynamic_cast<FastFedCablingAnalysis*>(ianal->second);
if (!anal) {
edm::LogError(mlDqmClient_) << "[FastFedCablingHistograms::" << __func__ << "]"
<< " NULL pointer to analysis object!";
continue;
}
std::stringstream ss;
anal->print(ss);
if (anal->isValid() && !(anal->isDirty()) && !(anal->badTrimDac())) {
LogTrace(mlDqmClient_) << ss.str();
} else {
edm::LogWarning(mlDqmClient_) << ss.str();
}
}
}
// -----------------------------------------------------------------------------
/** */
void FastFedCablingHistograms::printSummary() {
std::stringstream good;
std::stringstream bad;
Analyses::iterator ianal = data().begin();
Analyses::iterator janal = data().end();
for (; ianal != janal; ++ianal) {
FastFedCablingAnalysis* anal = dynamic_cast<FastFedCablingAnalysis*>(ianal->second);
if (!anal) {
edm::LogError(mlDqmClient_) << "[FastFedCablingHistograms::" << __func__ << "]"
<< " NULL pointer to analysis object!";
continue;
}
if (anal->isValid() && !(anal->isDirty()) && !(anal->badTrimDac())) {
anal->summary(good);
} else {
anal->summary(bad);
}
}
if (good.str().empty()) {
good << "None found!";
}
LogTrace(mlDqmClient_) << "[FastFedCablingHistograms::" << __func__ << "]"
<< " Printing summary of good analyses:"
<< "\n"
<< good.str();
if (bad.str().empty()) {
return;
} //@@ bad << "None found!"; }
LogTrace(mlDqmClient_) << "[FastFedCablingHistograms::" << __func__ << "]"
<< " Printing summary of bad analyses:"
<< "\n"
<< bad.str();
}
|