Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:10:42

0001 #include "FWCore/Framework/interface/Frameworkfwd.h"
0002 #include "FWCore/Framework/interface/global/EDProducer.h"
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/Framework/interface/MakerMacros.h"
0005 #include "FWCore/Framework/interface/ESHandle.h"
0006 
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 #include "FWCore/ServiceRegistry/interface/Service.h"
0009 #include "FWCore/Utilities/interface/ESGetToken.h"
0010 
0011 #include "EventFilter/HcalRawToDigi/interface/HcalUHTRData.h"
0012 #include "DataFormats/HcalDigi/interface/HcalQIESample.h"
0013 
0014 #include "CondFormats/HcalObjects/interface/HcalElectronicsMap.h"
0015 
0016 #include "DataFormats/Common/interface/Handle.h"
0017 #include "DataFormats/FEDRawData/interface/FEDRawDataCollection.h"
0018 #include "DataFormats/FEDRawData/interface/FEDHeader.h"
0019 #include "DataFormats/FEDRawData/interface/FEDTrailer.h"
0020 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
0021 #include "DataFormats/FEDRawData/interface/FEDRawData.h"
0022 #include "FWCore/Utilities/interface/CRC16.h"
0023 
0024 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
0025 #include "CalibFormats/HcalObjects/interface/HcalDbRecord.h"
0026 
0027 #include "PackerHelp.h"
0028 
0029 #include <fstream>
0030 #include <iostream>
0031 #include <memory>
0032 
0033 #include <sstream>
0034 #include <string>
0035 
0036 /* QUESTION: what about dual FED readout? */
0037 /* QUESTION: what do I do if the number of 16-bit words
0038    are not divisible by 4? -- these need to
0039    fit into the 64-bit words of the FEDRawDataFormat */
0040 
0041 using namespace std;
0042 
0043 class HcalDigiToRawuHTR : public edm::global::EDProducer<> {
0044 public:
0045   explicit HcalDigiToRawuHTR(const edm::ParameterSet&);
0046   ~HcalDigiToRawuHTR() override;
0047 
0048   void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0049 
0050   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0051 
0052 private:
0053   const int _verbosity;
0054   const vector<int> tdc1_;
0055   const vector<int> tdc2_;
0056   const bool packHBTDC_;
0057   static constexpr int tdcmax_ = 49;
0058 
0059   const edm::EDGetTokenT<HcalDataFrameContainer<QIE10DataFrame>> tok_QIE10DigiCollection_;
0060   const edm::EDGetTokenT<HcalDataFrameContainer<QIE11DataFrame>> tok_QIE11DigiCollection_;
0061   const edm::EDGetTokenT<HBHEDigiCollection> tok_HBHEDigiCollection_;
0062   const edm::EDGetTokenT<HFDigiCollection> tok_HFDigiCollection_;
0063   const edm::EDGetTokenT<HcalTrigPrimDigiCollection> tok_TPDigiCollection_;
0064   const edm::ESGetToken<HcalElectronicsMap, HcalElectronicsMapRcd> tok_electronicsMap_;
0065 
0066   const bool premix_;
0067 };
0068 
0069 HcalDigiToRawuHTR::HcalDigiToRawuHTR(const edm::ParameterSet& iConfig)
0070     : _verbosity(iConfig.getUntrackedParameter<int>("Verbosity", 0)),
0071       tdc1_(iConfig.getParameter<vector<int>>("tdc1")),
0072       tdc2_(iConfig.getParameter<vector<int>>("tdc2")),
0073       packHBTDC_(iConfig.getParameter<bool>("packHBTDC")),
0074       tok_QIE10DigiCollection_(
0075           consumes<HcalDataFrameContainer<QIE10DataFrame>>(iConfig.getParameter<edm::InputTag>("QIE10"))),
0076       tok_QIE11DigiCollection_(
0077           consumes<HcalDataFrameContainer<QIE11DataFrame>>(iConfig.getParameter<edm::InputTag>("QIE11"))),
0078       tok_HBHEDigiCollection_(consumes<HBHEDigiCollection>(iConfig.getParameter<edm::InputTag>("HBHEqie8"))),
0079       tok_HFDigiCollection_(consumes<HFDigiCollection>(iConfig.getParameter<edm::InputTag>("HFqie8"))),
0080       tok_TPDigiCollection_(consumes<HcalTrigPrimDigiCollection>(iConfig.getParameter<edm::InputTag>("TP"))),
0081       tok_electronicsMap_(esConsumes<HcalElectronicsMap, HcalElectronicsMapRcd>(
0082           edm::ESInputTag("", iConfig.getParameter<std::string>("ElectronicsMap")))),
0083       premix_(iConfig.getParameter<bool>("premix")) {
0084   produces<FEDRawDataCollection>("");
0085   for (size_t i = 0; i < tdc1_.size(); i++) {
0086     if (!(tdc1_.at(i) >= 0 && tdc1_.at(i) <= tdc2_.at(i) && tdc2_.at(i) <= tdcmax_))
0087       edm::LogWarning("HcalDigiToRawuHTR")
0088           << " incorrect TDC ranges " << i << "-th element: " << tdc1_.at(i) << ", " << tdc2_.at(i) << ", " << tdcmax_;
0089   }
0090 }
0091 
0092 HcalDigiToRawuHTR::~HcalDigiToRawuHTR() {}
0093 
0094 void HcalDigiToRawuHTR::produce(edm::StreamID id, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0095   using namespace edm;
0096 
0097   edm::ESHandle<HcalElectronicsMap> item = iSetup.getHandle(tok_electronicsMap_);
0098   const HcalElectronicsMap* readoutMap = item.product();
0099 
0100   //collection to be inserted into event
0101   std::unique_ptr<FEDRawDataCollection> fed_buffers(new FEDRawDataCollection());
0102 
0103   //
0104   //  Extracting All the Collections containing useful Info
0105   edm::Handle<QIE10DigiCollection> qie10DigiCollection;
0106   edm::Handle<QIE11DigiCollection> qie11DigiCollection;
0107   edm::Handle<HBHEDigiCollection> hbheDigiCollection;
0108   edm::Handle<HFDigiCollection> hfDigiCollection;
0109   edm::Handle<HcalTrigPrimDigiCollection> tpDigiCollection;
0110   iEvent.getByToken(tok_QIE10DigiCollection_, qie10DigiCollection);
0111   iEvent.getByToken(tok_QIE11DigiCollection_, qie11DigiCollection);
0112   iEvent.getByToken(tok_HBHEDigiCollection_, hbheDigiCollection);
0113   iEvent.getByToken(tok_HFDigiCollection_, hfDigiCollection);
0114   iEvent.getByToken(tok_TPDigiCollection_, tpDigiCollection);
0115 
0116   // first argument is the fedid (minFEDID+crateId)
0117   map<int, unique_ptr<HCalFED>> fedMap;
0118 
0119   // - - - - - - - - - - - - - - - - - - - - - - - - - - -
0120   // QIE10 precision data
0121   // - - - - - - - - - - - - - - - - - - - - - - - - - - -
0122   UHTRpacker uhtrs;
0123   // loop over each digi and allocate memory for each
0124   if (qie10DigiCollection.isValid()) {
0125     const QIE10DigiCollection& qie10dc = *(qie10DigiCollection);
0126     for (unsigned int j = 0; j < qie10dc.size(); j++) {
0127       QIE10DataFrame qiedf = static_cast<QIE10DataFrame>(qie10dc[j]);
0128       DetId detid = qiedf.detid();
0129       HcalElectronicsId eid(readoutMap->lookup(detid));
0130       int crateId = eid.crateId();
0131       int slotId = eid.slot();
0132       int uhtrIndex = ((slotId & 0xF) << 8) | (crateId & 0xFF);
0133       int presamples = qiedf.presamples();
0134 
0135       /* Defining a custom index that will encode only
0136      the information about the crate and slot of a
0137      given channel:   crate: bits 0-7
0138      slot:  bits 8-12 */
0139 
0140       if (!uhtrs.exist(uhtrIndex)) {
0141         uhtrs.newUHTR(uhtrIndex, presamples);
0142       }
0143       uhtrs.addChannel(uhtrIndex, qiedf, readoutMap, _verbosity);
0144     }
0145   }
0146   // - - - - - - - - - - - - - - - - - - - - - - - - - - -
0147   // QIE11 precision data
0148   // - - - - - - - - - - - - - - - - - - - - - - - - - - -
0149   //UHTRpacker uhtrs;
0150   // loop over each digi and allocate memory for each
0151   if (qie11DigiCollection.isValid()) {
0152     const QIE11DigiCollection& qie11dc = *(qie11DigiCollection);
0153     for (unsigned int j = 0; j < qie11dc.size(); j++) {
0154       QIE11DataFrame qiedf = static_cast<QIE11DataFrame>(qie11dc[j]);
0155       DetId detid = qiedf.detid();
0156       HcalElectronicsId eid(readoutMap->lookup(detid));
0157       int crateId = eid.crateId();
0158       int slotId = eid.slot();
0159       int uhtrIndex = ((slotId & 0xF) << 8) | (crateId & 0xFF);
0160       int presamples = qiedf.presamples();
0161 
0162       //   convert to hb qie data if hb
0163       if (packHBTDC_ && HcalDetId(detid.rawId()).subdet() == HcalSubdetector::HcalBarrel)
0164         qiedf = convertHB(qiedf, tdc1_, tdc2_, tdcmax_);
0165 
0166       if (!uhtrs.exist(uhtrIndex)) {
0167         uhtrs.newUHTR(uhtrIndex, presamples);
0168       }
0169       uhtrs.addChannel(uhtrIndex, qiedf, readoutMap, _verbosity);
0170     }
0171   }
0172   // - - - - - - - - - - - - - - - - - - - - - - - - - - -
0173   // HF (QIE8) precision data
0174   // - - - - - - - - - - - - - - - - - - - - - - - - - - -
0175   // loop over each digi and allocate memory for each
0176   if (hfDigiCollection.isValid()) {
0177     const HFDigiCollection& qie8hfdc = *(hfDigiCollection);
0178     for (HFDigiCollection::const_iterator qiedf = qie8hfdc.begin(); qiedf != qie8hfdc.end(); qiedf++) {
0179       DetId detid = qiedf->id();
0180 
0181       HcalElectronicsId eid(readoutMap->lookup(detid));
0182       int crateId = eid.crateId();
0183       int slotId = eid.slot();
0184       int uhtrIndex = (crateId & 0xFF) | ((slotId & 0xF) << 8);
0185       int presamples = qiedf->presamples();
0186 
0187       if (!uhtrs.exist(uhtrIndex)) {
0188         uhtrs.newUHTR(uhtrIndex, presamples);
0189       }
0190       uhtrs.addChannel(uhtrIndex, qiedf, readoutMap, premix_, _verbosity);
0191     }
0192   }
0193   // - - - - - - - - - - - - - - - - - - - - - - - - - - -
0194   // HBHE (QIE8) precision data
0195   // - - - - - - - - - - - - - - - - - - - - - - - - - - -
0196   // loop over each digi and allocate memory for each
0197   if (hbheDigiCollection.isValid()) {
0198     const HBHEDigiCollection& qie8hbhedc = *(hbheDigiCollection);
0199     for (HBHEDigiCollection::const_iterator qiedf = qie8hbhedc.begin(); qiedf != qie8hbhedc.end(); qiedf++) {
0200       DetId detid = qiedf->id();
0201 
0202       HcalElectronicsId eid(readoutMap->lookup(detid));
0203       int crateId = eid.crateId();
0204       int slotId = eid.slot();
0205       int uhtrIndex = (crateId & 0xFF) | ((slotId & 0xF) << 8);
0206       int presamples = qiedf->presamples();
0207 
0208       if (!uhtrs.exist(uhtrIndex)) {
0209         uhtrs.newUHTR(uhtrIndex, presamples);
0210       }
0211       uhtrs.addChannel(uhtrIndex, qiedf, readoutMap, premix_, _verbosity);
0212     }
0213   }
0214   // - - - - - - - - - - - - - - - - - - - - - - - - - - -
0215   // TP data
0216   // - - - - - - - - - - - - - - - - - - - - - - - - - - -
0217   // loop over each digi and allocate memory for each
0218   if (tpDigiCollection.isValid()) {
0219     const HcalTrigPrimDigiCollection& qietpdc = *(tpDigiCollection);
0220     for (HcalTrigPrimDigiCollection::const_iterator qiedf = qietpdc.begin(); qiedf != qietpdc.end(); qiedf++) {
0221       DetId detid = qiedf->id();
0222       HcalElectronicsId eid(readoutMap->lookupTrigger(detid));
0223 
0224       int crateId = eid.crateId();
0225       int slotId = eid.slot();
0226       int uhtrIndex = (crateId & 0xFF) | ((slotId & 0xF) << 8);
0227       int ilink = eid.fiberIndex();
0228       int itower = eid.fiberChanId();
0229       int channelid = (itower & 0xF) | ((ilink & 0xF) << 4);
0230       int presamples = qiedf->presamples();
0231 
0232       if (!uhtrs.exist(uhtrIndex)) {
0233         uhtrs.newUHTR(uhtrIndex, presamples);
0234       }
0235       uhtrs.addChannel(uhtrIndex, qiedf, channelid, _verbosity);
0236     }
0237   }
0238   // -----------------------------------------------------
0239   // -----------------------------------------------------
0240   // loop over each uHTR and format data
0241   // -----------------------------------------------------
0242   // -----------------------------------------------------
0243   // loop over each uHTR and format data
0244   for (UHTRpacker::UHTRMap::iterator uhtr = uhtrs.uhtrs.begin(); uhtr != uhtrs.uhtrs.end(); ++uhtr) {
0245     uint64_t crateId = (uhtr->first) & 0xFF;
0246     uint64_t slotId = (uhtr->first & 0xF00) >> 8;
0247 
0248     uhtrs.finalizeHeadTail(&(uhtr->second), _verbosity);
0249     int fedId = FEDNumbering::MINHCALuTCAFEDID + crateId;
0250     if (fedMap.find(fedId) == fedMap.end()) {
0251       /* QUESTION: where should the orbit number come from? */
0252       fedMap[fedId] =
0253           std::make_unique<HCalFED>(fedId, iEvent.id().event(), iEvent.orbitNumber(), iEvent.bunchCrossing());
0254     }
0255     fedMap[fedId]->addUHTR(uhtr->second, crateId, slotId);
0256   }  // end loop over uhtr containers
0257 
0258   /* ------------------------------------------------------
0259      ------------------------------------------------------
0260            putting together the FEDRawDataCollection
0261      ------------------------------------------------------
0262      ------------------------------------------------------ */
0263   for (map<int, unique_ptr<HCalFED>>::iterator fed = fedMap.begin(); fed != fedMap.end(); ++fed) {
0264     int fedId = fed->first;
0265 
0266     auto& rawData = fed_buffers->FEDData(fedId);
0267     fed->second->formatFEDdata(rawData);
0268 
0269     FEDHeader hcalFEDHeader(rawData.data());
0270     hcalFEDHeader.set(rawData.data(), 1, iEvent.id().event(), iEvent.bunchCrossing(), fedId);
0271     FEDTrailer hcalFEDTrailer(rawData.data() + (rawData.size() - 8));
0272     hcalFEDTrailer.set(rawData.data() + (rawData.size() - 8),
0273                        rawData.size() / 8,
0274                        evf::compute_crc(rawData.data(), rawData.size()),
0275                        0,
0276                        0);
0277 
0278   }  // end loop over FEDs with data
0279 
0280   iEvent.put(std::move(fed_buffers));
0281 }
0282 
0283 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0284 void HcalDigiToRawuHTR::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0285   //The following says we do not know what parameters are allowed so do no validation
0286   // Please change this to state exactly what you do use, even if it is no parameters
0287   edm::ParameterSetDescription desc;
0288   desc.addUntracked<int>("Verbosity", 0);
0289   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,
0290                                  12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12,
0291                                  12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12});
0292   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,
0293                                  14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14,
0294                                  14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14, 14});
0295   desc.add<bool>("packHBTDC", true);
0296   desc.add<std::string>("ElectronicsMap", "");
0297   desc.add<edm::InputTag>("QIE10", edm::InputTag("simHcalDigis", "HFQIE10DigiCollection"));
0298   desc.add<edm::InputTag>("QIE11", edm::InputTag("simHcalDigis", "HBHEQIE11DigiCollection"));
0299   desc.add<edm::InputTag>("HBHEqie8", edm::InputTag("simHcalDigis"));
0300   desc.add<edm::InputTag>("HFqie8", edm::InputTag("simHcalDigis"));
0301   desc.add<edm::InputTag>("TP", edm::InputTag("simHcalTriggerPrimitiveDigis"));
0302   desc.add<bool>("premix", false);
0303   descriptions.add("hcalDigiToRawuHTR", desc);
0304   descriptions.addDefault(desc);
0305 }
0306 
0307 //define this as a plug-in
0308 DEFINE_FWK_MODULE(HcalDigiToRawuHTR);