Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2022-08-25 01:53:27

0001 /****************************************************************************
0002 *
0003 * This is a part of TOTEM offline software.
0004 * Authors:
0005 *   Jan Kašpar (jan.kaspar@gmail.com)
0006 *   Seyed Mohsen Etesami (setesami@cern.ch)
0007 ****************************************************************************/
0008 
0009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0010 
0011 #include "EventFilter/CTPPSRawToDigi/interface/RawToDigiConverter.h"
0012 #include "EventFilter/CTPPSRawToDigi/interface/CounterChecker.h"
0013 #include "EventFilter/CTPPSRawToDigi/interface/DiamondVFATFrame.h"
0014 #include "EventFilter/CTPPSRawToDigi/interface/TotemSampicFrame.h"
0015 #include "EventFilter/CTPPSRawToDigi/interface/TotemT2VFATFrame.h"
0016 
0017 #include "DataFormats/CTPPSDetId/interface/TotemRPDetId.h"
0018 #include "DataFormats/CTPPSDetId/interface/CTPPSDiamondDetId.h"
0019 #include "DataFormats/CTPPSDetId/interface/TotemTimingDetId.h"
0020 #include "DataFormats/CTPPSDetId/interface/TotemT2DetId.h"
0021 
0022 using namespace std;
0023 using namespace edm;
0024 
0025 RawToDigiConverter::RawToDigiConverter(const edm::ParameterSet &conf)
0026     : verbosity(conf.getUntrackedParameter<unsigned int>("verbosity", 0)),
0027       printErrorSummary(conf.getUntrackedParameter<unsigned int>("printErrorSummary", 1)),
0028       printUnknownFrameSummary(conf.getUntrackedParameter<unsigned int>("printUnknownFrameSummary", 1)),
0029 
0030       testFootprint(conf.getParameter<unsigned int>("testFootprint")),
0031       testCRC(conf.getParameter<unsigned int>("testCRC")),
0032       testID(conf.getParameter<unsigned int>("testID")),
0033       testECMostFrequent(conf.getParameter<unsigned int>("testECMostFrequent")),
0034       testBCMostFrequent(conf.getParameter<unsigned int>("testBCMostFrequent")),
0035 
0036       EC_min(conf.getUntrackedParameter<unsigned int>("EC_min", 10)),
0037       BC_min(conf.getUntrackedParameter<unsigned int>("BC_min", 10)),
0038 
0039       EC_fraction(conf.getUntrackedParameter<double>("EC_fraction", 0.6)),
0040       BC_fraction(conf.getUntrackedParameter<double>("BC_fraction", 0.6)) {}
0041 
0042 void RawToDigiConverter::runCommon(const VFATFrameCollection &input,
0043                                    const TotemDAQMapping &mapping,
0044                                    map<TotemFramePosition, RawToDigiConverter::Record> &records) {
0045   // EC and BC checks (wrt. the most frequent value), BC checks per subsystem
0046   CounterChecker ECChecker(CounterChecker::ECChecker, "EC", EC_min, EC_fraction, verbosity);
0047   CounterChecker BCChecker(CounterChecker::BCChecker, "BC", BC_min, BC_fraction, verbosity);
0048 
0049   // initialise structure merging vfat frame data with the mapping
0050   for (auto &p : mapping.VFATMapping) {
0051     TotemVFATStatus st;
0052     st.setMissing(true);
0053     records[p.first] = {&p.second, nullptr, st};
0054   }
0055 
0056   // event error message buffer
0057   stringstream ees;
0058 
0059   // associate data frames with records
0060   for (VFATFrameCollection::Iterator fr(&input); !fr.IsEnd(); fr.Next()) {
0061     // frame error message buffer
0062     stringstream fes;
0063 
0064     bool problemsPresent = false;
0065     bool stopProcessing = false;
0066 
0067     // skip data frames not listed in the DAQ mapping
0068     auto records_it = records.find(fr.Position());
0069     if (records_it == records.end()) {
0070       unknownSummary[fr.Position()]++;
0071       continue;
0072     }
0073 
0074     // update record
0075     Record &record = records_it->second;
0076     record.frame = fr.Data();
0077     record.status.setMissing(false);
0078     record.status.setNumberOfClustersSpecified(record.frame->isNumberOfClustersPresent());
0079     record.status.setNumberOfClusters(record.frame->getNumberOfClusters());
0080 
0081     // check footprint
0082     if (testFootprint != tfNoTest && !record.frame->checkFootprint()) {
0083       problemsPresent = true;
0084 
0085       if (verbosity > 0)
0086         fes << "    invalid footprint" << endl;
0087 
0088       if (testFootprint == tfErr) {
0089         record.status.setFootprintError();
0090         stopProcessing = true;
0091       }
0092     }
0093 
0094     // check CRC
0095     if (testCRC != tfNoTest && !record.frame->checkCRC()) {
0096       problemsPresent = true;
0097 
0098       if (verbosity > 0)
0099         fes << "    CRC failure" << endl;
0100 
0101       if (testCRC == tfErr) {
0102         record.status.setCRCError();
0103         stopProcessing = true;
0104       }
0105     }
0106 
0107     // check the id mismatch
0108     if (testID != tfNoTest && record.frame->isIDPresent() &&
0109         (record.frame->getChipID() & 0xFFF) != (record.info->hwID & 0xFFF)) {
0110       if (verbosity > 0)
0111         fes << "    ID mismatch (data: 0x" << hex << record.frame->getChipID() << ", mapping: 0x" << record.info->hwID
0112             << dec << ", symbId: " << record.info->symbolicID.symbolicID << ")" << endl;
0113 
0114       if (testID == tfErr) {
0115         record.status.setIDMismatch();
0116         stopProcessing = true;
0117       }
0118     }
0119 
0120     // if there were errors, put the information to ees buffer
0121     if (verbosity > 0 && problemsPresent) {
0122       string message = (stopProcessing) ? "(and will be dropped)" : "(but will be used though)";
0123       if (verbosity > 2) {
0124         ees << "  Frame at " << fr.Position() << " seems corrupted " << message << ":" << endl;
0125         ees << fes.rdbuf();
0126       } else
0127         ees << "  Frame at " << fr.Position() << " seems corrupted " << message << "." << endl;
0128     }
0129 
0130     // if there were serious errors, do not process this frame
0131     if (stopProcessing)
0132       continue;
0133 
0134     // fill EC and BC values to the statistics
0135     if (fr.Data()->isECPresent())
0136       ECChecker.Fill(fr.Data()->getEC(), fr.Position());
0137 
0138     if (fr.Data()->isBCPresent())
0139       BCChecker.Fill(fr.Data()->getBC(), fr.Position());
0140   }
0141 
0142   // analyze EC and BC statistics
0143   if (testECMostFrequent != tfNoTest)
0144     ECChecker.Analyze(records, (testECMostFrequent == tfErr), ees);
0145 
0146   if (testBCMostFrequent != tfNoTest)
0147     BCChecker.Analyze(records, (testBCMostFrequent == tfErr), ees);
0148 
0149   // add error message for missing frames
0150   if (verbosity > 1) {
0151     for (const auto &p : records) {
0152       if (p.second.status.isMissing())
0153         ees << "Frame for VFAT " << p.first << " is not present in the data." << endl;
0154     }
0155   }
0156 
0157   // print error message
0158   if (verbosity > 0 && !ees.rdbuf()->str().empty()) {
0159     if (verbosity > 1)
0160       LogWarning("Totem") << "Error in RawToDigiConverter::runCommon > "
0161                           << "event contains the following problems:" << endl
0162                           << ees.rdbuf() << endl;
0163     else
0164       LogWarning("Totem") << "Error in RawToDigiConverter::runCommon > "
0165                           << "event contains problems." << endl;
0166   }
0167 
0168   // increase error counters
0169   if (printErrorSummary) {
0170     for (const auto &it : records) {
0171       if (!it.second.status.isOK()) {
0172         auto &m = errorSummary[it.first];
0173         m[it.second.status]++;
0174       }
0175     }
0176   }
0177 }
0178 
0179 void RawToDigiConverter::run(const VFATFrameCollection &input,
0180                              const TotemDAQMapping &mapping,
0181                              const TotemAnalysisMask &analysisMask,
0182                              DetSetVector<TotemRPDigi> &rpData,
0183                              DetSetVector<TotemVFATStatus> &finalStatus) {
0184   // structure merging vfat frame data with the mapping
0185   map<TotemFramePosition, Record> records;
0186 
0187   // common processing - frame validation
0188   runCommon(input, mapping, records);
0189 
0190   // second loop over data
0191   for (auto &p : records) {
0192     Record &record = p.second;
0193 
0194     // calculate ids
0195     TotemRPDetId chipId(record.info->symbolicID.symbolicID);
0196     uint8_t chipPosition = chipId.chip();
0197     TotemRPDetId detId = chipId.planeId();
0198 
0199     // update chipPosition in status
0200     record.status.setChipPosition(chipPosition);
0201 
0202     // produce digi only for good frames
0203     if (record.status.isOK()) {
0204       // find analysis mask (needs a default=no mask, if not in present the mapping)
0205       TotemVFATAnalysisMask anMa;
0206       anMa.fullMask = false;
0207 
0208       auto analysisIter = analysisMask.analysisMask.find(record.info->symbolicID);
0209       if (analysisIter != analysisMask.analysisMask.end()) {
0210         // if there is some information about masked channels - save it into conversionStatus
0211         anMa = analysisIter->second;
0212         if (anMa.fullMask)
0213           record.status.setFullyMaskedOut();
0214         else
0215           record.status.setPartiallyMaskedOut();
0216       }
0217 
0218       // create the digi
0219       unsigned short offset = chipPosition * 128;
0220       const vector<unsigned char> &activeChannels = record.frame->getActiveChannels();
0221 
0222       for (auto ch : activeChannels) {
0223         // skip masked channels
0224         if (!anMa.fullMask && anMa.maskedChannels.find(ch) == anMa.maskedChannels.end()) {
0225           DetSet<TotemRPDigi> &digiDetSet = rpData.find_or_insert(detId);
0226           digiDetSet.emplace_back(offset + ch);
0227         }
0228       }
0229     }
0230 
0231     // save status
0232     DetSet<TotemVFATStatus> &statusDetSet = finalStatus.find_or_insert(detId);
0233     statusDetSet.push_back(record.status);
0234   }
0235 }
0236 
0237 void RawToDigiConverter::run(const VFATFrameCollection &coll,
0238                              const TotemDAQMapping &mapping,
0239                              const TotemAnalysisMask &mask,
0240                              edm::DetSetVector<CTPPSDiamondDigi> &digi,
0241                              edm::DetSetVector<TotemVFATStatus> &status) {
0242   // structure merging vfat frame data with the mapping
0243   map<TotemFramePosition, Record> records;
0244 
0245   // common processing - frame validation
0246   runCommon(coll, mapping, records);
0247 
0248   // second loop over data
0249   for (auto &p : records) {
0250     Record &record = p.second;
0251 
0252     // calculate ids
0253     CTPPSDiamondDetId detId(record.info->symbolicID.symbolicID);
0254 
0255     if (record.status.isOK()) {
0256       // update Event Counter in status
0257       record.status.setEC(record.frame->getEC() & 0xFF);
0258 
0259       // create the digi
0260       DetSet<CTPPSDiamondDigi> &digiDetSet = digi.find_or_insert(detId);
0261       digiDetSet.emplace_back(pps::diamond::vfat::getLeadingEdgeTime(*record.frame),
0262                               pps::diamond::vfat::getTrailingEdgeTime(*record.frame),
0263                               pps::diamond::vfat::getThresholdVoltage(*record.frame),
0264                               pps::diamond::vfat::getMultihit(*record.frame),
0265                               pps::diamond::vfat::getHptdcErrorFlag(*record.frame));
0266     }
0267 
0268     // save status
0269     DetSet<TotemVFATStatus> &statusDetSet = status.find_or_insert(detId);
0270     statusDetSet.push_back(record.status);
0271   }
0272 }
0273 
0274 void RawToDigiConverter::run(const VFATFrameCollection &coll,
0275                              const TotemDAQMapping &mapping,
0276                              const TotemAnalysisMask &mask,
0277                              edm::DetSetVector<TotemTimingDigi> &digi,
0278                              edm::DetSetVector<TotemVFATStatus> &status) {
0279   // structure merging vfat frame data with the mapping
0280   map<TotemFramePosition, Record> records;
0281 
0282   // common processing - frame validation
0283   runCommon(coll, mapping, records);
0284 
0285   // second loop over data
0286   for (auto &p : records) {
0287     Record &record = p.second;
0288     if (!record.status.isOK())
0289       continue;
0290 
0291     const TotemFramePosition *framepos = &p.first;
0292 
0293     if (((framepos->getIdxInFiber() % 2) == 0) && (framepos->getIdxInFiber() < 14)) {
0294       //corresponding channel data are always in the neighbouring idx in fiber
0295 
0296       TotemFramePosition frameposdata(framepos->getSubSystemId(),
0297                                       framepos->getTOTFEDId(),
0298                                       framepos->getOptoRxId(),
0299                                       framepos->getGOHId(),
0300                                       (framepos->getIdxInFiber() + 1));
0301       TotemFramePosition frameposEvtInfo(
0302           framepos->getSubSystemId(), framepos->getTOTFEDId(), framepos->getOptoRxId(), framepos->getGOHId(), 0xe);
0303 
0304       auto channelwaveformPtr = records.find(frameposdata);
0305       auto eventInfoPtr = records.find(frameposEvtInfo);
0306 
0307       if (channelwaveformPtr != records.end() && eventInfoPtr != records.end()) {
0308         Record &channelwaveform = records[frameposdata];
0309         Record &eventInfo = records[frameposEvtInfo];
0310 
0311         // Extract all the waveform information from the raw data
0312         TotemSampicFrame totemSampicFrame((const uint8_t *)record.frame->getData(),
0313                                           (const uint8_t *)channelwaveform.frame->getData(),
0314                                           (const uint8_t *)eventInfo.frame->getData());
0315 
0316         if (totemSampicFrame.valid()) {
0317           // create the digi
0318           TotemTimingEventInfo eventInfoTmp(totemSampicFrame.getEventHardwareId(),
0319                                             totemSampicFrame.getL1ATimestamp(),
0320                                             totemSampicFrame.getBunchNumber(),
0321                                             totemSampicFrame.getOrbitNumber(),
0322                                             totemSampicFrame.getEventNumber(),
0323                                             totemSampicFrame.getChannelMap(),
0324                                             totemSampicFrame.getL1ALatency(),
0325                                             totemSampicFrame.getNumberOfSentSamples(),
0326                                             totemSampicFrame.getOffsetOfSamples(),
0327                                             totemSampicFrame.getPLLInfo());
0328           TotemTimingDigi digiTmp(totemSampicFrame.getHardwareId(),
0329                                   totemSampicFrame.getFPGATimestamp(),
0330                                   totemSampicFrame.getTimestampA(),
0331                                   totemSampicFrame.getTimestampB(),
0332                                   totemSampicFrame.getCellInfo(),
0333                                   totemSampicFrame.getSamples(),
0334                                   eventInfoTmp);
0335           // calculate ids
0336           TotemTimingDetId detId(record.info->symbolicID.symbolicID);
0337           const TotemDAQMapping::TotemTimingPlaneChannelPair SWpair =
0338               mapping.getTimingChannel(totemSampicFrame.getHardwareId());
0339           // for FW Version > 0 plane and channel are encoded in the dataframe
0340           if (totemSampicFrame.getFWVersion() == 0)  // Mapping not present in HW, read from SW for FW versions == 0
0341           {
0342             if (SWpair.plane == -1 || SWpair.channel == -1) {
0343               if (verbosity > 0)
0344                 LogWarning("Totem") << "Error in RawToDigiConverter::TotemTiming > "
0345                                     << "HwId not recognized!  hwId: " << std::hex
0346                                     << (unsigned int)totemSampicFrame.getHardwareId() << endl;
0347             } else {
0348               detId.setPlane(SWpair.plane % 4);
0349               detId.setChannel(SWpair.channel);
0350             }
0351           } else  // Mapping read from HW, checked by SW
0352           {
0353             const int HWplane = totemSampicFrame.getDetPlane() % 16;
0354             const int HWchannel = totemSampicFrame.getDetChannel() % 16;
0355 
0356             if (SWpair.plane == -1 || SWpair.channel == -1) {
0357               if (verbosity > 0)
0358                 LogWarning("Totem") << "Warning in RawToDigiConverter::TotemTiming > "
0359                                     << "HwId not recognized!  hwId: " << std::hex
0360                                     << (unsigned int)totemSampicFrame.getHardwareId()
0361                                     << "\tUsing plane and ch from HW without check!" << endl;
0362             } else {
0363               if (verbosity > 0 && (SWpair.plane != HWplane || SWpair.channel != HWchannel))
0364                 LogWarning("Totem") << "Warning in RawToDigiConverter::TotemTiming > "
0365                                     << "Hw mapping different from SW mapping. hwId: " << std::hex
0366                                     << (unsigned int)totemSampicFrame.getHardwareId() << "HW: " << std::dec << HWplane
0367                                     << ":" << HWchannel << "\tSW " << SWpair.plane << ":" << SWpair.channel
0368                                     << "\tUsing plane and ch from HW!" << endl;
0369             }
0370             detId.setPlane(HWplane % 4);
0371             detId.setChannel(HWchannel);
0372           }
0373 
0374           DetSet<TotemTimingDigi> &digiDetSet = digi.find_or_insert(detId);
0375           digiDetSet.push_back(digiTmp);
0376         }
0377       }
0378     }
0379   }
0380 }
0381 
0382 void RawToDigiConverter::run(const VFATFrameCollection &coll,
0383                              const TotemDAQMapping &mapping,
0384                              const TotemAnalysisMask &mask,
0385                              edmNew::DetSetVector<TotemT2Digi> &digi,
0386                              edm::DetSetVector<TotemVFATStatus> &status) {
0387   // structure merging vfat frame data with the mapping
0388   map<TotemFramePosition, Record> records;
0389 
0390   // common processing - frame validation
0391   runCommon(coll, mapping, records);
0392 
0393   // second loop over data
0394   for (auto &p : records) {
0395     Record &record = p.second;
0396 
0397     // calculate ids
0398     TotemT2DetId detId(record.info->symbolicID.symbolicID);
0399 
0400     if (record.status.isOK()) {
0401       // update Event Counter in status
0402       record.status.setEC(record.frame->getEC() & 0xFF);
0403 
0404       // create the digi
0405       edmNew::DetSetVector<TotemT2Digi>::FastFiller(digi, detId)
0406           .emplace_back(totem::nt2::vfat::geoId(*record.frame),
0407                         totem::nt2::vfat::channelId(*record.frame),
0408                         totem::nt2::vfat::channelMarker(*record.frame),
0409                         totem::nt2::vfat::leadingEdgeTime(*record.frame),
0410                         totem::nt2::vfat::trailingEdgeTime(*record.frame));
0411     }
0412 
0413     // save status
0414     DetSet<TotemVFATStatus> &statusDetSet = status.find_or_insert(detId);
0415     statusDetSet.push_back(record.status);
0416   }
0417 }
0418 
0419 void RawToDigiConverter::printSummaries() const {
0420   // print error summary
0421   if (printErrorSummary) {
0422     if (!errorSummary.empty()) {
0423       stringstream ees;
0424       for (const auto &vit : errorSummary) {
0425         ees << vit.first << endl;
0426 
0427         for (const auto &it : vit.second)
0428           ees << "    " << it.first << " : " << it.second << endl;
0429       }
0430 
0431       LogWarning("Totem") << "RawToDigiConverter: error summary (error signature : number of such events)\n"
0432                           << endl
0433                           << ees.rdbuf();
0434     } else {
0435       LogInfo("Totem") << "RawToDigiConverter: no errors to be reported.";
0436     }
0437   }
0438 
0439   // print summary of unknown frames (found in data but not in the mapping)
0440   if (printUnknownFrameSummary) {
0441     if (!unknownSummary.empty()) {
0442       stringstream ees;
0443       for (const auto &it : unknownSummary)
0444         ees << "  " << it.first << " : " << it.second << endl;
0445 
0446       LogWarning("Totem")
0447           << "RawToDigiConverter: frames found in data, but not in the mapping (frame position : number of events)\n"
0448           << endl
0449           << ees.rdbuf();
0450     } else {
0451       LogInfo("Totem") << "RawToDigiConverter: no unknown frames to be reported.";
0452     }
0453   }
0454 }