CMS 3D CMS Logo

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