File indexing completed on 2024-04-06 12:04:18
0001 #ifndef DIGIHCAL_HFDATAFRAME_H
0002 #define DIGIHCAL_HFDATAFRAME_H
0003
0004 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
0005 #include "DataFormats/HcalDetId/interface/HcalElectronicsId.h"
0006 #include "DataFormats/HcalDigi/interface/HcalQIESample.h"
0007 #include <ostream>
0008
0009
0010
0011
0012
0013
0014 class HFDataFrame {
0015 public:
0016 typedef HcalDetId key_type;
0017
0018 constexpr HFDataFrame() : id_(0), size_(0), hcalPresamples_(0) {}
0019 constexpr explicit HFDataFrame(const HcalDetId& id)
0020 : id_(id), size_(0), hcalPresamples_(0) {
0021 }
0022
0023 constexpr HcalDetId const& id() const { return id_; }
0024 constexpr HcalElectronicsId const& elecId() const { return electronicsId_; }
0025
0026
0027 constexpr int size() const { return size_ & 0xF; }
0028
0029 constexpr int presamples() const { return hcalPresamples_ & 0xF; }
0030
0031 constexpr bool zsMarkAndPass() const { return (hcalPresamples_ & 0x10); }
0032
0033 constexpr bool zsUnsuppressed() const { return (hcalPresamples_ & 0x20); }
0034
0035 constexpr uint32_t zsCrossingMask() const { return (hcalPresamples_ & 0x3FF000) >> 12; }
0036
0037
0038 constexpr HcalQIESample const& operator[](int i) const { return data_[i]; }
0039
0040 constexpr HcalQIESample const& sample(int i) const { return data_[i]; }
0041
0042
0043 constexpr int fiberIdleOffset() const {
0044 int val = (hcalPresamples_ & 0xF00) >> 8;
0045 return (val == 0) ? (-1000) : (((val & 0x8) == 0) ? (-(val & 0x7)) : (val & 0x7));
0046 }
0047
0048
0049 constexpr bool validate(int firstSample = 0, int nSamples = 100) const {
0050 int capid = -1;
0051 bool ok = true;
0052 for (int i = 0; ok && i < nSamples && i + firstSample < size_; i++) {
0053 if (data_[i + firstSample].er() || !data_[i + firstSample].dv())
0054 ok = false;
0055 if (i == 0)
0056 capid = data_[i + firstSample].capid();
0057 if (capid != data_[i + firstSample].capid())
0058 ok = false;
0059 capid = (capid + 1) % 4;
0060 }
0061 return ok;
0062 }
0063
0064 constexpr void setSize(int size) {
0065 if (size > MAXSAMPLES)
0066 size_ = MAXSAMPLES;
0067 else if (size <= 0)
0068 size_ = 0;
0069 else
0070 size_ = size;
0071 }
0072 constexpr void setPresamples(int ps) { hcalPresamples_ |= ps & 0xF; }
0073 constexpr void setZSInfo(bool unsuppressed, bool markAndPass, uint32_t crossingMask = 0) {
0074 hcalPresamples_ &= 0x7FC00F0F;
0075 if (markAndPass)
0076 hcalPresamples_ |= 0x10;
0077 if (unsuppressed)
0078 hcalPresamples_ |= 0x20;
0079 hcalPresamples_ |= (crossingMask & 0x3FF) << 12;
0080 }
0081 constexpr void setSample(int i, const HcalQIESample& sam) { data_[i] = sam; }
0082 constexpr void setReadoutIds(const HcalElectronicsId& eid) { electronicsId_ = eid; }
0083 constexpr void setFiberIdleOffset(int offset) {
0084 hcalPresamples_ &= 0x7FFFF0FF;
0085 if (offset >= 7)
0086 hcalPresamples_ |= 0xF00;
0087 else if (offset >= 0)
0088 hcalPresamples_ |= (0x800) | (offset << 8);
0089 else if (offset >= -7)
0090 hcalPresamples_ |= ((-offset) << 8);
0091 else
0092 hcalPresamples_ |= 0x700;
0093 }
0094
0095 static const int MAXSAMPLES = 10;
0096
0097 private:
0098 HcalDetId id_;
0099 HcalElectronicsId electronicsId_;
0100 int size_;
0101 int hcalPresamples_;
0102 HcalQIESample data_[MAXSAMPLES];
0103 };
0104
0105 std::ostream& operator<<(std::ostream&, const HFDataFrame&);
0106
0107 #endif