File indexing completed on 2023-10-25 10:00:39
0001 #include <cuda_runtime.h>
0002
0003 #include "CondFormats/SiPixelTransient/interface/SiPixelTemplate.h"
0004 #include "DataFormats/DetId/interface/DetId.h"
0005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0006 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0007 #include "Geometry/TrackerGeometryBuilder/interface/RectangularPixelTopology.h"
0008 #include "Geometry/CommonTopologies/interface/SimplePixelTopology.h"
0009 #include "HeterogeneousCore/CUDAUtilities/interface/cudaCheck.h"
0010 #include "MagneticField/Engine/interface/MagneticField.h"
0011 #include "RecoLocalTracker/SiPixelRecHits/interface/PixelCPEFast.h"
0012
0013
0014
0015
0016 namespace {
0017 constexpr float micronsToCm = 1.0e-4;
0018 }
0019
0020
0021
0022
0023 template <typename TrackerTraits>
0024 PixelCPEFast<TrackerTraits>::PixelCPEFast(edm::ParameterSet const& conf,
0025 const MagneticField* mag,
0026 const TrackerGeometry& geom,
0027 const TrackerTopology& ttopo,
0028 const SiPixelLorentzAngle* lorentzAngle,
0029 const SiPixelGenErrorDBObject* genErrorDBObject,
0030 const SiPixelLorentzAngle* lorentzAngleWidth)
0031 : PixelCPEGenericBase(conf, mag, geom, ttopo, lorentzAngle, genErrorDBObject, lorentzAngleWidth) {
0032
0033 if (useErrorsFromTemplates_) {
0034 if (!SiPixelGenError::pushfile(*genErrorDBObject_, thePixelGenError_))
0035 throw cms::Exception("InvalidCalibrationLoaded")
0036 << "ERROR: GenErrors not filled correctly. Check the sqlite file. Using SiPixelTemplateDBObject version "
0037 << (*genErrorDBObject_).version();
0038 }
0039
0040 fillParamsForGpu();
0041
0042 cpuData_ = {
0043 &commonParamsGPU_,
0044 detParamsGPU_.data(),
0045 &layerGeometry_,
0046 &averageGeometry_,
0047 };
0048 }
0049
0050 template <typename TrackerTraits>
0051 const pixelCPEforGPU::ParamsOnGPUT<TrackerTraits>* PixelCPEFast<TrackerTraits>::getGPUProductAsync(
0052 cudaStream_t cudaStream) const {
0053 using ParamsOnGPU = pixelCPEforGPU::ParamsOnGPUT<TrackerTraits>;
0054 using LayerGeometry = pixelCPEforGPU::LayerGeometryT<TrackerTraits>;
0055 using AverageGeometry = pixelTopology::AverageGeometryT<TrackerTraits>;
0056
0057 const auto& data = gpuData_.dataForCurrentDeviceAsync(cudaStream, [this](GPUData& data, cudaStream_t stream) {
0058
0059
0060 cudaCheck(cudaMalloc((void**)&data.paramsOnGPU_h.m_commonParams, sizeof(pixelCPEforGPU::CommonParams)));
0061 cudaCheck(cudaMalloc((void**)&data.paramsOnGPU_h.m_detParams,
0062 this->detParamsGPU_.size() * sizeof(pixelCPEforGPU::DetParams)));
0063 cudaCheck(cudaMalloc((void**)&data.paramsOnGPU_h.m_averageGeometry, sizeof(AverageGeometry)));
0064 cudaCheck(cudaMalloc((void**)&data.paramsOnGPU_h.m_layerGeometry, sizeof(LayerGeometry)));
0065 cudaCheck(cudaMalloc((void**)&data.paramsOnGPU_d, sizeof(ParamsOnGPU)));
0066 cudaCheck(cudaMemcpyAsync(data.paramsOnGPU_d, &data.paramsOnGPU_h, sizeof(ParamsOnGPU), cudaMemcpyDefault, stream));
0067 cudaCheck(cudaMemcpyAsync((void*)data.paramsOnGPU_h.m_commonParams,
0068 &this->commonParamsGPU_,
0069 sizeof(pixelCPEforGPU::CommonParams),
0070 cudaMemcpyDefault,
0071 stream));
0072 cudaCheck(cudaMemcpyAsync((void*)data.paramsOnGPU_h.m_averageGeometry,
0073 &this->averageGeometry_,
0074 sizeof(AverageGeometry),
0075 cudaMemcpyDefault,
0076 stream));
0077 cudaCheck(cudaMemcpyAsync((void*)data.paramsOnGPU_h.m_layerGeometry,
0078 &this->layerGeometry_,
0079 sizeof(LayerGeometry),
0080 cudaMemcpyDefault,
0081 stream));
0082 cudaCheck(cudaMemcpyAsync((void*)data.paramsOnGPU_h.m_detParams,
0083 this->detParamsGPU_.data(),
0084 this->detParamsGPU_.size() * sizeof(pixelCPEforGPU::DetParams),
0085 cudaMemcpyDefault,
0086 stream));
0087 });
0088 return data.paramsOnGPU_d;
0089 }
0090
0091 template <typename TrackerTraits>
0092 void PixelCPEFast<TrackerTraits>::fillParamsForGpu() {
0093
0094
0095
0096
0097
0098 commonParamsGPU_.theThicknessB = m_DetParams.front().theThickness;
0099 commonParamsGPU_.theThicknessE = m_DetParams.back().theThickness;
0100 commonParamsGPU_.thePitchX = m_DetParams[0].thePitchX;
0101 commonParamsGPU_.thePitchY = m_DetParams[0].thePitchY;
0102
0103 commonParamsGPU_.numberOfLaddersInBarrel = TrackerTraits::numberOfLaddersInBarrel;
0104
0105 LogDebug("PixelCPEFast") << "pitch & thickness " << commonParamsGPU_.thePitchX << ' ' << commonParamsGPU_.thePitchY
0106 << " " << commonParamsGPU_.theThicknessB << ' ' << commonParamsGPU_.theThicknessE;
0107
0108
0109 memset(&averageGeometry_, 0, sizeof(pixelTopology::AverageGeometryT<TrackerTraits>));
0110
0111 uint32_t oldLayer = 0;
0112 uint32_t oldLadder = 0;
0113 float rl = 0;
0114 float zl = 0;
0115 float miz = 500, mxz = 0;
0116 float pl = 0;
0117 int nl = 0;
0118 detParamsGPU_.resize(m_DetParams.size());
0119
0120 for (auto i = 0U; i < m_DetParams.size(); ++i) {
0121 auto& p = m_DetParams[i];
0122 auto& g = detParamsGPU_[i];
0123
0124 g.nRowsRoc = p.theDet->specificTopology().rowsperroc();
0125 g.nColsRoc = p.theDet->specificTopology().colsperroc();
0126 g.nRows = p.theDet->specificTopology().rocsX() * g.nRowsRoc;
0127 g.nCols = p.theDet->specificTopology().rocsY() * g.nColsRoc;
0128
0129 g.numPixsInModule = g.nRows * g.nCols;
0130
0131 assert(p.theDet->index() == int(i));
0132 assert(commonParamsGPU_.thePitchY == p.thePitchY);
0133 assert(commonParamsGPU_.thePitchX == p.thePitchX);
0134
0135 g.isBarrel = GeomDetEnumerators::isBarrel(p.thePart);
0136 g.isPosZ = p.theDet->surface().position().z() > 0;
0137 g.layer = ttopo_.layer(p.theDet->geographicalId());
0138 g.index = i;
0139 g.rawId = p.theDet->geographicalId();
0140 auto thickness = g.isBarrel ? commonParamsGPU_.theThicknessB : commonParamsGPU_.theThicknessE;
0141 assert(thickness == p.theThickness);
0142
0143 auto ladder = ttopo_.pxbLadder(p.theDet->geographicalId());
0144 if (oldLayer != g.layer) {
0145 oldLayer = g.layer;
0146 LogDebug("PixelCPEFast") << "new layer at " << i << (g.isBarrel ? " B " : (g.isPosZ ? " E+ " : " E- "))
0147 << g.layer << " starting at " << g.rawId << '\n'
0148 << "old layer had " << nl << " ladders";
0149 nl = 0;
0150 }
0151 if (oldLadder != ladder) {
0152 oldLadder = ladder;
0153 LogDebug("PixelCPEFast") << "new ladder at " << i << (g.isBarrel ? " B " : (g.isPosZ ? " E+ " : " E- "))
0154 << ladder << " starting at " << g.rawId << '\n'
0155 << "old ladder ave z,r,p mz " << zl / 8.f << " " << rl / 8.f << " " << pl / 8.f << ' '
0156 << miz << ' ' << mxz;
0157 rl = 0;
0158 zl = 0;
0159 pl = 0;
0160 miz = 500;
0161 mxz = 0;
0162 nl++;
0163 }
0164
0165 g.shiftX = 0.5f * p.lorentzShiftInCmX;
0166 g.shiftY = 0.5f * p.lorentzShiftInCmY;
0167 g.chargeWidthX = p.lorentzShiftInCmX * p.widthLAFractionX;
0168 g.chargeWidthY = p.lorentzShiftInCmY * p.widthLAFractionY;
0169
0170 g.x0 = p.theOrigin.x();
0171 g.y0 = p.theOrigin.y();
0172 g.z0 = p.theOrigin.z();
0173
0174 auto vv = p.theDet->surface().position();
0175 auto rr = pixelCPEforGPU::Rotation(p.theDet->surface().rotation());
0176 g.frame = pixelCPEforGPU::Frame(vv.x(), vv.y(), vv.z(), rr);
0177
0178 zl += vv.z();
0179 miz = std::min(miz, std::abs(vv.z()));
0180 mxz = std::max(mxz, std::abs(vv.z()));
0181 rl += vv.perp();
0182 pl += vv.phi();
0183
0184
0185 ClusterParamGeneric cp;
0186
0187 cp.with_track_angle = false;
0188
0189 auto lape = p.theDet->localAlignmentError();
0190 if (lape.invalid())
0191 lape = LocalError();
0192
0193 g.apeXX = lape.xx();
0194 g.apeYY = lape.yy();
0195
0196 auto toMicron = [&](float x) { return std::min(511, int(x * 1.e4f + 0.5f)); };
0197
0198
0199 auto gvx = p.theOrigin.x() + 40.f * commonParamsGPU_.thePitchX;
0200 auto gvy = p.theOrigin.y();
0201 auto gvz = 1.f / p.theOrigin.z();
0202
0203
0204 {
0205
0206 cp.cotalpha = gvx * gvz;
0207 cp.cotbeta = gvy * gvz;
0208
0209 errorFromTemplates(p, cp, 20000.);
0210 }
0211
0212 #ifdef EDM_ML_DEBUG
0213 auto m = 10000.f;
0214 for (float qclus = 15000; qclus < 35000; qclus += 15000) {
0215 errorFromTemplates(p, cp, qclus);
0216 LogDebug("PixelCPEFast") << i << ' ' << qclus << ' ' << cp.pixmx << ' ' << m * cp.sigmax << ' ' << m * cp.sx1
0217 << ' ' << m * cp.sx2 << ' ' << m * cp.sigmay << ' ' << m * cp.sy1 << ' ' << m * cp.sy2;
0218 }
0219 LogDebug("PixelCPEFast") << i << ' ' << m * std::sqrt(lape.xx()) << ' ' << m * std::sqrt(lape.yy());
0220 #endif
0221
0222 g.pixmx = std::max(0, cp.pixmx);
0223 g.sx2 = toMicron(cp.sx2);
0224 g.sy1 = std::max(21, toMicron(cp.sy1));
0225 g.sy2 = std::max(55, toMicron(cp.sy2));
0226
0227
0228
0229
0230 float moduleOffsetX = -(0.5f * float(g.nRows) + TrackerTraits::bigPixXCorrection);
0231 auto const xoff = moduleOffsetX * commonParamsGPU_.thePitchX;
0232
0233 for (int ix = 0; ix < CPEFastParametrisation::kNumErrorBins; ++ix) {
0234 auto x = xoff * (1.f - (0.5f + float(ix)) / 8.f);
0235 auto gvx = p.theOrigin.x() - x;
0236 auto gvy = p.theOrigin.y();
0237 auto gvz = 1.f / p.theOrigin.z();
0238 cp.cotbeta = gvy * gvz;
0239 cp.cotalpha = gvx * gvz;
0240 errorFromTemplates(p, cp, 20000.f);
0241 g.sigmax[ix] = toMicron(cp.sigmax);
0242 g.sigmax1[ix] = toMicron(cp.sx1);
0243 LogDebug("PixelCPEFast") << "sigmax vs x " << i << ' ' << x << ' ' << cp.cotalpha << ' ' << int(g.sigmax[ix])
0244 << ' ' << int(g.sigmax1[ix]) << ' ' << 10000.f * cp.sigmay << std::endl;
0245 }
0246 #ifdef EDM_ML_DEBUG
0247
0248
0249 float moduleOffsetY = 0.5f * float(g.nCols) + TrackerTraits::bigPixYCorrection;
0250 auto const yoff = -moduleOffsetY * commonParamsGPU_.thePitchY;
0251
0252 for (int ix = 0; ix < CPEFastParametrisation::kNumErrorBins; ++ix) {
0253 auto y = yoff * (1.f - (0.5f + float(ix)) / 8.f);
0254 auto gvx = p.theOrigin.x() + 40.f * commonParamsGPU_.thePitchY;
0255 auto gvy = p.theOrigin.y() - y;
0256 auto gvz = 1.f / p.theOrigin.z();
0257 cp.cotbeta = gvy * gvz;
0258 cp.cotalpha = gvx * gvz;
0259 errorFromTemplates(p, cp, 20000.f);
0260 LogDebug("PixelCPEFast") << "sigmay vs y " << i << ' ' << y << ' ' << cp.cotbeta << ' ' << 10000.f * cp.sigmay
0261 << std::endl;
0262 }
0263 #endif
0264
0265
0266 cp.cotalpha = gvx * gvz;
0267 cp.cotbeta = gvy * gvz;
0268 auto aveCB = cp.cotbeta;
0269
0270
0271 int qbin = CPEFastParametrisation::kGenErrorQBins;
0272 int k = 0;
0273 for (int qclus = 1000; qclus < 200000; qclus += 1000) {
0274 errorFromTemplates(p, cp, qclus);
0275 if (cp.qBin_ == qbin)
0276 continue;
0277 qbin = cp.qBin_;
0278 g.xfact[k] = cp.sigmax;
0279 g.yfact[k] = cp.sigmay;
0280 g.minCh[k++] = qclus;
0281 #ifdef EDM_ML_DEBUG
0282 LogDebug("PixelCPEFast") << i << ' ' << g.rawId << ' ' << cp.cotalpha << ' ' << qclus << ' ' << cp.qBin_ << ' '
0283 << cp.pixmx << ' ' << m * cp.sigmax << ' ' << m * cp.sx1 << ' ' << m * cp.sx2 << ' '
0284 << m * cp.sigmay << ' ' << m * cp.sy1 << ' ' << m * cp.sy2 << std::endl;
0285 #endif
0286 }
0287
0288 assert(k <= CPEFastParametrisation::kGenErrorQBins);
0289
0290
0291 for (int kk = k; kk < CPEFastParametrisation::kGenErrorQBins; ++kk) {
0292 g.xfact[kk] = g.xfact[k - 1];
0293 g.yfact[kk] = g.yfact[k - 1];
0294 g.minCh[kk] = g.minCh[k - 1];
0295 }
0296 auto detx = 1.f / g.xfact[0];
0297 auto dety = 1.f / g.yfact[0];
0298 for (int kk = 0; kk < CPEFastParametrisation::kGenErrorQBins; ++kk) {
0299 g.xfact[kk] *= detx;
0300 g.yfact[kk] *= dety;
0301 }
0302
0303 float ys = 8.f - 4.f;
0304
0305
0306 for (int iy = 0; iy < CPEFastParametrisation::kNumErrorBins; ++iy) {
0307 ys += 1.f;
0308 if (CPEFastParametrisation::kNumErrorBins - 1 == iy)
0309 ys += 8.f;
0310
0311 cp.cotbeta = std::copysign(ys * (commonParamsGPU_.thePitchY / (8.f * thickness)), aveCB);
0312 errorFromTemplates(p, cp, 20000.f);
0313 g.sigmay[iy] = toMicron(cp.sigmay);
0314 LogDebug("PixelCPEFast") << "sigmax/sigmay " << i << ' ' << (ys + 4.f) / 8.f << ' ' << cp.cotalpha << '/'
0315 << cp.cotbeta << ' ' << 10000.f * cp.sigmax << '/' << int(g.sigmay[iy]) << std::endl;
0316 }
0317 }
0318
0319 constexpr int numberOfModulesInLadder = TrackerTraits::numberOfModulesInLadder;
0320 constexpr int numberOfLaddersInBarrel = TrackerTraits::numberOfLaddersInBarrel;
0321 constexpr int numberOfModulesInBarrel = TrackerTraits::numberOfModulesInBarrel;
0322
0323 constexpr float ladderFactor = 1.f / float(numberOfModulesInLadder);
0324
0325 constexpr int firstEndcapPos = TrackerTraits::firstEndcapPos;
0326 constexpr int firstEndcapNeg = TrackerTraits::firstEndcapNeg;
0327
0328
0329
0330 auto& aveGeom = averageGeometry_;
0331 int il = 0;
0332 for (int im = 0, nm = numberOfModulesInBarrel; im < nm; ++im) {
0333 auto const& g = detParamsGPU_[im];
0334 il = im / numberOfModulesInLadder;
0335 assert(il < int(numberOfLaddersInBarrel));
0336 auto z = g.frame.z();
0337 aveGeom.ladderZ[il] += ladderFactor * z;
0338 aveGeom.ladderMinZ[il] = std::min(aveGeom.ladderMinZ[il], z);
0339 aveGeom.ladderMaxZ[il] = std::max(aveGeom.ladderMaxZ[il], z);
0340 aveGeom.ladderX[il] += ladderFactor * g.frame.x();
0341 aveGeom.ladderY[il] += ladderFactor * g.frame.y();
0342 aveGeom.ladderR[il] += ladderFactor * sqrt(g.frame.x() * g.frame.x() + g.frame.y() * g.frame.y());
0343 }
0344 assert(il + 1 == int(numberOfLaddersInBarrel));
0345
0346 constexpr float moduleLength = TrackerTraits::moduleLength;
0347 constexpr float module_tolerance = 0.2f;
0348 for (int il = 0, nl = numberOfLaddersInBarrel; il < nl; ++il) {
0349 aveGeom.ladderMinZ[il] -= (0.5f * moduleLength - module_tolerance);
0350 aveGeom.ladderMaxZ[il] += (0.5f * moduleLength - module_tolerance);
0351 }
0352
0353
0354 for (auto im = TrackerTraits::layerStart[firstEndcapPos]; im < TrackerTraits::layerStart[firstEndcapPos + 1]; ++im) {
0355 auto const& g = detParamsGPU_[im];
0356 aveGeom.endCapZ[0] = std::max(aveGeom.endCapZ[0], g.frame.z());
0357 }
0358 for (auto im = TrackerTraits::layerStart[firstEndcapNeg]; im < TrackerTraits::layerStart[firstEndcapNeg + 1]; ++im) {
0359 auto const& g = detParamsGPU_[im];
0360 aveGeom.endCapZ[1] = std::min(aveGeom.endCapZ[1], g.frame.z());
0361 }
0362
0363 aveGeom.endCapZ[0] -= TrackerTraits::endcapCorrection;
0364 aveGeom.endCapZ[1] += TrackerTraits::endcapCorrection;
0365 #ifdef EDM_ML_DEBUG
0366 for (int jl = 0, nl = numberOfLaddersInBarrel; jl < nl; ++jl) {
0367 LogDebug("PixelCPEFast") << jl << ':' << aveGeom.ladderR[jl] << '/'
0368 << std::sqrt(aveGeom.ladderX[jl] * aveGeom.ladderX[jl] +
0369 aveGeom.ladderY[jl] * aveGeom.ladderY[jl])
0370 << ',' << aveGeom.ladderZ[jl] << ',' << aveGeom.ladderMinZ[jl] << ','
0371 << aveGeom.ladderMaxZ[jl] << '\n';
0372 }
0373 LogDebug("PixelCPEFast") << aveGeom.endCapZ[0] << ' ' << aveGeom.endCapZ[1];
0374 #endif
0375
0376
0377 memset(&layerGeometry_, 0, sizeof(pixelCPEforGPU::LayerGeometryT<TrackerTraits>));
0378 memcpy(layerGeometry_.layerStart,
0379 TrackerTraits::layerStart,
0380 sizeof(pixelCPEforGPU::LayerGeometryT<TrackerTraits>::layerStart));
0381 memcpy(layerGeometry_.layer, pixelTopology::layer<TrackerTraits>.data(), pixelTopology::layer<TrackerTraits>.size());
0382 layerGeometry_.maxModuleStride = pixelTopology::maxModuleStride<TrackerTraits>;
0383 }
0384
0385 template <typename TrackerTraits>
0386 PixelCPEFast<TrackerTraits>::GPUData::~GPUData() {
0387 if (paramsOnGPU_d != nullptr) {
0388 cudaFree((void*)paramsOnGPU_h.m_commonParams);
0389 cudaFree((void*)paramsOnGPU_h.m_detParams);
0390 cudaFree((void*)paramsOnGPU_h.m_averageGeometry);
0391 cudaFree((void*)paramsOnGPU_h.m_layerGeometry);
0392 cudaFree(paramsOnGPU_d);
0393 }
0394 }
0395
0396 template <typename TrackerTraits>
0397 void PixelCPEFast<TrackerTraits>::errorFromTemplates(DetParam const& theDetParam,
0398 ClusterParamGeneric& theClusterParam,
0399 float qclus) const {
0400 float locBz = theDetParam.bz;
0401 float locBx = theDetParam.bx;
0402 LogDebug("PixelCPEFast") << "PixelCPEFast::localPosition(...) : locBz = " << locBz;
0403
0404 theClusterParam.pixmx = std::numeric_limits<int>::max();
0405
0406 theClusterParam.sigmay = -999.9;
0407 theClusterParam.sigmax = -999.9;
0408 theClusterParam.sy1 = -999.9;
0409 theClusterParam.sy2 = -999.9;
0410 theClusterParam.sx1 = -999.9;
0411 theClusterParam.sx2 = -999.9;
0412
0413 float dummy;
0414
0415 SiPixelGenError gtempl(thePixelGenError_);
0416 int gtemplID = theDetParam.detTemplateId;
0417
0418 theClusterParam.qBin_ = gtempl.qbin(gtemplID,
0419 theClusterParam.cotalpha,
0420 theClusterParam.cotbeta,
0421 locBz,
0422 locBx,
0423 qclus,
0424 false,
0425 theClusterParam.pixmx,
0426 theClusterParam.sigmay,
0427 dummy,
0428 theClusterParam.sigmax,
0429 dummy,
0430 theClusterParam.sy1,
0431 dummy,
0432 theClusterParam.sy2,
0433 dummy,
0434 theClusterParam.sx1,
0435 dummy,
0436 theClusterParam.sx2,
0437 dummy);
0438
0439 theClusterParam.sigmax = theClusterParam.sigmax * micronsToCm;
0440 theClusterParam.sx1 = theClusterParam.sx1 * micronsToCm;
0441 theClusterParam.sx2 = theClusterParam.sx2 * micronsToCm;
0442
0443 theClusterParam.sigmay = theClusterParam.sigmay * micronsToCm;
0444 theClusterParam.sy1 = theClusterParam.sy1 * micronsToCm;
0445 theClusterParam.sy2 = theClusterParam.sy2 * micronsToCm;
0446 }
0447
0448 template <>
0449 void PixelCPEFast<pixelTopology::Phase2>::errorFromTemplates(DetParam const& theDetParam,
0450 ClusterParamGeneric& theClusterParam,
0451 float qclus) const {
0452 theClusterParam.qBin_ = 0.0f;
0453 }
0454
0455
0456
0457
0458
0459
0460 template <typename TrackerTraits>
0461 LocalPoint PixelCPEFast<TrackerTraits>::localPosition(DetParam const& theDetParam,
0462 ClusterParam& theClusterParamBase) const {
0463 ClusterParamGeneric& theClusterParam = static_cast<ClusterParamGeneric&>(theClusterParamBase);
0464
0465 if (useErrorsFromTemplates_) {
0466 errorFromTemplates(theDetParam, theClusterParam, theClusterParam.theCluster->charge());
0467 } else {
0468 theClusterParam.qBin_ = 0;
0469 }
0470
0471 int q_f_X;
0472 int q_l_X;
0473 int q_f_Y;
0474 int q_l_Y;
0475 collect_edge_charges(theClusterParam, q_f_X, q_l_X, q_f_Y, q_l_Y, useErrorsFromTemplates_ && truncatePixelCharge_);
0476
0477
0478 pixelCPEforGPU::ClusParams cp;
0479
0480 cp.minRow[0] = theClusterParam.theCluster->minPixelRow();
0481 cp.maxRow[0] = theClusterParam.theCluster->maxPixelRow();
0482 cp.minCol[0] = theClusterParam.theCluster->minPixelCol();
0483 cp.maxCol[0] = theClusterParam.theCluster->maxPixelCol();
0484
0485 cp.q_f_X[0] = q_f_X;
0486 cp.q_l_X[0] = q_l_X;
0487 cp.q_f_Y[0] = q_f_Y;
0488 cp.q_l_Y[0] = q_l_Y;
0489
0490 cp.charge[0] = theClusterParam.theCluster->charge();
0491
0492 auto ind = theDetParam.theDet->index();
0493 pixelCPEforGPU::position<TrackerTraits>(commonParamsGPU_, detParamsGPU_[ind], cp, 0);
0494 auto xPos = cp.xpos[0];
0495 auto yPos = cp.ypos[0];
0496
0497
0498 pixelCPEforGPU::errorFromDB<TrackerTraits>(commonParamsGPU_, detParamsGPU_[ind], cp, 0);
0499 theClusterParam.sigmax = cp.xerr[0];
0500 theClusterParam.sigmay = cp.yerr[0];
0501
0502 LogDebug("PixelCPEFast") << " in PixelCPEFast:localPosition - pos = " << xPos << " " << yPos << " size "
0503 << cp.maxRow[0] - cp.minRow[0] << ' ' << cp.maxCol[0] - cp.minCol[0];
0504
0505
0506 LocalPoint pos_in_local(xPos, yPos);
0507 return pos_in_local;
0508 }
0509
0510
0511
0512
0513
0514
0515 template <typename TrackerTraits>
0516 LocalError PixelCPEFast<TrackerTraits>::localError(DetParam const& theDetParam,
0517 ClusterParam& theClusterParamBase) const {
0518 ClusterParamGeneric& theClusterParam = static_cast<ClusterParamGeneric&>(theClusterParamBase);
0519
0520 auto xerr = theClusterParam.sigmax;
0521 auto yerr = theClusterParam.sigmay;
0522
0523 LogDebug("PixelCPEFast") << " errors " << xerr << " " << yerr;
0524
0525 auto xerr_sq = xerr * xerr;
0526 auto yerr_sq = yerr * yerr;
0527
0528 return LocalError(xerr_sq, 0, yerr_sq);
0529 }
0530
0531 template <typename TrackerTraits>
0532 void PixelCPEFast<TrackerTraits>::fillPSetDescription(edm::ParameterSetDescription& desc) {
0533
0534 PixelCPEGenericBase::fillPSetDescription(desc);
0535 }
0536
0537 template class PixelCPEFast<pixelTopology::Phase1>;
0538 template class PixelCPEFast<pixelTopology::Phase2>;
0539 template class PixelCPEFast<pixelTopology::HIonPhase1>;