File indexing completed on 2025-01-30 06:08:07
0001 #ifndef RecoTracker_LSTCore_interface_Circle_h
0002 #define RecoTracker_LSTCore_interface_Circle_h
0003
0004 #include <tuple>
0005
0006 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0007
0008 namespace lst {
0009
0010 template <typename TAcc>
0011 ALPAKA_FN_ACC ALPAKA_FN_INLINE std::tuple<float, float, float> computeRadiusFromThreeAnchorHits(
0012 TAcc const& acc, float x1, float y1, float x2, float y2, float x3, float y3) {
0013 float radius = 0.f;
0014
0015
0016
0017 float denomInv = 1.0f / ((y1 - y3) * (x2 - x3) - (x1 - x3) * (y2 - y3));
0018
0019 float xy1sqr = x1 * x1 + y1 * y1;
0020
0021 float xy2sqr = x2 * x2 + y2 * y2;
0022
0023 float xy3sqr = x3 * x3 + y3 * y3;
0024
0025 float regressionCenterX = 0.5f * ((y3 - y2) * xy1sqr + (y1 - y3) * xy2sqr + (y2 - y1) * xy3sqr) * denomInv;
0026
0027 float regressionCenterY = 0.5f * ((x2 - x3) * xy1sqr + (x3 - x1) * xy2sqr + (x1 - x2) * xy3sqr) * denomInv;
0028
0029 float c = ((x2 * y3 - x3 * y2) * xy1sqr + (x3 * y1 - x1 * y3) * xy2sqr + (x1 * y2 - x2 * y1) * xy3sqr) * denomInv;
0030
0031 if (((y1 - y3) * (x2 - x3) - (x1 - x3) * (y2 - y3) == 0) ||
0032 (regressionCenterX * regressionCenterX + regressionCenterY * regressionCenterY - c < 0)) {
0033 #ifdef WARNINGS
0034 printf("three collinear points or FATAL! r^2 < 0!\n");
0035 #endif
0036 radius = -1.f;
0037 } else
0038 radius =
0039 alpaka::math::sqrt(acc, regressionCenterX * regressionCenterX + regressionCenterY * regressionCenterY - c);
0040
0041 return std::make_tuple(radius, regressionCenterX, regressionCenterY);
0042 }
0043
0044 }
0045
0046 #endif