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