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
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
|
#include <array>
#include <string>
#include <iostream>
#include <map>
#include "FWCore/Framework/interface/one/EDAnalyzer.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "CondFormats/PCLConfig/interface/AlignPCLThresholds.h"
#include "CondFormats/PCLConfig/interface/AlignPCLThresholdsHG.h"
#include "CondFormats/DataRecord/interface/AlignPCLThresholdsRcd.h"
#include "CondFormats/DataRecord/interface/AlignPCLThresholdsHGRcd.h"
#include "HLTrigger/HLTcore/interface/defaultModuleLabel.h"
namespace edmtest {
template <typename T, typename R>
class AlignPCLThresholdsReader : public edm::one::EDAnalyzer<> {
public:
explicit AlignPCLThresholdsReader(edm::ParameterSet const& p);
~AlignPCLThresholdsReader() override;
static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
private:
void analyze(const edm::Event& e, const edm::EventSetup& c) override;
// ----------member data ---------------------------
const edm::ESGetToken<T, R> thresholdToken_;
const bool printdebug_;
const std::string formatedOutput_;
};
template <typename T, typename R>
AlignPCLThresholdsReader<T, R>::AlignPCLThresholdsReader(edm::ParameterSet const& p)
: thresholdToken_(esConsumes()),
printdebug_(p.getUntrackedParameter<bool>("printDebug", true)),
formatedOutput_(p.getUntrackedParameter<std::string>("outputFile", "")) {
edm::LogInfo("AlignPCLThresholdsReader") << "AlignPCLThresholdsReader" << std::endl;
}
template <typename T, typename R>
AlignPCLThresholdsReader<T, R>::~AlignPCLThresholdsReader() {
edm::LogInfo("AlignPCLThresholdsReader") << "~AlignPCLThresholdsReader " << std::endl;
}
template <typename T, typename R>
void AlignPCLThresholdsReader<T, R>::analyze(const edm::Event& e, const edm::EventSetup& context) {
edm::LogInfo("AlignPCLThresholdsReader") << "### AlignPCLThresholdsReader::analyze ###" << std::endl;
edm::LogInfo("AlignPCLThresholdsReader") << " I AM IN RUN NUMBER " << e.id().run() << std::endl;
edm::LogInfo("AlignPCLThresholdsReader") << " ---EVENT NUMBER " << e.id().event() << std::endl;
edm::eventsetup::EventSetupRecordKey inputKey = edm::eventsetup::EventSetupRecordKey::makeKey<R>();
edm::eventsetup::EventSetupRecordKey recordKey(
edm::eventsetup::EventSetupRecordKey::TypeTag::findType(inputKey.type().name()));
if (recordKey.type() == edm::eventsetup::EventSetupRecordKey::TypeTag()) {
//record not found
edm::LogInfo("AlignPCLThresholdsReader")
<< "Record \"" << inputKey.type().name() << "\" does not exist " << std::endl;
}
//this part gets the handle of the event source and the record (i.e. the Database)
edm::ESHandle<T> thresholdHandle = context.getHandle(thresholdToken_);
edm::LogInfo("AlignPCLThresholdsReader") << "got eshandle" << std::endl;
if (!thresholdHandle.isValid()) {
edm::LogError("AlignPCLThresholdsReader") << " Could not get Handle" << std::endl;
return;
}
const T* thresholds = thresholdHandle.product();
edm::LogInfo("AlignPCLThresholdsReader") << "got AlignPCLThresholds* " << std::endl;
edm::LogInfo("AlignPCLThresholdsReader") << "print pointer address : ";
edm::LogInfo("AlignPCLThresholdsReader") << thresholds << std::endl;
edm::LogInfo("AlignPCLThresholdsReader") << "Size " << thresholds->size() << std::endl;
edm::LogInfo("AlignPCLThresholdsReader") << "Content of myThresholds " << std::endl;
// use built-in method in the CondFormat to print the content
if (thresholds && printdebug_) {
thresholds->printAll();
}
FILE* pFile = nullptr;
if (!formatedOutput_.empty())
pFile = fopen(formatedOutput_.c_str(), "w");
if (pFile) {
fprintf(pFile, "AlignPCLThresholds::printAll() \n");
fprintf(pFile,
" ======================================================================================================="
"============\n");
fprintf(pFile, "N records cut: %i \n", thresholds->getNrecords());
AlignPCLThresholds::threshold_map m_thresholds = thresholds->getThreshold_Map();
AlignPCLThresholdsHG::param_map m_floatMap{};
if constexpr (std::is_same_v<T, AlignPCLThresholdsHG>) {
m_floatMap = thresholds->getFloatMap();
}
for (auto it = m_thresholds.begin(); it != m_thresholds.end(); ++it) {
bool hasFractionCut = (m_floatMap.find(it->first) != m_floatMap.end());
fprintf(pFile,
" ====================================================================================================="
"==============\n");
fprintf(pFile, "key : %s \n", (it->first).c_str());
fprintf(pFile, "- Xcut : %8.3f um ", (it->second).getXcut());
fprintf(pFile, "| sigXcut : %8.3f ", (it->second).getSigXcut());
fprintf(pFile, "| maxMoveXcut : %8.3f um ", (it->second).getMaxMoveXcut());
fprintf(pFile, "| ErrorXcut : %8.3f um ", (it->second).getErrorXcut());
if constexpr (std::is_same_v<T, AlignPCLThresholdsHG>) {
if (hasFractionCut) {
fprintf(pFile,
"| X_fractionCut : %8.3f \n",
thresholds->getFractionCut(it->first, AlignPCLThresholds::coordType::X));
} else {
fprintf(pFile, "\n");
}
} else {
fprintf(pFile, "\n");
}
fprintf(pFile, "- thetaXcut : %8.3f urad ", (it->second).getThetaXcut());
fprintf(pFile, "| sigThetaXcut : %8.3f ", (it->second).getSigThetaXcut());
fprintf(pFile, "| maxMoveThetaXcut : %8.3f urad ", (it->second).getMaxMoveThetaXcut());
fprintf(pFile, "| ErrorThetaXcut : %8.3f urad ", (it->second).getErrorThetaXcut());
if constexpr (std::is_same_v<T, AlignPCLThresholdsHG>) {
if (hasFractionCut) {
fprintf(pFile,
"| thetaX_fractionCut : %8.3f \n",
thresholds->getFractionCut(it->first, AlignPCLThresholds::coordType::theta_X));
} else {
fprintf(pFile, "\n");
}
} else {
fprintf(pFile, "\n");
}
fprintf(pFile, "- Ycut : %8.3f um ", (it->second).getYcut());
fprintf(pFile, "| sigYcut : %8.3f ", (it->second).getSigXcut());
fprintf(pFile, "| maxMoveYcut : %8.3f um ", (it->second).getMaxMoveYcut());
fprintf(pFile, "| ErrorYcut : %8.3f um ", (it->second).getErrorYcut());
if constexpr (std::is_same_v<T, AlignPCLThresholdsHG>) {
if (hasFractionCut) {
fprintf(pFile,
"| Y_fractionCut : %8.3f \n",
thresholds->getFractionCut(it->first, AlignPCLThresholds::coordType::Y));
} else {
fprintf(pFile, "\n");
}
} else {
fprintf(pFile, "\n");
}
fprintf(pFile, "- thetaYcut : %8.3f urad ", (it->second).getThetaYcut());
fprintf(pFile, "| sigThetaYcut : %8.3f ", (it->second).getSigThetaYcut());
fprintf(pFile, "| maxMoveThetaYcut : %8.3f urad ", (it->second).getMaxMoveThetaYcut());
fprintf(pFile, "| ErrorThetaYcut : %8.3f urad ", (it->second).getErrorThetaYcut());
if constexpr (std::is_same_v<T, AlignPCLThresholdsHG>) {
if (hasFractionCut) {
fprintf(pFile,
"| thetaY_fractionCut : %8.3f \n",
thresholds->getFractionCut(it->first, AlignPCLThresholds::coordType::theta_Y));
} else {
fprintf(pFile, "\n");
}
} else {
fprintf(pFile, "\n");
}
fprintf(pFile, "- Zcut : %8.3f um ", (it->second).getZcut());
fprintf(pFile, "| sigZcut : %8.3f ", (it->second).getSigZcut());
fprintf(pFile, "| maxMoveZcut : %8.3f um ", (it->second).getMaxMoveZcut());
fprintf(pFile, "| ErrorZcut : %8.3f um ", (it->second).getErrorZcut());
if constexpr (std::is_same_v<T, AlignPCLThresholdsHG>) {
if (hasFractionCut) {
fprintf(pFile,
"| Z_fractionCut : %8.3f \n",
thresholds->getFractionCut(it->first, AlignPCLThresholds::coordType::Z));
} else {
fprintf(pFile, "\n");
}
} else {
fprintf(pFile, "\n");
}
fprintf(pFile, "- thetaZcut : %8.3f urad ", (it->second).getThetaZcut());
fprintf(pFile, "| sigThetaZcut : %8.3f ", (it->second).getSigThetaZcut());
fprintf(pFile, "| maxMoveThetaZcut : %8.3f urad ", (it->second).getMaxMoveThetaZcut());
fprintf(pFile, "| ErrorThetaZcut : %8.3f urad ", (it->second).getErrorThetaZcut());
if constexpr (std::is_same_v<T, AlignPCLThresholdsHG>) {
if (hasFractionCut) {
fprintf(pFile,
"| thetaZ_fractionCut : %8.3f \n",
thresholds->getFractionCut(it->first, AlignPCLThresholds::coordType::theta_Z));
} else {
fprintf(pFile, "\n");
}
} else {
fprintf(pFile, "\n");
}
if ((it->second).hasExtraDOF()) {
for (unsigned int j = 0; j < (it->second).extraDOFSize(); j++) {
std::array<float, 4> extraDOFCuts = thresholds->getExtraDOFCutsForAlignable(it->first, j);
fprintf(pFile,
"Extra DOF: %i with label %s \n ",
j,
thresholds->getExtraDOFLabelForAlignable(it->first, j).c_str());
fprintf(pFile, "- cut : %8.3f ", extraDOFCuts.at(0));
fprintf(pFile, "| sigCut : %8.3f ", extraDOFCuts.at(1));
fprintf(pFile, "| maxMoveCut : %8.3f ", extraDOFCuts.at(2));
fprintf(pFile, "| maxErrorCut : %8.3f \n ", extraDOFCuts.at(3));
}
}
}
}
}
template <typename T, typename R>
void AlignPCLThresholdsReader<T, R>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
desc.setComment("Reads payloads of type AlignPCLThresholds");
desc.addUntracked<bool>("printDebug", true);
desc.addUntracked<std::string>("outputFile", "");
descriptions.add(defaultModuleLabel<AlignPCLThresholdsReader<T, R>>(), desc);
}
typedef AlignPCLThresholdsReader<AlignPCLThresholds, AlignPCLThresholdsRcd> AlignPCLThresholdsLGReader;
typedef AlignPCLThresholdsReader<AlignPCLThresholdsHG, AlignPCLThresholdsHGRcd> AlignPCLThresholdsHGReader;
DEFINE_FWK_MODULE(AlignPCLThresholdsLGReader);
DEFINE_FWK_MODULE(AlignPCLThresholdsHGReader);
} // namespace edmtest
|