Back to home page

Project CMSSW displayed by LXR

 
 

    


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

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