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 //-----------------------------------------------------------------------------
24  const MagneticField* mag,
25  const TrackerGeometry& geom,
26  const TrackerTopology& ttopo,
27  const SiPixelLorentzAngle* lorentzAngle,
28  const SiPixelGenErrorDBObject* genErrorDBObject,
29  const SiPixelLorentzAngle* lorentzAngleWidth)
30  : PixelCPEGenericBase(conf, mag, geom, ttopo, lorentzAngle, genErrorDBObject, lorentzAngleWidth) {
31  // Use errors from templates or from GenError
34  throw cms::Exception("InvalidCalibrationLoaded")
35  << "ERROR: GenErrors not filled correctly. Check the sqlite file. Using SiPixelTemplateDBObject version "
36  << (*genErrorDBObject_).version();
37  }
38 
40 
41  cpuData_ = {
43  detParamsGPU_.data(),
46  };
47 }
48 
49 const pixelCPEforGPU::ParamsOnGPU* PixelCPEFast::getGPUProductAsync(cudaStream_t cudaStream) const {
50  const auto& data = gpuData_.dataForCurrentDeviceAsync(cudaStream, [this](GPUData& data, cudaStream_t stream) {
51  // and now copy to device...
52  cudaCheck(cudaMalloc((void**)&data.paramsOnGPU_h.m_commonParams, sizeof(pixelCPEforGPU::CommonParams)));
53  cudaCheck(cudaMalloc((void**)&data.paramsOnGPU_h.m_detParams,
54  this->detParamsGPU_.size() * sizeof(pixelCPEforGPU::DetParams)));
55  cudaCheck(cudaMalloc((void**)&data.paramsOnGPU_h.m_averageGeometry, sizeof(pixelCPEforGPU::AverageGeometry)));
56  cudaCheck(cudaMalloc((void**)&data.paramsOnGPU_h.m_layerGeometry, sizeof(pixelCPEforGPU::LayerGeometry)));
57  cudaCheck(cudaMalloc((void**)&data.paramsOnGPU_d, sizeof(pixelCPEforGPU::ParamsOnGPU)));
58 
59  cudaCheck(cudaMemcpyAsync(
60  data.paramsOnGPU_d, &data.paramsOnGPU_h, sizeof(pixelCPEforGPU::ParamsOnGPU), cudaMemcpyDefault, stream));
61  cudaCheck(cudaMemcpyAsync((void*)data.paramsOnGPU_h.m_commonParams,
62  &this->commonParamsGPU_,
64  cudaMemcpyDefault,
65  stream));
66  cudaCheck(cudaMemcpyAsync((void*)data.paramsOnGPU_h.m_averageGeometry,
67  &this->averageGeometry_,
69  cudaMemcpyDefault,
70  stream));
71  cudaCheck(cudaMemcpyAsync((void*)data.paramsOnGPU_h.m_layerGeometry,
72  &this->layerGeometry_,
74  cudaMemcpyDefault,
75  stream));
76  cudaCheck(cudaMemcpyAsync((void*)data.paramsOnGPU_h.m_detParams,
77  this->detParamsGPU_.data(),
78  this->detParamsGPU_.size() * sizeof(pixelCPEforGPU::DetParams),
79  cudaMemcpyDefault,
80  stream));
81  });
82  return data.paramsOnGPU_d;
83 }
84 
86  //
87  // this code executes only once per job, computation inefficiency is not an issue
88  // many code blocks are repeated: better keep the computation local and self oconsistent as blocks may in future move around, be deleted ...
89  // It is valid only for Phase1 and the version of GenError in DB used in late 2018 and in 2021
90 
91  commonParamsGPU_.theThicknessB = m_DetParams.front().theThickness;
92  commonParamsGPU_.theThicknessE = m_DetParams.back().theThickness;
93  commonParamsGPU_.thePitchX = m_DetParams[0].thePitchX;
94  commonParamsGPU_.thePitchY = m_DetParams[0].thePitchY;
95 
96  LogDebug("PixelCPEFast") << "pitch & thickness " << commonParamsGPU_.thePitchX << ' ' << commonParamsGPU_.thePitchY
98 
99  // zero average geometry
101 
102  uint32_t oldLayer = 0;
103  uint32_t oldLadder = 0;
104  float rl = 0;
105  float zl = 0;
106  float miz = 90, mxz = 0;
107  float pl = 0;
108  int nl = 0;
109  detParamsGPU_.resize(m_DetParams.size());
110  for (auto i = 0U; i < m_DetParams.size(); ++i) {
111  auto& p = m_DetParams[i];
112  auto& g = detParamsGPU_[i];
113 
114  assert(p.theDet->index() == int(i));
115  assert(commonParamsGPU_.thePitchY == p.thePitchY);
116  assert(commonParamsGPU_.thePitchX == p.thePitchX);
117 
118  g.isBarrel = GeomDetEnumerators::isBarrel(p.thePart);
119  g.isPosZ = p.theDet->surface().position().z() > 0;
120  g.layer = ttopo_.layer(p.theDet->geographicalId());
121  g.index = i; // better be!
122  g.rawId = p.theDet->geographicalId();
124  assert(thickness == p.theThickness);
125 
126  auto ladder = ttopo_.pxbLadder(p.theDet->geographicalId());
127  if (oldLayer != g.layer) {
128  oldLayer = g.layer;
129  LogDebug("PixelCPEFast") << "new layer at " << i << (g.isBarrel ? " B " : (g.isPosZ ? " E+ " : " E- "))
130  << g.layer << " starting at " << g.rawId << '\n'
131  << "old layer had " << nl << " ladders";
132  nl = 0;
133  }
134  if (oldLadder != ladder) {
135  oldLadder = ladder;
136  LogDebug("PixelCPEFast") << "new ladder at " << i << (g.isBarrel ? " B " : (g.isPosZ ? " E+ " : " E- "))
137  << ladder << " starting at " << g.rawId << '\n'
138  << "old ladder ave z,r,p mz " << zl / 8.f << " " << rl / 8.f << " " << pl / 8.f << ' '
139  << miz << ' ' << mxz;
140  rl = 0;
141  zl = 0;
142  pl = 0;
143  miz = 90;
144  mxz = 0;
145  nl++;
146  }
147 
148  g.shiftX = 0.5f * p.lorentzShiftInCmX;
149  g.shiftY = 0.5f * p.lorentzShiftInCmY;
150  g.chargeWidthX = p.lorentzShiftInCmX * p.widthLAFractionX;
151  g.chargeWidthY = p.lorentzShiftInCmY * p.widthLAFractionY;
152 
153  g.x0 = p.theOrigin.x();
154  g.y0 = p.theOrigin.y();
155  g.z0 = p.theOrigin.z();
156 
157  auto vv = p.theDet->surface().position();
158  auto rr = pixelCPEforGPU::Rotation(p.theDet->surface().rotation());
159  g.frame = pixelCPEforGPU::Frame(vv.x(), vv.y(), vv.z(), rr);
160 
161  zl += vv.z();
162  miz = std::min(miz, std::abs(vv.z()));
163  mxz = std::max(mxz, std::abs(vv.z()));
164  rl += vv.perp();
165  pl += vv.phi(); // (not obvious)
166 
167  // errors .....
169 
170  cp.with_track_angle = false;
171 
172  auto lape = p.theDet->localAlignmentError();
173  if (lape.invalid())
174  lape = LocalError(); // zero....
175 
176  g.apeXX = lape.xx();
177  g.apeYY = lape.yy();
178 
179  auto toMicron = [&](float x) { return std::min(511, int(x * 1.e4f + 0.5f)); };
180 
181  {
182  // average angle
183  auto gvx = p.theOrigin.x() + 40.f * commonParamsGPU_.thePitchX;
184  auto gvy = p.theOrigin.y();
185  auto gvz = 1.f / p.theOrigin.z();
186  cp.cotalpha = gvx * gvz;
187  cp.cotbeta = gvy * gvz;
188  errorFromTemplates(p, cp, 20000.);
189  }
190 
191 #ifdef EDM_ML_DEBUG
192  auto m = 10000.f;
193  for (float qclus = 15000; qclus < 35000; qclus += 15000) {
194  errorFromTemplates(p, cp, qclus);
195  LogDebug("PixelCPEFast") << i << ' ' << qclus << ' ' << cp.pixmx << ' ' << m * cp.sigmax << ' ' << m * cp.sx1
196  << ' ' << m * cp.sx2 << ' ' << m * cp.sigmay << ' ' << m * cp.sy1 << ' ' << m * cp.sy2;
197  }
198  LogDebug("PixelCPEFast") << i << ' ' << m * std::sqrt(lape.xx()) << ' ' << m * std::sqrt(lape.yy());
199 #endif // EDM_ML_DEBUG
200 
201  g.pixmx = std::max(0, cp.pixmx);
202  g.sx2 = toMicron(cp.sx2);
203  g.sy1 = std::max(21, toMicron(cp.sy1)); // for some angles sy1 is very small
204  g.sy2 = std::max(55, toMicron(cp.sy2)); // sometimes sy2 is smaller than others (due to angle?)
205 
206  // sample xerr as function of position
207  auto const xoff = -81.f * commonParamsGPU_.thePitchX;
208 
209  for (int ix = 0; ix < 16; ++ix) {
210  auto x = xoff * (1.f - (0.5f + float(ix)) / 8.f);
211  auto gvx = p.theOrigin.x() - x;
212  auto gvy = p.theOrigin.y();
213  auto gvz = 1.f / p.theOrigin.z();
214  cp.cotbeta = gvy * gvz;
215  cp.cotalpha = gvx * gvz;
216  errorFromTemplates(p, cp, 20000.f);
217  g.sigmax[ix] = toMicron(cp.sigmax);
218  g.sigmax1[ix] = toMicron(cp.sx1);
219  LogDebug("PixelCPEFast") << "sigmax vs x " << i << ' ' << x << ' ' << cp.cotalpha << ' ' << int(g.sigmax[ix])
220  << ' ' << int(g.sigmax1[ix]) << ' ' << 10000.f * cp.sigmay << std::endl;
221  }
222 #ifdef EDM_ML_DEBUG
223  // sample yerr as function of position
224  auto const yoff = -54.f * 4.f * commonParamsGPU_.thePitchY;
225  for (int ix = 0; ix < 16; ++ix) {
226  auto y = yoff * (1.f - (0.5f + float(ix)) / 8.f);
227  auto gvx = p.theOrigin.x() + 40.f * commonParamsGPU_.thePitchY;
228  auto gvy = p.theOrigin.y() - y;
229  auto gvz = 1.f / p.theOrigin.z();
230  cp.cotbeta = gvy * gvz;
231  cp.cotalpha = gvx * gvz;
232  errorFromTemplates(p, cp, 20000.f);
233  LogDebug("PixelCPEFast") << "sigmay vs y " << i << ' ' << y << ' ' << cp.cotbeta << ' ' << 10000.f * cp.sigmay
234  << std::endl;
235  }
236 #endif // EDM_ML_DEBUG
237 
238  // average angle
239  auto gvx = p.theOrigin.x() + 40.f * commonParamsGPU_.thePitchX;
240  auto gvy = p.theOrigin.y();
241  auto gvz = 1.f / p.theOrigin.z();
242  //--- Note that the normalization is not required as only the ratio used
243 
244  // calculate angles (fed into errorFromTemplates)
245  cp.cotalpha = gvx * gvz;
246  cp.cotbeta = gvy * gvz;
247  auto aveCB = cp.cotbeta;
248 
249  // sample x by charge
250  int qbin = 5; // low charge
251  int k = 0;
252  for (int qclus = 1000; qclus < 200000; qclus += 1000) {
253  errorFromTemplates(p, cp, qclus);
254  if (cp.qBin_ == qbin)
255  continue;
256  qbin = cp.qBin_;
257  g.xfact[k] = cp.sigmax;
258  g.yfact[k] = cp.sigmay;
259  g.minCh[k++] = qclus;
260 #ifdef EDM_ML_DEBUG
261  LogDebug("PixelCPEFast") << i << ' ' << g.rawId << ' ' << cp.cotalpha << ' ' << qclus << ' ' << cp.qBin_ << ' '
262  << cp.pixmx << ' ' << m * cp.sigmax << ' ' << m * cp.sx1 << ' ' << m * cp.sx2 << ' '
263  << m * cp.sigmay << ' ' << m * cp.sy1 << ' ' << m * cp.sy2 << std::endl;
264 #endif // EDM_ML_DEBUG
265  }
266  assert(k <= 5);
267  // fill the rest (sometimes bin 4 is missing)
268  for (int kk = k; kk < 5; ++kk) {
269  g.xfact[kk] = g.xfact[k - 1];
270  g.yfact[kk] = g.yfact[k - 1];
271  g.minCh[kk] = g.minCh[k - 1];
272  }
273  auto detx = 1.f / g.xfact[0];
274  auto dety = 1.f / g.yfact[0];
275  for (int kk = 0; kk < 5; ++kk) {
276  g.xfact[kk] *= detx;
277  g.yfact[kk] *= dety;
278  }
279  // sample y in "angle" (estimated from cluster size)
280  float ys = 8.f - 4.f; // apperent bias of half pixel (see plot)
281  // sample yerr as function of "size"
282  for (int iy = 0; iy < 16; ++iy) {
283  ys += 1.f; // first bin 0 is for size 9 (and size is in fixed point 2^3)
284  if (15 == iy)
285  ys += 8.f; // last bin for "overflow"
286  // cp.cotalpha = ys*(commonParamsGPU_.thePitchX/(8.f*thickness)); // use this to print sampling in "x" (and comment the line below)
287  cp.cotbeta = std::copysign(ys * (commonParamsGPU_.thePitchY / (8.f * thickness)), aveCB);
288  errorFromTemplates(p, cp, 20000.f);
289  g.sigmay[iy] = toMicron(cp.sigmay);
290  LogDebug("PixelCPEFast") << "sigmax/sigmay " << i << ' ' << (ys + 4.f) / 8.f << ' ' << cp.cotalpha << '/'
291  << cp.cotbeta << ' ' << 10000.f * cp.sigmax << '/' << int(g.sigmay[iy]) << std::endl;
292  }
293  } // loop over det
294 
295  //
296  // compute ladder baricenter (only in global z) for the barrel
297  //
298  auto& aveGeom = averageGeometry_;
299  int il = 0;
300  for (int im = 0, nm = phase1PixelTopology::numberOfModulesInBarrel; im < nm; ++im) {
301  auto const& g = detParamsGPU_[im];
302  il = im / 8;
304  auto z = g.frame.z();
305  aveGeom.ladderZ[il] += 0.125f * z;
306  aveGeom.ladderMinZ[il] = std::min(aveGeom.ladderMinZ[il], z);
307  aveGeom.ladderMaxZ[il] = std::max(aveGeom.ladderMaxZ[il], z);
308  aveGeom.ladderX[il] += 0.125f * g.frame.x();
309  aveGeom.ladderY[il] += 0.125f * g.frame.y();
310  aveGeom.ladderR[il] += 0.125f * sqrt(g.frame.x() * g.frame.x() + g.frame.y() * g.frame.y());
311  }
313  // add half_module and tollerance
314  constexpr float module_length = 6.7f;
315  constexpr float module_tolerance = 0.2f;
316  for (int il = 0, nl = phase1PixelTopology::numberOfLaddersInBarrel; il < nl; ++il) {
317  aveGeom.ladderMinZ[il] -= (0.5f * module_length - module_tolerance);
318  aveGeom.ladderMaxZ[il] += (0.5f * module_length - module_tolerance);
319  }
320 
321  // compute "max z" for first layer in endcap (should we restrict to the outermost ring?)
322  for (auto im = phase1PixelTopology::layerStart[4]; im < phase1PixelTopology::layerStart[5]; ++im) {
323  auto const& g = detParamsGPU_[im];
324  aveGeom.endCapZ[0] = std::max(aveGeom.endCapZ[0], g.frame.z());
325  }
326  for (auto im = phase1PixelTopology::layerStart[7]; im < phase1PixelTopology::layerStart[8]; ++im) {
327  auto const& g = detParamsGPU_[im];
328  aveGeom.endCapZ[1] = std::min(aveGeom.endCapZ[1], g.frame.z());
329  }
330  // correct for outer ring being closer
331  aveGeom.endCapZ[0] -= 1.5f;
332  aveGeom.endCapZ[1] += 1.5f;
333 
334  for (int jl = 0, nl = phase1PixelTopology::numberOfLaddersInBarrel; jl < nl; ++jl) {
335  LogDebug("PixelCPEFast") << jl << ':' << aveGeom.ladderR[jl] << '/'
336  << std::sqrt(aveGeom.ladderX[jl] * aveGeom.ladderX[jl] +
337  aveGeom.ladderY[jl] * aveGeom.ladderY[jl])
338  << ',' << aveGeom.ladderZ[jl] << ',' << aveGeom.ladderMinZ[jl] << ','
339  << aveGeom.ladderMaxZ[jl] << '\n';
340  }
341  LogDebug("PixelCPEFast") << aveGeom.endCapZ[0] << ' ' << aveGeom.endCapZ[1];
342 
343  // fill Layer and ladders geometry
346 }
347 
349  if (paramsOnGPU_d != nullptr) {
350  cudaFree((void*)paramsOnGPU_h.m_commonParams);
351  cudaFree((void*)paramsOnGPU_h.m_detParams);
352  cudaFree((void*)paramsOnGPU_h.m_averageGeometry);
353  cudaFree((void*)paramsOnGPU_h.m_layerGeometry);
354  cudaFree(paramsOnGPU_d);
355  }
356 }
357 
359  ClusterParamGeneric& theClusterParam,
360  float qclus) const {
361  float locBz = theDetParam.bz;
362  float locBx = theDetParam.bx;
363  LogDebug("PixelCPEFast") << "PixelCPEFast::localPosition(...) : locBz = " << locBz;
364 
365  theClusterParam.pixmx = std::numeric_limits<int>::max(); // max pixel charge for truncation of 2-D cluster
366 
367  theClusterParam.sigmay = -999.9; // CPE Generic y-error for multi-pixel cluster
368  theClusterParam.sigmax = -999.9; // CPE Generic x-error for multi-pixel cluster
369  theClusterParam.sy1 = -999.9; // CPE Generic y-error for single single-pixel
370  theClusterParam.sy2 = -999.9; // CPE Generic y-error for single double-pixel cluster
371  theClusterParam.sx1 = -999.9; // CPE Generic x-error for single single-pixel cluster
372  theClusterParam.sx2 = -999.9; // CPE Generic x-error for single double-pixel cluster
373 
374  float dummy;
375 
377  int gtemplID = theDetParam.detTemplateId;
378 
379  theClusterParam.qBin_ = gtempl.qbin(gtemplID,
380  theClusterParam.cotalpha,
381  theClusterParam.cotbeta,
382  locBz,
383  locBx,
384  qclus,
385  false,
386  theClusterParam.pixmx,
387  theClusterParam.sigmay,
388  dummy,
389  theClusterParam.sigmax,
390  dummy,
391  theClusterParam.sy1,
392  dummy,
393  theClusterParam.sy2,
394  dummy,
395  theClusterParam.sx1,
396  dummy,
397  theClusterParam.sx2,
398  dummy);
399 
400  theClusterParam.sigmax = theClusterParam.sigmax * micronsToCm;
401  theClusterParam.sx1 = theClusterParam.sx1 * micronsToCm;
402  theClusterParam.sx2 = theClusterParam.sx2 * micronsToCm;
403 
404  theClusterParam.sigmay = theClusterParam.sigmay * micronsToCm;
405  theClusterParam.sy1 = theClusterParam.sy1 * micronsToCm;
406  theClusterParam.sy2 = theClusterParam.sy2 * micronsToCm;
407 }
408 
409 //-----------------------------------------------------------------------------
413 //-----------------------------------------------------------------------------
414 LocalPoint PixelCPEFast::localPosition(DetParam const& theDetParam, ClusterParam& theClusterParamBase) const {
415  ClusterParamGeneric& theClusterParam = static_cast<ClusterParamGeneric&>(theClusterParamBase);
416 
417  assert(!theClusterParam.with_track_angle);
418 
420  errorFromTemplates(theDetParam, theClusterParam, theClusterParam.theCluster->charge());
421  } else {
422  theClusterParam.qBin_ = 0;
423  }
424 
425  int q_f_X;
426  int q_l_X;
427  int q_f_Y;
428  int q_l_Y;
429  collect_edge_charges(theClusterParam, q_f_X, q_l_X, q_f_Y, q_l_Y, useErrorsFromTemplates_ && truncatePixelCharge_);
430 
431  // do GPU like ...
433 
434  cp.minRow[0] = theClusterParam.theCluster->minPixelRow();
435  cp.maxRow[0] = theClusterParam.theCluster->maxPixelRow();
436  cp.minCol[0] = theClusterParam.theCluster->minPixelCol();
437  cp.maxCol[0] = theClusterParam.theCluster->maxPixelCol();
438 
439  cp.q_f_X[0] = q_f_X;
440  cp.q_l_X[0] = q_l_X;
441  cp.q_f_Y[0] = q_f_Y;
442  cp.q_l_Y[0] = q_l_Y;
443 
444  cp.charge[0] = theClusterParam.theCluster->charge();
445 
446  auto ind = theDetParam.theDet->index();
448  auto xPos = cp.xpos[0];
449  auto yPos = cp.ypos[0];
450 
451  // set the error (mind ape....)
453  theClusterParam.sigmax = cp.xerr[0];
454  theClusterParam.sigmay = cp.yerr[0];
455 
456  LogDebug("PixelCPEFast") << " in PixelCPEFast:localPosition - pos = " << xPos << " " << yPos << " size "
457  << cp.maxRow[0] - cp.minRow[0] << ' ' << cp.maxCol[0] - cp.minCol[0];
458 
459  //--- Now put the two together
460  LocalPoint pos_in_local(xPos, yPos);
461  return pos_in_local;
462 }
463 
464 //============== INFLATED ERROR AND ERRORS FROM DB BELOW ================
465 
466 //-------------------------------------------------------------------------
467 // Hit error in the local frame
468 //-------------------------------------------------------------------------
469 LocalError PixelCPEFast::localError(DetParam const& theDetParam, ClusterParam& theClusterParamBase) const {
470  ClusterParamGeneric& theClusterParam = static_cast<ClusterParamGeneric&>(theClusterParamBase);
471 
472  auto xerr = theClusterParam.sigmax;
473  auto yerr = theClusterParam.sigmay;
474 
475  LogDebug("PixelCPEFast") << " errors " << xerr << " " << yerr;
476 
477  auto xerr_sq = xerr * xerr;
478  auto yerr_sq = yerr * yerr;
479 
480  return LocalError(xerr_sq, 0, yerr_sq);
481 }
482 
484  // call PixelCPEGenericBase fillPSetDescription to add common rechit errors
486 }
PixelCPEFast::GPUData::~GPUData
~GPUData()
Definition: PixelCPEFast.cc:348
SiPixelCluster::minPixelCol
int minPixelCol() const
Definition: SiPixelCluster.h:152
PixelCPEBase::ClusterParam
Definition: PixelCPEBase.h:69
PixelCPEBase::DetParam::bz
float bz
Definition: PixelCPEBase.h:58
PixelCPEGenericBase::useErrorsFromTemplates_
bool useErrorsFromTemplates_
Definition: PixelCPEGenericBase.h:86
DDAxes::y
PixelCPEBase::DetParam::bx
float bx
Definition: PixelCPEBase.h:59
pixelCPEforGPU::Rotation
SOARotation< float > Rotation
Definition: pixelCPEforGPU.h:19
mps_fire.i
i
Definition: mps_fire.py:428
PixelCPEGenericBase::ClusterParamGeneric::sy2
float sy2
Definition: PixelCPEGenericBase.h:22
MessageLogger.h
dqmMemoryStats.float
float
Definition: dqmMemoryStats.py:127
pixelCPEforGPU::LayerGeometry::layer
uint8_t layer[phase1PixelTopology::layerIndexSize]
Definition: pixelCPEforGPU.h:57
PixelCPEFast::cpuData_
pixelCPEforGPU::ParamsOnGPU cpuData_
Definition: PixelCPEFast.h:49
RectangularPixelTopology.h
f
double f[11][100]
Definition: MuScleFitUtils.cc:78
PixelCPEFast::localError
LocalError localError(DetParam const &theDetParam, ClusterParam &theClusterParam) const override
Definition: PixelCPEFast.cc:469
pixelCPEforGPU::CommonParams::thePitchX
float thePitchX
Definition: pixelCPEforGPU.h:25
min
T min(T a, T b)
Definition: MathUtil.h:58
findQualityFiles.rr
string rr
Definition: findQualityFiles.py:185
TrackerTopology
Definition: TrackerTopology.h:16
PixelCPEBase::genErrorDBObject_
const SiPixelGenErrorDBObject * genErrorDBObject_
Definition: PixelCPEBase.h:235
PixelCPEGenericBase::ClusterParamGeneric::sx1
float sx1
Definition: PixelCPEGenericBase.h:23
TrackerTopology::pxbLadder
unsigned int pxbLadder(const DetId &id) const
Definition: TrackerTopology.h:155
PixelCPEBase::ttopo_
const TrackerTopology & ttopo_
Definition: PixelCPEBase.h:230
SiPixelCluster::maxPixelCol
int maxPixelCol() const
Definition: SiPixelCluster.h:153
PixelCPEFast::fillPSetDescription
static void fillPSetDescription(edm::ParameterSetDescription &desc)
Definition: PixelCPEFast.cc:483
PixelCPEGenericBase::ClusterParamGeneric::pixmx
int pixmx
Definition: PixelCPEGenericBase.h:16
phase1PixelTopology::numberOfModulesInBarrel
constexpr uint32_t numberOfModulesInBarrel
Definition: phase1PixelTopology.h:50
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
PixelCPEGenericBase::ClusterParamGeneric::sy1
float sy1
Definition: PixelCPEGenericBase.h:21
cms::cuda::stream
uint32_t const T *__restrict__ const uint32_t *__restrict__ int32_t int Histo::index_type cudaStream_t stream
Definition: HistoContainer.h:51
PixelCPEGenericBase::fillPSetDescription
static void fillPSetDescription(edm::ParameterSetDescription &desc)
Definition: PixelCPEGenericBase.cc:217
cms::cuda::assert
assert(be >=bs)
SiPixelCluster::maxPixelRow
int maxPixelRow() const
Definition: SiPixelCluster.h:151
SiPixelGenError::qbin
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)
Definition: SiPixelGenError.cc:538
TrackerTopology::layer
unsigned int layer(const DetId &id) const
Definition: TrackerTopology.cc:47
GeomDet::index
int index() const
Definition: GeomDet.h:83
hgcal_conditions::parameters
Definition: HGCConditions.h:86
DDAxes::x
PixelCPEFast.h
pixelCPEforGPU::ParamsOnGPU::m_commonParams
CommonParams const * m_commonParams
Definition: pixelCPEforGPU.h:61
PixelCPEGenericBase::ClusterParamGeneric::sigmax
float sigmax
Definition: PixelCPEGenericBase.h:20
PixelCPEBase::DetParam::theDet
const PixelGeomDetUnit * theDet
Definition: PixelCPEBase.h:47
PixelCPEGenericBase::ClusterParamGeneric
Definition: PixelCPEGenericBase.h:9
PixelCPEBase::ClusterParam::theCluster
const SiPixelCluster * theCluster
Definition: PixelCPEBase.h:75
PixelCPEFast::GPUData::paramsOnGPU_d
pixelCPEforGPU::ParamsOnGPU * paramsOnGPU_d
Definition: PixelCPEFast.h:55
SiPixelTemplate.h
PixelCPEGenericBase::collect_edge_charges
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)
Definition: PixelCPEGenericBase.cc:45
PixelCPEFast::thePixelGenError_
std::vector< SiPixelGenErrorStore > thePixelGenError_
Definition: PixelCPEFast.h:42
PixelCPEBase::m_DetParams
DetParams m_DetParams
Definition: PixelCPEBase.h:283
Calorimetry_cff.thickness
thickness
Definition: Calorimetry_cff.py:115
visualization-live-secondInstance_cfg.m
m
Definition: visualization-live-secondInstance_cfg.py:78
SiPixelCluster::minPixelRow
int minPixelRow() const
Definition: SiPixelCluster.h:150
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
DDAxes::z
SiPixelLorentzAngle
Definition: SiPixelLorentzAngle.h:11
relativeConstraints.geom
geom
Definition: relativeConstraints.py:72
GetRecoTauVFromDQM_MC_cff.kk
kk
Definition: GetRecoTauVFromDQM_MC_cff.py:84
dqmdumpme.k
k
Definition: dqmdumpme.py:60
Point3DBase< float, LocalTag >
pixelCPEforGPU::LayerGeometry::layerStart
uint32_t layerStart[phase1PixelTopology::numberOfLayers+1]
Definition: pixelCPEforGPU.h:56
mitigatedMETSequence_cff.U
U
Definition: mitigatedMETSequence_cff.py:36
phase1PixelTopology::layer
constexpr std::array< uint8_t, layerIndexSize > layer
Definition: phase1PixelTopology.h:99
PixelCPEGenericBase::ClusterParamGeneric::sx2
float sx2
Definition: PixelCPEGenericBase.h:24
pixelCPEforGPU::Frame
SOAFrame< float > Frame
Definition: pixelCPEforGPU.h:18
PixelCPEFast::errorFromTemplates
void errorFromTemplates(DetParam const &theDetParam, ClusterParamGeneric &theClusterParam, float qclus) const
Definition: PixelCPEFast.cc:358
PixelCPEFast::GPUData
Definition: PixelCPEFast.h:51
phase1PixelTopology::numberOfLaddersInBarrel
constexpr uint32_t numberOfLaddersInBarrel
Definition: phase1PixelTopology.h:51
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:233
SiPixelGenError
Definition: SiPixelGenError.h:113
edm::ParameterSet
Definition: ParameterSet.h:47
AlCaHLTBitMon_ParallelJobs.p
def p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
pixelCPEforGPU::LayerGeometry
Definition: pixelCPEforGPU.h:55
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
PixelCPEGenericBase::truncatePixelCharge_
const bool truncatePixelCharge_
Definition: PixelCPEGenericBase.h:87
LocalError
Definition: LocalError.h:12
createfilelist.int
int
Definition: createfilelist.py:10
GeomDetEnumerators::isBarrel
bool isBarrel(GeomDetEnumerators::SubDetector m)
Definition: GeomDetEnumerators.cc:57
cudaCheck.h
pixelCPEforGPU::ParamsOnGPU::m_averageGeometry
AverageGeometry const * m_averageGeometry
Definition: pixelCPEforGPU.h:64
MagneticField.h
pixelCPEforGPU::position
constexpr void position(CommonParams const &__restrict__ comParams, DetParams const &__restrict__ detParams, ClusParams &cp, uint32_t ic)
Definition: pixelCPEforGPU.h:183
PixelCPEFast::fillParamsForGpu
void fillParamsForGpu()
Definition: PixelCPEFast.cc:85
pixelCPEforGPU::CommonParams::theThicknessE
float theThicknessE
Definition: pixelCPEforGPU.h:24
pixelCPEforGPU::CommonParams
Definition: pixelCPEforGPU.h:22
PixelCPEGenericBase
Definition: PixelCPEGenericBase.h:7
phase1PixelTopology::AverageGeometry
Definition: phase1PixelTopology.h:161
SiPixelGenError::pushfile
static bool pushfile(int filenum, std::vector< SiPixelGenErrorStore > &pixelTemp, std::string dir="")
Definition: SiPixelGenError.cc:55
PixelCPEFast::GPUData::paramsOnGPU_h
pixelCPEforGPU::ParamsOnGPU paramsOnGPU_h
Definition: PixelCPEFast.h:54
PixelCPEFast::layerGeometry_
pixelCPEforGPU::LayerGeometry layerGeometry_
Definition: PixelCPEFast.h:47
PixelCPEBase::DetParam::detTemplateId
int detTemplateId
Definition: PixelCPEBase.h:65
PixelCPEFast::detParamsGPU_
std::vector< pixelCPEforGPU::DetParams > detParamsGPU_
Definition: PixelCPEFast.h:45
pixelCPEforGPU::ParamsOnGPU
Definition: pixelCPEforGPU.h:60
mag
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
Definition: Basic3DVectorLD.h:127
PixelCPEBase::DetParam
Definition: PixelCPEBase.h:45
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
PixelCPEFast::averageGeometry_
pixelCPEforGPU::AverageGeometry averageGeometry_
Definition: PixelCPEFast.h:48
cudaCheck
#define cudaCheck(ARG,...)
Definition: cudaCheck.h:69
pixelCPEforGPU::ParamsOnGPU::m_layerGeometry
LayerGeometry const * m_layerGeometry
Definition: pixelCPEforGPU.h:63
DetId.h
PixelCPEGenericBase::ClusterParamGeneric::sigmay
float sigmay
Definition: PixelCPEGenericBase.h:19
Exception
Definition: hltDiff.cc:245
PixelCPEFast::gpuData_
cms::cuda::ESProduct< GPUData > gpuData_
Definition: PixelCPEFast.h:57
pixelCPEforGPU::errorFromDB
constexpr void errorFromDB(CommonParams const &__restrict__ comParams, DetParams const &__restrict__ detParams, ClusParams &cp, uint32_t ic)
Definition: pixelCPEforGPU.h:323
PixelGeomDetUnit.h
pixelCPEforGPU::ParamsOnGPU::m_detParams
DetParams const * m_detParams
Definition: pixelCPEforGPU.h:62
PVValHelper::ladder
Definition: PVValidationHelpers.h:73
SiPixelCluster::charge
int charge() const
Definition: SiPixelCluster.h:142
phase1PixelTopology::layerStart
constexpr uint32_t layerStart[numberOfLayers+1]
Definition: phase1PixelTopology.h:26
data
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:79
PixelCPEBase::ClusterParam::cotbeta
float cotbeta
Definition: PixelCPEBase.h:79
SiPixelGenErrorDBObject
Definition: SiPixelGenErrorDBObject.h:16
PixelCPEFast::getGPUProductAsync
const pixelCPEforGPU::ParamsOnGPU * getGPUProductAsync(cudaStream_t cudaStream) const
Definition: PixelCPEFast.cc:49
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
dummy
Definition: DummySelector.h:38
PixelCPEBase::ClusterParam::cotalpha
float cotalpha
Definition: PixelCPEBase.h:78
pixelCPEforGPU::ClusParamsT
Definition: pixelCPEforGPU.h:84
pixelCPEforGPU::DetParams
Definition: pixelCPEforGPU.h:29
CommonMethods.cp
def cp(fromDir, toDir, listOfFiles, overwrite=False, smallList=False)
Definition: CommonMethods.py:191
SOAFrame::z
constexpr T z() const
Definition: SOARotation.h:133
MagneticField
Definition: MagneticField.h:19
PixelCPEFast::localPosition
LocalPoint localPosition(DetParam const &theDetParam, ClusterParam &theClusterParam) const override
Definition: PixelCPEFast.cc:414
phase1PixelTopology.h
pixelCPEforGPU::CommonParams::theThicknessB
float theThicknessB
Definition: pixelCPEforGPU.h:23
pixelCPEforGPU::CommonParams::thePitchY
float thePitchY
Definition: pixelCPEforGPU.h:26
PixelCPEBase::ClusterParam::with_track_angle
bool with_track_angle
Definition: PixelCPEBase.h:102
g
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
PixelCPEFast::PixelCPEFast
PixelCPEFast(edm::ParameterSet const &conf, const MagneticField *, const TrackerGeometry &, const TrackerTopology &, const SiPixelLorentzAngle *, const SiPixelGenErrorDBObject *, const SiPixelLorentzAngle *)
The constructor.
Definition: PixelCPEFast.cc:23
TrackerGeometry
Definition: TrackerGeometry.h:14
PixelCPEBase::ClusterParam::qBin_
int qBin_
Definition: PixelCPEBase.h:95
PixelCPEFast::commonParamsGPU_
pixelCPEforGPU::CommonParams commonParamsGPU_
Definition: PixelCPEFast.h:46