25#ifndef LSST_AFW_MATH_LeastSquares_h_INCLUDED
26#define LSST_AFW_MATH_LeastSquares_h_INCLUDED
29#include "ndarray/eigen.h"
99 template <
typename T1,
typename T2,
int C1,
int C2>
101 ndarray::Array<T2, 1, C2>
const&
data,
103 LeastSquares r(factorization, design.template getSize<1>());
104 r.setDesignMatrix(design,
data);
109 template <
typename D1,
typename D2>
111 Eigen::MatrixBase<D2>
const&
data,
114 r.setDesignMatrix(design,
data);
119 template <
typename T1,
typename T2,
int C1,
int C2>
125 _getDesignMatrix() = ndarray::asEigenMatrix(design).template cast<double>();
126 _getDataVector() = ndarray::asEigenMatrix(
data).template cast<double>();
131 template <
typename D1,
typename D2>
133 _getDesignMatrix() = design.template cast<double>();
134 _getDataVector() =
data.template cast<double>();
139 template <
typename T1,
int C1>
141 _getDesignMatrix() = ndarray::asEigenMatrix(design).template cast<double>();
146 template <
typename D1,
typename D2>
148 _getDesignMatrix() = design.template cast<double>();
153 template <
typename T1,
typename T2,
int C1,
int C2>
155 ndarray::Array<T2, 1, C2>
const& rhs,
157 LeastSquares r(factorization, fisher.template getSize<0>());
158 r.setNormalEquations(fisher, rhs);
163 template <
typename D1,
typename D2>
165 Eigen::MatrixBase<D2>
const& rhs,
168 r.setNormalEquations(fisher, rhs);
173 template <
typename T1,
typename T2,
int C1,
int C2>
174 void setNormalEquations(ndarray::Array<T1, 2, C1>
const& fisher, ndarray::Array<T2, 1, C2>
const& rhs) {
175 if ((C1 > 0) ==
bool(Eigen::MatrixXd::IsRowMajor)) {
176 _getFisherMatrix() = ndarray::asEigenMatrix(fisher).template cast<double>();
178 _getFisherMatrix() = ndarray::asEigenMatrix(fisher).transpose().template cast<double>();
180 _getRhsVector() = ndarray::asEigenMatrix(rhs).template cast<double>();
185 template <
typename D1,
typename D2>
187 if (
bool(Eigen::MatrixBase<D1>::IsRowMajor) ==
bool(Eigen::MatrixXd::IsRowMajor)) {
188 _getFisherMatrix() = fisher.template cast<double>();
190 _getFisherMatrix() = fisher.transpose().template cast<double>();
192 _getRhsVector() = rhs.template cast<double>();
324 Eigen::MatrixXd& _getDesignMatrix();
325 Eigen::VectorXd& _getDataVector();
330 Eigen::MatrixXd& _getFisherMatrix();
331 Eigen::VectorXd& _getRhsVector();
336 void _factor(
bool haveNormalEquations);
Solver for linear least-squares problems.
LeastSquares & operator=(LeastSquares &&)
void setDesignMatrix(ndarray::Array< T1, 2, C1 > const &design)
Reset the design matrix given as an ndarray; dimension and data are not changed.
ndarray::Array< double const, 1, 1 > getSolution()
Return the vector solution to the least squares problem.
static LeastSquares fromDesignMatrix(ndarray::Array< T1, 2, C1 > const &design, ndarray::Array< T2, 1, C2 > const &data, Factorization factorization=NORMAL_EIGENSYSTEM)
Initialize from the design matrix and data vector given as ndarrays.
ndarray::Array< double const, 2, 2 > getFisherMatrix()
Return the Fisher matrix (inverse of the covariance) of the parameters.
int getRank() const
Return the rank of the problem (number of nonzero Eigenvalues).
ndarray::Array< double const, 2, 2 > getCovariance()
Return the covariance matrix of the least squares problem.
LeastSquares & operator=(LeastSquares const &)
Factorization
Private implementation; forward-declared publicly so we can inherit from it in .cc.
@ NORMAL_EIGENSYSTEM
Use the normal equations with a symmetric Eigensystem decomposition.
@ NORMAL_CHOLESKY
Use the normal equations with a Cholesky decomposition.
@ DIRECT_SVD
Use a thin singular value decomposition of the design matrix.
LeastSquares(Factorization factorization, int dimension)
Construct a least-squares object for the given factorization and dimensionality.
int getDimension() const
Return the number of parameters.
void setDesignMatrix(Eigen::MatrixBase< D1 > const &design, Eigen::MatrixBase< D2 > const &data)
Reset the design matrix and data vector given as Eigen objects; dimension must not change.
static LeastSquares fromNormalEquations(Eigen::MatrixBase< D1 > const &fisher, Eigen::MatrixBase< D2 > const &rhs, Factorization factorization=NORMAL_EIGENSYSTEM)
Initialize from the terms in the normal equations, given as Eigen objects.
Factorization getFactorization() const
Retun the type of factorization used by the solver.
ndarray::Array< double const, 1, 1 > getDiagnostic(Factorization factorization)
Return a factorization-dependent vector that can be used to characterize the stability of the solutio...
void setNormalEquations(Eigen::MatrixBase< D1 > const &fisher, Eigen::MatrixBase< D2 > const &rhs)
Reset the terms in the normal equations given as Eigen objects; dimension must not change.
void setDesignMatrix(ndarray::Array< T1, 2, C1 > const &design, ndarray::Array< T2, 1, C2 > const &data)
Reset the design matrix and data vector given as ndarrays; dimension must not change.
double getThreshold() const
Get the threshold used to determine when to truncate Eigenvalues.
LeastSquares(LeastSquares &&)
void setNormalEquations(ndarray::Array< T1, 2, C1 > const &fisher, ndarray::Array< T2, 1, C2 > const &rhs)
Reset the terms in the normal equations given as ndarrays; dimension must not change.
static LeastSquares fromDesignMatrix(Eigen::MatrixBase< D1 > const &design, Eigen::MatrixBase< D2 > const &data, Factorization factorization=NORMAL_EIGENSYSTEM)
Initialize from the design matrix and data vector given as an Eigen objects.
void setThreshold(double threshold)
Set the threshold used to determine when to truncate Eigenvalues.
void setDesignMatrix(Eigen::MatrixBase< D1 > const &design)
Reset the design matrix given as an Eigen object; dimension and data are not changed.
LeastSquares(LeastSquares const &)
static LeastSquares fromNormalEquations(ndarray::Array< T1, 2, C1 > const &fisher, ndarray::Array< T2, 1, C2 > const &rhs, Factorization factorization=NORMAL_EIGENSYSTEM)
Initialize from the terms in the normal equations, given as ndarrays.