CosmoBolognaLib
Free Software C++/Python libraries for cosmological calculations
CovarianceMatrix.h
Go to the documentation of this file.
1 /********************************************************************
2  * Copyright (C) 2014 by Federico Marulli and Alfonso Veropalumbo *
3  * federico.marulli3@unibo.it *
4  * *
5  * This program is free software; you can redistribute it and/or *
6  * modify it under the terms of the GNU General Public License as *
7  * published by the Free Software Foundation; either version 2 of *
8  * the License, or (at your option) any later version. *
9  * *
10  * This program is distributed in the hope that it will be useful, *
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of *
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
13  * GNU General Public License for more details. *
14  * *
15  * You should have received a copy of the GNU General Public *
16  * License along with this program; if not, write to the Free *
17  * Software Foundation, Inc., *
18  * 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. *
19  ********************************************************************/
20 
33 #ifndef __COVMAT__
34 #define __COVMAT__
35 
36 #include "Func.h"
37 
38 namespace cbl {
39 
47  namespace data {
48 
49 
60  {
61 
62  protected:
63 
65  size_t m_order;
66 
68  Eigen::MatrixXd m_matrix;
69 
71  Eigen::MatrixXd m_precision;
72 
74  Eigen::MatrixXd m_correlation;
75 
77  Eigen::VectorXd m_variance;
78 
80  Eigen::VectorXd m_std;
81 
83  double m_determinant;
84 
86  // the covariance is measured from multiple dataset
88 
93  void m_set_default ();
94 
105  virtual void m_set (const std::vector<double> matrix, const double nmeasures=-1, const double prec=1.e-10);
106 
122  inline double hartlap_factor (const double order, const double nmeasures=-1) { return (nmeasures>1) ? 1-(order+1)/(nmeasures-1) : 1; }
123 
124  public:
125 
130 
137 
151  CovarianceMatrix (std::vector<std::vector<double>> covariance_matrix, const double nmeasures=-1, const double prec=1.e-10)
152  { set_from_matrix(covariance_matrix, nmeasures, prec); }
153 
163  CovarianceMatrix (std::vector<double> standard_deviation, const double nmeasures=-1)
165 
186  CovarianceMatrix (const std::string filename, const int cov_col=2, const int skipped_lines=0, const double nmeasures=-1, const double prec=1.e-10)
187  {read(filename, cov_col, skipped_lines, prec, nmeasures); }
188 
189 
193  virtual ~CovarianceMatrix () = default;
194 
196 
201 
211  double operator() (const int i, const int j) const {return m_matrix(i, j);}
212 
217  std::vector<std::vector<double>> operator() () const;
218 
229  double correlation (const int i, const int j) const {return m_correlation(i, j);}
230 
237  std::vector<std::vector<double>> correlation () const;
238 
248  double precision (const int i, const int j) const {return m_precision(i, j);}
249 
255  std::vector<std::vector<double>> precision () const;
256 
262  double determinant () const {return m_determinant;};
263 
284  double precision_hartlap (const int i, const int j) const {return m_hartlap_factor*m_precision(i, j);}
285 
302  std::vector<std::vector<double>> precision_hartlap () const;
303 
311  double standard_deviation (const int i) const { return m_std(i); }
312 
318  std::vector<double> standard_deviation () const;
319 
327  double variance (const int i) const { return m_variance(i); }
328 
334  std::vector<double> variance () const;
335 
341  size_t order () const { return m_order; }
342 
344 
345 
350 
362  void set_from_matrix (const std::vector<std::vector<double>> covariance, const double nmeasures=-1, const double prec=1.e-10)
363  { m_set(cbl::flatten(covariance), nmeasures, prec); }
364 
372  void set_from_standard_deviation (const std::vector<double> standard_deviation, const double nmeasures=-1);
373 
395  void measure (const std::vector<std::shared_ptr<Data>> dataset, const double normalization=1, const double prec=1.e-10);
396 
420  void measure (const std::vector<std::vector<std::shared_ptr<Data>>> dataset, const double normalization=1, const double prec=1.e-10);
421 
423 
428 
447  void read (const std::string filename, const int cov_col=2, const int skipped_lines=0, const double nmeasures=-1, const double prec=1.e-10);
448 
462  void write (const std::string dir, const std::string file, const int precision=4, const int rank=0) const;
463 
465 
471 
479  CovarianceMatrix cut (const std::vector<bool> mask) const;
480 
482 
488 
499  CovarianceMatrix operator += (const CovarianceMatrix covariance) const;
500 
512  CovarianceMatrix operator += (const std::shared_ptr<CovarianceMatrix> covariance) const;
513 
526  CovarianceMatrix operator += (const std::vector<CovarianceMatrix> covariance) const;
527 
540  CovarianceMatrix operator += (const std::vector<std::shared_ptr<CovarianceMatrix>> covariance) const;
541 
543 
544 
545  };
546 
547 
548  }
549 }
550 
551 #endif
Useful generic functions.
The class CovarianceMatrix.
std::vector< std::vector< double > > operator()() const
get the covariance matrix
void measure(const std::vector< std::shared_ptr< Data >> dataset, const double normalization=1, const double prec=1.e-10)
measure the covariance from a collection of dataset
Eigen::VectorXd m_std
standard deviation
std::vector< double > standard_deviation() const
get the standard deviation
void set_from_standard_deviation(const std::vector< double > standard_deviation, const double nmeasures=-1)
set the covariance by passing the diagonal square root
double determinant() const
get the determinant
CovarianceMatrix(std::vector< std::vector< double >> covariance_matrix, const double nmeasures=-1, const double prec=1.e-10)
constructor which sets the covariance matrix
CovarianceMatrix(std::vector< double > standard_deviation, const double nmeasures=-1)
constructor which gets the data from an input vector
void read(const std::string filename, const int cov_col=2, const int skipped_lines=0, const double nmeasures=-1, const double prec=1.e-10)
set the covariance matrix, reading from an input file
std::vector< std::vector< double > > precision_hartlap() const
get the value of the precision matrix at index i,j times the Hartlap factor
size_t m_order
number of data
virtual ~CovarianceMatrix()=default
default destructor
Eigen::VectorXd m_variance
diagonal of the covariance matrix
std::vector< std::vector< double > > correlation() const
get the value of the correlation matrix at index i,j
std::vector< double > variance() const
get the variance
double hartlap_factor(const double order, const double nmeasures=-1)
compute the hartlap factor. This is used to de-bias precision matrix measured from covariance measure...
double precision(const int i, const int j) const
get the value of the precision matrix at index i,j
double m_hartlap_factor
The hartlap factor, only set when.
void m_set_default()
set internal attributes to default values
double m_determinant
determinant
double precision_hartlap(const int i, const int j) const
get the value of the precision matrix at index i,j times the Hartlap factor
Eigen::MatrixXd m_correlation
correlation matrix
virtual void m_set(const std::vector< double > matrix, const double nmeasures=-1, const double prec=1.e-10)
set internal attributes
CovarianceMatrix()
default constructor
double correlation(const int i, const int j) const
get the value of the correlation matrix at index i,j
double standard_deviation(const int i) const
get value of the standard deviation at index i
double variance(const int i) const
get value of the variance at index i
void write(const std::string dir, const std::string file, const int precision=4, const int rank=0) const
write the covariance matrix
Eigen::MatrixXd m_matrix
covariance matrix
std::vector< std::vector< double > > precision() const
get the precision matrix
size_t order() const
return the covariance matrix order
CovarianceMatrix(const std::string filename, const int cov_col=2, const int skipped_lines=0, const double nmeasures=-1, const double prec=1.e-10)
constructs with sets the covariance matrix reading from an input file
void set_from_matrix(const std::vector< std::vector< double >> covariance, const double nmeasures=-1, const double prec=1.e-10)
set the covariance matrix by passing a std::vector<std::vector<double>> object;
CovarianceMatrix cut(const std::vector< bool > mask) const
cut the data, for Data1D
CovarianceMatrix operator+=(const CovarianceMatrix covariance) const
overloading of the += operator, to sum two catalogues
Eigen::MatrixXd m_precision
precision matrix
The global namespace of the CosmoBolognaLib
Definition: CAMB.h:38
std::vector< T > flatten(std::vector< std::vector< T >> matrix)
flatten a matrix in a vector of size
Definition: Kernel.h:1686
void covariance_matrix(const std::vector< std::vector< double >> mat, std::vector< std::vector< double >> &cov, const bool JK=false)
compute the covariance matrix from an input dataset
Definition: Func.cpp:761