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 namespace {
27  constexpr int LINK_bits = 6;
28  constexpr int ROC_bits = 5;
29  constexpr int DCOL_bits = 5;
30  constexpr int PXID_bits = 8;
31  constexpr int ADC_bits = 8;
32 
33  // Add phase1 constants
34  // For phase1
35  //GO BACK TO OLD VALUES. THE 48-CHAN FED DOES NOT NEED A NEW FORMAT
36  // 28/9/16 d.k.
37  constexpr int LINK_bits1 = 6; // 7;
38  constexpr int ROC_bits1 = 5; // 4;
39  // Special for layer 1 bpix rocs 6/9/16 d.k. THIS STAYS.
40  constexpr int COL_bits1_l1 = 6;
41  constexpr int ROW_bits1_l1 = 7;
42 
43  // Moved to the header file, keep commented out unti the final version is done/
44  // constexpr int ADC_shift = 0;
45  // constexpr int PXID_shift = ADC_shift + ADC_bits;
46  // constexpr int DCOL_shift = PXID_shift + PXID_bits;
47  // constexpr int ROC_shift = DCOL_shift + DCOL_bits;
48  // constexpr int LINK_shift = ROC_shift + ROC_bits;
49  // constexpr PixelDataFormatter::Word32 LINK_mask = ~(~PixelDataFormatter::Word32(0) << LINK_bits);
50  // constexpr PixelDataFormatter::Word32 ROC_mask = ~(~PixelDataFormatter::Word32(0) << ROC_bits);
51  // constexpr PixelDataFormatter::Word32 DCOL_mask = ~(~PixelDataFormatter::Word32(0) << DCOL_bits);
52  // constexpr PixelDataFormatter::Word32 PXID_mask = ~(~PixelDataFormatter::Word32(0) << PXID_bits);
53  // constexpr PixelDataFormatter::Word32 ADC_mask = ~(~PixelDataFormatter::Word32(0) << ADC_bits);
54  //const bool DANEK = false;
55 }
56 
58  : theDigiCounter(0), theWordCounter(0), theCablingTree(map), badPixelInfo(nullptr), modulesToUnpack(nullptr), phase1(phase)
59 {
60  int s32 = sizeof(Word32);
61  int s64 = sizeof(Word64);
62  int s8 = sizeof(char);
63  if ( s8 != 1 || s32 != 4*s8 || s64 != 2*s32) {
64  LogError("UnexpectedSizes")
65  <<" unexpected sizes: "
66  <<" size of char is: " << s8
67  <<", size of Word32 is: " << s32
68  <<", size of Word64 is: " << s64
69  <<", send exception" ;
70  }
71  includeErrors = false;
72  useQualityInfo = false;
73  allDetDigis = 0;
74  hasDetDigis = 0;
75 
76  ADC_shift = 0;
77  PXID_shift = ADC_shift + ADC_bits;
78  DCOL_shift = PXID_shift + PXID_bits;
79  ROC_shift = DCOL_shift + DCOL_bits;
80 
81  if(phase1) { // for phase 1
82  LINK_shift = ROC_shift + ROC_bits1;
83  LINK_mask = ~(~PixelDataFormatter::Word32(0) << LINK_bits1);
84  ROC_mask = ~(~PixelDataFormatter::Word32(0) << ROC_bits1);
85  // special for layer 1 ROC
86  ROW_shift = ADC_shift + ADC_bits;
87  COL_shift = ROW_shift + ROW_bits1_l1;
88  COL_mask = ~(~PixelDataFormatter::Word32(0) << COL_bits1_l1);
89  ROW_mask = ~(~PixelDataFormatter::Word32(0) << ROW_bits1_l1);
90  maxROCIndex=8;
91 
92  } else { // for phase 0
93  LINK_shift = ROC_shift + ROC_bits;
94  LINK_mask = ~(~PixelDataFormatter::Word32(0) << LINK_bits);
95  ROC_mask = ~(~PixelDataFormatter::Word32(0) << ROC_bits);
96  maxROCIndex=25;
97  }
98 
99  DCOL_mask = ~(~PixelDataFormatter::Word32(0) << DCOL_bits);
100  PXID_mask = ~(~PixelDataFormatter::Word32(0) << PXID_bits);
101  ADC_mask = ~(~PixelDataFormatter::Word32(0) << ADC_bits);
102 
103 }
104 
106 {
107  includeErrors = ErrorStatus;
109 }
110 
111 void PixelDataFormatter::setQualityStatus(bool QualityStatus, const SiPixelQuality* QualityInfo)
112 {
113  useQualityInfo = QualityStatus;
114  badPixelInfo = QualityInfo;
115 }
116 
117 void PixelDataFormatter::setModulesToUnpack(const std::set<unsigned int> * moduleIds)
118 {
119  modulesToUnpack = moduleIds;
120 }
121 
123 {
124  theFrameReverter = reverter;
125 }
126 
127 void PixelDataFormatter::interpretRawData(bool& errorsInEvent, int fedId, const FEDRawData& rawData, Collection & digis, Errors& errors)
128 {
129  using namespace sipixelobjects;
130 
131  int nWords = rawData.size()/sizeof(Word64);
132  if (nWords==0) return;
133 
135 
136  // check CRC bit
137  const Word64* trailer = reinterpret_cast<const Word64* >(rawData.data())+(nWords-1);
138  if(!errorcheck.checkCRC(errorsInEvent, fedId, trailer, errors)) return;
139 
140  // check headers
141  const Word64* header = reinterpret_cast<const Word64* >(rawData.data()); header--;
142  bool moreHeaders = true;
143  while (moreHeaders) {
144  header++;
145  LogTrace("")<<"HEADER: " << print(*header);
146  bool headerStatus = errorcheck.checkHeader(errorsInEvent, fedId, header, errors);
147  moreHeaders = headerStatus;
148  }
149 
150  // check trailers
151  bool moreTrailers = true;
152  trailer++;
153  while (moreTrailers) {
154  trailer--;
155  LogTrace("")<<"TRAILER: " << print(*trailer);
156  bool trailerStatus = errorcheck.checkTrailer(errorsInEvent, fedId, nWords, trailer, errors);
157  moreTrailers = trailerStatus;
158  }
159 
160  // data words
161  theWordCounter += 2*(nWords-2);
162  LogTrace("")<<"data words: "<< (trailer-header-1);
163 
164  int link = -1;
165  int roc = -1;
166  int layer = 0;
167  PixelROC const * rocp=nullptr;
168  bool skipROC=false;
169  edm::DetSet<PixelDigi> * detDigis=nullptr;
170 
171  const Word32 * bw =(const Word32 *)(header+1);
172  const Word32 * ew =(const Word32 *)(trailer);
173  if ( *(ew-1) == 0 ) { ew--; theWordCounter--;}
174  for (auto word = bw; word < ew; ++word) {
175  LogTrace("")<<"DATA: " << print(*word);
176 
177  auto ww = *word;
178  if unlikely(ww==0) { theWordCounter--; continue;}
179  int nlink = (ww >> LINK_shift) & LINK_mask;
180  int nroc = (ww >> ROC_shift) & ROC_mask;
181 
182  //if(DANEK) cout<<" fed, link, roc "<<fedId<<" "<<nlink<<" "<<nroc<<endl;
183 
184  if ( (nlink!=link) | (nroc!=roc) ) { // new roc
185  link = nlink; roc=nroc;
186  skipROC = likely(roc<maxROCIndex) ? false : !errorcheck.checkROC(errorsInEvent, fedId, &converter, theCablingTree, ww, errors);
187  if (skipROC) continue;
188  rocp = converter.toRoc(link,roc);
189  if unlikely(!rocp) {
190  errorsInEvent = true;
191  errorcheck.conversionError(fedId, &converter, 2, ww, errors);
192  skipROC=true;
193  continue;
194  }
195  auto rawId = rocp->rawId();
196  bool barrel = PixelModuleName::isBarrel(rawId);
197  if(barrel) layer = PixelROC::bpixLayerPhase1(rawId);
198  else layer=0;
199 
200  //if(DANEK) cout<<" rocp "<<rocp->print()<<" layer "<<rocp->bpixLayerPhase1(rawId)<<" "
201  // <<layer<<" phase1 "<<phase1<<" rawid "<<rawId<<endl;
202 
203  if (useQualityInfo&(nullptr!=badPixelInfo)) {
204  short rocInDet = (short) rocp->idInDetUnit();
205  skipROC = badPixelInfo->IsRocBad(rawId, rocInDet);
206  if (skipROC) continue;
207  }
208  skipROC= modulesToUnpack && ( modulesToUnpack->find(rawId) == modulesToUnpack->end());
209  if (skipROC) continue;
210 
211  detDigis = &digis.find_or_insert(rawId);
212  if ( (*detDigis).empty() ) (*detDigis).data.reserve(32); // avoid the first relocations
213  }
214 
215  // skip is roc to be skipped ot invalid
216  if unlikely(skipROC || !rocp) continue;
217 
218  int adc = (ww >> ADC_shift) & ADC_mask;
219  std::unique_ptr<LocalPixel> local;
220 
221  if(phase1 && layer==1) { // special case for layer 1ROC
222  // for l1 roc use the roc column and row index instead of dcol and pixel index.
223  int col = (ww >> COL_shift) & COL_mask;
224  int row = (ww >> ROW_shift) & ROW_mask;
225  //if(DANEK) cout<<" layer 1: raw2digi "<<link<<" "<<roc<<" "
226  // <<col<<" "<<row<<" "<<adc<<endl;
227 
228  LocalPixel::RocRowCol localCR = { row, col }; // build pixel
229  //if(DANEK)cout<<localCR.rocCol<<" "<<localCR.rocRow<<endl;
230  if unlikely(!localCR.valid()) {
231  LogDebug("PixelDataFormatter::interpretRawData")
232  << "status #3";
233  errorsInEvent = true;
234  errorcheck.conversionError(fedId, &converter, 3, ww, errors);
235  continue;
236  }
237  local = std::make_unique<LocalPixel>(localCR); // local pixel coordinate
238  //if(DANEK) cout<<local->dcol()<<" "<<local->pxid()<<" "<<local->rocCol()<<" "<<local->rocRow()<<endl;
239 
240  } else { // phase0 and phase1 except bpix layer 1
241  int dcol = (ww >> DCOL_shift) & DCOL_mask;
242  int pxid = (ww >> PXID_shift) & PXID_mask;
243  //if(DANEK) cout<<" raw2digi "<<link<<" "<<roc<<" "
244  //<<dcol<<" "<<pxid<<" "<<adc<<" "<<layer<<endl;
245 
246  LocalPixel::DcolPxid localDP = { dcol, pxid };
247  //if(DANEK) cout<<localDP.dcol<<" "<<localDP.pxid<<endl;
248 
249  if unlikely(!localDP.valid()) {
250  LogDebug("PixelDataFormatter::interpretRawData")
251  << "status #3";
252  errorsInEvent = true;
253  errorcheck.conversionError(fedId, &converter, 3, ww, errors);
254  continue;
255  }
256  local = std::make_unique<LocalPixel>(localDP); // local pixel coordinate
257  //if(DANEK) cout<<local->dcol()<<" "<<local->pxid()<<" "<<local->rocCol()<<" "<<local->rocRow()<<endl;
258  }
259 
260  GlobalPixel global = rocp->toGlobal( *local ); // global pixel coordinate (in module)
261  (*detDigis).data.emplace_back(global.row, global.col, adc);
262  //if(DANEK) cout<<global.row<<" "<<global.col<<" "<<adc<<endl;
263  LogTrace("") << (*detDigis).data.back();
264  }
265 
266 }
267 
268 // I do not know what this was for or if it is needed? d.k. 10.14
269 // Keep it commented out until we are sure that it is not needed.
270 // void doVectorize(int const * __restrict__ w, int * __restrict__ row, int * __restrict__ col, int * __restrict__ valid, int N, PixelROC const * rocp) {
271 // for (int i=0; i<N; ++i) {
272 // auto ww = w[i];
273 // int dcol = (ww >> DCOL_shift) & DCOL_mask;
274 // int pxid = (ww >> PXID_shift) & PXID_mask;
275 // // int adc = (ww >> ADC_shift) & ADC_mask;
276 // LocalPixel::DcolPxid local = { dcol, pxid };
277 // valid[i] = local.valid();
278 // GlobalPixel global = rocp->toGlobal( LocalPixel(local) );
279 // row[i]=global.row; col[i]=global.col;
280 // }
281 // }
282 
283 
284 void PixelDataFormatter::formatRawData(unsigned int lvl1_ID, RawData & fedRawData, const Digis & digis)
285 {
286  std::map<int, vector<Word32> > words;
287 
288  // translate digis into 32-bit raw words and store in map indexed by Fed
289  for (Digis::const_iterator im = digis.begin(); im != digis.end(); im++) {
290  allDetDigis++;
291  cms_uint32_t rawId = im->first;
292  int layer=0;
293  bool barrel = PixelModuleName::isBarrel(rawId);
294  if(barrel) layer = PixelROC::bpixLayerPhase1(rawId);
295  //if(DANEK) cout<<" layer "<<layer<<" "<<phase1<<endl;
296 
297  hasDetDigis++;
298  const DetDigis & detDigis = im->second;
299  for (DetDigis::const_iterator it = detDigis.begin(); it != detDigis.end(); it++) {
300  theDigiCounter++;
301  const PixelDigi & digi = (*it);
302  int status=0;
303 
304  if(layer==1 && phase1) status = digi2wordPhase1Layer1( rawId, digi, words);
305  else status = digi2word( rawId, digi, words);
306 
307  if (status) {
308  LogError("FormatDataException")
309  <<" digi2word returns error #"<<status
310  <<" Ndigis: "<<theDigiCounter << endl
311  <<" detector: "<<rawId<< endl
312  << print(digi) <<endl;
313  } // if (status)
314  } // for (DetDigis
315  } // for (Digis
316  LogTrace(" allDetDigis/hasDetDigis : ") << allDetDigis<<"/"<<hasDetDigis;
317 
318  typedef std::map<int, vector<Word32> >::const_iterator RI;
319  for (RI feddata = words.begin(); feddata != words.end(); feddata++) {
320  int fedId = feddata->first;
321  // since raw words are written in the form of 64-bit packets
322  // add extra 32-bit word to make number of words even if necessary
323  if (words.find(fedId)->second.size() %2 != 0) words[fedId].push_back( Word32(0) );
324 
325  // size in Bytes; create output structure
326  int dataSize = words.find(fedId)->second.size() * sizeof(Word32);
327  int nHeaders = 1;
328  int nTrailers = 1;
329  dataSize += (nHeaders+nTrailers)*sizeof(Word64);
330  FEDRawData * rawData = new FEDRawData(dataSize);
331 
332  // get begining of data;
333  Word64 * word = reinterpret_cast<Word64* >(rawData->data());
334 
335  // write one header
336  FEDHeader::set( reinterpret_cast<unsigned char*>(word), 0, lvl1_ID, 0, fedId);
337  word++;
338 
339  // write data
340  unsigned int nWord32InFed = words.find(fedId)->second.size();
341  for (unsigned int i=0; i < nWord32InFed; i+=2) {
342  *word = (Word64(words.find(fedId)->second[i+1]) << 32 ) | words.find(fedId)->second[i];
343  LogDebug("PixelDataFormatter") << print(*word);
344  word++;
345  }
346 
347  // write one trailer
348  FEDTrailer::set( reinterpret_cast<unsigned char*>(word), dataSize/sizeof(Word64), 0,0,0);
349  word++;
350 
351  // check memory
352  if (word != reinterpret_cast<Word64* >(rawData->data()+dataSize)) {
353  string s = "** PROBLEM in PixelDataFormatter !!!";
354  throw cms::Exception(s);
355  } // if (word !=
356  fedRawData[fedId] = *rawData;
357  delete rawData;
358  } // for (RI feddata
359 }
360 
362  std::map<int, vector<Word32> > & words) const
363 {
364  LogDebug("PixelDataFormatter")
365 // <<" detId: " << detId
366  <<print(digi);
367 
368  DetectorIndex detector = {detId, digi.row(), digi.column()};
369  ElectronicIndex cabling;
370  int fedId = theFrameReverter->toCabling(cabling, detector);
371  if (fedId<0) return fedId;
372 
373  //if(DANEK) cout<<" digi2raw "<<detId<<" "<<digi.column()<<" "<<digi.row()<<" "<<digi.adc()<<" "
374  // <<cabling.link<<" "<<cabling.roc<<" "<<cabling.dcol<<" "<<cabling.pxid<<endl;
375 
376  Word32 word =
377  (cabling.link << LINK_shift)
378  | (cabling.roc << ROC_shift)
379  | (cabling.dcol << DCOL_shift)
380  | (cabling.pxid << PXID_shift)
381  | (digi.adc() << ADC_shift);
382  words[fedId].push_back(word);
383  theWordCounter++;
384 
385  return 0;
386 }
388  std::map<int, vector<Word32> > & words) const
389 {
390  LogDebug("PixelDataFormatter")
391 // <<" detId: " << detId
392  <<print(digi);
393 
394  DetectorIndex detector = {detId, digi.row(), digi.column()};
395  ElectronicIndex cabling;
396  int fedId = theFrameReverter->toCabling(cabling, detector);
397  if (fedId<0) return fedId;
398 
399  int col = ((cabling.dcol)*2) + ((cabling.pxid)%2);
400  int row = LocalPixel::numRowsInRoc - ((cabling.pxid)/2);
401 
402  //if(DANEK) cout<<" layer 1: digi2raw "<<detId<<" "<<digi.column()<<" "<<digi.row()<<" "<<digi.adc()<<" "
403  // <<cabling.link<<" "<<cabling.roc<<" "<<cabling.dcol<<" "<<cabling.pxid<<" "
404  // <<col<<" "<<row<<endl;
405 
406  Word32 word =
407  (cabling.link << LINK_shift)
408  | (cabling.roc << ROC_shift)
409  | (col << COL_shift)
410  | (row << ROW_shift)
411  | (digi.adc() << ADC_shift);
412  words[fedId].push_back(word);
413  theWordCounter++;
414 
415  return 0;
416 }
417 
418 
419 // obsolete...
421  const bool includeErrors, const bool useQuality, const Word32 & word, Digis & digis) const
422 {
423  // do not interpret false digis
424  if (word == 0 ) return 0;
425 
426  ElectronicIndex cabling;
427  cabling.dcol = (word >> DCOL_shift) & DCOL_mask;
428  cabling.pxid = (word >> PXID_shift) & PXID_mask;
429  cabling.link = (word >> LINK_shift) & LINK_mask;
430  cabling.roc = (word >> ROC_shift) & ROC_mask;
431  int adc = (word >> ADC_shift) & ADC_mask;
432 
433  if (debug) {
434  LocalPixel::DcolPxid pixel = {cabling.dcol,cabling.pxid};
435  LocalPixel local(pixel);
436  LogTrace("")<<" link: "<<cabling.link<<", roc: "<<cabling.roc
437  <<" rocRow: "<<local.rocRow()<<", rocCol:"<<local.rocCol()
438  <<" (dcol: "<<cabling.dcol<<", pxid:"<<cabling.pxid<<"), adc:"<<adc;
439  }
440 
441  if (!converter) return 0;
442 
443  DetectorIndex detIdx;
444  int status = converter->toDetector(cabling, detIdx);
445  if (status) return status;
446 
447  // exclude ROC(raw) based on bad ROC list bad in SiPixelQuality
448  // enable: process.siPixelDigis.UseQualityInfo = True
449  // 20-10-2010 A.Y.
450  if (useQuality&&badPixelInfo) {
451  CablingPathToDetUnit path = {static_cast<unsigned int>(fedId),
452  static_cast<unsigned int>(cabling.link),
453  static_cast<unsigned int>(cabling.roc)};
454  const PixelROC * roc = theCablingTree->findItem(path);
455  short rocInDet = (short) roc->idInDetUnit();
456  bool badROC = badPixelInfo->IsRocBad(detIdx.rawId, rocInDet);
457  if (badROC) return 0;
458  }
459 
460  if (modulesToUnpack && modulesToUnpack->find(detIdx.rawId) == modulesToUnpack->end()) return 0;
461 
462  digis[detIdx.rawId].emplace_back(detIdx.row, detIdx.col, adc);
463 
464  theDigiCounter++;
465 
466  if (debug) LogTrace("") << digis[detIdx.rawId].back();
467  return 0;
468 }
469 
471 {
472  ostringstream str;
473  str << " DIGI: row: " << digi.row() <<", col: " << digi.column() <<", adc: " << digi.adc();
474  return str.str();
475 }
476 
478 {
479  ostringstream str;
480  str <<"word64: " << reinterpret_cast<const bitset<64>&> (word);
481  return str.str();
482 }
483 
#define LogDebug(id)
int adc(sample_type sample)
get the ADC sample (12 bits)
int row() const
Definition: PixelDigi.h:56
std::map< cms_uint32_t, DetDigis > Digis
void passFrameReverter(const SiPixelFrameReverter *reverter)
void setErrorStatus(bool ErrorStatus)
bool IsRocBad(const uint32_t &detid, const short &rocNb) const
int digi2wordPhase1Layer1(cms_uint32_t detId, const PixelDigi &digi, std::map< int, std::vector< Word32 > > &words) const
#define nullptr
#define constexpr
size_t size() const
Lenght of the data buffer in bytes.
Definition: FEDRawData.h:47
identify pixel inside single ROC
Definition: LocalPixel.h:7
bool checkROC(bool &errorsInEvent, int fedId, const SiPixelFrameConverter *converter, const SiPixelFedCabling *theCablingTree, Word32 &errorWord, Errors &errors)
int toCabling(sipixelobjects::ElectronicIndex &cabling, const sipixelobjects::DetectorIndex &detector) const
bool checkTrailer(bool &errorsInEvent, int fedId, unsigned int nWords, const Word64 *trailer, Errors &errors)
Definition: ErrorChecker.cc:89
reference find_or_insert(det_id_type id)
Definition: DetSetVector.h:254
global coordinates (row and column in DetUnit, as in PixelDigi)
Definition: GlobalPixel.h:6
#define unlikely(x)
#define likely(x)
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:60
bool checkCRC(bool &errorsInEvent, int fedId, const Word64 *trailer, Errors &errors)
Definition: ErrorChecker.cc:58
unsigned int idInDetUnit() const
id of this ROC in DetUnit etermined by token path
Definition: PixelROC.h:40
uint32_t rawId() const
return the DetUnit to which this ROC belongs to.
Definition: PixelROC.h:37
std::string print(const PixelDigi &digi) const
unsigned short adc() const
Definition: PixelDigi.h:59
void formatRawData(unsigned int lvl1_ID, RawData &fedRawData, const Digis &digis)
void setQualityStatus(bool QualityStatus, const SiPixelQuality *QualityInfo)
const std::set< unsigned int > * modulesToUnpack
const SiPixelQuality * badPixelInfo
void setErrorStatus(bool ErrorStatus)
Definition: ErrorChecker.cc:53
unsigned int cms_uint32_t
Definition: typedefs.h:15
#define LogTrace(id)
double collumn and pixel ID in double collumn representation
Definition: LocalPixel.h:22
std::map< cms_uint32_t, DetErrors > Errors
virtual bool isBarrel() const
true for barrel modules
void interpretRawData(bool &errorsInEvent, int fedId, const FEDRawData &data, Collection &digis, Errors &errors)
void setModulesToUnpack(const std::set< unsigned int > *moduleIds)
virtual const sipixelobjects::PixelROC * findItem(const sipixelobjects::CablingPathToDetUnit &) const =0
const SiPixelFrameReverter * theFrameReverter
void conversionError(int fedId, const SiPixelFrameConverter *converter, int status, Word32 &errorWord, Errors &errors)
row and collumn in ROC representation
Definition: LocalPixel.h:15
std::map< int, FEDRawData > RawData
std::vector< PixelDigi > DetDigis
int digi2word(cms_uint32_t detId, const PixelDigi &digi, std::map< int, std::vector< Word32 > > &words) const
collection_type data
Definition: DetSet.h:78
HLT enums.
sipixelobjects::PixelROC const * toRoc(int link, int roc) const
const unsigned char * data() const
Return a const pointer to the beginning of the data buffer.
Definition: FEDRawData.cc:28
col
Definition: cuy.py:1008
int toDetector(const sipixelobjects::ElectronicIndex &cabling, sipixelobjects::DetectorIndex &detector) const
int word2digi(const int fedId, const SiPixelFrameConverter *converter, const bool includeError, const bool useQuality, const Word32 &word, Digis &digis) const
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:47
SiPixelFedCabling const * theCablingTree
int column() const
Definition: PixelDigi.h:57
PixelDataFormatter(const SiPixelFedCabling *map, bool phase1=false)
GlobalPixel toGlobal(const LocalPixel &loc) const
Definition: PixelROC.h:59
bool checkHeader(bool &errorsInEvent, int fedId, const Word64 *header, Errors &errors)
Definition: ErrorChecker.cc:71