CMS 3D CMS Logo

PixelCPEFast.cc
Go to the documentation of this file.
1 #include <cuda_runtime.h>
2 
12 
13 // Services
14 // this is needed to get errors from templates
15 
16 namespace {
17  constexpr float micronsToCm = 1.0e-4;
18 }
19 
20 //-----------------------------------------------------------------------------
22 //-----------------------------------------------------------------------------
23 template <typename TrackerTraits>
25  const MagneticField* mag,
26  const TrackerGeometry& geom,
27  const TrackerTopology& ttopo,
28  const SiPixelLorentzAngle* lorentzAngle,
29  const SiPixelGenErrorDBObject* genErrorDBObject,
30  const SiPixelLorentzAngle* lorentzAngleWidth)
31  : PixelCPEGenericBase(conf, mag, geom, ttopo, lorentzAngle, genErrorDBObject, lorentzAngleWidth) {
32  // Use errors from templates or from GenError
35  throw cms::Exception("InvalidCalibrationLoaded")
36  << "ERROR: GenErrors not filled correctly. Check the sqlite file. Using SiPixelTemplateDBObject version "
37  << (*genErrorDBObject_).version();
38  }
39 
41 
42  cpuData_ = {
44  detParamsGPU_.data(),
47  };
48 }
49 
50 template <typename TrackerTraits>
52  cudaStream_t cudaStream) const {
56 
57  const auto& data = gpuData_.dataForCurrentDeviceAsync(cudaStream, [this](GPUData& data, cudaStream_t stream) {
58  // and now copy to device...
59 
60  cudaCheck(cudaMalloc((void**)&data.paramsOnGPU_h.m_commonParams, sizeof(pixelCPEforGPU::CommonParams)));
61  cudaCheck(cudaMalloc((void**)&data.paramsOnGPU_h.m_detParams,
62  this->detParamsGPU_.size() * sizeof(pixelCPEforGPU::DetParams)));
63  cudaCheck(cudaMalloc((void**)&data.paramsOnGPU_h.m_averageGeometry, sizeof(AverageGeometry)));
64  cudaCheck(cudaMalloc((void**)&data.paramsOnGPU_h.m_layerGeometry, sizeof(LayerGeometry)));
65  cudaCheck(cudaMalloc((void**)&data.paramsOnGPU_d, sizeof(ParamsOnGPU)));
66  cudaCheck(cudaMemcpyAsync(data.paramsOnGPU_d, &data.paramsOnGPU_h, sizeof(ParamsOnGPU), cudaMemcpyDefault, stream));
67  cudaCheck(cudaMemcpyAsync((void*)data.paramsOnGPU_h.m_commonParams,
68  &this->commonParamsGPU_,
70  cudaMemcpyDefault,
71  stream));
72  cudaCheck(cudaMemcpyAsync((void*)data.paramsOnGPU_h.m_averageGeometry,
73  &this->averageGeometry_,
74  sizeof(AverageGeometry),
75  cudaMemcpyDefault,
76  stream));
77  cudaCheck(cudaMemcpyAsync((void*)data.paramsOnGPU_h.m_layerGeometry,
78  &this->layerGeometry_,
79  sizeof(LayerGeometry),
80  cudaMemcpyDefault,
81  stream));
82  cudaCheck(cudaMemcpyAsync((void*)data.paramsOnGPU_h.m_detParams,
83  this->detParamsGPU_.data(),
84  this->detParamsGPU_.size() * sizeof(pixelCPEforGPU::DetParams),
85  cudaMemcpyDefault,
86  stream));
87  });
88  return data.paramsOnGPU_d;
89 }
90 
91 template <typename TrackerTraits>
93  //
94  // this code executes only once per job, computation inefficiency is not an issue
95  // many code blocks are repeated: better keep the computation local and self oconsistent as blocks may in future move around, be deleted ...
96  // It is valid only for Phase1 and the version of GenError in DB used in late 2018 and in 2021
97 
98  commonParamsGPU_.theThicknessB = m_DetParams.front().theThickness;
99  commonParamsGPU_.theThicknessE = m_DetParams.back().theThickness;
100  commonParamsGPU_.thePitchX = m_DetParams[0].thePitchX;
101  commonParamsGPU_.thePitchY = m_DetParams[0].thePitchY;
102 
103  commonParamsGPU_.numberOfLaddersInBarrel = TrackerTraits::numberOfLaddersInBarrel;
104 
105  LogDebug("PixelCPEFast") << "pitch & thickness " << commonParamsGPU_.thePitchX << ' ' << commonParamsGPU_.thePitchY
106  << " " << commonParamsGPU_.theThicknessB << ' ' << commonParamsGPU_.theThicknessE;
107 
108  // zero average geometry
109  memset(&averageGeometry_, 0, sizeof(pixelTopology::AverageGeometryT<TrackerTraits>));
110 
111  uint32_t oldLayer = 0;
112  uint32_t oldLadder = 0;
113  float rl = 0;
114  float zl = 0;
115  float miz = 500, mxz = 0;
116  float pl = 0;
117  int nl = 0;
118  detParamsGPU_.resize(m_DetParams.size());
119 
120  for (auto i = 0U; i < m_DetParams.size(); ++i) {
121  auto& p = m_DetParams[i];
122  auto& g = detParamsGPU_[i];
123 
124  g.nRowsRoc = p.theDet->specificTopology().rowsperroc();
125  g.nColsRoc = p.theDet->specificTopology().colsperroc();
126  g.nRows = p.theDet->specificTopology().rocsX() * g.nRowsRoc;
127  g.nCols = p.theDet->specificTopology().rocsY() * g.nColsRoc;
128 
129  g.numPixsInModule = g.nRows * g.nCols;
130 
131  assert(p.theDet->index() == int(i));
132  assert(commonParamsGPU_.thePitchY == p.thePitchY);
133  assert(commonParamsGPU_.thePitchX == p.thePitchX);
134 
135  g.isBarrel = GeomDetEnumerators::isBarrel(p.thePart);
136  g.isPosZ = p.theDet->surface().position().z() > 0;
137  g.layer = ttopo_.layer(p.theDet->geographicalId());
138  g.index = i; // better be!
139  g.rawId = p.theDet->geographicalId();
140  auto thickness = g.isBarrel ? commonParamsGPU_.theThicknessB : commonParamsGPU_.theThicknessE;
141  assert(thickness == p.theThickness);
142 
143  auto ladder = ttopo_.pxbLadder(p.theDet->geographicalId());
144  if (oldLayer != g.layer) {
145  oldLayer = g.layer;
146  LogDebug("PixelCPEFast") << "new layer at " << i << (g.isBarrel ? " B " : (g.isPosZ ? " E+ " : " E- "))
147  << g.layer << " starting at " << g.rawId << '\n'
148  << "old layer had " << nl << " ladders";
149  nl = 0;
150  }
151  if (oldLadder != ladder) {
152  oldLadder = ladder;
153  LogDebug("PixelCPEFast") << "new ladder at " << i << (g.isBarrel ? " B " : (g.isPosZ ? " E+ " : " E- "))
154  << ladder << " starting at " << g.rawId << '\n'
155  << "old ladder ave z,r,p mz " << zl / 8.f << " " << rl / 8.f << " " << pl / 8.f << ' '
156  << miz << ' ' << mxz;
157  rl = 0;
158  zl = 0;
159  pl = 0;
160  miz = 500;
161  mxz = 0;
162  nl++;
163  }
164 
165  g.shiftX = 0.5f * p.lorentzShiftInCmX;
166  g.shiftY = 0.5f * p.lorentzShiftInCmY;
167  g.chargeWidthX = p.lorentzShiftInCmX * p.widthLAFractionX;
168  g.chargeWidthY = p.lorentzShiftInCmY * p.widthLAFractionY;
169 
170  g.x0 = p.theOrigin.x();
171  g.y0 = p.theOrigin.y();
172  g.z0 = p.theOrigin.z();
173 
174  auto vv = p.theDet->surface().position();
175  auto rr = pixelCPEforGPU::Rotation(p.theDet->surface().rotation());
176  g.frame = pixelCPEforGPU::Frame(vv.x(), vv.y(), vv.z(), rr);
177 
178  zl += vv.z();
179  miz = std::min(miz, std::abs(vv.z()));
180  mxz = std::max(mxz, std::abs(vv.z()));
181  rl += vv.perp();
182  pl += vv.phi(); // (not obvious)
183 
184  // errors .....
186 
187  cp.with_track_angle = false;
188 
189  auto lape = p.theDet->localAlignmentError();
190  if (lape.invalid())
191  lape = LocalError(); // zero....
192 
193  g.apeXX = lape.xx();
194  g.apeYY = lape.yy();
195 
196  auto toMicron = [&](float x) { return std::min(511, int(x * 1.e4f + 0.5f)); };
197 
198  // average angle
199  auto gvx = p.theOrigin.x() + 40.f * commonParamsGPU_.thePitchX;
200  auto gvy = p.theOrigin.y();
201  auto gvz = 1.f / p.theOrigin.z();
202  //--- Note that the normalization is not required as only the ratio used
203 
204  {
205  // calculate angles (fed into errorFromTemplates)
206  cp.cotalpha = gvx * gvz;
207  cp.cotbeta = gvy * gvz;
208 
209  errorFromTemplates(p, cp, 20000.);
210  }
211 
212 #ifdef EDM_ML_DEBUG
213  auto m = 10000.f;
214  for (float qclus = 15000; qclus < 35000; qclus += 15000) {
215  errorFromTemplates(p, cp, qclus);
216  LogDebug("PixelCPEFast") << i << ' ' << qclus << ' ' << cp.pixmx << ' ' << m * cp.sigmax << ' ' << m * cp.sx1
217  << ' ' << m * cp.sx2 << ' ' << m * cp.sigmay << ' ' << m * cp.sy1 << ' ' << m * cp.sy2;
218  }
219  LogDebug("PixelCPEFast") << i << ' ' << m * std::sqrt(lape.xx()) << ' ' << m * std::sqrt(lape.yy());
220 #endif // EDM_ML_DEBUG
221 
222  g.pixmx = std::max(0, cp.pixmx);
223  g.sx2 = toMicron(cp.sx2);
224  g.sy1 = std::max(21, toMicron(cp.sy1)); // for some angles sy1 is very small
225  g.sy2 = std::max(55, toMicron(cp.sy2)); // sometimes sy2 is smaller than others (due to angle?)
226 
227  //sample xerr as function of position
228  // moduleOffsetX is the definition of TrackerTraits::xOffset,
229  // needs to be calculated because for Phase2 the modules are not uniform
230  float moduleOffsetX = -(0.5f * float(g.nRows) + TrackerTraits::bigPixXCorrection);
231  auto const xoff = moduleOffsetX * commonParamsGPU_.thePitchX;
232 
233  for (int ix = 0; ix < CPEFastParametrisation::kNumErrorBins; ++ix) {
234  auto x = xoff * (1.f - (0.5f + float(ix)) / 8.f);
235  auto gvx = p.theOrigin.x() - x;
236  auto gvy = p.theOrigin.y();
237  auto gvz = 1.f / p.theOrigin.z();
238  cp.cotbeta = gvy * gvz;
239  cp.cotalpha = gvx * gvz;
240  errorFromTemplates(p, cp, 20000.f);
241  g.sigmax[ix] = toMicron(cp.sigmax);
242  g.sigmax1[ix] = toMicron(cp.sx1);
243  LogDebug("PixelCPEFast") << "sigmax vs x " << i << ' ' << x << ' ' << cp.cotalpha << ' ' << int(g.sigmax[ix])
244  << ' ' << int(g.sigmax1[ix]) << ' ' << 10000.f * cp.sigmay << std::endl;
245  }
246 #ifdef EDM_ML_DEBUG
247  // sample yerr as function of position
248  // moduleOffsetY is the definition of TrackerTraits::yOffset (removed)
249  float moduleOffsetY = 0.5f * float(g.nCols) + TrackerTraits::bigPixYCorrection;
250  auto const yoff = -moduleOffsetY * commonParamsGPU_.thePitchY;
251 
252  for (int ix = 0; ix < CPEFastParametrisation::kNumErrorBins; ++ix) {
253  auto y = yoff * (1.f - (0.5f + float(ix)) / 8.f);
254  auto gvx = p.theOrigin.x() + 40.f * commonParamsGPU_.thePitchY;
255  auto gvy = p.theOrigin.y() - y;
256  auto gvz = 1.f / p.theOrigin.z();
257  cp.cotbeta = gvy * gvz;
258  cp.cotalpha = gvx * gvz;
259  errorFromTemplates(p, cp, 20000.f);
260  LogDebug("PixelCPEFast") << "sigmay vs y " << i << ' ' << y << ' ' << cp.cotbeta << ' ' << 10000.f * cp.sigmay
261  << std::endl;
262  }
263 #endif // EDM_ML_DEBUG
264 
265  // calculate angles (repeated)
266  cp.cotalpha = gvx * gvz;
267  cp.cotbeta = gvy * gvz;
268  auto aveCB = cp.cotbeta;
269 
270  // sample x by charge
271  int qbin = CPEFastParametrisation::kGenErrorQBins; // low charge
272  int k = 0;
273  for (int qclus = 1000; qclus < 200000; qclus += 1000) {
274  errorFromTemplates(p, cp, qclus);
275  if (cp.qBin_ == qbin)
276  continue;
277  qbin = cp.qBin_;
278  g.xfact[k] = cp.sigmax;
279  g.yfact[k] = cp.sigmay;
280  g.minCh[k++] = qclus;
281 #ifdef EDM_ML_DEBUG
282  LogDebug("PixelCPEFast") << i << ' ' << g.rawId << ' ' << cp.cotalpha << ' ' << qclus << ' ' << cp.qBin_ << ' '
283  << cp.pixmx << ' ' << m * cp.sigmax << ' ' << m * cp.sx1 << ' ' << m * cp.sx2 << ' '
284  << m * cp.sigmay << ' ' << m * cp.sy1 << ' ' << m * cp.sy2 << std::endl;
285 #endif // EDM_ML_DEBUG
286  }
287 
289 
290  // fill the rest (sometimes bin 4 is missing)
291  for (int kk = k; kk < CPEFastParametrisation::kGenErrorQBins; ++kk) {
292  g.xfact[kk] = g.xfact[k - 1];
293  g.yfact[kk] = g.yfact[k - 1];
294  g.minCh[kk] = g.minCh[k - 1];
295  }
296  auto detx = 1.f / g.xfact[0];
297  auto dety = 1.f / g.yfact[0];
298  for (int kk = 0; kk < CPEFastParametrisation::kGenErrorQBins; ++kk) {
299  g.xfact[kk] *= detx;
300  g.yfact[kk] *= dety;
301  }
302  // sample y in "angle" (estimated from cluster size)
303  float ys = 8.f - 4.f; // apperent bias of half pixel (see plot)
304  // plot: https://indico.cern.ch/event/934821/contributions/3974619/attachments/2091853/3515041/DigilessReco.pdf page 25
305  // sample yerr as function of "size"
306  for (int iy = 0; iy < CPEFastParametrisation::kNumErrorBins; ++iy) {
307  ys += 1.f; // first bin 0 is for size 9 (and size is in fixed point 2^3)
309  ys += 8.f; // last bin for "overflow"
310  // cp.cotalpha = ys*(commonParamsGPU_.thePitchX/(8.f*thickness)); // use this to print sampling in "x" (and comment the line below)
311  cp.cotbeta = std::copysign(ys * (commonParamsGPU_.thePitchY / (8.f * thickness)), aveCB);
312  errorFromTemplates(p, cp, 20000.f);
313  g.sigmay[iy] = toMicron(cp.sigmay);
314  LogDebug("PixelCPEFast") << "sigmax/sigmay " << i << ' ' << (ys + 4.f) / 8.f << ' ' << cp.cotalpha << '/'
315  << cp.cotbeta << ' ' << 10000.f * cp.sigmax << '/' << int(g.sigmay[iy]) << std::endl;
316  }
317  } // loop over det
318 
319  constexpr int numberOfModulesInLadder = TrackerTraits::numberOfModulesInLadder;
320  constexpr int numberOfLaddersInBarrel = TrackerTraits::numberOfLaddersInBarrel;
321  constexpr int numberOfModulesInBarrel = TrackerTraits::numberOfModulesInBarrel;
322 
323  constexpr float ladderFactor = 1.f / float(numberOfModulesInLadder);
324 
325  constexpr int firstEndcapPos = TrackerTraits::firstEndcapPos;
326  constexpr int firstEndcapNeg = TrackerTraits::firstEndcapNeg;
327 
328  // compute ladder baricenter (only in global z) for the barrel
329  //
330  auto& aveGeom = averageGeometry_;
331  int il = 0;
332  for (int im = 0, nm = numberOfModulesInBarrel; im < nm; ++im) {
333  auto const& g = detParamsGPU_[im];
334  il = im / numberOfModulesInLadder;
335  assert(il < int(numberOfLaddersInBarrel));
336  auto z = g.frame.z();
337  aveGeom.ladderZ[il] += ladderFactor * z;
338  aveGeom.ladderMinZ[il] = std::min(aveGeom.ladderMinZ[il], z);
339  aveGeom.ladderMaxZ[il] = std::max(aveGeom.ladderMaxZ[il], z);
340  aveGeom.ladderX[il] += ladderFactor * g.frame.x();
341  aveGeom.ladderY[il] += ladderFactor * g.frame.y();
342  aveGeom.ladderR[il] += ladderFactor * sqrt(g.frame.x() * g.frame.x() + g.frame.y() * g.frame.y());
343  }
344  assert(il + 1 == int(numberOfLaddersInBarrel));
345  // add half_module and tollerance
346  constexpr float moduleLength = TrackerTraits::moduleLength;
347  constexpr float module_tolerance = 0.2f;
348  for (int il = 0, nl = numberOfLaddersInBarrel; il < nl; ++il) {
349  aveGeom.ladderMinZ[il] -= (0.5f * moduleLength - module_tolerance);
350  aveGeom.ladderMaxZ[il] += (0.5f * moduleLength - module_tolerance);
351  }
352 
353  // compute "max z" for first layer in endcap (should we restrict to the outermost ring?)
354  for (auto im = TrackerTraits::layerStart[firstEndcapPos]; im < TrackerTraits::layerStart[firstEndcapPos + 1]; ++im) {
355  auto const& g = detParamsGPU_[im];
356  aveGeom.endCapZ[0] = std::max(aveGeom.endCapZ[0], g.frame.z());
357  }
358  for (auto im = TrackerTraits::layerStart[firstEndcapNeg]; im < TrackerTraits::layerStart[firstEndcapNeg + 1]; ++im) {
359  auto const& g = detParamsGPU_[im];
360  aveGeom.endCapZ[1] = std::min(aveGeom.endCapZ[1], g.frame.z());
361  }
362  // correct for outer ring being closer
363  aveGeom.endCapZ[0] -= TrackerTraits::endcapCorrection;
364  aveGeom.endCapZ[1] += TrackerTraits::endcapCorrection;
365 #ifdef EDM_ML_DEBUG
366  for (int jl = 0, nl = numberOfLaddersInBarrel; jl < nl; ++jl) {
367  LogDebug("PixelCPEFast") << jl << ':' << aveGeom.ladderR[jl] << '/'
368  << std::sqrt(aveGeom.ladderX[jl] * aveGeom.ladderX[jl] +
369  aveGeom.ladderY[jl] * aveGeom.ladderY[jl])
370  << ',' << aveGeom.ladderZ[jl] << ',' << aveGeom.ladderMinZ[jl] << ','
371  << aveGeom.ladderMaxZ[jl] << '\n';
372  }
373  LogDebug("PixelCPEFast") << aveGeom.endCapZ[0] << ' ' << aveGeom.endCapZ[1];
374 #endif // EDM_ML_DEBUG
375 
376  // fill Layer and ladders geometry
377  memset(&layerGeometry_, 0, sizeof(pixelCPEforGPU::LayerGeometryT<TrackerTraits>));
378  memcpy(layerGeometry_.layerStart,
381  memcpy(layerGeometry_.layer, pixelTopology::layer<TrackerTraits>.data(), pixelTopology::layer<TrackerTraits>.size());
382  layerGeometry_.maxModuleStride = pixelTopology::maxModuleStride<TrackerTraits>;
383 }
384 
385 template <typename TrackerTraits>
387  if (paramsOnGPU_d != nullptr) {
388  cudaFree((void*)paramsOnGPU_h.m_commonParams);
389  cudaFree((void*)paramsOnGPU_h.m_detParams);
390  cudaFree((void*)paramsOnGPU_h.m_averageGeometry);
391  cudaFree((void*)paramsOnGPU_h.m_layerGeometry);
392  cudaFree(paramsOnGPU_d);
393  }
394 }
395 
396 template <typename TrackerTraits>
398  ClusterParamGeneric& theClusterParam,
399  float qclus) const {
400  float locBz = theDetParam.bz;
401  float locBx = theDetParam.bx;
402  LogDebug("PixelCPEFast") << "PixelCPEFast::localPosition(...) : locBz = " << locBz;
403 
404  theClusterParam.pixmx = std::numeric_limits<int>::max(); // max pixel charge for truncation of 2-D cluster
405 
406  theClusterParam.sigmay = -999.9; // CPE Generic y-error for multi-pixel cluster
407  theClusterParam.sigmax = -999.9; // CPE Generic x-error for multi-pixel cluster
408  theClusterParam.sy1 = -999.9; // CPE Generic y-error for single single-pixel
409  theClusterParam.sy2 = -999.9; // CPE Generic y-error for single double-pixel cluster
410  theClusterParam.sx1 = -999.9; // CPE Generic x-error for single single-pixel cluster
411  theClusterParam.sx2 = -999.9; // CPE Generic x-error for single double-pixel cluster
412 
413  float dummy;
414 
416  int gtemplID = theDetParam.detTemplateId;
417 
418  theClusterParam.qBin_ = gtempl.qbin(gtemplID,
419  theClusterParam.cotalpha,
420  theClusterParam.cotbeta,
421  locBz,
422  locBx,
423  qclus,
424  false,
425  theClusterParam.pixmx,
426  theClusterParam.sigmay,
427  dummy,
428  theClusterParam.sigmax,
429  dummy,
430  theClusterParam.sy1,
431  dummy,
432  theClusterParam.sy2,
433  dummy,
434  theClusterParam.sx1,
435  dummy,
436  theClusterParam.sx2,
437  dummy);
438 
439  theClusterParam.sigmax = theClusterParam.sigmax * micronsToCm;
440  theClusterParam.sx1 = theClusterParam.sx1 * micronsToCm;
441  theClusterParam.sx2 = theClusterParam.sx2 * micronsToCm;
442 
443  theClusterParam.sigmay = theClusterParam.sigmay * micronsToCm;
444  theClusterParam.sy1 = theClusterParam.sy1 * micronsToCm;
445  theClusterParam.sy2 = theClusterParam.sy2 * micronsToCm;
446 }
447 
448 template <>
450  ClusterParamGeneric& theClusterParam,
451  float qclus) const {
452  theClusterParam.qBin_ = 0.0f;
453 }
454 
455 //-----------------------------------------------------------------------------
459 //-----------------------------------------------------------------------------
460 template <typename TrackerTraits>
462  ClusterParam& theClusterParamBase) const {
463  ClusterParamGeneric& theClusterParam = static_cast<ClusterParamGeneric&>(theClusterParamBase);
464 
466  errorFromTemplates(theDetParam, theClusterParam, theClusterParam.theCluster->charge());
467  } else {
468  theClusterParam.qBin_ = 0;
469  }
470 
471  int q_f_X;
472  int q_l_X;
473  int q_f_Y;
474  int q_l_Y;
475  collect_edge_charges(theClusterParam, q_f_X, q_l_X, q_f_Y, q_l_Y, useErrorsFromTemplates_ && truncatePixelCharge_);
476 
477  // do GPU like ...
479 
480  cp.minRow[0] = theClusterParam.theCluster->minPixelRow();
481  cp.maxRow[0] = theClusterParam.theCluster->maxPixelRow();
482  cp.minCol[0] = theClusterParam.theCluster->minPixelCol();
483  cp.maxCol[0] = theClusterParam.theCluster->maxPixelCol();
484 
485  cp.q_f_X[0] = q_f_X;
486  cp.q_l_X[0] = q_l_X;
487  cp.q_f_Y[0] = q_f_Y;
488  cp.q_l_Y[0] = q_l_Y;
489 
490  cp.charge[0] = theClusterParam.theCluster->charge();
491 
492  auto ind = theDetParam.theDet->index();
493  pixelCPEforGPU::position<TrackerTraits>(commonParamsGPU_, detParamsGPU_[ind], cp, 0);
494  auto xPos = cp.xpos[0];
495  auto yPos = cp.ypos[0];
496 
497  // set the error (mind ape....)
498  pixelCPEforGPU::errorFromDB<TrackerTraits>(commonParamsGPU_, detParamsGPU_[ind], cp, 0);
499  theClusterParam.sigmax = cp.xerr[0];
500  theClusterParam.sigmay = cp.yerr[0];
501 
502  LogDebug("PixelCPEFast") << " in PixelCPEFast:localPosition - pos = " << xPos << " " << yPos << " size "
503  << cp.maxRow[0] - cp.minRow[0] << ' ' << cp.maxCol[0] - cp.minCol[0];
504 
505  //--- Now put the two together
506  LocalPoint pos_in_local(xPos, yPos);
507  return pos_in_local;
508 }
509 
510 //============== INFLATED ERROR AND ERRORS FROM DB BELOW ================
511 
512 //-------------------------------------------------------------------------
513 // Hit error in the local frame
514 //-------------------------------------------------------------------------
515 template <typename TrackerTraits>
517  ClusterParam& theClusterParamBase) const {
518  ClusterParamGeneric& theClusterParam = static_cast<ClusterParamGeneric&>(theClusterParamBase);
519 
520  auto xerr = theClusterParam.sigmax;
521  auto yerr = theClusterParam.sigmay;
522 
523  LogDebug("PixelCPEFast") << " errors " << xerr << " " << yerr;
524 
525  auto xerr_sq = xerr * xerr;
526  auto yerr_sq = yerr * yerr;
527 
528  return LocalError(xerr_sq, 0, yerr_sq);
529 }
530 
531 template <typename TrackerTraits>
533  // call PixelCPEGenericBase fillPSetDescription to add common rechit errors
535 }
536 
static bool pushfile(int filenum, std::vector< SiPixelGenErrorStore > &pixelTemp, std::string dir="")
LocalError localError(DetParam const &theDetParam, ClusterParam &theClusterParam) const override
ParamsOnGPU * paramsOnGPU_d
Definition: PixelCPEFast.h:60
SOAFrame< float > Frame
int index() const
Definition: GeomDet.h:83
const SiPixelCluster * theCluster
Definition: PixelCPEBase.h:76
AverageGeometry const * m_averageGeometry
int maxPixelRow() const
bool isBarrel(GeomDetEnumerators::SubDetector m)
static void fillPSetDescription(edm::ParameterSetDescription &desc)
static void collect_edge_charges(ClusterParam &theClusterParam, int &q_f_X, int &q_l_X, int &q_f_Y, int &q_l_Y, bool truncate)
DetParams const * m_detParams
constexpr int kGenErrorQBins
uint32_t T const *__restrict__ uint32_t const *__restrict__ int32_t int Histo::index_type cudaStream_t stream
assert(be >=bs)
const PixelGeomDetUnit * theDet
Definition: PixelCPEBase.h:47
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
float float float z
int minPixelRow() const
int charge() const
std::vector< SiPixelGenErrorStore > thePixelGenError_
Definition: PixelCPEFast.h:47
int qbin(int id, float cotalpha, float cotbeta, float locBz, float locBx, float qclus, bool irradiationCorrections, int &pixmx, float &sigmay, float &deltay, float &sigmax, float &deltax, float &sy1, float &dy1, float &sy2, float &dy2, float &sx1, float &dx1, float &sx2, float &dx2)
PixelCPEFast(edm::ParameterSet const &conf, const MagneticField *, const TrackerGeometry &, const TrackerTopology &, const SiPixelLorentzAngle *, const SiPixelGenErrorDBObject *, const SiPixelLorentzAngle *)
The constructor.
Definition: PixelCPEFast.cc:24
CommonParams const * m_commonParams
std::vector< pixelCPEforGPU::DetParams > detParamsGPU_
Definition: PixelCPEFast.h:50
LayerGeometry const * m_layerGeometry
int minPixelCol() const
ParamsOnGPU cpuData_
Definition: PixelCPEFast.h:54
T sqrt(T t)
Definition: SSEVec.h:19
int maxPixelCol() const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
constexpr T z() const
Definition: SOARotation.h:133
LayerGeometry layerGeometry_
Definition: PixelCPEFast.h:52
double f[11][100]
void fillParamsForGpu()
Definition: PixelCPEFast.cc:92
static void fillPSetDescription(edm::ParameterSetDescription &desc)
LocalPoint localPosition(DetParam const &theDetParam, ClusterParam &theClusterParam) const override
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
SOARotation< float > Rotation
ParamsOnGPU paramsOnGPU_h
Definition: PixelCPEFast.h:59
const SiPixelGenErrorDBObject * genErrorDBObject_
Definition: PixelCPEBase.h:236
AverageGeometry averageGeometry_
Definition: PixelCPEFast.h:53
const ParamsOnGPU * getGPUProductAsync(cudaStream_t cudaStream) const
Definition: PixelCPEFast.cc:51
pixelCPEforGPU::CommonParams commonParamsGPU_
Definition: PixelCPEFast.h:51
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:80
static constexpr uint32_t layerStart[numberOfLayers+1]
float x
void errorFromTemplates(DetParam const &theDetParam, ClusterParamGeneric &theClusterParam, float qclus) const
#define cudaCheck(ARG,...)
Definition: cudaCheck.h:69
constexpr int kNumErrorBins
def cp(fromDir, toDir, listOfFiles, overwrite=False, smallList=False)
#define LogDebug(id)