CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
HcalFEDIntegrityTask Class Reference
Inheritance diagram for HcalFEDIntegrityTask:

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 
void dqmBeginRun (const edm::Run &, edm::EventSetup const &) override
 
 HcalFEDIntegrityTask (const edm::ParameterSet &ps)
 
 ~HcalFEDIntegrityTask () override
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Member Functions

void labelBins (MonitorElement *me)
 

Private Attributes

std::string dirName_
 
int maxFEDNum_
 
MonitorElementmeFEDEntries_
 
MonitorElementmeFEDFatal_
 
MonitorElementmeFEDNonFatal_
 
int minFEDNum_
 
int numOfFED_
 
edm::EDGetTokenT< FEDRawDataCollectiontokFEDs_
 
edm::EDGetTokenT< HcalUnpackerReporttokReport_
 

Detailed Description

Definition at line 24 of file HcalFEDIntegrityTask.cc.

Constructor & Destructor Documentation

◆ HcalFEDIntegrityTask()

HcalFEDIntegrityTask::HcalFEDIntegrityTask ( const edm::ParameterSet ps)

Definition at line 48 of file HcalFEDIntegrityTask.cc.

References maxFEDNum_, minFEDNum_, and numOfFED_.

49  : tokFEDs_(consumes<FEDRawDataCollection>(
50  ps.getUntrackedParameter<edm::InputTag>("tagFEDs", edm::InputTag("rawDataCollector")))),
51  tokReport_(consumes<HcalUnpackerReport>(
52  ps.getUntrackedParameter<edm::InputTag>("tagReport", edm::InputTag("hcalDigis")))),
55  dirName_(ps.getUntrackedParameter<std::string>("DirName", "Hcal/FEDIntegrity/")) {
56  LogInfo("HcalDQM") << "HcalFEDIntegrityTask::HcalFEDIntegrityTask: Constructor Initialization" << endl;
58 }
edm::EDGetTokenT< FEDRawDataCollection > tokFEDs_
T getUntrackedParameter(std::string const &, T const &) const
Log< level::Info, false > LogInfo
edm::EDGetTokenT< HcalUnpackerReport > tokReport_

◆ ~HcalFEDIntegrityTask()

HcalFEDIntegrityTask::~HcalFEDIntegrityTask ( )
override

Definition at line 60 of file HcalFEDIntegrityTask.cc.

60  {
61  LogInfo("HcalDQM") << "HcalFEDIntegrityTask::~HcalFEDIntegrityTask: Destructor" << endl;
62 }
Log< level::Info, false > LogInfo

Member Function Documentation

◆ analyze()

void HcalFEDIntegrityTask::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

Definition at line 80 of file HcalFEDIntegrityTask.cc.

References FEDRawDataCollection::FEDData(), l1t_dqm_sourceclient-live_cfg::fedRawData, dqm::impl::MonitorElement::Fill(), spr::find(), iEvent, HcalUHTRData::l1ANumber(), FEDNumbering::MAXHCALuTCAFEDID, meFEDEntries_, meFEDFatal_, FEDNumbering::MINHCALuTCAFEDID, edmIntegrityCheck::report, tokFEDs_, and tokReport_.

80  {
83  iEvent.getByToken(tokFEDs_, raw);
84  iEvent.getByToken(tokReport_, report);
85 
86  // FEDs with unpacking errors: https://github.com/cms-sw/cmssw/blob/master/EventFilter/HcalRawToDigi/plugins/HcalRawToDigi.cc#L235-L262
87  const std::vector<int> FedsError = (*report).getFedsError();
88  for (auto &fed : FedsError) {
89  if (fed < 1000)
90  LogWarning("HcalDQM") << "HcalFEDIntegrityTask::analyze: obsoleteVME FEDs from HcalUnpackerReport" << endl;
91  meFEDFatal_->Fill(fed);
92  }
93 
95  // Same checks as the RawTask Summary map
97  for (int fed = FEDNumbering::MINHCALuTCAFEDID; fed <= FEDNumbering::MAXHCALuTCAFEDID; fed++) {
98  const FEDRawData &fedRawData = raw->FEDData(fed);
99  if (fedRawData.size() != 0) {
100  meFEDEntries_->Fill(fed);
101  }
102 
103  hcal::AMC13Header const *amc13 = (hcal::AMC13Header const *)fedRawData.data();
104  if (!amc13) {
105  continue;
106  }
107 
108  uint32_t bcn = amc13->bunchId();
109  uint32_t orn = amc13->orbitNumber() & 0xFFFF; // LS 16bits only
110  uint32_t evn = amc13->l1aNumber();
111  int namc = amc13->NAMC();
112 
113  // looping over AMC in this packet
114  for (int iamc = 0; iamc < namc; iamc++) {
115  if (!amc13->AMCEnabled(iamc) || !amc13->AMCDataPresent(iamc) || !amc13->AMCCRCOk(iamc) ||
116  amc13->AMCSegmented(iamc)) {
117  LogWarning("HcalDQM") << "HcalFEDIntegrityTask::analyze: AMC issue on iamc" << iamc << endl;
118  continue;
119  }
120 
121  HcalUHTRData uhtr(amc13->AMCPayload(iamc), amc13->AMCSize(iamc));
122  uint32_t uhtr_evn = uhtr.l1ANumber();
123  uint32_t uhtr_bcn = uhtr.bunchNumber();
124  uint32_t uhtr_orn = uhtr.orbitNumber();
125 
126  if (uhtr_evn != evn || uhtr_bcn != bcn || uhtr_orn != orn) {
127  if (std::find(FedsError.begin(), FedsError.end(), fed) ==
128  FedsError.end()) // FED not already in the error list from unpacker report
129  meFEDFatal_->Fill(fed);
130  break; // one mismatch is sufficient enough to determine it's a bad data
131  }
132  }
133 
134  } // end of Hcal FED looping
135 }
edm::EDGetTokenT< FEDRawDataCollection > tokFEDs_
MonitorElement * meFEDEntries_
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:224
MonitorElement * meFEDFatal_
const FEDRawData & FEDData(int fedid) const
retrieve data for fed
uint32_t l1ANumber() const
Get the HTR event number.
Definition: HcalUHTRData.h:87
Log< level::Warning, false > LogWarning
edm::EDGetTokenT< HcalUnpackerReport > tokReport_

◆ bookHistograms()

void HcalFEDIntegrityTask::bookHistograms ( DQMStore::IBooker iBooker,
edm::Run const &  iRun,
edm::EventSetup const &  iSetup 
)
override

Definition at line 66 of file HcalFEDIntegrityTask.cc.

References dqm::implementation::IBooker::book1D(), dqm::implementation::NavigatorBase::cd(), dirName_, labelBins(), maxFEDNum_, meFEDEntries_, meFEDFatal_, meFEDNonFatal_, minFEDNum_, numOfFED_, and dqm::implementation::NavigatorBase::setCurrentFolder().

68  {
69  iBooker.cd();
70  iBooker.setCurrentFolder(dirName_);
71 
72  meFEDEntries_ = iBooker.book1D("FEDEntries", "FED Entries", numOfFED_, minFEDNum_, maxFEDNum_ + 1);
73  this->labelBins(meFEDEntries_);
74  meFEDFatal_ = iBooker.book1D("FEDFatal", "FED Fatal Errors", numOfFED_, minFEDNum_, maxFEDNum_ + 1);
75  this->labelBins(meFEDFatal_);
76  meFEDNonFatal_ = iBooker.book1D("FEDNonFatal", "FED NON Fatal Errors", numOfFED_, minFEDNum_, maxFEDNum_ + 1);
78 }
void labelBins(MonitorElement *me)
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
MonitorElement * meFEDEntries_
MonitorElement * meFEDFatal_
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
MonitorElement * meFEDNonFatal_

◆ dqmBeginRun()

void HcalFEDIntegrityTask::dqmBeginRun ( const edm::Run iRun,
edm::EventSetup const &  iSetup 
)
override

Definition at line 64 of file HcalFEDIntegrityTask.cc.

64 {}

◆ fillDescriptions()

void HcalFEDIntegrityTask::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 149 of file HcalFEDIntegrityTask.cc.

References edm::ConfigurationDescriptions::addWithDefaultLabel(), submitPVResolutionJobs::desc, HLT_2022v15_cff::InputTag, FEDNumbering::MAXHCALuTCAFEDID, FEDNumbering::MINHCALuTCAFEDID, and AlCaHLTBitMon_QueryRunRegistry::string.

149  {
151  desc.addUntracked<std::string>("name", "HcalFEDIntegrityTask");
152  desc.addUntracked<int>("debug", 0);
153  desc.addUntracked<edm::InputTag>("tagFEDs", edm::InputTag("rawDataCollector"));
154  desc.addUntracked<edm::InputTag>("tagReport", edm::InputTag("hcalDigis"));
155  desc.addUntracked<int>(
156  "MinHcalFEDID",
157  FEDNumbering::MINHCALuTCAFEDID); // Assuming no more VME FEDs after LS2, according to Hcal Phase1 upgrade.
158  desc.addUntracked<int>("MaxHcalFEDID", FEDNumbering::MAXHCALuTCAFEDID);
159  desc.addUntracked<std::string>("DirName", "Hcal/FEDIntegrity/");
160  descriptions.addWithDefaultLabel(desc);
161 }
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)

◆ labelBins()

void HcalFEDIntegrityTask::labelBins ( MonitorElement me)
private

Definition at line 137 of file HcalFEDIntegrityTask.cc.

References dqm-mbProfile::format, mps_fire::i, hlt_dqm_clientPB-live_cfg::me, minFEDNum_, numOfFED_, AlCaHLTBitMon_QueryRunRegistry::string, and fw3dlego::xbins.

Referenced by bookHistograms().

137  {
138  int xbins = me->getNbinsX();
139 
140  if (xbins != numOfFED_)
141  return;
142 
143  for (int i = 0; i < xbins; i++) {
144  const std::string xLabel = fmt::format("{}", minFEDNum_ + i);
145  me->setBinLabel(i + 1, xLabel, 1);
146  }
147 }
const double xbins[]

Member Data Documentation

◆ dirName_

std::string HcalFEDIntegrityTask::dirName_
private

Definition at line 42 of file HcalFEDIntegrityTask.cc.

Referenced by bookHistograms().

◆ maxFEDNum_

int HcalFEDIntegrityTask::maxFEDNum_
private

Definition at line 41 of file HcalFEDIntegrityTask.cc.

Referenced by bookHistograms(), and HcalFEDIntegrityTask().

◆ meFEDEntries_

MonitorElement* HcalFEDIntegrityTask::meFEDEntries_
private

Definition at line 43 of file HcalFEDIntegrityTask.cc.

Referenced by analyze(), and bookHistograms().

◆ meFEDFatal_

MonitorElement* HcalFEDIntegrityTask::meFEDFatal_
private

Definition at line 44 of file HcalFEDIntegrityTask.cc.

Referenced by analyze(), and bookHistograms().

◆ meFEDNonFatal_

MonitorElement* HcalFEDIntegrityTask::meFEDNonFatal_
private

Definition at line 45 of file HcalFEDIntegrityTask.cc.

Referenced by bookHistograms().

◆ minFEDNum_

int HcalFEDIntegrityTask::minFEDNum_
private

Definition at line 41 of file HcalFEDIntegrityTask.cc.

Referenced by bookHistograms(), HcalFEDIntegrityTask(), and labelBins().

◆ numOfFED_

int HcalFEDIntegrityTask::numOfFED_
private

Definition at line 41 of file HcalFEDIntegrityTask.cc.

Referenced by bookHistograms(), HcalFEDIntegrityTask(), and labelBins().

◆ tokFEDs_

edm::EDGetTokenT<FEDRawDataCollection> HcalFEDIntegrityTask::tokFEDs_
private

Definition at line 38 of file HcalFEDIntegrityTask.cc.

Referenced by analyze().

◆ tokReport_

edm::EDGetTokenT<HcalUnpackerReport> HcalFEDIntegrityTask::tokReport_
private

Definition at line 39 of file HcalFEDIntegrityTask.cc.

Referenced by analyze().