CMS 3D CMS Logo

SiPixelGainCalibrationForHLTGPU.cc
Go to the documentation of this file.
1 #include <cuda.h>
2 
10 
12  const TrackerGeometry& geom)
13  : gains_(&gains) {
14  // bizzarre logic (looking for fist strip-det) don't ask
15  auto const& dus = geom.detUnits();
16  unsigned int n_detectors = dus.size();
17  for (unsigned int i = 1; i < 7; ++i) {
18  const auto offset = geom.offsetDU(GeomDetEnumerators::tkDetEnum[i]);
19  if (offset != dus.size() && dus[offset]->type().isTrackerStrip()) {
20  if (n_detectors > offset)
21  n_detectors = offset;
22  }
23  }
24 
25  LogDebug("SiPixelGainCalibrationForHLTGPU")
26  << "caching calibs for " << n_detectors << " pixel detectors of size " << gains.data().size() << '\n'
27  << "sizes " << sizeof(char) << ' ' << sizeof(uint8_t) << ' ' << sizeof(SiPixelGainForHLTonGPU::DecodingStructure);
28 
29  cudaCheck(cudaMallocHost((void**)&gainForHLTonHost_, sizeof(SiPixelGainForHLTonGPU)));
31  (SiPixelGainForHLTonGPU_DecodingStructure*)this->gains_->data().data(); // so it can be used on CPU as well...
32 
33  // do not read back from the (possibly write-combined) memory buffer
34  auto minPed = gains.getPedLow();
35  auto maxPed = gains.getPedHigh();
36  auto minGain = gains.getGainLow();
37  auto maxGain = gains.getGainHigh();
38  auto nBinsToUseForEncoding = 253;
39 
40  // we will simplify later (not everything is needed....)
41  gainForHLTonHost_->minPed_ = minPed;
42  gainForHLTonHost_->maxPed_ = maxPed;
43  gainForHLTonHost_->minGain_ = minGain;
44  gainForHLTonHost_->maxGain_ = maxGain;
45 
47  gainForHLTonHost_->nBinsToUseForEncoding_ = nBinsToUseForEncoding;
50 
51  gainForHLTonHost_->pedPrecision_ = static_cast<float>(maxPed - minPed) / nBinsToUseForEncoding;
52  gainForHLTonHost_->gainPrecision_ = static_cast<float>(maxGain - minGain) / nBinsToUseForEncoding;
53 
54  LogDebug("SiPixelGainCalibrationForHLTGPU")
55  << "precisions g " << gainForHLTonHost_->pedPrecision_ << ' ' << gainForHLTonHost_->gainPrecision_;
56 
57  // fill the index map
58  auto const& ind = gains.getIndexes();
59  LogDebug("SiPixelGainCalibrationForHLTGPU") << ind.size() << " " << n_detectors;
60 
61  for (auto i = 0U; i < n_detectors; ++i) {
62  auto p = std::lower_bound(
63  ind.begin(), ind.end(), dus[i]->geographicalId().rawId(), SiPixelGainCalibrationForHLT::StrictWeakOrdering());
64  assert(p != ind.end() && p->detid == dus[i]->geographicalId());
65  assert(p->iend <= gains.data().size());
66  assert(p->iend >= p->ibegin);
67  assert(0 == p->ibegin % 2);
68  assert(0 == p->iend % 2);
69  assert(p->ibegin != p->iend);
70  assert(p->ncols > 0);
71  gainForHLTonHost_->rangeAndCols_[i] = std::make_pair(SiPixelGainForHLTonGPU::Range(p->ibegin, p->iend), p->ncols);
72  if (ind[i].detid != dus[i]->geographicalId())
73  LogDebug("SiPixelGainCalibrationForHLTGPU") << ind[i].detid << "!=" << dus[i]->geographicalId();
74  }
75 }
76 
78 
80  cudaCheck(cudaFree(gainForHLTonGPU));
81  cudaCheck(cudaFree(gainDataOnGPU));
82 }
83 
85  const auto& data = gpuData_.dataForCurrentDeviceAsync(cudaStream, [this](GPUData& data, cudaStream_t stream) {
86  cudaCheck(cudaMalloc((void**)&data.gainForHLTonGPU, sizeof(SiPixelGainForHLTonGPU)));
87  cudaCheck(cudaMalloc((void**)&data.gainDataOnGPU, this->gains_->data().size()));
88  // gains.data().data() is used also for non-GPU code, we cannot allocate it on aligned and write-combined memory
89  cudaCheck(cudaMemcpyAsync(
90  data.gainDataOnGPU, this->gains_->data().data(), this->gains_->data().size(), cudaMemcpyDefault, stream));
91 
92  cudaCheck(cudaMemcpyAsync(
93  data.gainForHLTonGPU, this->gainForHLTonHost_, sizeof(SiPixelGainForHLTonGPU), cudaMemcpyDefault, stream));
94  cudaCheck(cudaMemcpyAsync(&(data.gainForHLTonGPU->v_pedestals_),
95  &(data.gainDataOnGPU),
97  cudaMemcpyDefault,
98  stream));
99  });
100  return data.gainForHLTonGPU;
101 }
const SiPixelGainCalibrationForHLT * gains_
const SiPixelGainForHLTonGPU * getGPUProductAsync(cudaStream_t cudaStream) const
DecodingStructure * v_pedestals_
std::pair< uint32_t, uint32_t > Range
uint32_t T const *__restrict__ uint32_t const *__restrict__ int32_t int Histo::index_type cudaStream_t stream
assert(be >=bs)
SiPixelGainForHLTonGPU_DecodingStructure DecodingStructure
std::pair< Range, int > rangeAndCols_[phase1PixelTopology::numberOfModules]
cms::cuda::ESProduct< GPUData > gpuData_
constexpr float gains[NGAINS]
Definition: EcalConstants.h:20
std::vector< char > const & data() const
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:80
SiPixelGainForHLTonGPU_DecodingStructure * gainDataOnGPU
SiPixelGainCalibrationForHLTGPU(const SiPixelGainCalibrationForHLT &gains, const TrackerGeometry &geom)
#define cudaCheck(ARG,...)
Definition: cudaCheck.h:69
constexpr SubDetector tkDetEnum[8]
#define LogDebug(id)