CMS 3D CMS Logo

SiPixelGainCalibrationForHLTSoAESProducer.cc
Go to the documentation of this file.
14 
20 
21 #include <memory>
22 
24 
26  public:
28  std::unique_ptr<SiPixelGainCalibrationForHLTHost> produce(const SiPixelGainCalibrationForHLTSoARcd& iRecord);
29 
30  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
31 
32  private:
35  };
36 
38  : ESProducer(iConfig) {
39  auto cc = setWhatProduced(this);
40  gainsToken_ = cc.consumes();
41  geometryToken_ = cc.consumes();
42  }
43 
46  descriptions.addWithDefaultLabel(desc);
47  }
48 
49  std::unique_ptr<SiPixelGainCalibrationForHLTHost> SiPixelGainCalibrationForHLTSoAESProducer::produce(
50  const SiPixelGainCalibrationForHLTSoARcd& iRecord) {
51  auto const& gains = iRecord.get(gainsToken_);
52  auto const& geom = iRecord.get(geometryToken_);
53 
54  auto product = std::make_unique<SiPixelGainCalibrationForHLTHost>(gains.data().size(), cms::alpakatools::host());
55 
56  // bizzarre logic (looking for fist strip-det) don't ask
57  auto const& dus = geom.detUnits();
58  unsigned int n_detectors = dus.size();
59  for (unsigned int i = 1; i < 7; ++i) {
60  const auto offset = geom.offsetDU(GeomDetEnumerators::tkDetEnum[i]);
61  if (offset != dus.size() && dus[offset]->type().isTrackerStrip()) {
62  if (n_detectors > offset)
63  n_detectors = offset;
64  }
65  }
66 
67  LogDebug("SiPixelGainCalibrationForHLTSoA")
68  << "caching calibs for " << n_detectors << " pixel detectors of size " << gains.data().size() << '\n'
69  << "sizes " << sizeof(char) << ' ' << sizeof(uint8_t) << ' ' << sizeof(siPixelGainsSoA::DecodingStructure);
70 
71  for (size_t i = 0; i < gains.data().size(); i = i + 2) {
72  product->view().v_pedestals()[i / 2].gain = gains.data()[i];
73  product->view().v_pedestals()[i / 2].ped = gains.data()[i + 1];
74  }
75 
76  //std::copy here
77  // do not read back from the (possibly write-combined) memory buffer
78  auto minPed = gains.getPedLow();
79  auto maxPed = gains.getPedHigh();
80  auto minGain = gains.getGainLow();
81  auto maxGain = gains.getGainHigh();
82  auto nBinsToUseForEncoding = 253;
83 
84  // we will simplify later (not everything is needed....)
85  product->view().minPed() = minPed;
86  product->view().maxPed() = maxPed;
87  product->view().minGain() = minGain;
88  product->view().maxGain() = maxGain;
89 
90  product->view().numberOfRowsAveragedOver() = 80;
91  product->view().nBinsToUseForEncoding() = nBinsToUseForEncoding;
92  product->view().deadFlag() = 255;
93  product->view().noisyFlag() = 254;
94 
95  product->view().pedPrecision() = static_cast<float>(maxPed - minPed) / nBinsToUseForEncoding;
96  product->view().gainPrecision() = static_cast<float>(maxGain - minGain) / nBinsToUseForEncoding;
97 
98  LogDebug("SiPixelGainCalibrationForHLTSoA")
99  << "precisions g " << product->view().pedPrecision() << ' ' << product->view().gainPrecision();
100 
101  // fill the index map
102  auto const& ind = gains.getIndexes();
103  LogDebug("SiPixelGainCalibrationForHLTSoA") << ind.size() << " " << n_detectors;
104 
105  for (auto i = 0U; i < n_detectors; ++i) {
106  auto p = std::lower_bound(
107  ind.begin(), ind.end(), dus[i]->geographicalId().rawId(), SiPixelGainCalibrationForHLT::StrictWeakOrdering());
108  assert(p != ind.end() && p->detid == dus[i]->geographicalId());
109  assert(p->iend <= gains.data().size());
110  assert(p->iend >= p->ibegin);
111  assert(0 == p->ibegin % 2);
112  assert(0 == p->iend % 2);
113  assert(p->ibegin != p->iend);
114  assert(p->ncols > 0);
115 
116  product->view().modStarts()[i] = p->ibegin;
117  product->view().modEnds()[i] = p->iend;
118  product->view().modCols()[i] = p->ncols;
119 
120  if (ind[i].detid != dus[i]->geographicalId())
121  LogDebug("SiPixelGainCalibrationForHLTSoA") << ind[i].detid << "!=" << dus[i]->geographicalId();
122  }
123 
124  return product;
125  }
126 
127 } // namespace ALPAKA_ACCELERATOR_NAMESPACE
128 DEFINE_FWK_EVENTSETUP_ALPAKA_MODULE(SiPixelGainCalibrationForHLTSoAESProducer);
auto setWhatProduced(T *iThis, const es::Label &iLabel={})
Definition: ESProducer.h:166
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
assert(be >=bs)
std::unique_ptr< SiPixelGainCalibrationForHLTHost > produce(const SiPixelGainCalibrationForHLTSoARcd &iRecord)
constexpr float gains[NGAINS]
Definition: EcalConstants.h:20
alpaka::DevCpu const & host()
Definition: host.h:14
#define DEFINE_FWK_EVENTSETUP_ALPAKA_MODULE(type)
Definition: ModuleFactory.h:17
ProductT const & get(ESGetToken< ProductT, DepRecordT > const &iToken) const
constexpr SubDetector tkDetEnum[8]
edm::ESGetToken< SiPixelGainCalibrationForHLT, SiPixelGainCalibrationForHLTRcd > gainsToken_
#define LogDebug(id)