CMS 3D CMS Logo

VectorHitBuilderAlgorithm.cc
Go to the documentation of this file.
6 
8  VectorHitCollection& vhAcc,
9  VectorHitCollection& vhRej,
12  LogDebug("VectorHitBuilderAlgorithm") << "Run VectorHitBuilderAlgorithm ... \n";
13  const auto* clustersPhase2Collection = clusters.product();
14 
15  //loop over the DetSetVector
16  LogDebug("VectorHitBuilderAlgorithm") << "with #clusters : " << clustersPhase2Collection->size() << std::endl;
17  for (auto dSViter : *clustersPhase2Collection) {
18  unsigned int rawDetId1(dSViter.detId());
19  DetId detId1(rawDetId1);
20  DetId lowerDetId, upperDetId;
21  if (tkTopo_->isLower(detId1)) {
22  lowerDetId = detId1;
23  upperDetId = tkTopo_->partnerDetId(detId1);
24  } else
25  continue;
26 
27  DetId detIdStack = tkTopo_->stack(detId1);
28 
29  //debug
30  LogDebug("VectorHitBuilderAlgorithm") << " DetId stack : " << detIdStack.rawId() << std::endl;
31  LogDebug("VectorHitBuilderAlgorithm") << " DetId lower set of clusters : " << lowerDetId.rawId();
32  LogDebug("VectorHitBuilderAlgorithm") << " DetId upper set of clusters : " << upperDetId.rawId() << std::endl;
33 
34  const GeomDet* gd;
35  const StackGeomDet* stackDet;
36  const auto& it_detLower = dSViter;
37  const auto& it_detUpper = clustersPhase2Collection->find(upperDetId);
38 
39  if (it_detUpper != clustersPhase2Collection->end()) {
40  gd = tkGeom_->idToDet(detIdStack);
41  stackDet = dynamic_cast<const StackGeomDet*>(gd);
42  buildVectorHits(vhAcc, vhRej, detIdStack, stackDet, clusters, it_detLower, *it_detUpper);
43  }
44  }
45  LogDebug("VectorHitBuilderAlgorithm") << "End run VectorHitBuilderAlgorithm ... \n";
46 }
47 
50  const Detset& theLowerDetSet,
51  const Detset& theUpperDetSet) const {
52  if (theLowerDetSet.size() == 1 && theUpperDetSet.size() == 1)
53  return true;
54 
55  //order lower clusters in u
56  std::vector<Phase2TrackerCluster1D> lowerClusters;
57  lowerClusters.reserve(theLowerDetSet.size());
58  if (theLowerDetSet.size() > 1)
59  LogDebug("VectorHitBuilderAlgorithm") << " more than 1 lower cluster! " << std::endl;
60  if (theUpperDetSet.size() > 1)
61  LogDebug("VectorHitBuilderAlgorithm") << " more than 1 upper cluster! " << std::endl;
62  for (const_iterator cil = theLowerDetSet.begin(); cil != theLowerDetSet.end(); ++cil) {
64  lowerClusters.push_back(*clusterLower);
65  }
66  return true;
67 }
68 
70  Local3DPoint& posupper,
71  LocalError& errlower,
72  LocalError& errupper) const {
73  return true;
74 }
75 
76 //----------------------------------------------------------------------------
77 //ERICA::in the DT code the global position is used to compute the alpha angle and put a cut on that.
79  VectorHitCollection& vhRej,
80  DetId detIdStack,
81  const StackGeomDet* stack,
83  const Detset& theLowerDetSet,
84  const Detset& theUpperDetSet,
85  const std::vector<bool>& phase2OTClustersToSkip) const {
86  if (checkClustersCompatibilityBeforeBuilding(clusters, theLowerDetSet, theUpperDetSet)) {
87  LogDebug("VectorHitBuilderAlgorithm") << " compatible -> continue ... " << std::endl;
88  } else {
89  LogTrace("VectorHitBuilderAlgorithm") << " not compatible, going to the next cluster";
90  }
91 
92  edmNew::DetSetVector<VectorHit>::FastFiller vh_colAcc(vhAcc, detIdStack);
93  edmNew::DetSetVector<VectorHit>::FastFiller vh_colRej(vhRej, detIdStack);
94 
95  unsigned int layerStack = tkTopo_->layer(stack->geographicalId());
96  if (stack->subDetector() == GeomDetEnumerators::SubDetector::P2OTB)
97  LogDebug("VectorHitBuilderAlgorithm") << " \t is barrel. " << std::endl;
98  if (stack->subDetector() == GeomDetEnumerators::SubDetector::P2OTEC)
99  LogDebug("VectorHitBuilderAlgorithm") << " \t is endcap. " << std::endl;
100  LogDebug("VectorHitBuilderAlgorithm") << " \t layer is : " << layerStack << std::endl;
101 
102  float cut = 0.0;
103  if (stack->subDetector() == GeomDetEnumerators::SubDetector::P2OTB)
104  cut = barrelCut_.at(layerStack);
105  if (stack->subDetector() == GeomDetEnumerators::SubDetector::P2OTEC)
106  cut = endcapCut_.at(layerStack);
107  LogDebug("VectorHitBuilderAlgorithm") << " \t the cut is:" << cut << std::endl;
108 
109  //only cache local parameters for upper cluster as we loop over lower clusters only once anyway
110  std::vector<std::pair<LocalPoint, LocalError>> localParamsUpper;
111  std::vector<const PixelGeomDetUnit*> localGDUUpper;
112 
113  const PixelGeomDetUnit* gduUpp = dynamic_cast<const PixelGeomDetUnit*>(stack->upperDet());
114  for (auto const& clusterUpper : theUpperDetSet) {
115  localGDUUpper.push_back(gduUpp);
116  localParamsUpper.push_back(cpe_->localParameters(clusterUpper, *gduUpp));
117  }
118  int upperIterator = 0;
119  const PixelGeomDetUnit* gduLow = dynamic_cast<const PixelGeomDetUnit*>(stack->lowerDet());
120  for (const_iterator cil = theLowerDetSet.begin(); cil != theLowerDetSet.end(); ++cil) {
121  LogDebug("VectorHitBuilderAlgorithm") << " lower clusters " << std::endl;
123 #ifdef EDM_ML_DEBUG
124  printCluster(stack->lowerDet(), &*cluL);
125 #endif
126  auto&& lparamsLow = cpe_->localParameters(*cluL, *gduLow);
127  upperIterator = 0;
128  for (const_iterator ciu = theUpperDetSet.begin(); ciu != theUpperDetSet.end(); ++ciu) {
129  LogDebug("VectorHitBuilderAlgorithm") << "\t upper clusters " << std::endl;
131 #ifdef EDM_ML_DEBUG
132  printCluster(stack->upperDet(), &*cluU);
133 #endif
134  //applying the parallax correction
135  double pC = computeParallaxCorrection(
136  gduLow, lparamsLow.first, localGDUUpper[upperIterator], localParamsUpper[upperIterator].first);
137  LogDebug("VectorHitBuilderAlgorithm") << " \t parallax correction:" << pC << std::endl;
138  double lpos_upp_corr = 0.0;
139  double lpos_low_corr = 0.0;
140  auto const localUpperX = localParamsUpper[upperIterator].first.x();
141  if (localUpperX > lparamsLow.first.x()) {
142  if (localUpperX > 0) {
143  lpos_low_corr = lparamsLow.first.x();
144  lpos_upp_corr = localParamsUpper[upperIterator].first.x() - std::abs(pC);
145  } else if (localUpperX < 0) {
146  lpos_low_corr = lparamsLow.first.x() + std::abs(pC);
147  lpos_upp_corr = localUpperX;
148  }
149  } else if (localUpperX < lparamsLow.first.x()) {
150  if (localUpperX > 0) {
151  lpos_low_corr = lparamsLow.first.x() - std::abs(pC);
152  lpos_upp_corr = localUpperX;
153  } else if (localUpperX < 0) {
154  lpos_low_corr = lparamsLow.first.x();
155  lpos_upp_corr = localUpperX + std::abs(pC);
156  }
157  } else {
158  if (localUpperX > 0) {
159  lpos_low_corr = lparamsLow.first.x();
160  lpos_upp_corr = localUpperX - std::abs(pC);
161  } else if (localUpperX < 0) {
162  lpos_low_corr = lparamsLow.first.x();
163  lpos_upp_corr = localUpperX + std::abs(pC);
164  }
165  }
166 
167  LogDebug("VectorHitBuilderAlgorithm") << " \t local pos upper corrected (x):" << lpos_upp_corr << std::endl;
168  LogDebug("VectorHitBuilderAlgorithm") << " \t local pos lower corrected (x):" << lpos_low_corr << std::endl;
169 
170  double width = lpos_low_corr - lpos_upp_corr;
171  LogDebug("VectorHitBuilderAlgorithm") << " \t width: " << width << std::endl;
172 
173  //old cut: indipendent from layer
174  //building my tolerance : 10*sigma
175  //double delta = 10.0 * sqrt(lparamsLow.second.xx() + localParamsUpper[upperIterator].second.xx());
176  //LogDebug("VectorHitBuilderAlgorithm") << " \t delta: " << delta << std::endl;
177  //if( (lpos_upp_corr < lpos_low_corr + delta) &&
178  // (lpos_upp_corr > lpos_low_corr - delta) ){
179  //new cut: dependent on layers
180  if (std::abs(width) < cut) {
181  LogDebug("VectorHitBuilderAlgorithm") << " accepting VH! " << std::endl;
182  VectorHit vh = buildVectorHit(stack, cluL, cluU);
183  //protection: the VH can also be empty!!
184  if (vh.isValid()) {
185  vh_colAcc.push_back(vh);
186  }
187 
188  } else {
189  LogDebug("VectorHitBuilderAlgorithm") << " rejecting VH: " << std::endl;
190  //storing vh rejected for combinatiorial studies
191  VectorHit vh = buildVectorHit(stack, cluL, cluU);
192  if (vh.isValid()) {
193  vh_colRej.push_back(vh);
194  }
195  }
196  upperIterator += 1;
197  }
198  }
199 }
200 
204  LogTrace("VectorHitBuilderAlgorithm") << "Build VH with: ";
205 #ifdef EDM_ML_DEBUG
206  printCluster(stack->upperDet(), &*upper);
207 #endif
208  const PixelGeomDetUnit* geomDetLower = static_cast<const PixelGeomDetUnit*>(stack->lowerDet());
209  const PixelGeomDetUnit* geomDetUpper = static_cast<const PixelGeomDetUnit*>(stack->upperDet());
210 
211  auto&& lparamsLower = cpe_->localParameters(*lower, *geomDetLower); // x, y, z, e2_xx, e2_xy, e2_yy
212  Global3DPoint gparamsLower = geomDetLower->surface().toGlobal(lparamsLower.first);
213  LogTrace("VectorHitBuilderAlgorithm") << "\t lower global pos: " << gparamsLower;
214 
215  auto&& lparamsUpper = cpe_->localParameters(*upper, *geomDetUpper);
216  Global3DPoint gparamsUpper = geomDetUpper->surface().toGlobal(lparamsUpper.first);
217  LogTrace("VectorHitBuilderAlgorithm") << "\t upper global pos: " << gparamsUpper;
218 
219  //local parameters of upper cluster in lower system of reference
220  Local3DPoint lparamsUpperInLower = geomDetLower->surface().toLocal(gparamsUpper);
221 
222  LogTrace("VectorHitBuilderAlgorithm") << "\t lower global pos: " << gparamsLower;
223  LogTrace("VectorHitBuilderAlgorithm") << "\t upper global pos: " << gparamsUpper;
224 
225  LogTrace("VectorHitBuilderAlgorithm") << "A:\t lower local pos: " << lparamsLower.first
226  << " with error: " << lparamsLower.second << std::endl;
227  LogTrace("VectorHitBuilderAlgorithm") << "A:\t upper local pos in the lower sof " << lparamsUpperInLower
228  << " with error: " << lparamsUpper.second << std::endl;
229 
230  bool ok =
231  checkClustersCompatibility(lparamsLower.first, lparamsUpper.first, lparamsLower.second, lparamsUpper.second);
232 
233  if (ok) {
234  AlgebraicSymMatrix22 covMat2Dzx;
235  double chi22Dzx = 0.0;
236  Local3DPoint pos2Dzx;
237  Local3DVector dir2Dzx;
238  fit2Dzx(lparamsLower.first,
239  lparamsUpperInLower,
240  lparamsLower.second,
241  lparamsUpper.second,
242  pos2Dzx,
243  dir2Dzx,
244  covMat2Dzx,
245  chi22Dzx);
246  LogTrace("VectorHitBuilderAlgorithm") << "\t pos2Dzx: " << pos2Dzx;
247  LogTrace("VectorHitBuilderAlgorithm") << "\t dir2Dzx: " << dir2Dzx;
248  LogTrace("VectorHitBuilderAlgorithm") << "\t cov2Dzx: " << covMat2Dzx;
249  VectorHit2D vh2Dzx = VectorHit2D(pos2Dzx, dir2Dzx, covMat2Dzx, chi22Dzx);
250 
251  AlgebraicSymMatrix22 covMat2Dzy;
252  double chi22Dzy = 0.0;
253  Local3DPoint pos2Dzy;
254  Local3DVector dir2Dzy;
255  fit2Dzy(lparamsLower.first,
256  lparamsUpperInLower,
257  lparamsLower.second,
258  lparamsUpper.second,
259  pos2Dzy,
260  dir2Dzy,
261  covMat2Dzy,
262  chi22Dzy);
263  LogTrace("VectorHitBuilderAlgorithm") << "\t pos2Dzy: " << pos2Dzy;
264  LogTrace("VectorHitBuilderAlgorithm") << "\t dir2Dzy: " << dir2Dzy;
265  LogTrace("VectorHitBuilderAlgorithm") << "\t cov2Dzy: " << covMat2Dzy;
266  VectorHit2D vh2Dzy = VectorHit2D(pos2Dzy, dir2Dzy, covMat2Dzy, chi22Dzy);
267 
268  OmniClusterRef lowerOmni(lower);
269  OmniClusterRef upperOmni(upper);
270 
271  Global3DPoint gPositionLower = VectorHit::phase2clusterGlobalPos(geomDetLower, lower);
272  Global3DPoint gPositionUpper = VectorHit::phase2clusterGlobalPos(geomDetUpper, upper);
273  GlobalError gErrorLower = VectorHit::phase2clusterGlobalPosErr(geomDetLower);
274  GlobalError gErrorUpper = VectorHit::phase2clusterGlobalPosErr(geomDetUpper);
275 
276  if (gPositionLower.perp() > gPositionUpper.perp()) {
277  std::swap(gPositionLower, gPositionUpper);
278  std::swap(gErrorLower, gErrorUpper);
279  }
280 
281  const CurvatureAndPhi curvatureAndPhi = curvatureANDphi(gPositionLower, gPositionUpper, gErrorLower, gErrorUpper);
282  VectorHit vh = VectorHit(*stack,
283  vh2Dzx,
284  vh2Dzy,
285  lowerOmni,
286  upperOmni,
287  curvatureAndPhi.curvature,
288  curvatureAndPhi.curvatureError,
289  curvatureAndPhi.phi);
290  return vh;
291  }
292 
293  return VectorHit();
294 }
295 
297  const Local3DPoint lpCO,
298  const LocalError leCI,
299  const LocalError leCO,
300  Local3DPoint& pos,
302  AlgebraicSymMatrix22& covMatrix,
303  double& chi2) const {
304  float x[2] = {lpCI.z(), lpCO.z()};
305  float y[2] = {lpCI.x(), lpCO.x()};
306  float sqCI = sqrt(leCI.xx());
307  float sqCO = sqrt(leCO.xx());
308  float sigy2[2] = {sqCI * sqCI, sqCO * sqCO};
309 
310  fit(x, y, sigy2, pos, dir, covMatrix, chi2);
311 
312  return;
313 }
314 
316  const Local3DPoint lpCO,
317  const LocalError leCI,
318  const LocalError leCO,
319  Local3DPoint& pos,
321  AlgebraicSymMatrix22& covMatrix,
322  double& chi2) const {
323  float x[2] = {lpCI.z(), lpCO.z()};
324  float y[2] = {lpCI.y(), lpCO.y()};
325  float sqCI = sqrt(leCI.yy());
326  float sqCO = sqrt(leCO.yy());
327  float sigy2[2] = {sqCI * sqCI, sqCO * sqCO};
328 
329  fit(x, y, sigy2, pos, dir, covMatrix, chi2);
330 
331  return;
332 }
333 
335  float y[2],
336  float sigy2[2],
337  Local3DPoint& pos,
339  AlgebraicSymMatrix22& covMatrix,
340  double& chi2) const {
341  float slope = 0.;
342  float intercept = 0.;
343  float covss = 0.;
344  float covii = 0.;
345  float covsi = 0.;
346 
347  linearFit(x, y, 2, sigy2, slope, intercept, covss, covii, covsi);
348 
349  covMatrix[0][0] = covss; // this is var(dy/dz)
350  covMatrix[1][1] = covii; // this is var(y)
351  covMatrix[1][0] = covsi; // this is cov(dy/dz,y)
352 
353  for (unsigned int j = 0; j < 2; j++) {
354  const double ypred = intercept + slope * x[j];
355  const double dy = (y[j] - ypred) / sqrt(sigy2[j]);
356  chi2 += dy * dy;
357  }
358 
359  pos = Local3DPoint(intercept, 0., 0.);
360  //difference in z is the difference of the lowermost and the uppermost cluster z pos
361  float slopeZ = x[1] - x[0];
362  dir = LocalVector(slope, 0., slopeZ);
363 }
364 
366  Global3DPoint gPositionUpper,
367  GlobalError gErrorLower,
368  GlobalError gErrorUpper) const {
370 
371  float curvature = -999.;
372  float errorCurvature = -999.;
373  float phi = -999.;
374 
375  float h1 = gPositionLower.x() * gPositionUpper.y() - gPositionUpper.x() * gPositionLower.y();
376 
377  //determine sign of curvature
378  AlgebraicVector2 n1;
379  n1[0] = -gPositionLower.y();
380  n1[1] = gPositionLower.x();
381  AlgebraicVector2 n2;
382  n2[0] = gPositionUpper.x() - gPositionLower.x();
383  n2[1] = gPositionUpper.y() - gPositionLower.y();
384 
385  double n3 = n1[0] * n2[0] + n1[1] * n2[1];
386  double signCurv = -copysign(1.0, n3);
387  double phi1 = atan2(gPositionUpper.y() - gPositionLower.y(), gPositionUpper.x() - gPositionLower.x());
388 
389  double x2Low = pow(gPositionLower.x(), 2);
390  double y2Low = pow(gPositionLower.y(), 2);
391  double x2Up = pow(gPositionUpper.x(), 2);
392  double y2Up = pow(gPositionUpper.y(), 2);
393 
394  if (h1 != 0) {
395  double h2 = 2 * h1;
396  double h2Inf = 1. / (2 * h1);
397  double r12 = gPositionLower.perp2();
398  double r22 = gPositionUpper.perp2();
399  double h3 = pow(n2[0], 2) + pow(n2[1], 2);
400  double h4 = -x2Low * gPositionUpper.x() + gPositionLower.x() * x2Up + gPositionLower.x() * y2Up -
401  gPositionUpper.x() * y2Low;
402  double h5 =
403  x2Low * gPositionUpper.y() - x2Up * gPositionLower.y() + y2Low * gPositionUpper.y() - gPositionLower.y() * y2Up;
404 
405  //radius of circle
406  double invRho2 = (4. * h1 * h1) / (r12 * r22 * h3);
407  curvature = sqrt(invRho2);
408 
409  //center of circle
410  double xcentre = h5 / h2;
411  double ycentre = h4 / h2;
412 
413  //to compute phi at the cluster points
414  double xtg = gPositionLower.y() - ycentre;
415  double ytg = -(gPositionLower.x() - xcentre);
416 
417  //to compute phi at the origin
418  phi = atan2(ytg, xtg);
419 
421 
422  double denom1 = 1. / sqrt(r12 * r22 * h3);
423  double denom2 = 1. / (pow(r12 * r22 * h3, 1.5));
424  jacobian[0][0] = 1.0; // dx1/dx1 dx1/dy1 dx2/dx1 dy2/dx1
425  jacobian[1][1] = 1.0; //dy1/dx1 dy1/dy1 dy2/dx1 dy2/dx1
426  jacobian[2][0] =
427  -2. * ((h1 * (gPositionLower.x() * r22 * h3 + (gPositionLower.x() - gPositionUpper.x()) * r12 * r22)) * denom2 -
428  (gPositionUpper.y()) * denom1); // dkappa/dx1
429  jacobian[2][1] =
430  -2. * ((gPositionUpper.x()) * denom1 +
431  (h1 * (gPositionLower.y() * r22 * h3 + r12 * r22 * (gPositionLower.y() - gPositionUpper.y()))) *
432  denom2); // dkappa/dy1
433  jacobian[2][2] =
434  -2. * ((gPositionLower.y()) * denom1 +
435  (h1 * (gPositionUpper.x() * r12 * h3 - (gPositionLower.x() - gPositionUpper.x()) * r12 * r22)) *
436  denom2); // dkappa/dx2
437  jacobian[2][3] =
438  -2. * ((h1 * (gPositionUpper.y() * r12 * h3 - r12 * r22 * (gPositionLower.y() - gPositionUpper.y()))) * denom2 -
439  (gPositionLower.x()) * denom1); // dkappa/dy2
440  AlgebraicVector2 mVector;
441  //to compute phi at the cluster points
442  mVector[0] = (gPositionLower.y() - ycentre) * invRho2; // dphi/dxcentre
443  mVector[1] = -(gPositionLower.x() - xcentre) * invRho2; // dphi/dycentre
444  //to compute phi at the origin
445 
446  double h22Inv = 1. / pow(h2, 2);
447 
449  kMatrix[0][0] =
450  2. * ((gPositionLower.x() * gPositionUpper.y()) * h2Inf - (gPositionUpper.y() * h5) * h22Inv); // dxm/dx1
451  kMatrix[0][1] = (2. * gPositionUpper.x() * h5) * h22Inv -
452  (x2Up + y2Up - 2. * gPositionLower.y() * gPositionUpper.y()) * h2Inf; // dxm/dy1
453  kMatrix[0][2] =
454  2. * ((gPositionLower.y() * h5) * h22Inv - (gPositionUpper.x() * gPositionLower.y()) * h2Inf); // dxm/dx2
455  kMatrix[0][3] = (x2Low + y2Low - 2. * gPositionUpper.y() * gPositionLower.y()) * h2Inf -
456  (2. * gPositionLower.x() * h5) * h22Inv; // dxm/dy2
457  kMatrix[1][0] = (x2Up - 2. * gPositionLower.x() * gPositionUpper.x() + y2Up) * h2Inf -
458  (2. * gPositionUpper.y() * h4) * h22Inv; // dym/dx1
459  kMatrix[1][1] =
460  2. * ((gPositionUpper.x() * h4) * h22Inv - (gPositionUpper.x() * gPositionLower.y()) * h2Inf); // dym/dy1
461  kMatrix[1][2] = (2. * gPositionLower.y() * h4) * h22Inv -
462  (x2Low - 2. * gPositionUpper.x() * gPositionLower.x() + y2Low) * h2Inf; // dym/dx2
463  kMatrix[1][3] =
464  2. * (gPositionLower.x() * gPositionUpper.y()) * h2Inf - (gPositionLower.x() * h4) * h22Inv; // dym/dy2
465 
466  AlgebraicVector4 nMatrix = mVector * kMatrix;
467  jacobian[3][0] = nMatrix[0]; // dphi/(dx1,dy1,dx2,dy2)
468  jacobian[3][1] = nMatrix[1]; // dphi/(dx1,dy1,dx2,dy2)
469  jacobian[3][2] = nMatrix[2]; // dphi/(dx1,dy1,dx2,dy2)
470  jacobian[3][3] = nMatrix[3]; // dphi/(dx1,dy1,dx2,dy2)
471 
472  //assign correct sign to the curvature errors
473  if ((signCurv < 0 && curvature > 0) || (signCurv > 0 && curvature < 0)) {
474  curvature = -curvature;
475  for (int i = 0; i < 4; i++) {
476  jacobian[2][i] = -jacobian[2][i];
477  }
478  }
479 
480  // bring phi in the same quadrant as phi1
481  if (deltaPhi(phi, phi1) > M_PI / 2.) {
482  phi = phi + M_PI;
483  if (phi > M_PI)
484  phi = phi - 2. * M_PI;
485  }
486 
487  //computing the curvature error
488  AlgebraicVector4 curvatureJacobian;
489  for (int i = 0; i < 4; i++) {
490  curvatureJacobian[i] = jacobian[2][i];
491  }
492 
494 
495  gErrors[0][0] = gErrorLower.cxx();
496  gErrors[0][1] = gErrorLower.cyx();
497  gErrors[1][0] = gErrorLower.cyx();
498  gErrors[1][1] = gErrorLower.cyy();
499  gErrors[2][2] = gErrorUpper.cxx();
500  gErrors[2][3] = gErrorUpper.cyx();
501  gErrors[3][2] = gErrorUpper.cyx();
502  gErrors[3][3] = gErrorUpper.cyy();
503 
504  AlgebraicVector4 temp = curvatureJacobian;
505  temp = temp * gErrors;
506  errorCurvature = temp[0] * curvatureJacobian[0] + temp[1] * curvatureJacobian[1] + temp[2] * curvatureJacobian[2] +
507  temp[3] * curvatureJacobian[3];
508  }
509 
510  result.curvature = curvature;
511  result.curvatureError = errorCurvature;
512  result.phi = phi;
513  return result;
514 }
Vector3DBase< float, LocalTag >
TrackerGeometry::idToDet
const TrackerGeomDet * idToDet(DetId) const override
Definition: TrackerGeometry.cc:193
DDAxes::y
VectorHitBuilderAlgorithm::buildVectorHit
VectorHit buildVectorHit(const StackGeomDet *stack, Phase2TrackerCluster1DRef lower, Phase2TrackerCluster1DRef upper) const override
Definition: VectorHitBuilderAlgorithm.cc:201
ApeEstimator_cff.width
width
Definition: ApeEstimator_cff.py:24
VectorHitBuilderAlgorithm::checkClustersCompatibility
bool checkClustersCompatibility(Local3DPoint &posinner, Local3DPoint &posouter, LocalError &errinner, LocalError &errouter) const
Definition: VectorHitBuilderAlgorithm.cc:69
mps_fire.i
i
Definition: mps_fire.py:428
VectorHitBuilderAlgorithm::CurvatureAndPhi::curvatureError
float curvatureError
Definition: VectorHitBuilderAlgorithm.h:42
GeomDet
Definition: GeomDet.h:27
PixelTopology.h
TkAlMuonSelectors_cfi.cut
cut
Definition: TkAlMuonSelectors_cfi.py:5
VectorHitBuilderAlgorithmBase::tkTopo_
const TrackerTopology * tkTopo_
Definition: VectorHitBuilderAlgorithmBase.h:61
TrackerTopology::isLower
bool isLower(const DetId &id) const
Definition: TrackerTopology.cc:195
PV3DBase::x
T x() const
Definition: PV3DBase.h:59
VectorHitBuilderAlgorithm.h
VectorHitBuilderAlgorithm::fit2Dzy
void fit2Dzy(const Local3DPoint lpCI, const Local3DPoint lpCO, const LocalError leCI, const LocalError leCO, Local3DPoint &pos, Local3DVector &dir, AlgebraicSymMatrix22 &covMatrix, double &chi2) const
Definition: VectorHitBuilderAlgorithm.cc:315
pos
Definition: PixelAliasList.h:18
edmNew::makeRefTo
edm::Ref< typename HandleT::element_type, typename HandleT::element_type::value_type::value_type > makeRefTo(const HandleT &iHandle, typename HandleT::element_type::value_type::const_iterator itIter)
Definition: DetSetVectorNew.h:689
GeomDetEnumerators::P2OTB
Definition: GeomDetEnumerators.h:23
VectorHitBuilderAlgorithm::fit
void fit(float x[2], float y[2], float sigy[2], Local3DPoint &pos, Local3DVector &dir, AlgebraicSymMatrix22 &covMatrix, double &chi2) const
Definition: VectorHitBuilderAlgorithm.cc:334
VectorHitBuilderAlgorithm::checkClustersCompatibilityBeforeBuilding
bool checkClustersCompatibilityBeforeBuilding(edm::Handle< edmNew::DetSetVector< Phase2TrackerCluster1D >> clusters, const Detset &theLowerDetSet, const Detset &theUpperDetSet) const
Definition: VectorHitBuilderAlgorithm.cc:48
TrackerTopology::layer
unsigned int layer(const DetId &id) const
Definition: TrackerTopology.cc:47
TrackerTopology::stack
uint32_t stack(const DetId &id) const
Definition: TrackerTopology.cc:104
DDAxes::x
edmNew::DetSetVector::FastFiller::push_back
void push_back(data_type const &d)
Definition: DetSetVectorNew.h:274
VectorHitBuilderAlgorithmBase::cpe_
const ClusterParameterEstimator< Phase2TrackerCluster1D > * cpe_
Definition: VectorHitBuilderAlgorithmBase.h:62
hltPixelTracks_cff.chi2
chi2
Definition: hltPixelTracks_cff.py:25
PixelRecoUtilities::curvature
T curvature(T InversePt, const edm::EventSetup &iSetup)
Definition: PixelRecoUtilities.h:42
convertSQLiteXML.ok
bool ok
Definition: convertSQLiteXML.py:98
edm::Handle
Definition: AssociativeIterator.h:50
groupFilesInBlocks.temp
list temp
Definition: groupFilesInBlocks.py:142
OmniClusterRef
Definition: OmniClusterRef.h:12
VectorHitBuilderAlgorithm::curvatureANDphi
CurvatureAndPhi curvatureANDphi(Global3DPoint gPositionLower, Global3DPoint gPositionUpper, GlobalError gErrorLower, GlobalError gErrorUpper) const
Definition: VectorHitBuilderAlgorithm.cc:365
edm::Ref
Definition: AssociativeIterator.h:58
edmNew::DetSet::size
size_type size() const
Definition: DetSetNew.h:68
PV3DBase::z
T z() const
Definition: PV3DBase.h:61
DetId
Definition: DetId.h:17
edmNew::DetSet::end
iterator end()
Definition: DetSetNew.h:56
GlobalErrorBase::cyy
T cyy() const
Definition: GlobalErrorBase.h:101
VectorHitBuilderAlgorithm::CurvatureAndPhi::phi
float phi
Definition: VectorHitBuilderAlgorithm.h:43
std::swap
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
Definition: DataFrameContainer.h:209
SiPixelRawToDigiRegional_cfi.deltaPhi
deltaPhi
Definition: SiPixelRawToDigiRegional_cfi.py:9
LocalError::xx
float xx() const
Definition: LocalError.h:22
PixelGeomDetUnit
Definition: PixelGeomDetUnit.h:15
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
linearFit
void linearFit(T const *__restrict__ x, T const *__restrict__ y, int ndat, T const *__restrict__ sigy2, T &slope, T &intercept, T &covss, T &covii, T &covsi)
Definition: LinearFit.h:29
edmNew::DetSet
Definition: DetSetNew.h:22
TrackerTopology::partnerDetId
DetId partnerDetId(const DetId &id) const
Definition: TrackerTopology.cc:233
GeomDetEnumerators::P2OTEC
Definition: GeomDetEnumerators.h:24
Point3DBase< float, LocalTag >
GlobalErrorBase::cxx
T cxx() const
Definition: GlobalErrorBase.h:97
VectorHitBuilderAlgorithm::CurvatureAndPhi
Definition: VectorHitBuilderAlgorithm.h:40
VectorHitBuilderAlgorithmBase::printCluster
void printCluster(const GeomDet *geomDetUnit, const Phase2TrackerCluster1D *cluster) const
Definition: VectorHitBuilderAlgorithmBase.cc:65
bsc_activity_cfg.clusters
clusters
Definition: bsc_activity_cfg.py:36
svgfig.stack
stack
Definition: svgfig.py:559
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:223
AlgebraicROOTObject::Matrix
ROOT::Math::SMatrix< double, D1, D2, ROOT::Math::MatRepStd< double, D1, D2 > > Matrix
Definition: AlgebraicROOTObjects.h:69
VectorHitBuilderAlgorithmBase::const_iterator
Detset::const_iterator const_iterator
Definition: VectorHitBuilderAlgorithmBase.h:22
GlobalErrorBase::cyx
T cyx() const
Definition: GlobalErrorBase.h:99
LocalError
Definition: LocalError.h:12
AlgebraicVector2
ROOT::Math::SVector< double, 2 > AlgebraicVector2
Definition: AlgebraicROOTObjects.h:11
VectorHit2D
Definition: VectorHit2D.h:11
VectorHitBuilderAlgorithmBase::computeParallaxCorrection
double computeParallaxCorrection(const PixelGeomDetUnit *, const Point3DBase< float, LocalTag > &, const PixelGeomDetUnit *, const Point3DBase< float, LocalTag > &) const
Definition: VectorHitBuilderAlgorithmBase.cc:23
PV3DBase::y
T y() const
Definition: PV3DBase.h:60
LocalVector
Local3DVector LocalVector
Definition: LocalVector.h:12
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:49
ClusterParameterEstimator::localParameters
virtual LocalValues localParameters(const T &, const GeomDetUnit &) const =0
edmNew::DetSet::begin
iterator begin()
Definition: DetSetNew.h:54
GlobalErrorBase< double, ErrorMatrixTag >
PVValHelper::dy
Definition: PVValidationHelpers.h:49
VectorHit::phase2clusterGlobalPosErr
static GlobalError phase2clusterGlobalPosErr(const PixelGeomDetUnit *geomDet)
Definition: VectorHit.cc:135
VectorHitBuilderAlgorithm::CurvatureAndPhi::curvature
float curvature
Definition: VectorHitBuilderAlgorithm.h:41
VectorHitBuilderAlgorithmBase::tkGeom_
const TrackerGeometry * tkGeom_
Definition: VectorHitBuilderAlgorithmBase.h:60
DDAxes::phi
VectorHit
Definition: VectorHit.h:28
VectorHit::phase2clusterGlobalPos
static Global3DPoint phase2clusterGlobalPos(const PixelGeomDetUnit *geomDet, ClusterRef cluster)
Definition: VectorHit.cc:114
GeomDet.h
edmNew::DetSetVector
Definition: DetSetNew.h:13
StackGeomDet
Definition: StackGeomDet.h:7
AlgebraicVector4
ROOT::Math::SVector< double, 4 > AlgebraicVector4
Definition: AlgebraicROOTObjects.h:13
VectorHitBuilderAlgorithmBase::barrelCut_
std::vector< double > barrelCut_
Definition: VectorHitBuilderAlgorithmBase.h:64
DetId::rawId
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
VectorHitBuilderAlgorithm::buildVectorHits
void buildVectorHits(VectorHitCollection &vhAcc, VectorHitCollection &vhRej, DetId detIdStack, const StackGeomDet *stack, edm::Handle< edmNew::DetSetVector< Phase2TrackerCluster1D >> clusters, const Detset &DSVinner, const Detset &DSVouter, const std::vector< bool > &phase2OTClustersToSkip=std::vector< bool >()) const override
Definition: VectorHitBuilderAlgorithm.cc:78
genVertex_cff.x
x
Definition: genVertex_cff.py:12
detailsBasic3DVector::y
float float y
Definition: extBasic3DVector.h:14
pileupCalc.upper
upper
Definition: pileupCalc.py:214
VectorHitBuilderAlgorithm::run
void run(edm::Handle< edmNew::DetSetVector< Phase2TrackerCluster1D >> clusters, VectorHitCollection &vhAcc, VectorHitCollection &vhRej, edmNew::DetSetVector< Phase2TrackerCluster1D > &clustersAcc, edmNew::DetSetVector< Phase2TrackerCluster1D > &clustersRej) const override
Definition: VectorHitBuilderAlgorithm.cc:7
VectorHit2D.h
edmNew::DetSetVector::FastFiller
Definition: DetSetVectorNew.h:202
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
VectorHitBuilderAlgorithm::fit2Dzx
void fit2Dzx(const Local3DPoint lpCI, const Local3DPoint lpCO, const LocalError leCI, const LocalError leCO, Local3DPoint &pos, Local3DVector &dir, AlgebraicSymMatrix22 &covMatrix, double &chi2) const
Definition: VectorHitBuilderAlgorithm.cc:296
mps_fire.result
result
Definition: mps_fire.py:311
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ClusterParameterEstimator.h
LogTrace
#define LogTrace(id)
Definition: MessageLogger.h:224
VectorHitBuilderAlgorithmBase::endcapCut_
std::vector< double > endcapCut_
Definition: VectorHitBuilderAlgorithmBase.h:65
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
slope
static const double slope[3]
Definition: CastorTimeSlew.cc:6
TrackingRecHit::isValid
bool isValid() const
Definition: TrackingRecHit.h:141
LocalError::yy
float yy() const
Definition: LocalError.h:24
PV3DBase::perp
T perp() const
Definition: PV3DBase.h:69
AlgebraicSymMatrix22
ROOT::Math::SMatrix< double, 2, 2, ROOT::Math::MatRepSym< double, 2 > > AlgebraicSymMatrix22
Definition: AlgebraicROOTObjects.h:20
PV3DBase::perp2
T perp2() const
Definition: PV3DBase.h:68
Local3DPoint
Point3DBase< float, LocalTag > Local3DPoint
Definition: LocalPoint.h:9
DeadROC_duringRun.dir
dir
Definition: DeadROC_duringRun.py:23