ATLAS Offline Software
EventPrimitivesHelpers.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_EVENTPRIMITIVESHELPERS_H
10 #define EVENTPRIMITIVES_EVENTPRIMITIVESHELPERS_H
11 
13 //
14 #include <cmath>
15 #include <vector>
16 
24 namespace Amg {
27 template <int N>
28 inline bool saneVector(const AmgVector(N) & vec) {
29  for (int i = 0; i < N; ++i) {
30  if (std::isnan(vec[i]) || std::isinf(vec[i]))
31  return false;
32  }
33  constexpr double max_length2 = 1.e16;
34  return vec.dot(vec) < max_length2;
35 }
36 
40 inline double error(const Amg::MatrixX& mat, int index) {
41  return std::sqrt(mat(index, index));
42 }
43 template <int N>
44 inline double error(const AmgSymMatrix(N) & mat, int index) {
45  assert(index < N);
46  return std::sqrt(mat(index, index));
47 }
48 
49 // expression template to evaluate the required size of the compressed matrix at
50 // compile time
51 inline constexpr int CalculateCompressedSize(int n) {
52  return (n * (n + 1)) / 2;
53 }
54 
55 template <int N>
56 inline void compress(const AmgSymMatrix(N) & covMatrix,
57  std::vector<float>& vec) {
58  vec.reserve(CalculateCompressedSize(N));
59  for (unsigned int i = 0; i < N; ++i)
60  for (unsigned int j = 0; j <= i; ++j)
61  vec.emplace_back(covMatrix(i, j));
62 }
63 
64 inline void compress(const MatrixX& covMatrix, std::vector<float>& vec) {
65  int rows = covMatrix.rows();
67  for (int i = 0; i < rows; ++i) {
68  for (int j = 0; j <= i; ++j) {
69  vec.emplace_back(covMatrix(i, j));
70  }
71  }
72 }
73 
74 template <int N>
75 inline void expand(std::vector<float>::const_iterator it,
76  std::vector<float>::const_iterator,
78  for (unsigned int i = 0; i < N; ++i) {
79  for (unsigned int j = 0; j <= i; ++j) {
80  covMatrix.fillSymmetric(i, j, *it);
81  ++it;
82  }
83  }
84 }
85 
86 inline void expand(std::vector<float>::const_iterator it,
87  std::vector<float>::const_iterator it_end,
88  MatrixX& covMatrix) {
89  unsigned int dist = std::distance(it, it_end);
90  unsigned int n;
91  for (n = 1; dist > n; ++n) {
92  dist = dist - n;
93  }
94  covMatrix = MatrixX(n, n);
95  for (unsigned int i = 0; i < n; ++i) {
96  for (unsigned int j = 0; j <= i; ++j) {
97  covMatrix.fillSymmetric(i, j, *it);
98  ++it;
99  }
100  }
101 }
102 
108 template <int N>
109 std::pair<int, int> compare(const AmgSymMatrix(N) & m1,
110  const AmgSymMatrix(N) & m2, double precision = 1e-9,
111  bool relative = false) {
112  for (int i = 0; i < N; ++i) {
113  for (int j = 0; j < N; ++j) {
114  if (relative) {
115  if (0.5 * std::abs(m1(i, j) - m2(i, j)) / (m1(i, j) + m2(i, j)) >
116  precision) {
117  return std::make_pair(i, j);
118  }
119  } else {
120  if (std::abs(m1(i, j) - m2(i, j)) > precision) {
121  return std::make_pair(i, j);
122  }
123  }
124  }
125  }
126  return std::make_pair(-1, -1);
127 }
128 
134 template <int N>
135 int compare(const AmgVector(N) & m1, const AmgVector(N) & m2,
136  double precision = 1e-9, bool relative = false) {
137  for (int i = 0; i < N; ++i) {
138  if (relative) {
139  if (0.5 * std::abs(m1(i) - m2(i)) / (m1(i) + m2(i)) > precision)
140  return i;
141  } else {
142  if (std::abs(m1(i) - m2(i)) > precision)
143  return i;
144  }
145  }
146  return -1;
147 }
148 
149 } // namespace Amg
150 
151 #endif
AllowedVariables::e
e
Definition: AsgElectronSelectorTool.cxx:37
Amg::compare
std::pair< int, int > compare(const AmgSymMatrix(N) &m1, const AmgSymMatrix(N) &m2, double precision=1e-9, bool relative=false)
compare two matrices, returns the indices of the first element that fails the condition,...
Definition: EventPrimitivesHelpers.h:109
Amg::MatrixX
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Dynamic Matrix - dynamic allocation.
Definition: EventPrimitives.h:27
python.SystemOfUnits.m2
int m2
Definition: SystemOfUnits.py:92
Amg::compress
void compress(const AmgSymMatrix(N) &covMatrix, std::vector< float > &vec)
Definition: EventPrimitivesHelpers.h:56
index
Definition: index.py:1
mat
GeoMaterial * mat
Definition: LArDetectorConstructionTBEC.cxx:55
Amg::CalculateCompressedSize
constexpr int CalculateCompressedSize(int n)
Definition: EventPrimitivesHelpers.h:51
skel.it
it
Definition: skel.GENtoEVGEN.py:396
Amg::saneVector
bool saneVector(const AmgVector(N) &vec)
Definition: EventPrimitivesHelpers.h:28
JetTiledMap::N
@ N
Definition: TiledEtaPhiMap.h:44
vec
std::vector< size_t > vec
Definition: CombinationsGeneratorTest.cxx:12
AmgSymMatrix
#define AmgSymMatrix(dim)
Definition: EventPrimitives.h:50
Amg::expand
void expand(std::vector< float >::const_iterator it, std::vector< float >::const_iterator, AmgSymMatrix(N) &covMatrix)
Definition: EventPrimitivesHelpers.h:75
python.changerun.m1
m1
Definition: changerun.py:32
lumiFormat.i
int i
Definition: lumiFormat.py:85
beamspotman.n
n
Definition: beamspotman.py:731
xAOD::covMatrix
covMatrix
Definition: TrackMeasurement_v1.cxx:19
AmgVector
AmgVector(4) T2BSTrackFilterTool
Definition: T2BSTrackFilterTool.cxx:114
beamspotnt.rows
list rows
Definition: bin/beamspotnt.py:1112
EventPrimitives.h
Amg
Definition of ATLAS Math & Geometry primitives (Amg)
Definition: AmgStringHelpers.h:19
Amg::error
double error(const Amg::MatrixX &mat, int index)
return diagonal error of the matrix caller should ensure the matrix is symmetric and the index is in ...
Definition: EventPrimitivesHelpers.h:40
Amg::distance
float distance(const Amg::Vector3D &p1, const Amg::Vector3D &p2)
calculates the distance between two point in 3D space
Definition: GeoPrimitivesHelpers.h:54