CMS 3D CMS Logo

QIE11DataFrame.h
Go to the documentation of this file.
1 #ifndef DATAFORMATS_HCALDIGI_QIE11DATAFRAME_H
2 #define DATAFORMATS_HCALDIGI_QIE11DATAFRAME_H
3 
6 #include <ostream>
7 
12 public:
13 
14  static const int WORDS_PER_SAMPLE = 1;
15  static const int HEADER_WORDS = 1;
16  static const int FLAG_WORDS = 1;
17 
18  static const int OFFSET_FLAVOR = 12;
19  static const int MASK_FLAVOR = 0x7;
20  static const int FLAVOR_HB = 3;
21  static const int MASK_LINKERROR = 0x800;
22 
25 
26  class Sample {
27  public:
29  static const int MASK_ADC = 0xFF;
30  static const int MASK_TDC_HE = 0x3F;
31  static const int MASK_TDC_HB = 0x3;
32  static const int OFFSET_TDC = 8; // 8 bits
33  static const int MASK_SOI = 0x4000;
34  static const int MASK_LE_HB = 0x2000;
35  static const int MASK_CAPID = 0x3;
36  static const int MASK_CAPID_INV_HB = 0xF3FF;
37  static const int MASK_CAPID_KEEP_HB = 0x0C00;
38  static const int OFFSET_CAPID_HE = 8;
39  static const int OFFSET_CAPID_HB = 10;
40  constexpr int flavor() const { return ((frame_[0]>>OFFSET_FLAVOR)&MASK_FLAVOR); }
41  constexpr int adc() const { return frame_[i_]&MASK_ADC; }
42  constexpr int tdc() const { return (frame_[i_]>>OFFSET_TDC)&((flavor()==FLAVOR_HB)?(MASK_TDC_HB):(MASK_TDC_HE)); }
43  constexpr bool soi() const { return frame_[i_]&MASK_SOI; }
44  constexpr int capid() const { return (flavor()==FLAVOR_HB)? ((frame_[i_]>>OFFSET_CAPID_HB)&MASK_CAPID) :
45  ((((frame_[0]>>OFFSET_CAPID_HE)&MASK_CAPID)+i_-HEADER_WORDS)&MASK_CAPID); }
46  constexpr bool linkError() const { return (flavor()==FLAVOR_HB)? (frame_[i_]&MASK_LE_HB) : (frame_[0]&MASK_LINKERROR); }
47  private:
50  };
51 
52  constexpr void copyContent(const QIE11DataFrame& digi) {
53  for (edm::DataFrame::size_type i=0; i<size() && i<digi.size();i++){
54  Sample sam = digi[i];
55  setSample(i,sam.adc(),sam.tdc(),sam.soi());
56  }
57  }
58 
60  constexpr DetId detid() const { return DetId(m_data.id()); }
70  constexpr int samples() const { return (size()-HEADER_WORDS-FLAG_WORDS)/WORDS_PER_SAMPLE; }
72  constexpr int presamples() const {
73  for (int i=0; i<samples(); i++) {
74  if ((*this)[i].soi()) return i;
75  }
76  return -1;
77  }
79  constexpr int flavor() const { return ((m_data[0]>>OFFSET_FLAVOR)&MASK_FLAVOR); }
81  constexpr bool linkError() const { return m_data[0]&MASK_LINKERROR; }
83  static const int MASK_CAPIDERROR = 0x400;
84  constexpr bool capidError() const { return m_data[0]&MASK_CAPIDERROR; }
86  constexpr bool zsMarkAndPass() const {return (flavor()==1); }
89  if(markAndPass) m_data[0] |= (markAndPass&MASK_FLAVOR)<<OFFSET_FLAVOR;
90  }
92  constexpr inline Sample operator[](edm::DataFrame::size_type i) const { return Sample(m_data,i+HEADER_WORDS); }
93  constexpr void setCapid0(int cap0) {
94  if (flavor()==FLAVOR_HB) {
95  for (int i=0; i<samples(); i++) {
98  }
99  } else {
100  m_data[0]&=0xFCFF; // inversion of the capid0 mask
102  }
103  }
106  int tdc, bool soi=false) {
107  if (isample>=size()) return;
108  if (flavor()==FLAVOR_HB) m_data[isample+1]=(adc&Sample::MASK_ADC)|(soi?(Sample::MASK_SOI):(0))|((tdc&Sample::MASK_TDC_HB)<<Sample::OFFSET_TDC)|(m_data[isample+1]&Sample::MASK_CAPID_KEEP_HB);
109  else m_data[isample+1]=(adc&Sample::MASK_ADC)|(soi?(Sample::MASK_SOI):(0))|((tdc&Sample::MASK_TDC_HE)<<Sample::OFFSET_TDC);
110  }
112  constexpr uint16_t flags() const { return m_data[size()-1]; }
114  constexpr void setFlags(uint16_t v) {
115  m_data[size()-1]=v;
116  }
117 
118  private:
120 
121 };
122 
123 std::ostream& operator<<(std::ostream&, const QIE11DataFrame&);
124 
125 #endif // DATAFORMATS_HCALDIGI_QIE11DATAFRAME_H
DetId detid() const
Get the detector id.
unsigned int id_type
Definition: DataFrame.h:19
void setSample(edm::DataFrame::size_type isample, int adc, int tdc, bool soi=false)
set the sample contents
static const int FLAG_WORDS
Sample(const edm::DataFrame &frame, edm::DataFrame::size_type i)
edm::DataFrame::size_type size() const
more accessors
QIE11DataFrame(edm::DataFrame const &df)
edm::DataFrame::const_iterator end() const
edm::DataFrame::size_type i_
constexpr iterator end()
Definition: DataFrame.h:51
void setCapid0(int cap0)
static const int MASK_CAPID
static const int HEADER_WORDS
int flavor() const
get the flavor of the frame
static const int MASK_ADC
static const int MASK_FLAVOR
void setFlags(uint16_t v)
set the flag word
uint16_t flags() const
get the flag word
static const int WORDS_PER_SAMPLE
static const int MASK_LINKERROR
static const int MASK_CAPIDERROR
was there a capid rotation error?
unsigned int size_type
Definition: DataFrame.h:18
static const int OFFSET_CAPID_HB
constexpr id_type id() const
Definition: DataFrame.h:61
void setZSInfo(bool markAndPass)
set ZS params
bool linkError() const
was there a link error?
bool zsMarkAndPass() const
was this a mark-and-pass ZS event?
Sample operator[](edm::DataFrame::size_type i) const
get the sample
constexpr size_type size() const
Definition: DataFrame.h:64
edm::DataFrame::const_iterator begin() const
edm::DataFrame::id_type id() const
const edm::DataFrame & frame_
static const int OFFSET_FLAVOR
bool capidError() const
edm::DataFrame m_data
void copyContent(const QIE11DataFrame &digi)
Definition: DetId.h:18
static const int OFFSET_TDC
std::ostream & operator<<(std::ostream &, const QIE11DataFrame &)
constexpr iterator begin()
Definition: DataFrame.h:48
static const int MASK_CAPID_KEEP_HB
static const int MASK_TDC_HB
static const int MASK_SOI
int presamples() const
for backward compatibility
static const int MASK_CAPID_INV_HB
data_type * iterator
Definition: DataFrame.h:21
data_type const * const_iterator
Definition: DataFrame.h:22
bool linkError() const
edm::DataFrame::iterator begin()
iterators
static const int OFFSET_CAPID_HE
edm::DataFrame::iterator end()
int samples() const
total number of samples in the digi
static const int MASK_LE_HB
#define constexpr
static const int FLAVOR_HB
static const int MASK_TDC_HE