Belle II Software  release-05-02-19
DecorrelationMatrix.h
1 /*************************************************************************
2 * BASF2 (Belle Analysis Framework 2) *
3 * Copyright(C) 2015 - Belle II Collaboration *
4 * *
5 * Author: The Belle II Collaboration *
6 * Contributors: Thomas Madlener *
7 * *
8 * This software is provided "as is" without any warranty. *
9 **************************************************************************/
10 
11 #pragma once
12 
13 #include <tracking/trackFindingVXD/filterTools/DecorrelationMatrixHelper.h>
14 
15 #include <Eigen/Core>
16 #include <Eigen/Eigenvalues>
17 
18 #include <array>
19 #include <vector>
20 #include <string>
21 #include <sstream>
22 #include <limits> // needed for getting full precision of doubles for printing
23 #include <iomanip>
24 
25 namespace Belle2 {
37  template<size_t Ndims>
38  class DecorrelationMatrix {
39  public:
41  using MatrixT = Eigen::Matrix<double, Ndims, Ndims, Eigen::RowMajor>;
42 
44  explicit DecorrelationMatrix(const MatrixT& matrix = MatrixT::Identity()) : m_matrix(matrix) { }
45 
47  DecorrelationMatrix(const DecorrelationMatrix& matrix) = default;
48 
50  DecorrelationMatrix& operator=(const DecorrelationMatrix& rhs) = default;
51 
52  const MatrixT& getMatrix() const { return m_matrix; }
60  void calculateDecorrMatrix(std::array<std::vector<double>, Ndims> inputData, bool normalise = true);
61 
62 
64  void calculateDecorrMatrix(const MatrixT& covMatrix, bool normalise = true);
65 
67  std::vector<double> decorrelate(const std::vector<double>& inputVec) const;
68 
70  std::vector<double> decorrelate(const std::array<double, Ndims>& inputs) const;
71 
75  std::array<std::vector<double>, Ndims> decorrelate(const std::array<std::vector<double>, Ndims>& inputMat) const;
76 
80  std::string print() const;
81 
88  bool readFromStream(std::istream& is);
89 
90  private:
94  };
95 
96 
97  // ======================================== IMPLEMENTATION ========================================
98 
99  template<size_t Ndims>
100  void DecorrelationMatrix<Ndims>::calculateDecorrMatrix(std::array<std::vector<double>, Ndims> inputData, bool normalise)
101  {
102  calculateDecorrMatrix(calculateCovMatrix(inputData), normalise);
103  }
104 
105  template<size_t Ndims>
107  {
108  // make an eigen decomposition of the covariance matrix to define the transformation matrix
109  Eigen::SelfAdjointEigenSolver<MatrixT> eigenSolver(covMatrix);
110  const MatrixT U = eigenSolver.eigenvectors();
111 
112  if (normalise) {
113  // numerically better (faster hopefully) to not invert a full matrix but to calculate the diagonal elements of the
114  // (diagonal) matrix containing the eigenvalues used for normalization
115  const MatrixT D = eigenSolver.eigenvalues().cwiseSqrt().cwiseInverse().asDiagonal();
116  m_matrix = U * D;
117  } else {
118  m_matrix = U;
119  }
120  }
121 
122  template<size_t Ndims>
123  std::vector<double> DecorrelationMatrix<Ndims>::decorrelate(const std::vector<double>& inputVec) const
124  {
125  using RVector = Eigen::Matrix<double, 1, Ndims, Eigen::RowMajor>;
126  Eigen::Map<const RVector> input(inputVec.data(), Ndims);
127 
128  const RVector transform = input * m_matrix;
129 
130  return std::vector<double>(transform.data(), transform.data() + transform.size());
131  }
132 
133  template<size_t Ndims>
134  std::vector<double> DecorrelationMatrix<Ndims>::decorrelate(const std::array<double, Ndims>& input) const
135  {
136  return decorrelate(std::vector<double>(input.begin(), input.end()));
137  }
138 
139  template<size_t Ndims>
140  std::array<std::vector<double>, Ndims>
141  DecorrelationMatrix<Ndims>::decorrelate(const std::array<std::vector<double>, Ndims>& inputMat) const
142  {
143  using DVector = Eigen::Matrix<double, Eigen::Dynamic, 1, Eigen::ColMajor>;
144  using DMatrix = Eigen::Matrix<double, Eigen::Dynamic, Ndims, Eigen::ColMajor>;
145 
146  size_t nSamples = inputMat[0].size();
147  DMatrix dataMatrix(nSamples, Ndims);
148  for (size_t i = 0; i < Ndims; ++i) {
149  dataMatrix.col(i) = Eigen::Map<const DVector>(inputMat[i].data(), inputMat[i].size());
150  }
151 
152  DMatrix transform = dataMatrix * m_matrix;
153 
154  std::array<std::vector<double>, Ndims> output;
155  for (size_t i = 0; i < Ndims; ++i) {
156  output[i] = std::vector<double>(transform.col(i).data(), transform.col(i).data() + transform.col(i).size());
157  }
158 
159  return output;
160  }
161 
162  template<size_t Ndims>
163  std::string DecorrelationMatrix<Ndims>::print() const
164  {
165  std::stringstream sstr{};
166  sstr << std::setprecision(std::numeric_limits<double>::digits10 + 1) << m_matrix; // ensure precision at write-out
167  return sstr.str();
168  }
169 
170  template<size_t Ndims>
172  {
173  MatrixT inMat{};
174  using RVector = Eigen::Matrix<double, 1, Ndims, Eigen::RowMajor>; // typedef for row vector
175 
176  size_t iRow = 0;
177  while (iRow < Ndims) { // read only as many lines as needed
178  std::string line{};
179  if (is.eof()) return false;
180  std::getline(is, line);
181  if (line.empty()) continue; // skip empty lines
182  std::stringstream lstr(line);
183  std::array<double, Ndims> linevalues;
184  for (double& value : linevalues) { // read only as many values in the line as needed
185  if (lstr.eof()) return false;
186  lstr >> value;
187  }
188  inMat.row(iRow) = Eigen::Map<const RVector>(linevalues.data(), linevalues.size());
189  ++iRow;
190  }
191 
192  m_matrix = inMat;
193  return true;
194  }
195 
197 } // end namespace Belle 2
Belle2::DecorrelationMatrix::MatrixT
Eigen::Matrix< double, Ndims, Ndims, Eigen::RowMajor > MatrixT
typedef for consistent type across class
Definition: DecorrelationMatrix.h:49
Belle2::DecorrelationMatrix::print
std::string print() const
print the matrix to a string.
Definition: DecorrelationMatrix.h:171
Belle2::DecorrelationMatrix::operator=
DecorrelationMatrix & operator=(const DecorrelationMatrix &rhs)=default
assignment operator
Belle2::DecorrelationMatrix::m_matrix
MatrixT m_matrix
internal matrix holding the data
Definition: DecorrelationMatrix.h:99
Belle2::DecorrelationMatrix::DecorrelationMatrix
DecorrelationMatrix(const MatrixT &matrix=MatrixT::Identity())
default constructor: initializes to identity matrix or to passed matrix
Definition: DecorrelationMatrix.h:52
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::DecorrelationMatrix::decorrelate
std::vector< double > decorrelate(const std::vector< double > &inputVec) const
"decorrelate" one measurement (i.e.
Definition: DecorrelationMatrix.h:131
Belle2::DecorrelationMatrix::calculateDecorrMatrix
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 ...
Definition: DecorrelationMatrix.h:108
Belle2::calculateCovMatrix
const Eigen::Matrix< double, Ndims, Ndims, Eigen::RowMajor > calculateCovMatrix(std::array< std::vector< double >, Ndims > inputData)
calculates the empirical covariance matrix from the inputData.
Definition: DecorrelationMatrixHelper.h:47
Belle2::DecorrelationMatrix::readFromStream
bool readFromStream(std::istream &is)
read from stream.
Definition: DecorrelationMatrix.h:179