File indexing completed on 2024-04-06 12:04:18
0001 #ifndef DATAFORMATS_HCALDIGI_QIE11DATAFRAME_H
0002 #define DATAFORMATS_HCALDIGI_QIE11DATAFRAME_H
0003
0004 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
0005 #include "DataFormats/Common/interface/DataFrame.h"
0006 #include <ostream>
0007
0008
0009
0010
0011 class QIE11DataFrame {
0012 public:
0013 static const int WORDS_PER_SAMPLE = 1;
0014 static const int HEADER_WORDS = 1;
0015 static const int FLAG_WORDS = 1;
0016
0017 static const int OFFSET_FLAVOR = 12;
0018 static const int MASK_FLAVOR = 0x7;
0019 static const int FLAVOR_HB = 3;
0020 static const int MASK_LINKERROR = 0x800;
0021
0022 constexpr QIE11DataFrame() {}
0023 constexpr QIE11DataFrame(edm::DataFrame const& df) : m_data(df) {}
0024
0025 class Sample {
0026 public:
0027 constexpr Sample(const edm::DataFrame& frame, edm::DataFrame::size_type i) : frame_(frame), i_(i) {}
0028 static const int MASK_ADC = 0xFF;
0029 static const int MASK_TDC_HE = 0x3F;
0030 static const int MASK_TDC_HB = 0x3;
0031 static const int OFFSET_TDC = 8;
0032 static const int MASK_SOI = 0x4000;
0033 static const int MASK_LE_HB = 0x2000;
0034 static const int MASK_CAPID = 0x3;
0035 static const int MASK_CAPID_INV_HB = 0xF3FF;
0036 static const int MASK_CAPID_KEEP_HB = 0x0C00;
0037 static const int OFFSET_CAPID_HE = 8;
0038 static const int OFFSET_CAPID_HB = 10;
0039 constexpr int flavor() const { return ((frame_[0] >> OFFSET_FLAVOR) & MASK_FLAVOR); }
0040 constexpr int adc() const { return frame_[i_] & MASK_ADC; }
0041 constexpr int tdc() const {
0042 return (frame_[i_] >> OFFSET_TDC) & ((flavor() == FLAVOR_HB) ? (MASK_TDC_HB) : (MASK_TDC_HE));
0043 }
0044 constexpr bool soi() const { return frame_[i_] & MASK_SOI; }
0045 constexpr int capid() const {
0046 return (flavor() == FLAVOR_HB)
0047 ? ((frame_[i_] >> OFFSET_CAPID_HB) & MASK_CAPID)
0048 : ((((frame_[0] >> OFFSET_CAPID_HE) & MASK_CAPID) + i_ - HEADER_WORDS) & MASK_CAPID);
0049 }
0050 constexpr bool linkError() const {
0051 return (flavor() == FLAVOR_HB) ? (frame_[i_] & MASK_LE_HB) : (frame_[0] & MASK_LINKERROR);
0052 }
0053
0054 private:
0055 const edm::DataFrame& frame_;
0056 edm::DataFrame::size_type i_;
0057 };
0058
0059 constexpr void copyContent(const QIE11DataFrame& digi) {
0060 for (edm::DataFrame::size_type i = 0; i < size() && i < digi.size(); i++) {
0061 Sample sam = digi[i];
0062 setSample(i, sam.adc(), sam.tdc(), sam.soi());
0063 }
0064 }
0065
0066
0067 constexpr DetId detid() const { return DetId(m_data.id()); }
0068 constexpr edm::DataFrame::id_type id() const { return m_data.id(); }
0069
0070 constexpr edm::DataFrame::size_type size() const { return m_data.size(); }
0071
0072 constexpr edm::DataFrame::iterator begin() { return m_data.begin(); }
0073 constexpr edm::DataFrame::iterator end() { return m_data.end(); }
0074 constexpr edm::DataFrame::const_iterator begin() const { return m_data.begin(); }
0075 constexpr edm::DataFrame::const_iterator end() const { return m_data.end(); }
0076
0077 constexpr int samples() const { return (size() - HEADER_WORDS - FLAG_WORDS) / WORDS_PER_SAMPLE; }
0078
0079 constexpr int presamples() const {
0080 for (int i = 0; i < samples(); i++) {
0081 if ((*this)[i].soi())
0082 return i;
0083 }
0084 return -1;
0085 }
0086
0087 constexpr int flavor() const { return ((m_data[0] >> OFFSET_FLAVOR) & MASK_FLAVOR); }
0088
0089 constexpr bool linkError() const { return m_data[0] & MASK_LINKERROR; }
0090
0091 static const int MASK_CAPIDERROR = 0x400;
0092 constexpr bool capidError() const { return m_data[0] & MASK_CAPIDERROR; }
0093
0094 constexpr bool zsMarkAndPass() const { return (flavor() == 1); }
0095
0096 constexpr void setZSInfo(bool markAndPass) {
0097 if (markAndPass)
0098 m_data[0] |= (markAndPass & MASK_FLAVOR) << OFFSET_FLAVOR;
0099 }
0100
0101 constexpr inline Sample operator[](edm::DataFrame::size_type i) const { return Sample(m_data, i + HEADER_WORDS); }
0102
0103
0104 constexpr void setFlavor(int flavor) {
0105 m_data[0] &= 0x9FFF;
0106 m_data[0] |= ((flavor & MASK_FLAVOR) << OFFSET_FLAVOR);
0107 }
0108
0109 constexpr void setCapid0(int cap0) {
0110 if (flavor() == FLAVOR_HB) {
0111 for (int i = 0; i < samples(); i++) {
0112 m_data[i + 1] &= Sample::MASK_CAPID_INV_HB;
0113 m_data[i + 1] |= ((cap0 + i) & Sample::MASK_CAPID) << Sample::OFFSET_CAPID_HB;
0114 }
0115 } else {
0116 m_data[0] &= 0xFCFF;
0117 m_data[0] |= ((cap0 & Sample::MASK_CAPID) << Sample::OFFSET_CAPID_HE);
0118 }
0119 }
0120
0121 constexpr void setSample(edm::DataFrame::size_type isample, int adc, int tdc, bool soi = false) {
0122 if (isample >= size())
0123 return;
0124 if (flavor() == FLAVOR_HB)
0125 m_data[isample + 1] = (adc & Sample::MASK_ADC) | (soi ? (Sample::MASK_SOI) : (0)) |
0126 ((tdc & Sample::MASK_TDC_HB) << Sample::OFFSET_TDC) |
0127 (m_data[isample + 1] & Sample::MASK_CAPID_KEEP_HB);
0128 else
0129 m_data[isample + 1] = (adc & Sample::MASK_ADC) | (soi ? (Sample::MASK_SOI) : (0)) |
0130 ((tdc & Sample::MASK_TDC_HE) << Sample::OFFSET_TDC);
0131 }
0132
0133 constexpr uint16_t flags() const { return m_data[size() - 1]; }
0134
0135 constexpr void setFlags(uint16_t v) { m_data[size() - 1] = v; }
0136
0137 private:
0138 edm::DataFrame m_data;
0139 };
0140
0141 std::ostream& operator<<(std::ostream&, const QIE11DataFrame&);
0142
0143 #endif