CMS 3D CMS Logo

PixelDataFormatter.cc
Go to the documentation of this file.
2 
6 
8 
12 
15 
18 
19 #include <bitset>
20 #include <sstream>
21 #include <iostream>
22 
23 using namespace std;
24 using namespace edm;
25 using namespace sipixelobjects;
26 using namespace sipixelconstants;
27 
29  : theDigiCounter_(0),
30  theWordCounter_(0),
31  theCablingTree_(map),
32  badPixelInfo_(nullptr),
33  modulesToUnpack_(nullptr),
34  phase1_(phase) {
35  int s32 = sizeof(Word32);
36  int s64 = sizeof(Word64);
37  int s8 = sizeof(char);
38  if (s8 != 1 || s32 != 4 * s8 || s64 != 2 * s32) {
39  LogError("UnexpectedSizes") << " unexpected sizes: "
40  << " size of char is: " << s8 << ", size of Word32 is: " << s32
41  << ", size of Word64 is: " << s64 << ", send exception";
42  }
43  includeErrors_ = false;
44  useQualityInfo_ = false;
45  allDetDigis_ = 0;
46  hasDetDigis_ = 0;
47 
48  if (phase1_) {
49  maxROCIndex_ = 8;
50  errorcheck_ = std::unique_ptr<ErrorCheckerBase>(new ErrorChecker());
51  } else {
52  maxROCIndex_ = 25;
53  errorcheck_ = std::unique_ptr<ErrorCheckerBase>(new ErrorCheckerPhase0());
54  }
55 }
56 
57 void PixelDataFormatter::setErrorStatus(bool ErrorStatus) {
58  includeErrors_ = ErrorStatus;
59  errorcheck_->setErrorStatus(includeErrors_);
60 }
61 
62 void PixelDataFormatter::setQualityStatus(bool QualityStatus, const SiPixelQuality* QualityInfo) {
63  useQualityInfo_ = QualityStatus;
64  badPixelInfo_ = QualityInfo;
65 }
66 
67 void PixelDataFormatter::setModulesToUnpack(const std::set<unsigned int>* moduleIds) { modulesToUnpack_ = moduleIds; }
68 
70 
72  bool& errorsInEvent, int fedId, const FEDRawData& rawData, Collection& digis, Errors& errors) {
73  using namespace sipixelobjects;
74 
75  int nWords = rawData.size() / sizeof(Word64);
76  if (nWords == 0)
77  return;
78 
80 
81  // check CRC bit
82  const Word64* trailer = reinterpret_cast<const Word64*>(rawData.data()) + (nWords - 1);
83  if (!errorcheck_->checkCRC(errorsInEvent, fedId, trailer, errors))
84  return;
85 
86  // check headers
87  const Word64* header = reinterpret_cast<const Word64*>(rawData.data());
88  header--;
89  bool moreHeaders = true;
90  while (moreHeaders) {
91  header++;
92  LogTrace("") << "HEADER: " << print(*header);
93  bool headerStatus = errorcheck_->checkHeader(errorsInEvent, fedId, header, errors);
94  moreHeaders = headerStatus;
95  }
96 
97  // check trailers
98  bool moreTrailers = true;
99  trailer++;
100  while (moreTrailers) {
101  trailer--;
102  LogTrace("") << "TRAILER: " << print(*trailer);
103  bool trailerStatus = errorcheck_->checkTrailer(errorsInEvent, fedId, nWords, trailer, errors);
104  moreTrailers = trailerStatus;
105  }
106 
107  // data words
108  theWordCounter_ += 2 * (nWords - 2);
109  LogTrace("") << "data words: " << (trailer - header - 1);
110 
111  int link = -1;
112  int roc = -1;
113  int layer = 0;
114  unsigned int rawId = 0;
115  unsigned int nrawId = 0;
116  PixelROC const* rocp = nullptr;
117  bool skipROC = false;
118  edm::DetSet<PixelDigi>* detDigis = nullptr;
119 
120  const Word32* bw = (const Word32*)(header + 1);
121  const Word32* ew = (const Word32*)(trailer);
122  if (*(ew - 1) == 0) {
123  ew--;
124  theWordCounter_--;
125  }
126  for (auto word = bw; word < ew; ++word) {
127  LogTrace("") << "DATA: " << print(*word);
128 
129  auto ww = *word;
130  if UNLIKELY (ww == 0) {
131  theWordCounter_--;
132  continue;
133  }
134  int nlink = getLink(ww);
135  int nroc = getROC(ww);
136 
137  if ((nlink != link) | (nroc != roc)) { // new roc
138  link = nlink;
139  roc = nroc;
140  skipROC = LIKELY(roc < maxROCIndex_)
141  ? false
142  : !errorcheck_->checkROC(errorsInEvent, fedId, &converter, theCablingTree_, ww, errors);
143  if (skipROC)
144  continue;
145  rocp = converter.toRoc(link, roc);
146  if UNLIKELY (!rocp) {
147  errorsInEvent = true;
148  errorcheck_->conversionError(fedId, &converter, 2, ww, errors);
149  skipROC = true;
150  continue;
151  }
152  rawId = rocp->rawId();
154  if (barrel)
155  layer = PixelROC::bpixLayerPhase1(rawId);
156  else
157  layer = 0;
158 
159  if (useQualityInfo_ & (nullptr != badPixelInfo_)) {
160  short rocInDet = (short)rocp->idInDetUnit();
161  skipROC = badPixelInfo_->IsRocBad(rawId, rocInDet);
162  if (skipROC)
163  continue;
164  }
165  skipROC = modulesToUnpack_ && (modulesToUnpack_->find(rawId) == modulesToUnpack_->end());
166  if (skipROC)
167  continue;
168  }
169 
170  // skip is roc to be skipped ot invalid
171  if UNLIKELY (skipROC || !rocp)
172  continue;
173 
174  int adc = getADC(ww);
175  std::unique_ptr<LocalPixel> local;
176 
177  if (phase1_ && layer == 1) { // special case for layer 1ROC
178  // for l1 roc use the roc column and row index instead of dcol and pixel index.
179  int col = getCol(ww);
180  int row = getRow(ww);
181 
182  LocalPixel::RocRowCol localCR = {row, col}; // build pixel
183  if UNLIKELY (!localCR.valid()) {
184  LogDebug("PixelDataFormatter::interpretRawData") << "status #3";
185  errorsInEvent = true;
186  errorcheck_->conversionError(fedId, &converter, 3, ww, errors);
187  continue;
188  }
189  local = std::make_unique<LocalPixel>(localCR); // local pixel coordinate
190 
191  } else { // phase0 and phase1 except bpix layer 1
192  int dcol = getDCol(ww);
193  int pxid = getPxId(ww);
194  LocalPixel::DcolPxid localDP = {dcol, pxid};
195 
196  if UNLIKELY (!localDP.valid()) {
197  LogDebug("PixelDataFormatter::interpretRawData") << "status #3";
198  errorsInEvent = true;
199  errorcheck_->conversionError(fedId, &converter, 3, ww, errors);
200  continue;
201  }
202  local = std::make_unique<LocalPixel>(localDP); // local pixel coordinate
203  }
204 
205  if (nrawId != rawId) {
206  nrawId = rawId;
207  detDigis = &digis.find_or_insert(rawId);
208  if ((*detDigis).empty()) {
209  (*detDigis).data.reserve(32); // avoid the first relocations
210  }
211  }
212 
213  if (detDigis) {
214  GlobalPixel global = rocp->toGlobal(*local); // global pixel coordinate (in module)
215  (*detDigis).data.emplace_back(global.row, global.col, adc);
216  LogTrace("") << (*detDigis).data.back();
217  } else {
218  LogError("NullPointerException") << "@SUB=PixelDataFormatter::interpretRawData"
219  << "DetSet pointer not set. This is not supposed to happen.";
220  }
221  }
222 }
223 
224 void PixelDataFormatter::formatRawData(unsigned int lvl1_ID,
226  const Digis& digis,
227  const BadChannels& badChannels) {
228  std::map<int, vector<Word32> > words;
229 
230  // translate digis into 32-bit raw words and store in map indexed by Fed
231  for (Digis::const_iterator im = digis.begin(); im != digis.end(); im++) {
232  allDetDigis_++;
233  cms_uint32_t rawId = im->first;
234  int layer = 0;
236  if (barrel)
237  layer = PixelROC::bpixLayerPhase1(rawId);
238 
239  BadChannels::const_iterator detBadChannels = badChannels.find(rawId);
240 
241  hasDetDigis_++;
242  const DetDigis& detDigis = im->second;
243  for (DetDigis::const_iterator it = detDigis.begin(); it != detDigis.end(); it++) {
244  theDigiCounter_++;
245  const PixelDigi& digi = (*it);
246  int fedId = 0;
247 
248  if (layer == 1 && phase1_)
249  fedId = digi2wordPhase1Layer1(rawId, digi, words);
250  else
251  fedId = digi2word(rawId, digi, words);
252 
253  if (fedId < 0) {
254  LogError("FormatDataException") << " digi2word returns error #" << fedId << " Ndigis: " << theDigiCounter_
255  << endl
256  << " detector: " << rawId << endl
257  << print(digi) << endl;
258  } else if (detBadChannels != badChannels.end()) {
259  auto badChannel =
260  std::find_if(detBadChannels->second.begin(), detBadChannels->second.end(), [&](const PixelFEDChannel& ch) {
261  return (int(ch.fed) == fedId && ch.link == getLink(words[fedId].back()));
262  });
263  if (badChannel != detBadChannels->second.end()) {
264  LogError("FormatDataException") << " while marked bad, found digi for FED " << fedId << " Link "
265  << getLink(words[fedId].back()) << " on module " << rawId << endl
266  << print(digi) << endl;
267  }
268  } // if (fedId)
269  } // for (DetDigis
270  } // for (Digis
271  LogTrace(" allDetDigis_/hasDetDigis_ : ") << allDetDigis_ << "/" << hasDetDigis_;
272 
273  // fill FED error 25 words
274  for (const auto& detBadChannels : badChannels) {
275  for (const auto& badChannel : detBadChannels.second) {
276  unsigned int FEDError25 = 25;
277  Word32 word = (badChannel.link << LINK_shift) | (FEDError25 << ROC_shift);
278  words[badChannel.fed].push_back(word);
279  theWordCounter_++;
280  }
281  }
282 
283  typedef std::map<int, vector<Word32> >::const_iterator RI;
284  for (RI feddata = words.begin(); feddata != words.end(); feddata++) {
285  int fedId = feddata->first;
286  // since raw words are written in the form of 64-bit packets
287  // add extra 32-bit word to make number of words even if necessary
288  if (words.find(fedId)->second.size() % 2 != 0)
289  words[fedId].push_back(Word32(0));
290 
291  // size in Bytes; create output structure
292  int dataSize = words.find(fedId)->second.size() * sizeof(Word32);
293  int nHeaders = 1;
294  int nTrailers = 1;
295  dataSize += (nHeaders + nTrailers) * sizeof(Word64);
296  FEDRawData* rawData = new FEDRawData(dataSize);
297 
298  // get begining of data;
299  Word64* word = reinterpret_cast<Word64*>(rawData->data());
300 
301  // write one header
302  FEDHeader::set(reinterpret_cast<unsigned char*>(word), 0, lvl1_ID, 0, fedId);
303  word++;
304 
305  // write data
306  unsigned int nWord32InFed = words.find(fedId)->second.size();
307  for (unsigned int i = 0; i < nWord32InFed; i += 2) {
308  *word = (Word64(words.find(fedId)->second[i + 1]) << 32) | words.find(fedId)->second[i];
309  LogDebug("PixelDataFormatter") << print(*word);
310  word++;
311  }
312 
313  // write one trailer
314  FEDTrailer::set(reinterpret_cast<unsigned char*>(word), dataSize / sizeof(Word64), 0, 0, 0);
315  word++;
316 
317  // check memory
318  if (word != reinterpret_cast<Word64*>(rawData->data() + dataSize)) {
319  string s = "** PROBLEM in PixelDataFormatter !!!";
320  throw cms::Exception(s);
321  } // if (word !=
323  delete rawData;
324  } // for (RI feddata
325 }
326 
328  const PixelDigi& digi,
329  std::map<int, vector<Word32> >& words) const {
330  LogDebug("PixelDataFormatter") << print(digi);
331 
332  DetectorIndex detector = {detId, digi.row(), digi.column()};
333  ElectronicIndex cabling;
334  int fedId = theFrameReverter_->toCabling(cabling, detector);
335  if (fedId < 0)
336  return fedId;
337 
338  Word32 word = (cabling.link << LINK_shift) | (cabling.roc << ROC_shift) | (cabling.dcol << DCOL_shift) |
339  (cabling.pxid << PXID_shift) | (digi.adc() << ADC_shift);
340  words[fedId].push_back(word);
341  theWordCounter_++;
342 
343  return fedId;
344 }
346  const PixelDigi& digi,
347  std::map<int, vector<Word32> >& words) const {
348  LogDebug("PixelDataFormatter") << print(digi);
349 
350  DetectorIndex detector = {detId, digi.row(), digi.column()};
351  ElectronicIndex cabling;
352  int fedId = theFrameReverter_->toCabling(cabling, detector);
353  if (fedId < 0)
354  return fedId;
355 
356  int col = ((cabling.dcol) * 2) + ((cabling.pxid) % 2);
357  int row = LocalPixel::numRowsInRoc - ((cabling.pxid) / 2);
358 
359  Word32 word = (cabling.link << LINK_shift) | (cabling.roc << ROC_shift) | (col << COL_shift) | (row << ROW_shift) |
360  (digi.adc() << ADC_shift);
361  words[fedId].push_back(word);
362  theWordCounter_++;
363 
364  return fedId;
365 }
366 
368  ostringstream str;
369  str << " DIGI: row: " << digi.row() << ", col: " << digi.column() << ", adc: " << digi.adc();
370  return str.str();
371 }
372 
374  ostringstream str;
375  str << "word64: " << reinterpret_cast<const bitset<64>&>(word);
376  return str.str();
377 }
378 
380  std::vector<int> const& tkerrorlist,
381  std::vector<int> const& usererrorlist,
383  DetIdCollection& tkerror_detidcollection,
384  DetIdCollection& usererror_detidcollection,
385  edmNew::DetSetVector<PixelFEDChannel>& disabled_channelcollection,
386  DetErrors& nodeterrors) {
387  const uint32_t dummyDetId = 0xffffffff;
388  for (const auto& [errorDetId, rawErrorsVec] : errors) {
389  if (errorDetId == dummyDetId) { // errors given dummy detId must be sorted by Fed
390  nodeterrors.insert(nodeterrors.end(), rawErrorsVec.begin(), rawErrorsVec.end());
391  } else {
392  edm::DetSet<SiPixelRawDataError>& errorDetSet = errorcollection.find_or_insert(errorDetId);
393  errorDetSet.data.insert(errorDetSet.data.end(), rawErrorsVec.begin(), rawErrorsVec.end());
394  // Fill detid of the detectors where there is error AND the error number is listed
395  // in the configurable error list in the job option cfi.
396  // Code needs to be here, because there can be a set of errors for each
397  // entry in the for loop over PixelDataFormatter::Errors
398 
399  std::vector<PixelFEDChannel> disabledChannelsDetSet;
400 
401  for (auto const& aPixelError : errorDetSet) {
402  // For the time being, we extend the error handling functionality with ErrorType 25
403  // In the future, we should sort out how the usage of tkerrorlist can be generalized
404  if (phase1_ && aPixelError.getType() == 25) {
405  int fedId = aPixelError.getFedId();
407  if (fed) {
408  cms_uint32_t linkId = getLink(aPixelError.getWord32());
409  const sipixelobjects::PixelFEDLink* link = fed->link(linkId);
410  if (link) {
411  // The "offline" 0..15 numbering is fixed by definition, also, the FrameConversion depends on it
412  // in contrast, the ROC-in-channel numbering is determined by hardware --> better to use the "offline" scheme
413  PixelFEDChannel ch = {fed->id(), linkId, 25, 0};
414  for (unsigned int iRoc = 1; iRoc <= link->numberOfROCs(); iRoc++) {
415  const sipixelobjects::PixelROC* roc = link->roc(iRoc);
416  if (roc->idInDetUnit() < ch.roc_first)
417  ch.roc_first = roc->idInDetUnit();
418  if (roc->idInDetUnit() > ch.roc_last)
419  ch.roc_last = roc->idInDetUnit();
420  }
421  disabledChannelsDetSet.push_back(ch);
422  }
423  }
424  } else {
425  // fill list of detIds to be turned off by tracking
426  if (!tkerrorlist.empty()) {
427  auto it_find = std::find(tkerrorlist.begin(), tkerrorlist.end(), aPixelError.getType());
428  if (it_find != tkerrorlist.end()) {
429  tkerror_detidcollection.push_back(errorDetId);
430  }
431  }
432  }
433 
434  // fill list of detIds with errors to be studied
435  if (!usererrorlist.empty()) {
436  auto it_find = std::find(usererrorlist.begin(), usererrorlist.end(), aPixelError.getType());
437  if (it_find != usererrorlist.end()) {
438  usererror_detidcollection.push_back(errorDetId);
439  }
440  }
441 
442  } // loop on DetSet of errors
443 
444  if (!disabledChannelsDetSet.empty()) {
445  disabled_channelcollection.insert(errorDetId, disabledChannelsDetSet.data(), disabledChannelsDetSet.size());
446  }
447 
448  } // if error assigned to a real DetId
449  } // loop on errors in event for this FED
450 }
std::string print(const PixelDigi &digi) const
constexpr uint32_t getRow(uint32_t ww)
constexpr cms_uint32_t dummyDetId
const PixelFEDLink * link(unsigned int id) const
return link identified by id. Link id&#39;s are ranged [1, numberOfLinks]
int toCabling(sipixelobjects::ElectronicIndex &cabling, const sipixelobjects::DetectorIndex &detector) const
constexpr uint32_t DCOL_shift
constexpr uint32_t PXID_shift
bool IsRocBad(const uint32_t &detid, const short &rocNb) const
int digi2word(cms_uint32_t detId, const PixelDigi &digi, FEDWordsMap &words) const
void passFrameReverter(const SiPixelFrameReverter *reverter)
void formatRawData(unsigned int lvl1_ID, RawData &fedRawData, const Digis &digis, const BadChannels &badChannels)
int digi2wordPhase1Layer1(cms_uint32_t detId, const PixelDigi &digi, FEDWordsMap &words) const
constexpr uint32_t getCol(uint32_t ww)
cms_uint32_t errorDetId(const SiPixelFrameConverter *converter, int fedId, int errorType, const Word32 &word) const
constexpr uint32_t ROC_shift
#define LIKELY(x)
Definition: Likely.h:20
void push_back(T const &t)
Definition: EDCollection.h:60
constexpr uint32_t ADC_shift
void setErrorStatus(bool ErrorStatus)
Log< level::Error, false > LogError
std::vector< PixelDigi > DetDigis
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
unsigned short adc() const
Definition: PixelDigi.h:64
const ModuleIDSet * modulesToUnpack_
#define LogTrace(id)
unsigned int roc_last
reference find_or_insert(det_id_type id)
Definition: DetSetVector.h:234
void setModulesToUnpack(const ModuleIDSet *moduleIds)
global coordinates (row and column in DetUnit, as in PixelDigi)
Definition: GlobalPixel.h:6
std::map< cms_uint32_t, DetBadChannels > BadChannels
constexpr uint32_t getPxId(uint32_t ww)
uint64_t word
static void set(unsigned char *trailer, uint32_t lenght, uint16_t crc, uint8_t evt_stat, uint8_t tts, bool moreTrailers=false)
Set all fields in the trailer.
Definition: FEDTrailer.cc:31
const SiPixelFrameReverter * theFrameReverter_
GlobalPixel toGlobal(const LocalPixel &loc) const
Definition: PixelROC.h:55
std::map< cms_uint32_t, DetDigis > Digis
constexpr uint32_t getDCol(uint32_t ww)
cms_uint32_t Word32
constexpr uint32_t getADC(uint32_t ww)
constexpr uint32_t getLink(uint32_t ww)
void setQualityStatus(bool QualityStatus, const SiPixelQuality *QualityInfo)
unsigned int cms_uint32_t
Definition: typedefs.h:15
int column() const
Definition: PixelDigi.h:57
DetSet insert(id_type iid, data_type const *idata, size_type isize)
double collumn and pixel ID in double collumn representation
Definition: LocalPixel.h:19
std::map< int, FEDRawData > RawData
constexpr uint32_t getROC(uint32_t ww)
const PixelFEDCabling * fed(unsigned int idFed) const
get fed identified by its id
uint64_t Word64
void interpretRawData(bool &errorsInEvent, int fedId, const FEDRawData &data, Collection &digis, Errors &errors)
constexpr uint32_t LINK_shift
virtual bool isBarrel() const
true for barrel modules
std::unique_ptr< ErrorCheckerBase > errorcheck_
SiPixelFedCablingTree const * theCablingTree_
unsigned int roc_first
const SiPixelQuality * badPixelInfo_
row and collumn in ROC representation
Definition: LocalPixel.h:13
deadvectors [0] push_back({0.0175431, 0.538005, 6.80997, 13.29})
PixelDataFormatter(const SiPixelFedCablingTree *map, bool phase1_=false)
collection_type data
Definition: DetSet.h:80
HLT enums.
col
Definition: cuy.py:1009
void unpackFEDErrors(Errors const &errors, std::vector< int > const &tkerrorlist, std::vector< int > const &usererrorlist, edm::DetSetVector< SiPixelRawDataError > &errorcollection, DetIdCollection &tkerror_detidcollection, DetIdCollection &usererror_detidcollection, edmNew::DetSetVector< PixelFEDChannel > &disabled_channelcollection, DetErrors &nodeterrors)
Definition: errors.py:1
static void set(unsigned char *header, uint8_t triggerType, uint32_t lvl1ID, uint16_t bxID, uint16_t sourceID, uint8_t version=0, bool moreHeaders=false)
Set all fields in the header.
Definition: FEDHeader.cc:25
#define UNLIKELY(x)
Definition: Likely.h:21
int row() const
Definition: PixelDigi.h:54
#define str(s)
uint32_t rawId() const
return the DetUnit to which this ROC belongs to.
Definition: PixelROC.h:34
unsigned int idInDetUnit() const
id of this ROC in DetUnit etermined by token path
Definition: PixelROC.h:37
std::vector< SiPixelRawDataError > DetErrors
std::map< cms_uint32_t, DetErrors > Errors
uint16_t *__restrict__ uint16_t const *__restrict__ adc
cms_uint64_t Word64
#define LogDebug(id)