CMS 3D CMS Logo

HcalRecoParamsWithPulseShapesGPU.cc
Go to the documentation of this file.
2 
6 
9 
10 #include <unordered_map>
11 
12 // FIXME: add proper getters to conditions
14  : totalChannels_{recoParams.getAllContainers()[0].second.size() + recoParams.getAllContainers()[1].second.size()},
15  param1_(totalChannels_),
16  param2_(totalChannels_),
17  ids_(totalChannels_) {
18 #ifdef HCAL_MAHI_CPUDEBUG
19  printf("hello from a reco params with pulse shapes\n");
20 #endif
21 
22  auto const containers = recoParams.getAllContainers();
23 
24  HcalPulseShapes pulseShapes;
25  std::unordered_map<unsigned int, uint32_t> idCache;
26 
27  // fill in eb
28  auto const& barrelValues = containers[0].second;
29  for (uint64_t i = 0; i < barrelValues.size(); ++i) {
30  param1_[i] = barrelValues[i].param1();
31  param2_[i] = barrelValues[i].param2();
32 
33  auto const pulseShapeId = barrelValues[i].pulseShapeID();
34  // FIXME: 0 throws upon look up to HcalPulseShapes
35  // although comments state that 0 is reserved,
36  // HcalPulseShapes::getShape throws on 0!
37  if (pulseShapeId == 0) {
38  ids_[i] = 0;
39  continue;
40  }
41  if (auto const iter = idCache.find(pulseShapeId); iter == idCache.end()) {
42  // new guy
43  auto const newId = idCache.size();
44  idCache[pulseShapeId] = newId;
45  // this will be the id
46  ids_[i] = newId;
47 
48  // resize value arrays
49  acc25nsVec_.resize(acc25nsVec_.size() + hcal::constants::maxPSshapeBin);
50  diff25nsItvlVec_.resize(diff25nsItvlVec_.size() + hcal::constants::maxPSshapeBin);
51  accVarLenIdxMinusOneVec_.resize(accVarLenIdxMinusOneVec_.size() + hcal::constants::nsPerBX);
52  diffVarItvlIdxMinusOneVec_.resize(diffVarItvlIdxMinusOneVec_.size() + hcal::constants::nsPerBX);
53  accVarLenIdxZEROVec_.resize(accVarLenIdxZEROVec_.size() + hcal::constants::nsPerBX);
54  diffVarItvlIdxZEROVec_.resize(diffVarItvlIdxZEROVec_.size() + hcal::constants::nsPerBX);
55 
56  // precompute and get values from the functor
57  auto const& pulseShape = pulseShapes.getShape(pulseShapeId);
58  FitterFuncs::PulseShapeFunctor functor{pulseShape, false, false, false, 1, 0, 0, hcal::constants::maxSamples};
59  auto const offset256 = newId * hcal::constants::maxPSshapeBin;
60  auto const offset25 = newId * hcal::constants::nsPerBX;
61  auto const numShapes = newId;
62  for (int i = 0; i < hcal::constants::maxPSshapeBin; i++) {
63  acc25nsVec_[offset256 * numShapes + i] = functor.acc25nsVec()[i];
64  diff25nsItvlVec_[offset256 * numShapes + i] = functor.diff25nsItvlVec()[i];
65  }
66 
67  for (int i = 0; i < hcal::constants::nsPerBX; i++) {
68  accVarLenIdxMinusOneVec_[offset25 * numShapes + i] = functor.accVarLenIdxMinusOneVec()[i];
69  diffVarItvlIdxMinusOneVec_[offset25 * numShapes + i] = functor.diffVarItvlIdxMinusOneVec()[i];
70  accVarLenIdxZEROVec_[offset25 * numShapes + i] = functor.accVarLenIdxZEROVec()[i];
71  diffVarItvlIdxZEROVec_[offset25 * numShapes + i] = functor.diffVarItvlIdxZEROVec()[i];
72  }
73  } else {
74  // already recorded this pulse shape, just set id
75  ids_[i] = iter->second;
76  }
77 #ifdef HCAL_MAHI_CPUDEBUG
78  if (barrelValues[i].rawId() == DETID_TO_DEBUG) {
79  printf("recoShapeId = %u myid = %u\n", pulseShapeId, ids_[i]);
80  }
81 #endif
82  }
83 
84  // fill in ee
85  auto const& endcapValues = containers[1].second;
86  auto const offset = barrelValues.size();
87  for (uint64_t i = 0; i < endcapValues.size(); ++i) {
88  param1_[i + offset] = endcapValues[i].param1();
89  param2_[i + offset] = endcapValues[i].param2();
90 
91  auto const pulseShapeId = endcapValues[i].pulseShapeID();
92  // FIXME: 0 throws upon look up to HcalPulseShapes
93  // although comments state that 0 is reserved,
94  // HcalPulseShapes::getShape throws on 0!
95  if (pulseShapeId == 0) {
96  ids_[i + offset] = 0;
97  continue;
98  }
99  if (auto const iter = idCache.find(pulseShapeId); iter == idCache.end()) {
100  // new guy
101  auto const newId = idCache.size();
102  idCache[pulseShapeId] = newId;
103  // this will be the id
104  ids_[i + offset] = newId;
105 
106  // resize value arrays
107  acc25nsVec_.resize(acc25nsVec_.size() + hcal::constants::maxPSshapeBin);
108  diff25nsItvlVec_.resize(diff25nsItvlVec_.size() + hcal::constants::maxPSshapeBin);
109  accVarLenIdxMinusOneVec_.resize(accVarLenIdxMinusOneVec_.size() + hcal::constants::nsPerBX);
110  diffVarItvlIdxMinusOneVec_.resize(diffVarItvlIdxMinusOneVec_.size() + hcal::constants::nsPerBX);
111  accVarLenIdxZEROVec_.resize(accVarLenIdxZEROVec_.size() + hcal::constants::nsPerBX);
112  diffVarItvlIdxZEROVec_.resize(diffVarItvlIdxZEROVec_.size() + hcal::constants::nsPerBX);
113 
114  // precompute and get values from the functor
115  auto const& pulseShape = pulseShapes.getShape(pulseShapeId);
116  FitterFuncs::PulseShapeFunctor functor{pulseShape, false, false, false, 1, 0, 0, hcal::constants::maxSamples};
117  auto const offset256 = newId * hcal::constants::maxPSshapeBin;
118  auto const offset25 = newId * hcal::constants::nsPerBX;
119  auto const numShapes = newId;
120  for (int i = 0; i < hcal::constants::maxPSshapeBin; i++) {
121  acc25nsVec_[offset256 * numShapes + i] = functor.acc25nsVec()[i];
122  diff25nsItvlVec_[offset256 * numShapes + i] = functor.diff25nsItvlVec()[i];
123  }
124 
125  for (int i = 0; i < hcal::constants::nsPerBX; i++) {
126  accVarLenIdxMinusOneVec_[offset25 * numShapes + i] = functor.accVarLenIdxMinusOneVec()[i];
127  diffVarItvlIdxMinusOneVec_[offset25 * numShapes + i] = functor.diffVarItvlIdxMinusOneVec()[i];
128  accVarLenIdxZEROVec_[offset25 * numShapes + i] = functor.accVarLenIdxZEROVec()[i];
129  diffVarItvlIdxZEROVec_[offset25 * numShapes + i] = functor.diffVarItvlIdxZEROVec()[i];
130  }
131  } else {
132  // already recorded this pulse shape, just set id
133  ids_[i + offset] = iter->second;
134  }
135  }
136 
137 #ifdef HCAL_MAHI_CPUDEBUG
138  for (auto const& p : idCache)
139  printf("recoPulseShapeId = %u id = %u\n", p.first, p.second);
140 #endif
141 }
142 
144  // deallocation
145  cudaCheck(cudaFree(param1));
146  cudaCheck(cudaFree(param2));
147  cudaCheck(cudaFree(ids));
148  cudaCheck(cudaFree(acc25nsVec));
149  cudaCheck(cudaFree(diff25nsItvlVec));
152  cudaCheck(cudaFree(accVarLenIdxZEROVec));
153  cudaCheck(cudaFree(diffVarItvlIdxZEROVec));
154 }
155 
157  cudaStream_t cudaStream) const {
158  auto const& product = product_.dataForCurrentDeviceAsync(
159  cudaStream, [this](HcalRecoParamsWithPulseShapesGPU::Product& product, cudaStream_t cudaStream) {
160  // malloc
161  cudaCheck(cudaMalloc((void**)&product.param1, this->param1_.size() * sizeof(uint32_t)));
162  cudaCheck(cudaMalloc((void**)&product.param2, this->param2_.size() * sizeof(uint32_t)));
163  cudaCheck(cudaMalloc((void**)&product.ids, this->ids_.size() * sizeof(uint32_t)));
164  cudaCheck(cudaMalloc((void**)&product.acc25nsVec, this->acc25nsVec_.size() * sizeof(float)));
165  cudaCheck(cudaMalloc((void**)&product.diff25nsItvlVec, this->diff25nsItvlVec_.size() * sizeof(float)));
166  cudaCheck(cudaMalloc((void**)&product.accVarLenIdxMinusOneVec,
167  this->accVarLenIdxMinusOneVec_.size() * sizeof(float)));
168  cudaCheck(cudaMalloc((void**)&product.diffVarItvlIdxMinusOneVec,
169  this->diffVarItvlIdxMinusOneVec_.size() * sizeof(float)));
170  cudaCheck(cudaMalloc((void**)&product.accVarLenIdxZEROVec, this->accVarLenIdxZEROVec_.size() * sizeof(float)));
171  cudaCheck(
172  cudaMalloc((void**)&product.diffVarItvlIdxZEROVec, this->diffVarItvlIdxZEROVec_.size() * sizeof(float)));
173 
174  // transfer
175  cudaCheck(cudaMemcpyAsync(product.param1,
176  this->param1_.data(),
177  this->param1_.size() * sizeof(uint32_t),
178  cudaMemcpyHostToDevice,
179  cudaStream));
180  cudaCheck(cudaMemcpyAsync(product.param2,
181  this->param2_.data(),
182  this->param2_.size() * sizeof(uint32_t),
183  cudaMemcpyHostToDevice,
184  cudaStream));
185  cudaCheck(cudaMemcpyAsync(
186  product.ids, this->ids_.data(), this->ids_.size() * sizeof(uint32_t), cudaMemcpyHostToDevice, cudaStream));
187  cudaCheck(cudaMemcpyAsync(product.acc25nsVec,
188  this->acc25nsVec_.data(),
189  this->acc25nsVec_.size() * sizeof(float),
190  cudaMemcpyHostToDevice,
191  cudaStream));
192  cudaCheck(cudaMemcpyAsync(product.diff25nsItvlVec,
193  this->diff25nsItvlVec_.data(),
194  this->diff25nsItvlVec_.size() * sizeof(float),
195  cudaMemcpyHostToDevice,
196  cudaStream));
197  cudaCheck(cudaMemcpyAsync(product.accVarLenIdxMinusOneVec,
198  this->accVarLenIdxMinusOneVec_.data(),
199  this->accVarLenIdxMinusOneVec_.size() * sizeof(float),
200  cudaMemcpyHostToDevice,
201  cudaStream));
202  cudaCheck(cudaMemcpyAsync(product.diffVarItvlIdxMinusOneVec,
203  this->diffVarItvlIdxMinusOneVec_.data(),
204  this->diffVarItvlIdxMinusOneVec_.size() * sizeof(float),
205  cudaMemcpyHostToDevice,
206  cudaStream));
207  cudaCheck(cudaMemcpyAsync(product.accVarLenIdxZEROVec,
208  this->accVarLenIdxZEROVec_.data(),
209  this->accVarLenIdxZEROVec_.size() * sizeof(float),
210  cudaMemcpyHostToDevice,
211  cudaStream));
212  cudaCheck(cudaMemcpyAsync(product.diffVarItvlIdxZEROVec,
213  this->diffVarItvlIdxZEROVec_.data(),
214  this->diffVarItvlIdxZEROVec_.size() * sizeof(float),
215  cudaMemcpyHostToDevice,
216  cudaStream));
217  });
218 
219  return product;
220 }
221 
std::vector< float, cms::cuda::HostAllocator< float > > diffVarItvlIdxZEROVec_
std::vector< float, cms::cuda::HostAllocator< float > > accVarLenIdxMinusOneVec_
std::vector< uint32_t, cms::cuda::HostAllocator< uint32_t > > param1_
std::vector< float, cms::cuda::HostAllocator< float > > diffVarItvlIdxMinusOneVec_
std::vector< float, cms::cuda::HostAllocator< float > > accVarLenIdxZEROVec_
Product const & getProduct(cudaStream_t) const
constexpr int maxPSshapeBin
Definition: HcalConstants.h:7
constexpr int nsPerBX
Definition: HcalConstants.h:8
const tAllContWithNames getAllContainers() const
#define TYPELOOKUP_DATA_REG(_dataclass_)
Definition: typelookup.h:102
const Shape & getShape(int shapeType) const
std::vector< float, cms::cuda::HostAllocator< float > > acc25nsVec_
std::vector< uint32_t, cms::cuda::HostAllocator< uint32_t > > param2_
unsigned long long uint64_t
Definition: Time.h:13
constexpr int maxSamples
Definition: HcalConstants.h:6
def functor(code, kwds, debug=0)
Definition: utils.py:70
#define cudaCheck(ARG,...)
Definition: cudaCheck.h:69
std::vector< float, cms::cuda::HostAllocator< float > > diff25nsItvlVec_
std::vector< uint32_t, cms::cuda::HostAllocator< uint32_t > > ids_