25 #include "Eigen/Eigenvalues"
27 #include "Eigen/Cholesky"
28 #include "boost/format.hpp"
29 #include "boost/make_shared.hpp"
35 namespace lsst {
namespace afw {
namespace math {
69 void setRank(Eigen::MatrixBase<D>
const & values) {
81 int toAdd = ~
state & desired;
85 fisher.selfadjointView<Eigen::Lower>().rankUpdate(
design.adjoint());
87 if (toAdd & FULL_FISHER_MATRIX) {
88 fisher.triangularView<Eigen::StrictlyUpper>() =
fisher.adjoint();
109 virtual void factor() = 0;
118 Impl(
int dimension_,
double threshold_=std::numeric_limits<double>::epsilon()) :
120 log(
"afw.math.LeastSquares")
131 explicit EigensystemSolver(
int dimension) :
132 Impl(dimension),
_eig(dimension),
_svd(),
_tmp(dimension)
135 virtual void factor() {
136 ensure(LOWER_FISHER_MATRIX | RHS_VECTOR);
137 _eig.compute(fisher);
138 if (
_eig.info() == Eigen::Success) {
139 setRank(
_eig.eigenvalues().reverse());
140 log.debug<5>(
"SelfAdjointEigenSolver succeeded: dimension=%d, rank=%d", dimension, rank);
145 ensure(FULL_FISHER_MATRIX);
146 _svd.compute(fisher, Eigen::ComputeFullU);
147 setRank(
_svd.singularValues());
149 "SelfAdjointEigenSolver failed; falling back to equivalent SVD: dimension=%d, rank=%d",
155 virtual void updateRank() {
156 if (
_eig.info() == Eigen::Success) {
157 setRank(
_eig.eigenvalues().reverse());
159 setRank(
_svd.singularValues());
163 virtual void updateDiagnostic() {
166 pex::exceptions::LogicError,
167 "Cannot compute NORMAL_CHOLESKY diagnostic from NORMAL_EIGENSYSTEM factorization."
170 if (
_eig.info() == Eigen::Success) {
171 diagnostic.asEigen() =
_eig.eigenvalues().reverse();
173 diagnostic.asEigen() =
_svd.singularValues();
176 diagnostic.asEigen<Eigen::ArrayXpr>() = diagnostic.asEigen<Eigen::ArrayXpr>().sqrt();
180 virtual void updateSolution() {
182 solution.asEigen().setZero();
185 if (
_eig.info() == Eigen::Success) {
186 _tmp.head(rank) =
_eig.eigenvectors().rightCols(rank).adjoint() * rhs;
187 _tmp.head(rank).array() /=
_eig.eigenvalues().tail(rank).array();
188 solution.asEigen() =
_eig.eigenvectors().rightCols(rank) *
_tmp.head(rank);
190 _tmp.head(rank) =
_svd.matrixU().leftCols(rank).adjoint() * rhs;
191 _tmp.head(rank).array() /=
_svd.singularValues().head(rank).array();
192 solution.asEigen() =
_svd.matrixU().leftCols(rank) *
_tmp.head(rank);
196 virtual void updateCovariance() {
198 covariance.asEigen().setZero();
201 if (
_eig.info() == Eigen::Success) {
202 covariance.asEigen() =
203 _eig.eigenvectors().rightCols(rank)
204 *
_eig.eigenvalues().tail(rank).array().inverse().matrix().asDiagonal()
205 *
_eig.eigenvectors().rightCols(rank).adjoint();
207 covariance.asEigen() =
208 _svd.matrixU().leftCols(rank)
209 *
_svd.singularValues().head(rank).array().inverse().matrix().asDiagonal()
210 *
_svd.matrixU().leftCols(rank).adjoint();
215 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd>
_eig;
216 Eigen::JacobiSVD<Eigen::MatrixXd>
_svd;
220 class CholeskySolver :
public LeastSquares::Impl {
223 explicit CholeskySolver(
int dimension) : Impl(dimension, 0.0),
_ldlt(dimension) {}
225 virtual void factor() {
226 ensure(LOWER_FISHER_MATRIX | RHS_VECTOR);
227 _ldlt.compute(fisher);
230 virtual void updateRank() {}
232 virtual void updateDiagnostic() {
235 pex::exceptions::LogicError,
236 "Can only compute NORMAL_CHOLESKY diagnostic from NORMAL_CHOLESKY factorization."
239 diagnostic.asEigen() =
_ldlt.vectorD();
242 virtual void updateSolution() { solution.asEigen() =
_ldlt.solve(rhs); }
244 virtual void updateCovariance() {
247 cov =
_ldlt.solve(cov);
254 class SvdSolver :
public LeastSquares::Impl {
257 explicit SvdSolver(
int dimension) : Impl(dimension),
_svd(),
_tmp(dimension) {}
259 virtual void factor() {
260 if (!(state & DESIGN_AND_DATA)) {
262 pex::exceptions::InvalidParameterError,
263 "Cannot initialize DIRECT_SVD solver with normal equations."
266 _svd.compute(design, Eigen::ComputeThinU | Eigen::ComputeThinV);
267 setRank(
_svd.singularValues());
268 log.debug<5>(
"Using direct SVD method; dimension=%d, rank=%d", dimension, rank);
271 virtual void updateRank() { setRank(
_svd.singularValues()); }
273 virtual void updateDiagnostic() {
274 switch (whichDiagnostic) {
276 diagnostic.asEigen<Eigen::ArrayXpr>() =
_svd.singularValues().array().square();
280 pex::exceptions::LogicError,
281 "Can only compute NORMAL_CHOLESKY diagnostic from DIRECT_SVD factorization."
284 diagnostic.asEigen() =
_svd.singularValues();
289 virtual void updateSolution() {
291 solution.asEigen().setZero();
294 _tmp.head(rank) =
_svd.matrixU().leftCols(rank).adjoint() * data;
295 _tmp.head(rank).array() /=
_svd.singularValues().head(rank).array();
296 solution.asEigen() =
_svd.matrixV().leftCols(rank) *
_tmp.head(rank);
299 virtual void updateCovariance() {
301 covariance.asEigen().setZero();
304 covariance.asEigen() =
305 _svd.matrixV().leftCols(rank)
306 *
_svd.singularValues().head(rank).array().inverse().square().matrix().asDiagonal()
307 *
_svd.matrixV().leftCols(rank).adjoint();
311 Eigen::JacobiSVD<Eigen::MatrixXd>
_svd;
312 Eigen::VectorXd
_tmp;
318 _impl->threshold = threshold;
328 return _impl->solution;
333 return _impl->covariance;
341 _impl->fisher.data(),
349 if (
_impl->whichDiagnostic != factorization) {
351 _impl->whichDiagnostic = factorization;
354 return _impl->diagnostic;
364 switch (factorization) {
366 _impl = boost::make_shared<EigensystemSolver>(dimension);
369 _impl = boost::make_shared<CholeskySolver>(dimension);
372 _impl = boost::make_shared<SvdSolver>(dimension);
375 _impl->factorization = factorization;
387 if (haveNormalEquations) {
390 pex::exceptions::InvalidParameterError,
391 (
boost::format(
"Number of rows of Fisher matrix (%d) does not match"
392 " dimension of LeastSquares solver.")
398 pex::exceptions::InvalidParameterError,
399 (
boost::format(
"Number of columns of Fisher matrix (%d) does not match"
400 " dimension of LeastSquares solver.")
406 pex::exceptions::InvalidParameterError,
407 (
boost::format(
"Number of elements in RHS vector (%d) does not match"
408 " dimension of LeastSquares solver.")
416 pex::exceptions::InvalidParameterError,
417 "Number of columns of design matrix does not match dimension of LeastSquares solver."
422 pex::exceptions::InvalidParameterError,
423 (
boost::format(
"Number of rows of design matrix (%d) does not match number of "
430 pex::exceptions::InvalidParameterError,
431 (
boost::format(
"Number of columns of design matrix (%d) must be smaller than number of "
virtual void updateRank()=0
Eigen::LDLT< Eigen::MatrixXd > _ldlt
Impl(int dimension_, double threshold_=std::numeric_limits< double >::epsilon())
Eigen::SelfAdjointEigenSolver< Eigen::MatrixXd > _eig
Eigen3 view into an ndarray::Array.
void _factor(bool haveNormalEquations)
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.
Include files required for standard LSST Exception handling.
Eigen::VectorXd & _getDataVector()
detail::ExternalInitializer< T, N, Owner > external(T *data, Vector< int, N > const &shape, Vector< int, N > const &strides, Owner const &owner)
Create an expression that initializes an Array with externally allocated memory.
ndarray::Array< double, 1, 1 > diagnostic
LeastSquares(Factorization factorization, int dimension)
Construct a least-squares object for the given factorization and dimensionality.
virtual void updateDiagnostic()=0
Use a thin singular value decomposition of the design matrix.
Use the normal equations with a Cholesky decomposition.
Factorization factorization
int getRank() const
Return the rank of the problem (number of nonzero Eigenvalues).
Use the normal equations with a symmetric Eigensystem decomposition.
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...
boost::shared_ptr< Impl > _impl
Vector< T, N > makeVector(T v1, T v2,..., T vN)
Variadic constructor for Vector.
ndarray::Array< double, 1, 1 > solution
ndarray::Array< double const, 2, 2 > getCovariance()
Return the covariance matrix of the least squares problem.
virtual void updateSolution()=0
virtual void updateCovariance()=0
detail::SimpleInitializer< N > allocate(Vector< int, N > const &shape)
Create an expression that allocates uninitialized memory for an array.
bool isEmpty() const
Return true if the array has a null data point.
#define LSST_EXCEPT(type,...)
Factorization whichDiagnostic
Factorization
Private implementation; forward-declared publicly so we can inherit from it in .cc.
Factorization getFactorization() const
Retun the type of factorization used by the solver.
ndarray::Array< double const, 1, 1 > getSolution()
Return the vector solution to the least squares problem.
Eigen::MatrixXd & _getFisherMatrix()
void setRank(Eigen::MatrixBase< D > const &values)
Eigen::MatrixXd & _getDesignMatrix()
ndarray::Array< double, 2, 2 > covariance
int getDimension() const
Return the number of parameters.
ndarray::Array< double const, 2, 2 > getFisherMatrix()
Return the Fisher matrix (inverse of the covariance) of the parameters.
Eigen::JacobiSVD< Eigen::MatrixXd > _svd
Eigen::VectorXd & _getRhsVector()