25#include "Eigen/Eigenvalues"
27#include "Eigen/Cholesky"
28#include "boost/format.hpp"
72 void setRank(Eigen::MatrixBase<D>
const& values) {
85 int toAdd = ~state & desired;
89 fisher.selfadjointView<Eigen::Lower>().rankUpdate(
design.adjoint());
92 fisher.triangularView<Eigen::StrictlyUpper>() =
fisher.adjoint();
142 explicit EigensystemSolver(
int dimension)
145 void factor()
override {
146 ensure(LOWER_FISHER_MATRIX | RHS_VECTOR);
147 _eig.compute(fisher);
148 if (_eig.info() == Eigen::Success) {
149 setRank(_eig.eigenvalues().reverse());
150 LOGL_DEBUG(_log,
"SelfAdjointEigenSolver succeeded: dimension=%d, rank=%d", dimension, rank);
155 ensure(FULL_FISHER_MATRIX);
156 _svd.compute(fisher, Eigen::ComputeFullU);
157 setRank(_svd.singularValues());
159 "SelfAdjointEigenSolver failed; falling back to equivalent SVD: dimension=%d, rank=%d",
164 void updateRank()
override {
165 if (_eig.info() == Eigen::Success) {
166 setRank(_eig.eigenvalues().reverse());
168 setRank(_svd.singularValues());
172 void updateDiagnostic()
override {
176 "Cannot compute NORMAL_CHOLESKY diagnostic from NORMAL_EIGENSYSTEM factorization.");
178 if (_eig.info() == Eigen::Success) {
179 ndarray::asEigenMatrix(diagnostic) = _eig.eigenvalues().reverse();
181 ndarray::asEigenMatrix(diagnostic) = _svd.singularValues();
184 ndarray::asEigenArray(diagnostic) = ndarray::asEigenArray(diagnostic).sqrt();
188 void updateSolution()
override {
190 ndarray::asEigenMatrix(solution).setZero();
193 if (_eig.info() == Eigen::Success) {
194 _tmp.head(rank) = _eig.eigenvectors().rightCols(rank).adjoint() * rhs;
195 _tmp.head(rank).array() /= _eig.eigenvalues().tail(rank).array();
196 ndarray::asEigenMatrix(solution) = _eig.eigenvectors().rightCols(rank) * _tmp.head(rank);
198 _tmp.head(rank) = _svd.matrixU().leftCols(rank).adjoint() * rhs;
199 _tmp.head(rank).array() /= _svd.singularValues().head(rank).array();
200 ndarray::asEigenMatrix(solution) = _svd.matrixU().leftCols(rank) * _tmp.head(rank);
204 void updateCovariance()
override {
209 if (_eig.info() == Eigen::Success) {
211 _eig.eigenvectors().rightCols(rank) *
212 _eig.eigenvalues().tail(rank).array().inverse().matrix().asDiagonal() *
213 _eig.eigenvectors().rightCols(rank).adjoint();
216 _svd.matrixU().leftCols(rank) *
217 _svd.singularValues().head(rank).array().inverse().matrix().asDiagonal() *
218 _svd.matrixU().leftCols(rank).adjoint();
223 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> _eig;
224 Eigen::JacobiSVD<Eigen::MatrixXd> _svd;
225 Eigen::VectorXd _tmp;
230 explicit CholeskySolver(
int dimension)
231 : Impl(dimension,
LeastSquares::NORMAL_CHOLESKY, 0.0), _ldlt(dimension) {}
233 void factor()
override {
234 ensure(LOWER_FISHER_MATRIX | RHS_VECTOR);
235 _ldlt.compute(fisher);
238 void updateRank()
override {}
240 void updateDiagnostic()
override {
243 pex::exceptions::LogicError,
244 "Can only compute NORMAL_CHOLESKY diagnostic from NORMAL_CHOLESKY factorization.");
246 ndarray::asEigenMatrix(diagnostic) = _ldlt.vectorD();
249 void updateSolution()
override { ndarray::asEigenMatrix(solution) = _ldlt.solve(rhs); }
251 void updateCovariance()
override {
252 auto cov = ndarray::asEigenMatrix(
covariance);
254 cov = _ldlt.solve(cov);
258 Eigen::LDLT<Eigen::MatrixXd> _ldlt;
261class SvdSolver :
public LeastSquares::Impl {
263 explicit SvdSolver(
int dimension)
264 : Impl(dimension, LeastSquares::DIRECT_SVD), _svd(), _tmp(dimension) {}
266 void factor()
override {
267 if (!(state & DESIGN_AND_DATA)) {
268 throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
269 "Cannot initialize DIRECT_SVD solver with normal equations.");
271 _svd.compute(design, Eigen::ComputeThinU | Eigen::ComputeThinV);
272 setRank(_svd.singularValues());
273 LOGL_DEBUG(_log,
"Using direct SVD method; dimension=%d, rank=%d", dimension, rank);
276 void updateRank()
override { setRank(_svd.singularValues()); }
278 void updateDiagnostic()
override {
279 switch (whichDiagnostic) {
281 ndarray::asEigenArray(diagnostic) = _svd.singularValues().array().square();
285 pex::exceptions::LogicError,
286 "Can only compute NORMAL_CHOLESKY diagnostic from DIRECT_SVD factorization.");
288 ndarray::asEigenMatrix(diagnostic) = _svd.singularValues();
293 void updateSolution()
override {
295 ndarray::asEigenMatrix(solution).setZero();
298 _tmp.head(rank) = _svd.matrixU().leftCols(rank).adjoint() *
data;
299 _tmp.head(rank).array() /= _svd.singularValues().head(rank).array();
300 ndarray::asEigenMatrix(solution) = _svd.matrixV().leftCols(rank) * _tmp.head(rank);
303 void updateCovariance()
override {
309 _svd.matrixV().leftCols(rank) *
310 _svd.singularValues().head(rank).array().inverse().square().matrix().asDiagonal() *
311 _svd.matrixV().leftCols(rank).adjoint();
315 Eigen::JacobiSVD<Eigen::MatrixXd> _svd;
316 Eigen::VectorXd _tmp;
345 ndarray::makeVector(_impl->
dimension, 1), _impl);
364 switch (factorization) {
366 _impl = std::make_shared<EigensystemSolver>(dimension);
369 _impl = std::make_shared<CholeskySolver>(dimension);
372 _impl = std::make_shared<SvdSolver>(dimension);
384Eigen::MatrixXd& LeastSquares::_getDesignMatrix() {
return _impl->
design; }
385Eigen::VectorXd& LeastSquares::_getDataVector() {
return _impl->
data; }
387Eigen::MatrixXd& LeastSquares::_getFisherMatrix() {
return _impl->
fisher; }
388Eigen::VectorXd& LeastSquares::_getRhsVector() {
return _impl->
rhs; }
390void LeastSquares::_factor(
bool haveNormalEquations) {
391 if (haveNormalEquations) {
392 if (_getFisherMatrix().rows() != _impl->
dimension) {
393 throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
394 (boost::format(
"Number of rows of Fisher matrix (%d) does not match"
395 " dimension of LeastSquares solver.") %
396 _getFisherMatrix().rows() % _impl->
dimension)
399 if (_getFisherMatrix().cols() != _impl->
dimension) {
400 throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
401 (boost::format(
"Number of columns of Fisher matrix (%d) does not match"
402 " dimension of LeastSquares solver.") %
403 _getFisherMatrix().cols() % _impl->
dimension)
406 if (_getRhsVector().size() != _impl->
dimension) {
407 throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
408 (boost::format(
"Number of elements in RHS vector (%d) does not match"
409 " dimension of LeastSquares solver.") %
410 _getRhsVector().size() % _impl->
dimension)
415 if (_getDesignMatrix().cols() != _impl->
dimension) {
417 pex::exceptions::InvalidParameterError,
418 "Number of columns of design matrix does not match dimension of LeastSquares solver.");
420 if (_getDesignMatrix().rows() != _getDataVector().size()) {
421 throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
422 (boost::format(
"Number of rows of design matrix (%d) does not match number of "
423 "data points (%d)") %
424 _getDesignMatrix().rows() % _getDataVector().size())
427 if (_getDesignMatrix().cols() > _getDataVector().size()) {
429 pex::exceptions::InvalidParameterError,
430 (boost::format(
"Number of columns of design matrix (%d) must be smaller than number of "
431 "data points (%d)") %
432 _getDesignMatrix().cols() % _getDataVector().size())
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
LSST DM logging module built on log4cxx.
#define LOG_GET(logger)
Returns a Log object associated with logger.
#define LOGL_DEBUG(logger, message...)
Log a debug-level message using a varargs/printf style interface.
virtual void updateRank()=0
ndarray::Array< double, 1, 1 > solution
ndarray::Array< double, 2, 2 > covariance
virtual void updateSolution()=0
virtual void updateCovariance()=0
Factorization whichDiagnostic
Impl(int dimension_, Factorization factorization_, double threshold_=std::numeric_limits< double >::epsilon())
ndarray::Array< double, 1, 1 > diagnostic
virtual void updateDiagnostic()=0
Factorization factorization
void setRank(Eigen::MatrixBase< D > const &values)
Solver for linear least-squares problems.
ndarray::Array< double const, 1, 1 > getSolution()
Return the vector solution to the least squares problem.
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.
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...
double getThreshold() const
Get the threshold used to determine when to truncate Eigenvalues.
void setThreshold(double threshold)
Set the threshold used to determine when to truncate Eigenvalues.
Reports errors in the logical structure of the program.