CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HcalDigiToRawuHTR.cc
Go to the documentation of this file.
6 
10 
13 
15 
23 
26 
27 #include "PackerHelp.h"
28 
29 #include <fstream>
30 #include <iostream>
31 #include <memory>
32 
33 #include <sstream>
34 #include <string>
35 
36 /* QUESTION: what about dual FED readout? */
37 /* QUESTION: what do I do if the number of 16-bit words
38  are not divisible by 4? -- these need to
39  fit into the 64-bit words of the FEDRawDataFormat */
40 
41 using namespace std;
42 
44 public:
45  explicit HcalDigiToRawuHTR(const edm::ParameterSet&);
46  ~HcalDigiToRawuHTR() override;
47 
48  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
49 
50  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
51 
52 private:
53  const int _verbosity;
54  const vector<int> tdc1_;
55  const vector<int> tdc2_;
56  const bool packHBTDC_;
57  static constexpr int tdcmax_ = 49;
58 
65 
66  const bool premix_;
67 };
68 
70  : _verbosity(iConfig.getUntrackedParameter<int>("Verbosity", 0)),
71  tdc1_(iConfig.getParameter<vector<int>>("tdc1")),
72  tdc2_(iConfig.getParameter<vector<int>>("tdc2")),
73  packHBTDC_(iConfig.getParameter<bool>("packHBTDC")),
74  tok_QIE10DigiCollection_(
75  consumes<HcalDataFrameContainer<QIE10DataFrame>>(iConfig.getParameter<edm::InputTag>("QIE10"))),
76  tok_QIE11DigiCollection_(
77  consumes<HcalDataFrameContainer<QIE11DataFrame>>(iConfig.getParameter<edm::InputTag>("QIE11"))),
78  tok_HBHEDigiCollection_(consumes<HBHEDigiCollection>(iConfig.getParameter<edm::InputTag>("HBHEqie8"))),
79  tok_HFDigiCollection_(consumes<HFDigiCollection>(iConfig.getParameter<edm::InputTag>("HFqie8"))),
80  tok_TPDigiCollection_(consumes<HcalTrigPrimDigiCollection>(iConfig.getParameter<edm::InputTag>("TP"))),
82  edm::ESInputTag("", iConfig.getParameter<std::string>("ElectronicsMap")))),
83  premix_(iConfig.getParameter<bool>("premix")) {
84  produces<FEDRawDataCollection>("");
85  for (size_t i = 0; i < tdc1_.size(); i++) {
86  if (!(tdc1_.at(i) >= 0 && tdc1_.at(i) <= tdc2_.at(i) && tdc2_.at(i) <= tdcmax_))
87  edm::LogWarning("HcalDigiToRawuHTR")
88  << " incorrect TDC ranges " << i << "-th element: " << tdc1_.at(i) << ", " << tdc2_.at(i) << ", " << tdcmax_;
89  }
90 }
91 
93 
95  using namespace edm;
96 
98  const HcalElectronicsMap* readoutMap = item.product();
99 
100  //collection to be inserted into event
101  std::unique_ptr<FEDRawDataCollection> fed_buffers(new FEDRawDataCollection());
102 
103  //
104  // Extracting All the Collections containing useful Info
105  edm::Handle<QIE10DigiCollection> qie10DigiCollection;
106  edm::Handle<QIE11DigiCollection> qie11DigiCollection;
107  edm::Handle<HBHEDigiCollection> hbheDigiCollection;
108  edm::Handle<HFDigiCollection> hfDigiCollection;
110  iEvent.getByToken(tok_QIE10DigiCollection_, qie10DigiCollection);
111  iEvent.getByToken(tok_QIE11DigiCollection_, qie11DigiCollection);
112  iEvent.getByToken(tok_HBHEDigiCollection_, hbheDigiCollection);
113  iEvent.getByToken(tok_HFDigiCollection_, hfDigiCollection);
114  iEvent.getByToken(tok_TPDigiCollection_, tpDigiCollection);
115 
116  // first argument is the fedid (minFEDID+crateId)
117  map<int, unique_ptr<HCalFED>> fedMap;
118 
119  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
120  // QIE10 precision data
121  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
122  UHTRpacker uhtrs;
123  // loop over each digi and allocate memory for each
124  if (qie10DigiCollection.isValid()) {
125  const QIE10DigiCollection& qie10dc = *(qie10DigiCollection);
126  for (unsigned int j = 0; j < qie10dc.size(); j++) {
127  QIE10DataFrame qiedf = static_cast<QIE10DataFrame>(qie10dc[j]);
128  DetId detid = qiedf.detid();
129  HcalElectronicsId eid(readoutMap->lookup(detid));
130  int crateId = eid.crateId();
131  int slotId = eid.slot();
132  int uhtrIndex = ((slotId & 0xF) << 8) | (crateId & 0xFF);
133  int presamples = qiedf.presamples();
134 
135  /* Defining a custom index that will encode only
136  the information about the crate and slot of a
137  given channel: crate: bits 0-7
138  slot: bits 8-12 */
139 
140  if (!uhtrs.exist(uhtrIndex)) {
141  uhtrs.newUHTR(uhtrIndex, presamples);
142  }
143  uhtrs.addChannel(uhtrIndex, qiedf, readoutMap, _verbosity);
144  }
145  }
146  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
147  // QIE11 precision data
148  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
149  //UHTRpacker uhtrs;
150  // loop over each digi and allocate memory for each
151  if (qie11DigiCollection.isValid()) {
152  const QIE11DigiCollection& qie11dc = *(qie11DigiCollection);
153  for (unsigned int j = 0; j < qie11dc.size(); j++) {
154  QIE11DataFrame qiedf = static_cast<QIE11DataFrame>(qie11dc[j]);
155  DetId detid = qiedf.detid();
156  HcalElectronicsId eid(readoutMap->lookup(detid));
157  int crateId = eid.crateId();
158  int slotId = eid.slot();
159  int uhtrIndex = ((slotId & 0xF) << 8) | (crateId & 0xFF);
160  int presamples = qiedf.presamples();
161 
162  // convert to hb qie data if hb
163  if (packHBTDC_ && HcalDetId(detid.rawId()).subdet() == HcalSubdetector::HcalBarrel)
164  qiedf = convertHB(qiedf, tdc1_, tdc2_, tdcmax_);
165 
166  if (!uhtrs.exist(uhtrIndex)) {
167  uhtrs.newUHTR(uhtrIndex, presamples);
168  }
169  uhtrs.addChannel(uhtrIndex, qiedf, readoutMap, _verbosity);
170  }
171  }
172  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
173  // HF (QIE8) precision data
174  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
175  // loop over each digi and allocate memory for each
176  if (hfDigiCollection.isValid()) {
177  const HFDigiCollection& qie8hfdc = *(hfDigiCollection);
178  for (HFDigiCollection::const_iterator qiedf = qie8hfdc.begin(); qiedf != qie8hfdc.end(); qiedf++) {
179  DetId detid = qiedf->id();
180 
181  HcalElectronicsId eid(readoutMap->lookup(detid));
182  int crateId = eid.crateId();
183  int slotId = eid.slot();
184  int uhtrIndex = (crateId & 0xFF) | ((slotId & 0xF) << 8);
185  int presamples = qiedf->presamples();
186 
187  if (!uhtrs.exist(uhtrIndex)) {
188  uhtrs.newUHTR(uhtrIndex, presamples);
189  }
190  uhtrs.addChannel(uhtrIndex, qiedf, readoutMap, premix_, _verbosity);
191  }
192  }
193  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
194  // HBHE (QIE8) precision data
195  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
196  // loop over each digi and allocate memory for each
197  if (hbheDigiCollection.isValid()) {
198  const HBHEDigiCollection& qie8hbhedc = *(hbheDigiCollection);
199  for (HBHEDigiCollection::const_iterator qiedf = qie8hbhedc.begin(); qiedf != qie8hbhedc.end(); qiedf++) {
200  DetId detid = qiedf->id();
201 
202  HcalElectronicsId eid(readoutMap->lookup(detid));
203  int crateId = eid.crateId();
204  int slotId = eid.slot();
205  int uhtrIndex = (crateId & 0xFF) | ((slotId & 0xF) << 8);
206  int presamples = qiedf->presamples();
207 
208  if (!uhtrs.exist(uhtrIndex)) {
209  uhtrs.newUHTR(uhtrIndex, presamples);
210  }
211  uhtrs.addChannel(uhtrIndex, qiedf, readoutMap, premix_, _verbosity);
212  }
213  }
214  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
215  // TP data
216  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
217  // loop over each digi and allocate memory for each
218  if (tpDigiCollection.isValid()) {
219  const HcalTrigPrimDigiCollection& qietpdc = *(tpDigiCollection);
220  for (HcalTrigPrimDigiCollection::const_iterator qiedf = qietpdc.begin(); qiedf != qietpdc.end(); qiedf++) {
221  DetId detid = qiedf->id();
222  HcalElectronicsId eid(readoutMap->lookupTrigger(detid));
223 
224  int crateId = eid.crateId();
225  int slotId = eid.slot();
226  int uhtrIndex = (crateId & 0xFF) | ((slotId & 0xF) << 8);
227  int ilink = eid.fiberIndex();
228  int itower = eid.fiberChanId();
229  int channelid = (itower & 0xF) | ((ilink & 0xF) << 4);
230  int presamples = qiedf->presamples();
231 
232  if (!uhtrs.exist(uhtrIndex)) {
233  uhtrs.newUHTR(uhtrIndex, presamples);
234  }
235  uhtrs.addChannel(uhtrIndex, qiedf, channelid, _verbosity);
236  }
237  }
238  // -----------------------------------------------------
239  // -----------------------------------------------------
240  // loop over each uHTR and format data
241  // -----------------------------------------------------
242  // -----------------------------------------------------
243  // loop over each uHTR and format data
244  int idxuhtr = -1;
245  for (UHTRpacker::UHTRMap::iterator uhtr = uhtrs.uhtrs.begin(); uhtr != uhtrs.uhtrs.end(); ++uhtr) {
246  idxuhtr++;
247 
248  uint64_t crateId = (uhtr->first) & 0xFF;
249  uint64_t slotId = (uhtr->first & 0xF00) >> 8;
250 
251  uhtrs.finalizeHeadTail(&(uhtr->second), _verbosity);
252  int fedId = FEDNumbering::MINHCALuTCAFEDID + crateId;
253  if (fedMap.find(fedId) == fedMap.end()) {
254  /* QUESTION: where should the orbit number come from? */
255  fedMap[fedId] =
256  std::make_unique<HCalFED>(fedId, iEvent.id().event(), iEvent.orbitNumber(), iEvent.bunchCrossing());
257  }
258  fedMap[fedId]->addUHTR(uhtr->second, crateId, slotId);
259  } // end loop over uhtr containers
260 
261  /* ------------------------------------------------------
262  ------------------------------------------------------
263  putting together the FEDRawDataCollection
264  ------------------------------------------------------
265  ------------------------------------------------------ */
266  for (map<int, unique_ptr<HCalFED>>::iterator fed = fedMap.begin(); fed != fedMap.end(); ++fed) {
267  int fedId = fed->first;
268 
269  auto& rawData = fed_buffers->FEDData(fedId);
270  fed->second->formatFEDdata(rawData);
271 
272  FEDHeader hcalFEDHeader(rawData.data());
273  hcalFEDHeader.set(rawData.data(), 1, iEvent.id().event(), iEvent.bunchCrossing(), fedId);
274  FEDTrailer hcalFEDTrailer(rawData.data() + (rawData.size() - 8));
275  hcalFEDTrailer.set(rawData.data() + (rawData.size() - 8),
276  rawData.size() / 8,
277  evf::compute_crc(rawData.data(), rawData.size()),
278  0,
279  0);
280 
281  } // end loop over FEDs with data
282 
283  iEvent.put(std::move(fed_buffers));
284 }
285 
286 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
288  //The following says we do not know what parameters are allowed so do no validation
289  // Please change this to state exactly what you do use, even if it is no parameters
291  desc.addUntracked<int>("Verbosity", 0);
292  desc.add<vector<int>>("tdc1", {8, 14, 15, 17, 8, 14, 15, 17, 8, 14, 14, 17, 8, 14, 14, 17, 8, 13, 14, 16, 8, 13,
293  14, 16, 8, 12, 14, 15, 8, 12, 14, 15, 7, 12, 13, 15, 7, 12, 13, 15, 7, 12, 13, 15,
294  7, 12, 13, 15, 7, 11, 12, 14, 7, 11, 12, 14, 7, 11, 12, 14, 7, 11, 12, 7});
295  desc.add<vector<int>>("tdc2", {10, 16, 17, 19, 10, 16, 17, 19, 10, 16, 16, 19, 10, 16, 16, 19, 10, 15, 16, 18, 10, 15,
296  16, 18, 10, 14, 16, 17, 10, 14, 16, 17, 9, 14, 15, 17, 9, 14, 15, 17, 9, 14, 15, 17,
297  9, 14, 15, 17, 9, 13, 14, 16, 9, 13, 14, 16, 9, 13, 14, 16, 9, 13, 14, 9});
298  desc.add<bool>("packHBTDC", true);
299  desc.add<std::string>("ElectronicsMap", "");
300  desc.add<edm::InputTag>("QIE10", edm::InputTag("simHcalDigis", "HFQIE10DigiCollection"));
301  desc.add<edm::InputTag>("QIE11", edm::InputTag("simHcalDigis", "HBHEQIE11DigiCollection"));
302  desc.add<edm::InputTag>("HBHEqie8", edm::InputTag("simHcalDigis"));
303  desc.add<edm::InputTag>("HFqie8", edm::InputTag("simHcalDigis"));
304  desc.add<edm::InputTag>("TP", edm::InputTag("simHcalTriggerPrimitiveDigis"));
305  desc.add<bool>("premix", false);
306  descriptions.add("hcalDigiToRawuHTR", desc);
307  descriptions.addDefault(desc);
308 }
309 
310 //define this as a plug-in
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
EventNumber_t event() const
Definition: EventID.h:40
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
const edm::EDGetTokenT< HBHEDigiCollection > tok_HBHEDigiCollection_
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
QIE11DataFrame convertHB(QIE11DataFrame qiehe, std::vector< int > const &tdc1, std::vector< int > const &tdc2, const int tdcmax)
Definition: PackerHelp.h:620
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
std::vector< T >::const_iterator const_iterator
static constexpr int tdcmax_
int bunchCrossing() const
Definition: EventBase.h:64
constexpr DetId detid() const
Get the detector id.
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int iEvent
Definition: GenABIO.cc:224
void addDefault(ParameterSetDescription const &psetDescription)
static void set(unsigned char *trailer, uint32_t lenght, uint16_t crc, uint8_t evt_stat, uint8_t tts, bool moreTrailers=false)
Set all fields in the trailer.
Definition: FEDTrailer.cc:31
constexpr int presamples() const
for backward compatibility
unsigned short compute_crc(unsigned char *buffer, unsigned int bufSize)
Definition: CRC16.h:46
def move
Definition: eostools.py:511
const edm::ESGetToken< HcalElectronicsMap, HcalElectronicsMapRcd > tok_electronicsMap_
~HcalDigiToRawuHTR() override
int orbitNumber() const
Definition: EventBase.h:65
HcalDigiToRawuHTR(const edm::ParameterSet &)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
const_iterator end() const
Definition: DetId.h:17
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
unsigned long long uint64_t
Definition: Time.h:13
constexpr int presamples() const
for backward compatibility
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EventID id() const
Definition: EventBase.h:59
const vector< int > tdc1_
static void set(unsigned char *header, uint8_t triggerType, uint32_t lvl1ID, uint16_t bxID, uint16_t sourceID, uint8_t version=0, bool moreHeaders=false)
Set all fields in the header.
Definition: FEDHeader.cc:25
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:151
Log< level::Warning, false > LogWarning
Readout chain identification for Hcal.
const edm::EDGetTokenT< HcalDataFrameContainer< QIE11DataFrame > > tok_QIE11DigiCollection_
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
const edm::EDGetTokenT< HcalTrigPrimDigiCollection > tok_TPDigiCollection_
const edm::EDGetTokenT< HFDigiCollection > tok_HFDigiCollection_
const_iterator begin() const
const vector< int > tdc2_
constexpr DetId detid() const
Get the detector id.
const edm::EDGetTokenT< HcalDataFrameContainer< QIE10DataFrame > > tok_QIE10DigiCollection_