Belle II Software  release-06-00-14
DecorrelationMatrix.h
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 #pragma once
10 
11 #include <tracking/trackFindingVXD/filterTools/DecorrelationMatrixHelper.h>
12 
13 #include <Eigen/Core>
14 #include <Eigen/Eigenvalues>
15 
16 #include <array>
17 #include <vector>
18 #include <string>
19 #include <sstream>
20 #include <limits> // needed for getting full precision of doubles for printing
21 #include <iomanip>
22 
23 namespace Belle2 {
35  template<size_t Ndims>
37  public:
39  using MatrixT = Eigen::Matrix<double, Ndims, Ndims, Eigen::RowMajor>;
40 
42  explicit DecorrelationMatrix(const MatrixT& matrix = MatrixT::Identity()) : m_matrix(matrix) { }
43 
45  DecorrelationMatrix(const DecorrelationMatrix& matrix) = default;
46 
49 
50  const MatrixT& getMatrix() const { return m_matrix; }
58  void calculateDecorrMatrix(std::array<std::vector<double>, Ndims> inputData, bool normalise = true);
59 
60 
62  void calculateDecorrMatrix(const MatrixT& covMatrix, bool normalise = true);
63 
65  std::vector<double> decorrelate(const std::vector<double>& inputVec) const;
66 
68  std::vector<double> decorrelate(const std::array<double, Ndims>& inputs) const;
69 
73  std::array<std::vector<double>, Ndims> decorrelate(const std::array<std::vector<double>, Ndims>& inputMat) const;
74 
78  std::string print() const;
79 
86  bool readFromStream(std::istream& is);
87 
88  private:
92  };
93 
94 
95  // ======================================== IMPLEMENTATION ========================================
96 
97  template<size_t Ndims>
98  void DecorrelationMatrix<Ndims>::calculateDecorrMatrix(std::array<std::vector<double>, Ndims> inputData, bool normalise)
99  {
100  calculateDecorrMatrix(calculateCovMatrix(inputData), normalise);
101  }
102 
103  template<size_t Ndims>
105  {
106  // make an eigen decomposition of the covariance matrix to define the transformation matrix
107  Eigen::SelfAdjointEigenSolver<MatrixT> eigenSolver(covMatrix);
108  const MatrixT U = eigenSolver.eigenvectors();
109 
110  if (normalise) {
111  // numerically better (faster hopefully) to not invert a full matrix but to calculate the diagonal elements of the
112  // (diagonal) matrix containing the eigenvalues used for normalization
113  const MatrixT D = eigenSolver.eigenvalues().cwiseSqrt().cwiseInverse().asDiagonal();
114  m_matrix = U * D;
115  } else {
116  m_matrix = U;
117  }
118  }
119 
120  template<size_t Ndims>
121  std::vector<double> DecorrelationMatrix<Ndims>::decorrelate(const std::vector<double>& inputVec) const
122  {
123  using RVector = Eigen::Matrix<double, 1, Ndims, Eigen::RowMajor>;
124  Eigen::Map<const RVector> input(inputVec.data(), Ndims);
125 
126  const RVector transform = input * m_matrix;
127 
128  return std::vector<double>(transform.data(), transform.data() + transform.size());
129  }
130 
131  template<size_t Ndims>
132  std::vector<double> DecorrelationMatrix<Ndims>::decorrelate(const std::array<double, Ndims>& input) const
133  {
134  return decorrelate(std::vector<double>(input.begin(), input.end()));
135  }
136 
137  template<size_t Ndims>
138  std::array<std::vector<double>, Ndims>
139  DecorrelationMatrix<Ndims>::decorrelate(const std::array<std::vector<double>, Ndims>& inputMat) const
140  {
141  using DVector = Eigen::Matrix<double, Eigen::Dynamic, 1, Eigen::ColMajor>;
142  using DMatrix = Eigen::Matrix<double, Eigen::Dynamic, Ndims, Eigen::ColMajor>;
143 
144  size_t nSamples = inputMat[0].size();
145  DMatrix dataMatrix(nSamples, Ndims);
146  for (size_t i = 0; i < Ndims; ++i) {
147  dataMatrix.col(i) = Eigen::Map<const DVector>(inputMat[i].data(), inputMat[i].size());
148  }
149 
150  DMatrix transform = dataMatrix * m_matrix;
151 
152  std::array<std::vector<double>, Ndims> output;
153  for (size_t i = 0; i < Ndims; ++i) {
154  output[i] = std::vector<double>(transform.col(i).data(), transform.col(i).data() + transform.col(i).size());
155  }
156 
157  return output;
158  }
159 
160  template<size_t Ndims>
162  {
163  std::stringstream sstr{};
164  sstr << std::setprecision(std::numeric_limits<double>::digits10 + 1) << m_matrix; // ensure precision at write-out
165  return sstr.str();
166  }
167 
168  template<size_t Ndims>
170  {
171  MatrixT inMat{};
172  using RVector = Eigen::Matrix<double, 1, Ndims, Eigen::RowMajor>; // typedef for row vector
173 
174  size_t iRow = 0;
175  while (iRow < Ndims) { // read only as many lines as needed
176  std::string line{};
177  if (is.eof()) return false;
178  std::getline(is, line);
179  if (line.empty()) continue; // skip empty lines
180  std::stringstream lstr(line);
181  std::array<double, Ndims> linevalues;
182  for (double& value : linevalues) { // read only as many values in the line as needed
183  if (lstr.eof()) return false;
184  lstr >> value;
185  }
186  inMat.row(iRow) = Eigen::Map<const RVector>(linevalues.data(), linevalues.size());
187  ++iRow;
188  }
189 
190  m_matrix = inMat;
191  return true;
192  }
193 
195 } // end namespace Belle 2
Class holding a Matrix that can be used to decorrelate input data to Machine Learning classifiers.
DecorrelationMatrix(const MatrixT &matrix=MatrixT::Identity())
default constructor: initializes to identity matrix or to passed matrix
const MatrixT & getMatrix() const
get the currently stored matrix
MatrixT m_matrix
internal matrix holding the data
void calculateDecorrMatrix(const MatrixT &covMatrix, bool normalise=true)
calculate the matrix that can be used to decorrelate the data that yield the passed covariance matrix...
DecorrelationMatrix & operator=(const DecorrelationMatrix &rhs)=default
assignment operator
DecorrelationMatrix(const DecorrelationMatrix &matrix)=default
copy constructor
Eigen::Matrix< double, Ndims, Ndims, Eigen::RowMajor > MatrixT
typedef for consistent type across class
const Eigen::Matrix< double, Ndims, Ndims, Eigen::RowMajor > calculateCovMatrix(std::array< std::vector< double >, Ndims > inputData)
calculates the empirical covariance matrix from the inputData.
std::vector< double > decorrelate(const std::vector< double > &inputVec) const
"decorrelate" one measurement (i.e.
std::vector< double > decorrelate(const std::array< double, Ndims > &inputs) const
"decorrelate" one measurement (i.e.
void calculateDecorrMatrix(std::array< std::vector< double >, Ndims > inputData, bool normalise=true)
calculate the transformation matrix that when applied to the input data yields linearly uncorrelated ...
bool readFromStream(std::istream &is)
read from stream.
std::string print() const
print the matrix to a string.
std::array< std::vector< double >, Ndims > decorrelate(const std::array< std::vector< double >, Ndims > &inputMat) const
decorrelate M measurements (i.e.
Abstract base class for different kinds of events.