CMS 3D CMS Logo

HGCalMappingCellIndexer.h
Go to the documentation of this file.
1 #ifndef CondFormats_HGCalObjects_interface_HGCalMappingCellParameterIndex_h
2 #define CondFormats_HGCalObjects_interface_HGCalMappingCellParameterIndex_h
3 
4 #include <iostream>
5 #include <cstdint>
6 #include <vector>
7 #include <numeric>
11 
16 public:
17  // typedef HGCalDenseIndexerBase WaferDenseIndexerBase;
18 
19  HGCalMappingCellIndexer() = default;
20 
24  void processNewCell(std::string typecode, uint16_t chip, uint16_t half) {
25  //assign index to this typecode and resize the max e-Rx vector
26  if (typeCodeIndexer_.count(typecode) == 0) {
27  typeCodeIndexer_[typecode] = typeCodeIndexer_.size();
28  maxErx_.resize(typeCodeIndexer_.size(), 0);
29  }
30 
31  size_t idx = typeCodeIndexer_[typecode];
32  uint16_t erx = chip * 2 + half + 1; //use the number not the index here
33  maxErx_[idx] = std::max(maxErx_[idx], erx);
34  }
35 
39  void update() {
40  uint32_t n = typeCodeIndexer_.size();
41  offsets_ = std::vector<uint32_t>(n, 0);
42  di_ = std::vector<HGCalDenseIndexerBase>(n, HGCalDenseIndexerBase(2));
43  for (uint32_t idx = 0; idx < n; idx++) {
44  uint16_t nerx = maxErx_[idx];
45  di_[idx].updateRanges({{nerx, maxChPerErx_}});
46  if (idx < n - 1)
47  offsets_[idx + 1] = di_[idx].getMaxIndex();
48  }
49 
50  //accumulate the offsets in the array
51  std::partial_sum(offsets_.begin(), offsets_.end(), offsets_.begin());
52  }
53 
57  size_t getEnumFromTypecode(std::string typecode) const {
58  auto it = typeCodeIndexer_.find(typecode);
59  if (it == typeCodeIndexer_.end())
60  throw cms::Exception("ValueError") << " unable to find typecode=" << typecode << " in cell indexer";
61  return it->second;
62  }
63 
68  for (const auto& it : typeCodeIndexer_)
69  if (it.second == idx)
70  return it.first;
71  throw cms::Exception("ValueError") << " unable to find typecode corresponding to idx=" << idx;
72  }
73 
78  return getDenseIndexerFor(getEnumFromTypecode(typecode));
79  }
80 
85  if (idx >= di_.size())
86  throw cms::Exception("ValueError") << " index requested for cell dense indexer (i=" << idx
87  << ") is larger than allocated";
88  return di_[idx];
89  }
90 
94  uint32_t denseIndex(std::string typecode, uint32_t chip, uint32_t half, uint32_t seq) const {
95  return denseIndex(getEnumFromTypecode(typecode), chip, half, seq);
96  }
97  uint32_t denseIndex(std::string typecode, uint32_t erx, uint32_t seq) const {
98  return denseIndex(getEnumFromTypecode(typecode), erx, seq);
99  }
100  uint32_t denseIndex(size_t idx, uint32_t chip, uint32_t half, uint32_t seq) const {
101  uint16_t erx = chip * maxHalfPerROC_ + half;
102  return denseIndex(idx, erx, seq);
103  }
104  uint32_t denseIndex(size_t idx, uint32_t erx, uint32_t seq) const {
105  return di_[idx].denseIndex({{erx, seq}}) + offsets_[idx];
106  }
107 
111  uint32_t elecIdFromIndex(uint32_t rtn, std::string typecode) const {
112  return elecIdFromIndex(rtn, getEnumFromTypecode(typecode));
113  }
114  uint32_t elecIdFromIndex(uint32_t rtn, size_t idx) const {
115  if (idx >= di_.size())
116  throw cms::Exception("ValueError") << " index requested for cell dense indexer (i=" << idx
117  << ") is larger than allocated";
118  rtn -= offsets_[idx];
119  auto rtn_codes = di_[idx].unpackDenseIndex(rtn);
120  return HGCalElectronicsId(false, 0, 0, 0, rtn_codes[0], rtn_codes[1]).raw();
121  }
122 
126  uint32_t maxDenseIndex() const {
127  size_t i = maxErx_.size();
128  if (i == 0)
129  return 0;
130  return offsets_.back() + maxErx_.back() * maxChPerErx_;
131  }
132 
136  size_t getNWordsExpectedFor(std::string typecode) const {
137  auto it = getEnumFromTypecode(typecode);
138  return getNWordsExpectedFor(it);
139  }
140  size_t getNWordsExpectedFor(size_t typecodeidx) const { return maxErx_[typecodeidx] * maxChPerErx_; }
141 
145  size_t getNErxExpectedFor(std::string typecode) const {
146  auto it = getEnumFromTypecode(typecode);
147  return getNErxExpectedFor(it);
148  }
149  size_t getNErxExpectedFor(size_t typecodeidx) const { return maxErx_[typecodeidx]; }
150 
151  constexpr static char maxHalfPerROC_ = 2;
152  constexpr static uint16_t maxChPerErx_ = 37; //36 channels + 1 calib
153 
154  std::map<std::string, size_t> typeCodeIndexer_;
155  std::vector<uint16_t> maxErx_;
156  std::vector<uint32_t> offsets_;
157  std::vector<HGCalDenseIndexerBase> di_;
158 
160 
162 };
163 
164 #endif
static constexpr uint16_t maxChPerErx_
uint32_t raw() const
uint32_t elecIdFromIndex(uint32_t rtn, size_t idx) const
uint32_t denseIndex(size_t idx, uint32_t chip, uint32_t half, uint32_t seq) const
uint32_t maxDenseIndex() const
returns the max. dense index expected
void processNewCell(std::string typecode, uint16_t chip, uint16_t half)
size_t getNErxExpectedFor(size_t typecodeidx) const
uint32_t denseIndex(std::string typecode, uint32_t chip, uint32_t half, uint32_t seq) const
builders for the dense index
size_t getEnumFromTypecode(std::string typecode) const
gets index given typecode string
HGCalDenseIndexerBase getDenseIndexerFor(size_t idx) const
returns the dense indexer for a given internal index
HGCalMappingCellIndexer()=default
this is a simple class that takes care of building a dense index for a set of categories the maximum ...
wrapper for a 32b data word identifying a readout channel in the raw data The format is the following...
uint32_t denseIndex(size_t idx, uint32_t erx, uint32_t seq) const
HGCalDenseIndexerBase getDenseIndexFor(std::string typecode) const
returns the dense indexer for a typecode
void update()
process the current list of type codes handled and updates the dense indexers
utility class to assign dense readout cell indexing
uint32_t denseIndex(std::string typecode, uint32_t erx, uint32_t seq) const
#define COND_SERIALIZABLE
Definition: Serializable.h:39
std::map< std::string, size_t > typeCodeIndexer_
size_t getNWordsExpectedFor(std::string typecode) const
gets the number of words for a given typecode
uint32_t elecIdFromIndex(uint32_t rtn, std::string typecode) const
decodes the dense index code
static constexpr char maxHalfPerROC_
std::string getTypecodeFromEnum(size_t idx) const
checks if there is a typecode corresponding to an index
std::vector< uint32_t > offsets_
std::vector< uint16_t > maxErx_
std::vector< HGCalDenseIndexerBase > di_
size_t getNWordsExpectedFor(size_t typecodeidx) const
size_t getNErxExpectedFor(std::string typecode) const
gets the number of e-Rx for a given typecode