ATLAS Offline Software
EventPrimitivesCovarianceHelpers.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 // EventPrimitivesHelpers.h, (c) ATLAS Detector software
8 
9 #ifndef EVENTPRIMITIVES_EVENTPRIMITIVESCOVARIANCEHELPERS_H
10 #define EVENTPRIMITIVES_EVENTPRIMITIVESCOVARIANCEHELPERS_H
11 
13 //
14 #include <Eigen/Cholesky>
15 //
16 #include <limits>
22 namespace Amg {
61 
63 inline bool saneCovarianceElement(double ele) {
64  constexpr double upper_covariance_cutoff = std::numeric_limits<float>::max();
65  return !(std::isnan(ele) || std::isinf(ele) ||
66  std::abs(ele) > upper_covariance_cutoff);
67 }
68 
72 template <int N>
74  constexpr int dim = N;
75  for (int i = 0; i < dim; ++i) {
76  if (mat(i, i) < 0.0 || !saneCovarianceElement(mat(i, i)))
77  return false;
78  }
79  return true;
80 }
81 
83  int dim = mat.rows();
84  for (int i = 0; i < dim; ++i) {
85  if (mat(i, i) < 0.0 || !saneCovarianceElement(mat(i, i)))
86  return false;
87  }
88  return true;
89 }
90 
95 template <int N>
96 inline bool hasPositiveDiagElems(const AmgSymMatrix(N) & mat) {
97  constexpr double MIN_COV_EPSILON = std::numeric_limits<float>::min();
98  constexpr int dim = N;
99  for (int i = 0; i < dim; ++i) {
100  if (mat(i, i) < MIN_COV_EPSILON || !saneCovarianceElement(mat(i, i)))
101  return false;
102  }
103  return true;
104 }
105 inline bool hasPositiveDiagElems(const Amg::MatrixX& mat) {
106  constexpr double MIN_COV_EPSILON = std::numeric_limits<float>::min();
107  int dim = mat.rows();
108  for (int i = 0; i < dim; ++i) {
109  if (mat(i, i) < MIN_COV_EPSILON || !saneCovarianceElement(mat(i, i)))
110  return false;
111  }
112  return true;
113 }
114 
117 template <int N>
118 inline bool isPositiveSemiDefinite(const AmgSymMatrix(N) & mat) {
120  // fast check for necessary condition
121  return false;
122  }
123  Eigen::LDLT<AmgSymMatrix(N)> ldltCov(mat);
124  return (ldltCov.info() == Eigen::Success && ldltCov.isPositive());
125 }
128  // fast check for necessary condition
129  return false;
130  }
131  Eigen::LDLT<Amg::MatrixX> ldltCov(mat);
132  return (ldltCov.info() == Eigen::Success && ldltCov.isPositive());
133 }
134 
137 template <int N>
138 inline bool isPositiveDefinite(const AmgSymMatrix(N) & mat) {
139  if (!hasPositiveDiagElems(mat)) {
140  // fast check for necessary condition
141  return false;
142  }
143  Eigen::LLT<AmgSymMatrix(N)> lltCov(mat);
144  return (lltCov.info() == Eigen::Success);
145 }
146 inline bool isPositiveDefinite(const Amg::MatrixX& mat) {
147  if (!hasPositiveDiagElems(mat)) {
148  // fast check for necessary condition
149  return false;
150  }
151  Eigen::LLT<Amg::MatrixX> lltCov(mat);
152  return (lltCov.info() == Eigen::Success);
153 }
154 
157 template <int N>
160  // fast check for necessary condition
161  return false;
162  }
163  Eigen::SelfAdjointEigenSolver<AmgSymMatrix(5)> eigensolver(mat);
164  auto res = eigensolver.eigenvalues();
165  for (size_t i = 0; i < N; ++i) {
166  if (res[i] < 0) {
167  return false;
168  }
169  }
170  return true;
171 }
174  // fast check for necessary condition
175  return false;
176  }
177  Eigen::SelfAdjointEigenSolver<Amg::MatrixX> eigensolver(mat);
178  auto res = eigensolver.eigenvalues();
179  int dim = mat.rows();
180  for (int i = 0; i < dim; ++i) {
181  if (res[i] < 0) {
182  return false;
183  }
184  }
185  return true;
186 }
187 template <int N>
188 inline bool isPositiveDefiniteSlow(const AmgSymMatrix(N) & mat) {
189  if (!hasPositiveDiagElems(mat)) {
190  // fast check for necessary condition
191  return false;
192  }
193  Eigen::SelfAdjointEigenSolver<AmgSymMatrix(5)> eigensolver(mat);
194  auto res = eigensolver.eigenvalues();
195  for (size_t i = 0; i < N; ++i) {
196  if (res[i] <= 0) {
197  return false;
198  }
199  }
200  return true;
201 }
203  if (!hasPositiveDiagElems(mat)) {
204  // fast check for necessary condition
205  return false;
206  }
207  Eigen::SelfAdjointEigenSolver<Amg::MatrixX> eigensolver(mat);
208  auto res = eigensolver.eigenvalues();
209  int dim = mat.rows();
210  for (int i = 0; i < dim; ++i) {
211  if (res[i] <= 0) {
212  return false;
213  }
214  }
215  return true;
216 }
217 // Chi squared statistic of a linear regression as defined in Frühwirth,
218 // section 3.2.1, equation 3.19. Precision is defined as the inverse of
219 // of the covariance.
220 template <typename T, typename U>
221 inline double chi2(const T& precision, const U& residual, const int sign = 1) {
222  return sign * residual.transpose() * precision * residual;
223 }
224 
225 } // namespace Amg
226 
227 #endif
Amg::hasPositiveDiagElems
bool hasPositiveDiagElems(const AmgSymMatrix(N) &mat)
Returns true if all diagonal elements of the covariance matrix are finite aka sane in the above defin...
Definition: EventPrimitivesCovarianceHelpers.h:96
Amg::MatrixX
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Dynamic Matrix - dynamic allocation.
Definition: EventPrimitives.h:29
yodamerge_tmp.dim
dim
Definition: yodamerge_tmp.py:239
max
#define max(a, b)
Definition: cfImp.cxx:41
ClusterSeg::residual
@ residual
Definition: ClusterNtuple.h:20
Amg::isPositiveDefinite
bool isPositiveDefinite(const AmgSymMatrix(N) &mat)
Check if is positive semidefinit using that fact that is needed for Cholesky decomposition.
Definition: EventPrimitivesCovarianceHelpers.h:138
mat
GeoMaterial * mat
Definition: LArDetectorConstructionTBEC.cxx:53
Amg::saneCovarianceElement
bool saneCovarianceElement(double ele)
A covariance matrix formally needs to be positive semi definite.
Definition: EventPrimitivesCovarianceHelpers.h:63
JetTiledMap::N
@ N
Definition: TiledEtaPhiMap.h:44
AmgSymMatrix
#define AmgSymMatrix(dim)
Definition: EventPrimitives.h:52
lumiFormat.i
int i
Definition: lumiFormat.py:92
Amg::isPositiveSemiDefiniteSlow
bool isPositiveSemiDefiniteSlow(const AmgSymMatrix(N) &mat)
These are the slow test following the definition.
Definition: EventPrimitivesCovarianceHelpers.h:158
res
std::pair< std::vector< unsigned int >, bool > res
Definition: JetGroupProductTest.cxx:14
sign
int sign(int a)
Definition: TRT_StrawNeighbourSvc.h:127
Amg::isPositiveSemiDefinite
bool isPositiveSemiDefinite(const AmgSymMatrix(N) &mat)
Check if is positive semidefinit using that fact that is needed for Cholesky decomposition.
Definition: EventPrimitivesCovarianceHelpers.h:118
min
#define min(a, b)
Definition: cfImp.cxx:40
EventPrimitives.h
Amg
Definition of ATLAS Math & Geometry primitives (Amg)
Definition: AmgStringHelpers.h:19
Amg::isPositiveDefiniteSlow
bool isPositiveDefiniteSlow(const AmgSymMatrix(N) &mat)
Definition: EventPrimitivesCovarianceHelpers.h:188
Amg::hasPositiveOrZeroDiagElems
bool hasPositiveOrZeroDiagElems(const AmgSymMatrix(N) &mat)
Returns true if all diagonal elements of the covariance matrix are finite aka sane in the above defin...
Definition: EventPrimitivesCovarianceHelpers.h:73
Amg::chi2
double chi2(const T &precision, const U &residual, const int sign=1)
Definition: EventPrimitivesCovarianceHelpers.h:221