CMS 3D CMS Logo

SiPixelPhase1RawDataErrorComparator.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 // Package: SiPixelPhase1RawDataErrorComparator
3 // Class: SiPixelPhase1RawDataErrorComparator
4 //
7 //
8 // Author: Marco Musich
9 //
29 // for string manipulations
30 #include <fmt/printf.h>
31 
32 namespace {
33  // same logic used for the MTV:
34  // cf https://github.com/cms-sw/cmssw/blob/master/Validation/RecoTrack/src/MTVHistoProducerAlgoForTracker.cc
36 
37  void setBinLog(TAxis* axis) {
38  int bins = axis->GetNbins();
39  float from = axis->GetXmin();
40  float to = axis->GetXmax();
41  float width = (to - from) / bins;
42  std::vector<float> new_bins(bins + 1, 0);
43  for (int i = 0; i <= bins; i++) {
44  new_bins[i] = TMath::Power(10, from + i * width);
45  }
46  axis->Set(bins, new_bins.data());
47  }
48 
49  void setBinLogX(TH1* h) {
50  TAxis* axis = h->GetXaxis();
51  setBinLog(axis);
52  }
53  void setBinLogY(TH1* h) {
54  TAxis* axis = h->GetYaxis();
55  setBinLog(axis);
56  }
57 
58  template <typename... Args>
59  dqm::reco::MonitorElement* make2DIfLog(DQMStore::IBooker& ibook, bool logx, bool logy, Args&&... args) {
60  auto h = std::make_unique<TH2I>(std::forward<Args>(args)...);
61  if (logx)
62  setBinLogX(h.get());
63  if (logy)
64  setBinLogY(h.get());
65  const auto& name = h->GetName();
66  return ibook.book2I(name, h.release());
67  }
68 
69  //errorType - a number (25-38) indicating the type of error recorded.
71  k_FED25 = 25, // 25 indicates an invalid ROC of 25
72  k_FED26 = 26, // 26 indicates a gap word
73  k_FED27 = 27, // 27 indicates a dummy word
74  k_FED28 = 28, // 28 indicates a FIFO full error
75  k_FED29 = 29, // 29 indicates a timeout error
76  k_FED30 = 30, // 30 indicates a TBM error trailer
77  k_FED31 = 31, // 31 indicates an event number error (TBM and FED event number mismatch)
78  k_FED32 = 32, // 32 indicates an incorrectly formatted Slink Header
79  k_FED33 = 33, // 33 indicates an incorrectly formatted Slink Trailer
80  k_FED34 = 34, // 34 indicates the evt size encoded in Slink Trailer is different than size found at raw2digi
81  k_FED35 = 35, // 35 indicates an invalid FED channel number
82  k_FED36 = 36, // 36 indicates an invalid ROC value
83  k_FED37 = 37, // 37 indicates an invalid dcol or pixel value
84  k_FED38 = 38 // 38 indicates the pixels on a ROC weren't read out from lowest to highest row and dcol value
85  };
86 
87  using MapToCodes = std::map<SiPixelFEDErrorCodes, std::string>;
88 
89  const MapToCodes errorCodeToStringMap = {{k_FED25, "FED25 error"},
90  {k_FED26, "FED26 error"},
91  {k_FED27, "FED27 error"},
92  {k_FED28, "FED28 error"},
93  {k_FED29, "FED29 error"},
94  {k_FED30, "FED30 error"},
95  {k_FED31, "FED31 error"}};
96 
97  const MapToCodes errorCodeToTypeMap = {{k_FED25, "ROC of 25"},
98  {k_FED26, "Gap word"},
99  {k_FED27, "Dummy word"},
100  {k_FED28, "FIFO full"},
101  {k_FED29, "Timeout"},
102  {k_FED30, "TBM error trailer"},
103  {k_FED31, "Event number"},
104  {k_FED32, "Slink header"},
105  {k_FED33, "Slink trailer"},
106  {k_FED34, "Event size"},
107  {k_FED35, "Invalid channel#"},
108  {k_FED36, "ROC value"},
109  {k_FED37, "Dcol or pixel value"},
110  {k_FED38, "Readout order"}};
111 } // namespace
112 
114 public:
116  ~SiPixelPhase1RawDataErrorComparator() override = default;
117  void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const& iRun, edm::EventSetup const& iSetup) override;
118  void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
119  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
120 
121 private:
127 
130  std::unordered_map<SiPixelFEDErrorCodes, MonitorElement*> h_nFEDErrors_;
131 
132  // name of the plugins
133  static constexpr const char* kName = "SiPixelPhase1RawDataErrorComparator";
134 
135  // Define the dimensions of the 2D array
137  static constexpr int nErrors = k_FED38 - k_FED25;
138 };
139 
140 //
141 // constructors
142 //
144  : pixelErrorSrcGPU_(iConfig.getParameter<edm::InputTag>("pixelErrorSrcGPU")),
145  pixelErrorSrcCPU_(iConfig.getParameter<edm::InputTag>("pixelErrorSrcCPU")),
146  tokenErrorsGPU_(consumes<edm::DetSetVector<SiPixelRawDataError>>(pixelErrorSrcGPU_)),
147  tokenErrorsCPU_(consumes<edm::DetSetVector<SiPixelRawDataError>>(pixelErrorSrcCPU_)),
148  topFolderName_(iConfig.getParameter<std::string>("topFolderName")) {}
149 
150 //
151 // -- Analyze
152 //
154  std::map<int, int> countsOnCPU;
155  std::map<int, int> countsOnGPU;
156 
157  std::array<std::array<int, nErrors>, nFEDs> countsMatrixOnCPU;
158  std::array<std::array<int, nErrors>, nFEDs> countsMatrixOnGPU;
159 
160  // initialize the counts for FED/error matrix
161  for (int i = 0; i < nFEDs; i++) {
162  for (int j = 0; j < nErrors; j++) {
163  countsMatrixOnCPU[i][j] = 0;
164  countsMatrixOnGPU[i][j] = 0;
165  }
166  }
167 
168  // initialize the counts for errors per type scatter plots
169  for (unsigned int j = k_FED25; j <= k_FED31; j++) {
170  countsOnCPU[j] = 0.;
171  countsOnGPU[j] = 0.;
172  }
173 
174  // check upfront if the error collection is present
176  iEvent.getByToken(tokenErrorsCPU_, inputFromCPU);
177  if (!inputFromCPU.isValid()) {
178  edm::LogError(kName) << "reference (cpu) SiPixelRawDataErrors collection (" << pixelErrorSrcCPU_.encode()
179  << ") not found; \n"
180  << "the comparison will not run.";
181  return;
182  }
183 
185  iEvent.getByToken(tokenErrorsGPU_, inputFromGPU);
186  if (!inputFromGPU.isValid()) {
187  edm::LogError(kName) << "target (gpu) SiPixelRawDataErrors collection (" << pixelErrorSrcGPU_.encode()
188  << ") not found; \n"
189  << "the comparison will not run.";
190  return;
191  }
192 
193  // fill the counters on host
194  uint errorsOnCPU{0};
195  for (auto it = inputFromCPU->begin(); it != inputFromCPU->end(); ++it) {
196  for (auto& siPixelRawDataError : *it) {
197  int fed = siPixelRawDataError.getFedId();
198  int type = siPixelRawDataError.getType();
199  DetId id = it->detId();
200 
201  // fill the error matrices for CPU
202  countsOnCPU[type] += 1;
203  countsMatrixOnCPU[fed - FEDNumbering::MINSiPixeluTCAFEDID][type - k_FED25] += 1;
204 
205  LogDebug(kName) << " (on cpu) FED: " << fed << " detid: " << id.rawId() << " type:" << type;
206  errorsOnCPU++;
207  }
208  }
209 
210  // fill the counters on device
211  uint errorsOnGPU{0};
212  for (auto it = inputFromGPU->begin(); it != inputFromGPU->end(); ++it) {
213  for (auto& siPixelRawDataError : *it) {
214  int fed = siPixelRawDataError.getFedId();
215  int type = siPixelRawDataError.getType();
216  DetId id = it->detId();
217 
218  // fill the error matrices for GPU
219  countsOnGPU[type] += 1;
220  countsMatrixOnGPU[fed - FEDNumbering::MINSiPixeluTCAFEDID][type - k_FED25] += 1;
221 
222  LogDebug(kName) << " (on gpu) FED: " << fed << " detid: " << id.rawId() << " type:" << type;
223  errorsOnGPU++;
224  }
225  }
226 
227  LogDebug(kName) << " on gpu found: " << errorsOnGPU << " on cpu found: " << errorsOnCPU;
228 
229  h_totFEDErrors_->Fill(errorsOnCPU, errorsOnGPU);
230 
231  // fill the correlations per error type
232  for (unsigned int j = k_FED25; j <= k_FED31; j++) {
233  const SiPixelFEDErrorCodes code = static_cast<SiPixelFEDErrorCodes>(j);
234  h_nFEDErrors_[code]->Fill(std::min(1000, countsOnCPU[j]), std::min(1000, countsOnGPU[j]));
235  }
236 
237  // fill the error unbalance per FEDid per error type
238  for (int i = 0; i < nFEDs; i++) {
239  for (int j = 0; j < nErrors; j++) {
240  if (countsMatrixOnGPU[i][j] != 0 || countsMatrixOnCPU[i][j] != 0) {
241  edm::LogVerbatim(kName) << "FED: " << i + FEDNumbering::MINSiPixeluTCAFEDID << " error: " << j + k_FED25
242  << " | GPU counts: " << countsMatrixOnGPU[i][j]
243  << " CPU counts:" << countsMatrixOnCPU[i][j];
245  j, i + FEDNumbering::MINSiPixeluTCAFEDID, countsMatrixOnGPU[i][j] - countsMatrixOnCPU[i][j]);
246  }
247  }
248  }
249 }
250 
251 //
252 // -- Book Histograms
253 //
255  edm::Run const& iRun,
256  edm::EventSetup const& iSetup) {
257  iBook.cd();
259 
261  iBook.book2I("FEErrorVsFEDIdUnbalance",
262  "difference (GPU-CPU) of FED errors per FEDid per error type;;FED Id number;GPU counts - CPU counts",
263  nErrors,
264  -0.5,
265  nErrors - 0.5,
266  nFEDs,
267  static_cast<double>(FEDNumbering::MINSiPixeluTCAFEDID) - 0.5,
268  static_cast<double>(FEDNumbering::MAXSiPixeluTCAFEDID) - 0.5);
269  for (int j = 0; j < nErrors; j++) {
270  const auto& errorCode = static_cast<SiPixelFEDErrorCodes>(j + k_FED25);
271  h_FEDerrorVsFEDIdUnbalance_->setBinLabel(j + 1, errorCodeToTypeMap.at(errorCode));
272  }
273 
274  h_totFEDErrors_ = make2DIfLog(iBook,
275  true,
276  true,
277  "nTotalFEDErrors",
278  "n. of total Pixel FEDErrors per event; CPU; GPU",
279  200,
280  log10(0.5),
281  log10(1000.),
282  200,
283  log10(0.5),
284  log10(1000.));
285 
286  for (const auto& element : errorCodeToStringMap) {
287  h_nFEDErrors_[element.first] = iBook.book2I(fmt::sprintf("nFED%i_Errors", static_cast<int>(element.first)),
288  fmt::sprintf("n. of %ss per event; CPU; GPU", element.second),
289  1000,
290  -0.5,
291  1000.5,
292  1000,
293  -0.5,
294  1000.5);
295  }
296 }
297 
300  desc.add<edm::InputTag>("pixelErrorSrcGPU", edm::InputTag("siPixelDigis@cuda"))
301  ->setComment("input GPU SiPixel FED errors");
302  desc.add<edm::InputTag>("pixelErrorSrcCPU", edm::InputTag("siPixelDigis@cpu"))
303  ->setComment("input CPU SiPixel FED errors");
304  desc.add<std::string>("topFolderName", "SiPixelHeterogeneous/PixelErrorCompareGPUvsCPU");
305  descriptions.addWithDefaultLabel(desc);
306 }
307 
Log< level::Info, true > LogVerbatim
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
std::unordered_map< SiPixelFEDErrorCodes, MonitorElement * > h_nFEDErrors_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
std::string encode() const
Definition: InputTag.cc:159
const edm::EDGetTokenT< edm::DetSetVector< SiPixelRawDataError > > tokenErrorsGPU_
Log< level::Error, false > LogError
~SiPixelPhase1RawDataErrorComparator() override=default
dqm::reco::DQMStore DQMStore
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
virtual void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
Definition: DetId.h:17
bool isValid() const
Definition: HandleBase.h:70
HLT enums.
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
MonitorElement * book2I(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:305
const edm::EDGetTokenT< edm::DetSetVector< SiPixelRawDataError > > tokenErrorsCPU_
void bookHistograms(DQMStore::IBooker &ibooker, edm::Run const &iRun, edm::EventSetup const &iSetup) override
Pixel error – collection of errors and error information.
Definition: Run.h:45
void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
#define LogDebug(id)