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
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
|
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/global/EDProducer.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/Framework/interface/ESHandle.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/ServiceRegistry/interface/Service.h"
#include "FWCore/Utilities/interface/ESGetToken.h"
#include "EventFilter/HcalRawToDigi/interface/HcalUHTRData.h"
#include "DataFormats/HcalDigi/interface/HcalQIESample.h"
#include "CondFormats/HcalObjects/interface/HcalElectronicsMap.h"
#include "DataFormats/Common/interface/Handle.h"
#include "DataFormats/FEDRawData/interface/FEDRawDataCollection.h"
#include "DataFormats/FEDRawData/interface/FEDHeader.h"
#include "DataFormats/FEDRawData/interface/FEDTrailer.h"
#include "DataFormats/FEDRawData/interface/FEDNumbering.h"
#include "DataFormats/FEDRawData/interface/FEDRawData.h"
#include "FWCore/Utilities/interface/CRC16.h"
#include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
#include "CalibFormats/HcalObjects/interface/HcalDbRecord.h"
#include "PackerHelp.h"
#include <fstream>
#include <iostream>
#include <memory>
#include <sstream>
#include <string>
/* QUESTION: what about dual FED readout? */
/* QUESTION: what do I do if the number of 16-bit words
are not divisible by 4? -- these need to
fit into the 64-bit words of the FEDRawDataFormat */
using namespace std;
class HcalDigiToRawuHTR : public edm::global::EDProducer<> {
public:
explicit HcalDigiToRawuHTR(const edm::ParameterSet&);
~HcalDigiToRawuHTR() override;
void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
private:
const int _verbosity;
const vector<int> tdc1_;
const vector<int> tdc2_;
const bool packHBTDC_;
static constexpr int tdcmax_ = 49;
const edm::EDGetTokenT<HcalDataFrameContainer<QIE10DataFrame>> tok_QIE10DigiCollection_;
const edm::EDGetTokenT<HcalDataFrameContainer<QIE11DataFrame>> tok_QIE11DigiCollection_;
const edm::EDGetTokenT<HBHEDigiCollection> tok_HBHEDigiCollection_;
const edm::EDGetTokenT<HFDigiCollection> tok_HFDigiCollection_;
const edm::EDGetTokenT<HcalTrigPrimDigiCollection> tok_TPDigiCollection_;
const edm::ESGetToken<HcalElectronicsMap, HcalElectronicsMapRcd> tok_electronicsMap_;
const bool premix_;
};
HcalDigiToRawuHTR::HcalDigiToRawuHTR(const edm::ParameterSet& iConfig)
: _verbosity(iConfig.getUntrackedParameter<int>("Verbosity", 0)),
tdc1_(iConfig.getParameter<vector<int>>("tdc1")),
tdc2_(iConfig.getParameter<vector<int>>("tdc2")),
packHBTDC_(iConfig.getParameter<bool>("packHBTDC")),
tok_QIE10DigiCollection_(
consumes<HcalDataFrameContainer<QIE10DataFrame>>(iConfig.getParameter<edm::InputTag>("QIE10"))),
tok_QIE11DigiCollection_(
consumes<HcalDataFrameContainer<QIE11DataFrame>>(iConfig.getParameter<edm::InputTag>("QIE11"))),
tok_HBHEDigiCollection_(consumes<HBHEDigiCollection>(iConfig.getParameter<edm::InputTag>("HBHEqie8"))),
tok_HFDigiCollection_(consumes<HFDigiCollection>(iConfig.getParameter<edm::InputTag>("HFqie8"))),
tok_TPDigiCollection_(consumes<HcalTrigPrimDigiCollection>(iConfig.getParameter<edm::InputTag>("TP"))),
tok_electronicsMap_(esConsumes<HcalElectronicsMap, HcalElectronicsMapRcd>(
edm::ESInputTag("", iConfig.getParameter<std::string>("ElectronicsMap")))),
premix_(iConfig.getParameter<bool>("premix")) {
produces<FEDRawDataCollection>("");
for (size_t i = 0; i < tdc1_.size(); i++) {
if (!(tdc1_.at(i) >= 0 && tdc1_.at(i) <= tdc2_.at(i) && tdc2_.at(i) <= tdcmax_))
edm::LogWarning("HcalDigiToRawuHTR")
<< " incorrect TDC ranges " << i << "-th element: " << tdc1_.at(i) << ", " << tdc2_.at(i) << ", " << tdcmax_;
}
}
HcalDigiToRawuHTR::~HcalDigiToRawuHTR() {}
void HcalDigiToRawuHTR::produce(edm::StreamID id, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
using namespace edm;
edm::ESHandle<HcalElectronicsMap> item = iSetup.getHandle(tok_electronicsMap_);
const HcalElectronicsMap* readoutMap = item.product();
//collection to be inserted into event
std::unique_ptr<FEDRawDataCollection> fed_buffers(new FEDRawDataCollection());
//
// Extracting All the Collections containing useful Info
edm::Handle<QIE10DigiCollection> qie10DigiCollection;
edm::Handle<QIE11DigiCollection> qie11DigiCollection;
edm::Handle<HBHEDigiCollection> hbheDigiCollection;
edm::Handle<HFDigiCollection> hfDigiCollection;
edm::Handle<HcalTrigPrimDigiCollection> tpDigiCollection;
iEvent.getByToken(tok_QIE10DigiCollection_, qie10DigiCollection);
iEvent.getByToken(tok_QIE11DigiCollection_, qie11DigiCollection);
iEvent.getByToken(tok_HBHEDigiCollection_, hbheDigiCollection);
iEvent.getByToken(tok_HFDigiCollection_, hfDigiCollection);
iEvent.getByToken(tok_TPDigiCollection_, tpDigiCollection);
// first argument is the fedid (minFEDID+crateId)
map<int, unique_ptr<HCalFED>> fedMap;
// - - - - - - - - - - - - - - - - - - - - - - - - - - -
// QIE10 precision data
// - - - - - - - - - - - - - - - - - - - - - - - - - - -
UHTRpacker uhtrs;
// loop over each digi and allocate memory for each
if (qie10DigiCollection.isValid()) {
const QIE10DigiCollection& qie10dc = *(qie10DigiCollection);
for (unsigned int j = 0; j < qie10dc.size(); j++) {
QIE10DataFrame qiedf = static_cast<QIE10DataFrame>(qie10dc[j]);
DetId detid = qiedf.detid();
HcalElectronicsId eid(readoutMap->lookup(detid));
int crateId = eid.crateId();
int slotId = eid.slot();
int uhtrIndex = ((slotId & 0xF) << 8) | (crateId & 0xFF);
int presamples = qiedf.presamples();
/* Defining a custom index that will encode only
the information about the crate and slot of a
given channel: crate: bits 0-7
slot: bits 8-12 */
if (!uhtrs.exist(uhtrIndex)) {
uhtrs.newUHTR(uhtrIndex, presamples);
}
uhtrs.addChannel(uhtrIndex, qiedf, readoutMap, _verbosity);
}
}
// - - - - - - - - - - - - - - - - - - - - - - - - - - -
// QIE11 precision data
// - - - - - - - - - - - - - - - - - - - - - - - - - - -
//UHTRpacker uhtrs;
// loop over each digi and allocate memory for each
if (qie11DigiCollection.isValid()) {
const QIE11DigiCollection& qie11dc = *(qie11DigiCollection);
for (unsigned int j = 0; j < qie11dc.size(); j++) {
QIE11DataFrame qiedf = static_cast<QIE11DataFrame>(qie11dc[j]);
DetId detid = qiedf.detid();
HcalElectronicsId eid(readoutMap->lookup(detid));
int crateId = eid.crateId();
int slotId = eid.slot();
int uhtrIndex = ((slotId & 0xF) << 8) | (crateId & 0xFF);
int presamples = qiedf.presamples();
// convert to hb qie data if hb
if (packHBTDC_ && HcalDetId(detid.rawId()).subdet() == HcalSubdetector::HcalBarrel)
qiedf = convertHB(qiedf, tdc1_, tdc2_, tdcmax_);
if (!uhtrs.exist(uhtrIndex)) {
uhtrs.newUHTR(uhtrIndex, presamples);
}
uhtrs.addChannel(uhtrIndex, qiedf, readoutMap, _verbosity);
}
}
// - - - - - - - - - - - - - - - - - - - - - - - - - - -
// HF (QIE8) precision data
// - - - - - - - - - - - - - - - - - - - - - - - - - - -
// loop over each digi and allocate memory for each
if (hfDigiCollection.isValid()) {
const HFDigiCollection& qie8hfdc = *(hfDigiCollection);
for (HFDigiCollection::const_iterator qiedf = qie8hfdc.begin(); qiedf != qie8hfdc.end(); qiedf++) {
DetId detid = qiedf->id();
HcalElectronicsId eid(readoutMap->lookup(detid));
int crateId = eid.crateId();
int slotId = eid.slot();
int uhtrIndex = (crateId & 0xFF) | ((slotId & 0xF) << 8);
int presamples = qiedf->presamples();
if (!uhtrs.exist(uhtrIndex)) {
uhtrs.newUHTR(uhtrIndex, presamples);
}
uhtrs.addChannel(uhtrIndex, qiedf, readoutMap, premix_, _verbosity);
}
}
// - - - - - - - - - - - - - - - - - - - - - - - - - - -
// HBHE (QIE8) precision data
// - - - - - - - - - - - - - - - - - - - - - - - - - - -
// loop over each digi and allocate memory for each
if (hbheDigiCollection.isValid()) {
const HBHEDigiCollection& qie8hbhedc = *(hbheDigiCollection);
for (HBHEDigiCollection::const_iterator qiedf = qie8hbhedc.begin(); qiedf != qie8hbhedc.end(); qiedf++) {
DetId detid = qiedf->id();
HcalElectronicsId eid(readoutMap->lookup(detid));
int crateId = eid.crateId();
int slotId = eid.slot();
int uhtrIndex = (crateId & 0xFF) | ((slotId & 0xF) << 8);
int presamples = qiedf->presamples();
if (!uhtrs.exist(uhtrIndex)) {
uhtrs.newUHTR(uhtrIndex, presamples);
}
uhtrs.addChannel(uhtrIndex, qiedf, readoutMap, premix_, _verbosity);
}
}
// - - - - - - - - - - - - - - - - - - - - - - - - - - -
// TP data
// - - - - - - - - - - - - - - - - - - - - - - - - - - -
// loop over each digi and allocate memory for each
if (tpDigiCollection.isValid()) {
const HcalTrigPrimDigiCollection& qietpdc = *(tpDigiCollection);
for (HcalTrigPrimDigiCollection::const_iterator qiedf = qietpdc.begin(); qiedf != qietpdc.end(); qiedf++) {
DetId detid = qiedf->id();
HcalElectronicsId eid(readoutMap->lookupTrigger(detid));
int crateId = eid.crateId();
int slotId = eid.slot();
int uhtrIndex = (crateId & 0xFF) | ((slotId & 0xF) << 8);
int ilink = eid.fiberIndex();
int itower = eid.fiberChanId();
int channelid = (itower & 0xF) | ((ilink & 0xF) << 4);
int presamples = qiedf->presamples();
if (!uhtrs.exist(uhtrIndex)) {
uhtrs.newUHTR(uhtrIndex, presamples);
}
uhtrs.addChannel(uhtrIndex, qiedf, channelid, _verbosity);
}
}
// -----------------------------------------------------
// -----------------------------------------------------
// loop over each uHTR and format data
// -----------------------------------------------------
// -----------------------------------------------------
// loop over each uHTR and format data
for (UHTRpacker::UHTRMap::iterator uhtr = uhtrs.uhtrs.begin(); uhtr != uhtrs.uhtrs.end(); ++uhtr) {
uint64_t crateId = (uhtr->first) & 0xFF;
uint64_t slotId = (uhtr->first & 0xF00) >> 8;
uhtrs.finalizeHeadTail(&(uhtr->second), _verbosity);
int fedId = FEDNumbering::MINHCALuTCAFEDID + crateId;
if (fedMap.find(fedId) == fedMap.end()) {
/* QUESTION: where should the orbit number come from? */
fedMap[fedId] =
std::make_unique<HCalFED>(fedId, iEvent.id().event(), iEvent.orbitNumber(), iEvent.bunchCrossing());
}
fedMap[fedId]->addUHTR(uhtr->second, crateId, slotId);
} // end loop over uhtr containers
/* ------------------------------------------------------
------------------------------------------------------
putting together the FEDRawDataCollection
------------------------------------------------------
------------------------------------------------------ */
for (map<int, unique_ptr<HCalFED>>::iterator fed = fedMap.begin(); fed != fedMap.end(); ++fed) {
int fedId = fed->first;
auto& rawData = fed_buffers->FEDData(fedId);
fed->second->formatFEDdata(rawData);
FEDHeader hcalFEDHeader(rawData.data());
hcalFEDHeader.set(rawData.data(), 1, iEvent.id().event(), iEvent.bunchCrossing(), fedId);
FEDTrailer hcalFEDTrailer(rawData.data() + (rawData.size() - 8));
hcalFEDTrailer.set(rawData.data() + (rawData.size() - 8),
rawData.size() / 8,
evf::compute_crc(rawData.data(), rawData.size()),
0,
0);
} // end loop over FEDs with data
iEvent.put(std::move(fed_buffers));
}
// ------------ method fills 'descriptions' with the allowed parameters for the module ------------
void HcalDigiToRawuHTR::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
//The following says we do not know what parameters are allowed so do no validation
// Please change this to state exactly what you do use, even if it is no parameters
edm::ParameterSetDescription desc;
desc.addUntracked<int>("Verbosity", 0);
desc.add<vector<int>>("tdc1", {12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12,
12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12,
12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12});
desc.add<vector<int>>("tdc2", {14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14,
14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14,
14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14});
desc.add<bool>("packHBTDC", true);
desc.add<std::string>("ElectronicsMap", "");
desc.add<edm::InputTag>("QIE10", edm::InputTag("simHcalDigis", "HFQIE10DigiCollection"));
desc.add<edm::InputTag>("QIE11", edm::InputTag("simHcalDigis", "HBHEQIE11DigiCollection"));
desc.add<edm::InputTag>("HBHEqie8", edm::InputTag("simHcalDigis"));
desc.add<edm::InputTag>("HFqie8", edm::InputTag("simHcalDigis"));
desc.add<edm::InputTag>("TP", edm::InputTag("simHcalTriggerPrimitiveDigis"));
desc.add<bool>("premix", false);
descriptions.add("hcalDigiToRawuHTR", desc);
descriptions.addDefault(desc);
}
//define this as a plug-in
DEFINE_FWK_MODULE(HcalDigiToRawuHTR);
|