LSSTApplications
19.0.0-14-gb0260a2+72efe9b372,20.0.0+7927753e06,20.0.0+8829bf0056,20.0.0+995114c5d2,20.0.0+b6f4b2abd1,20.0.0+bddc4f4cbe,20.0.0-1-g253301a+8829bf0056,20.0.0-1-g2b7511a+0d71a2d77f,20.0.0-1-g5b95a8c+7461dd0434,20.0.0-12-g321c96ea+23efe4bbff,20.0.0-16-gfab17e72e+fdf35455f6,20.0.0-2-g0070d88+ba3ffc8f0b,20.0.0-2-g4dae9ad+ee58a624b3,20.0.0-2-g61b8584+5d3db074ba,20.0.0-2-gb780d76+d529cf1a41,20.0.0-2-ged6426c+226a441f5f,20.0.0-2-gf072044+8829bf0056,20.0.0-2-gf1f7952+ee58a624b3,20.0.0-20-geae50cf+e37fec0aee,20.0.0-25-g3dcad98+544a109665,20.0.0-25-g5eafb0f+ee58a624b3,20.0.0-27-g64178ef+f1f297b00a,20.0.0-3-g4cc78c6+e0676b0dc8,20.0.0-3-g8f21e14+4fd2c12c9a,20.0.0-3-gbd60e8c+187b78b4b8,20.0.0-3-gbecbe05+48431fa087,20.0.0-38-ge4adf513+a12e1f8e37,20.0.0-4-g97dc21a+544a109665,20.0.0-4-gb4befbc+087873070b,20.0.0-4-gf910f65+5d3db074ba,20.0.0-5-gdfe0fee+199202a608,20.0.0-5-gfbfe500+d529cf1a41,20.0.0-6-g64f541c+d529cf1a41,20.0.0-6-g9a5b7a1+a1cd37312e,20.0.0-68-ga3f3dda+5fca18c6a4,20.0.0-9-g4aef684+e18322736b,w.2020.45
LSSTDataManagementBasePackage
|
Go to the documentation of this file.
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>());
109 template <
typename D1,
typename D2>
111 Eigen::MatrixBase<D2>
const&
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>());
163 template <
typename D1,
typename D2>
165 Eigen::MatrixBase<D2>
const& 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);
344 #endif // !LSST_AFW_MATH_LeastSquares_h_INCLUDED
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.
@ DIRECT_SVD
Use a thin singular value decomposition of the design matrix.
LeastSquares & operator=(LeastSquares &&)
LeastSquares(LeastSquares const &)
void setThreshold(double threshold)
Set the threshold used to determine when to truncate Eigenvalues.
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...
@ NORMAL_CHOLESKY
Use the normal equations with a Cholesky decomposition.
Factorization getFactorization() const
Retun the type of factorization used by the solver.
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.
Solver for linear least-squares problems.
@ NORMAL_EIGENSYSTEM
Use the normal equations with a symmetric Eigensystem decomposition.
LeastSquares & operator=(LeastSquares const &)
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.
double getThreshold() const
Get the threshold used to determine when to truncate Eigenvalues.
LeastSquares(Factorization factorization, int dimension)
Construct a least-squares object for the given factorization and dimensionality.
void setDesignMatrix(Eigen::MatrixBase< D1 > const &design)
Reset the design matrix given as an Eigen object; dimension and data are not changed.
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 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, 2, 2 > getFisherMatrix()
Return the Fisher matrix (inverse of the covariance) of the parameters.
A base class for image defects.
int getDimension() const
Return the number of parameters.
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.
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.
Factorization
Private implementation; forward-declared publicly so we can inherit from it in .cc.
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.
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.
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.
LeastSquares(LeastSquares &&)