Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "EventFilter/RPCRawToDigi/interface/RPCPackingModule.h"
0002 
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/Framework/interface/EventSetup.h"
0005 
0006 #include "DataFormats/FEDRawData/interface/FEDRawDataCollection.h"
0007 #include "DataFormats/FEDRawData/interface/FEDHeader.h"
0008 #include "DataFormats/FEDRawData/interface/FEDTrailer.h"
0009 
0010 #include "DataFormats/Common/interface/Handle.h"
0011 #include "FWCore/Framework/interface/ESHandle.h"
0012 
0013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0015 
0016 #include "CondFormats/RPCObjects/interface/RPCReadOutMapping.h"
0017 
0018 #include "EventFilter/RPCRawToDigi/interface/RPCRecordFormatter.h"
0019 #include "EventFilter/RPCRawToDigi/interface/DebugDigisPrintout.h"
0020 #include "DataFormats/RPCDigi/interface/EmptyWord.h"
0021 
0022 #include <string>
0023 #include <sstream>
0024 
0025 using namespace std;
0026 using namespace edm;
0027 using namespace rpcrawtodigi;
0028 
0029 typedef uint64_t Word64;
0030 
0031 RPCPackingModule::RPCPackingModule(const ParameterSet& pset) {
0032   dataLabel_ = consumes<RPCDigiCollection>(pset.getParameter<edm::InputTag>("InputLabel"));
0033   readoutMappingToken_ = esConsumes<RPCEMap, RPCEMapRcd>();
0034   theCabling = new RPCReadOutMapping("");
0035   produces<FEDRawDataCollection>();
0036 }
0037 
0038 RPCPackingModule::~RPCPackingModule() { delete theCabling; }
0039 
0040 void RPCPackingModule::produce(edm::Event& ev, const edm::EventSetup& es) {
0041   Handle<RPCDigiCollection> digiCollection;
0042   ev.getByToken(dataLabel_, digiCollection);
0043   LogDebug("") << DebugDigisPrintout()(digiCollection.product());
0044 
0045   if (recordWatcher_.check(es)) {
0046     delete theCabling;
0047     LogTrace("") << "record has CHANGED!!, initialise readout map!";
0048     ESHandle<RPCEMap> readoutMapping = es.getHandle(readoutMappingToken_);
0049     theCabling = readoutMapping->convert();
0050     LogTrace("") << " READOUT MAP VERSION: " << theCabling->version() << endl;
0051   }
0052 
0053   auto buffers = std::make_unique<FEDRawDataCollection>();
0054 
0055   //  pair<int,int> rpcFEDS=FEDNumbering::getRPCFEDIds();
0056   pair<int, int> rpcFEDS(790, 792);
0057   for (int id = rpcFEDS.first; id <= rpcFEDS.second; ++id) {
0058     RPCRecordFormatter formatter(id, theCabling);
0059     unsigned int lvl1_ID = ev.id().event();
0060     FEDRawData* rawData = RPCPackingModule::rawData(id, lvl1_ID, digiCollection.product(), formatter);
0061     FEDRawData& fedRawData = buffers->FEDData(id);
0062 
0063     fedRawData = *rawData;
0064     delete rawData;
0065   }
0066   ev.put(std::move(buffers));
0067 }
0068 
0069 FEDRawData* RPCPackingModule::rawData(int fedId,
0070                                       unsigned int lvl1_ID,
0071                                       const RPCDigiCollection* digis,
0072                                       const RPCRecordFormatter& formatter) const {
0073   //
0074   // get merged records
0075   //
0076   int trigger_BX = 200;  // FIXME - set event by event but correct bx assigment in digi
0077   vector<EventRecords> merged = RPCPackingModule::eventRecords(fedId, trigger_BX, digis, formatter);
0078 
0079   //
0080   // create data words
0081   //
0082   vector<Word64> dataWords;
0083   EmptyWord empty;
0084   typedef vector<EventRecords>::const_iterator IR;
0085   for (IR ir = merged.begin(), irEnd = merged.end(); ir != irEnd; ++ir) {
0086     Word64 w = (((Word64(ir->recordBX().data()) << 16) | ir->recordSLD().data()) << 16 | ir->recordCD().data()) << 16 |
0087                empty.data();
0088     dataWords.push_back(w);
0089   }
0090 
0091   //
0092   // create raw data
0093   //
0094   int nHeaders = 1;
0095   int nTrailers = 1;
0096   int dataSize = (nHeaders + nTrailers + dataWords.size()) * sizeof(Word64);
0097   FEDRawData* raw = new FEDRawData(dataSize);
0098 
0099   //
0100   // add header
0101   //
0102   unsigned char* pHeader = raw->data();
0103   int evt_ty = 3;
0104   int source_ID = fedId;
0105   FEDHeader::set(pHeader, evt_ty, lvl1_ID, trigger_BX, source_ID);
0106 
0107   //
0108   // add datawords
0109   //
0110   for (unsigned int idata = 0; idata < dataWords.size(); idata++) {
0111     Word64* word = reinterpret_cast<Word64*>(pHeader + (idata + 1) * sizeof(Word64));
0112     *word = dataWords[idata];
0113   }
0114 
0115   //
0116   // add trailer
0117   //
0118   unsigned char* pTrailer = pHeader + raw->size() - sizeof(Word64);
0119   int crc = 0;
0120   int evt_stat = 15;
0121   int tts = 0;
0122   int datasize = raw->size() / sizeof(Word64);
0123   FEDTrailer::set(pTrailer, datasize, crc, evt_stat, tts);
0124 
0125   return raw;
0126 }
0127 
0128 vector<EventRecords> RPCPackingModule::eventRecords(int fedId,
0129                                                     int trigger_BX,
0130                                                     const RPCDigiCollection* digis,
0131                                                     const RPCRecordFormatter& formatter) {
0132   typedef DigiContainerIterator<RPCDetId, RPCDigi> DigiRangeIterator;
0133   vector<EventRecords> dataRecords;
0134 
0135   LogDebug("RPCRawDataPacker") << "Packing Fed id=" << fedId;
0136   for (DigiRangeIterator it = digis->begin(); it != digis->end(); it++) {
0137     RPCDetId rpcDetId = (*it).first;
0138     uint32_t rawDetId = rpcDetId.rawId();
0139     RPCDigiCollection::Range range = digis->get(rpcDetId);
0140     for (vector<RPCDigi>::const_iterator id = range.first; id != range.second; id++) {
0141       const RPCDigi& digi = (*id);
0142       vector<EventRecords> rawFromDigi = formatter.recordPack(rawDetId, digi, trigger_BX);
0143       dataRecords.insert(dataRecords.end(), rawFromDigi.begin(), rawFromDigi.end());
0144     }
0145   }
0146 
0147   //
0148   // merge data words
0149   //
0150   LogTrace("RPCRawDataPacker") << " size of   data: " << dataRecords.size();
0151   vector<EventRecords> merged = EventRecords::mergeRecords(dataRecords);
0152   LogTrace("") << " size of megred: " << merged.size();
0153 
0154   return merged;
0155 }