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();
132 explicit EigensystemSolver(
int dimension) : Impl(dimension), _eig(dimension), _svd(), _tmp(dimension) {}
134 void factor()
override {
135 ensure(LOWER_FISHER_MATRIX | RHS_VECTOR);
136 _eig.compute(fisher);
137 if (_eig.info() == Eigen::Success) {
138 setRank(_eig.eigenvalues().reverse());
139 LOGL_DEBUG(_log,
"SelfAdjointEigenSolver succeeded: dimension=%d, rank=%d", dimension, rank);
144 ensure(FULL_FISHER_MATRIX);
145 _svd.compute(fisher, Eigen::ComputeFullU);
146 setRank(_svd.singularValues());
148 "SelfAdjointEigenSolver failed; falling back to equivalent SVD: dimension=%d, rank=%d",
153 void updateRank()
override {
154 if (_eig.info() == Eigen::Success) {
155 setRank(_eig.eigenvalues().reverse());
157 setRank(_svd.singularValues());
161 void updateDiagnostic()
override {
164 pex::exceptions::LogicError,
165 "Cannot compute NORMAL_CHOLESKY diagnostic from NORMAL_EIGENSYSTEM factorization.");
167 if (_eig.info() == Eigen::Success) {
168 ndarray::asEigenMatrix(diagnostic) = _eig.eigenvalues().reverse();
170 ndarray::asEigenMatrix(diagnostic) = _svd.singularValues();
173 ndarray::asEigenArray(diagnostic) = ndarray::asEigenArray(diagnostic).sqrt();
177 void updateSolution()
override {
179 ndarray::asEigenMatrix(solution).setZero();
182 if (_eig.info() == Eigen::Success) {
183 _tmp.head(rank) = _eig.eigenvectors().rightCols(rank).adjoint() * rhs;
184 _tmp.head(rank).array() /= _eig.eigenvalues().tail(rank).array();
185 ndarray::asEigenMatrix(solution) = _eig.eigenvectors().rightCols(rank) * _tmp.head(rank);
187 _tmp.head(rank) = _svd.matrixU().leftCols(rank).adjoint() * rhs;
188 _tmp.head(rank).array() /= _svd.singularValues().head(rank).array();
189 ndarray::asEigenMatrix(solution) = _svd.matrixU().leftCols(rank) * _tmp.head(rank);
193 void updateCovariance()
override {
198 if (_eig.info() == Eigen::Success) {
200 _eig.eigenvectors().rightCols(rank) *
201 _eig.eigenvalues().tail(rank).array().inverse().matrix().asDiagonal() *
202 _eig.eigenvectors().rightCols(rank).adjoint();
205 _svd.matrixU().leftCols(rank) *
206 _svd.singularValues().head(rank).array().inverse().matrix().asDiagonal() *
207 _svd.matrixU().leftCols(rank).adjoint();
212 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> _eig;
213 Eigen::JacobiSVD<Eigen::MatrixXd> _svd;
214 Eigen::VectorXd _tmp;
217 class CholeskySolver :
public LeastSquares::Impl {
219 explicit CholeskySolver(
int dimension) : Impl(dimension, 0.0), _ldlt(dimension) {}
221 void factor()
override {
222 ensure(LOWER_FISHER_MATRIX | RHS_VECTOR);
223 _ldlt.compute(fisher);
226 void updateRank()
override {}
228 void updateDiagnostic()
override {
231 pex::exceptions::LogicError,
232 "Can only compute NORMAL_CHOLESKY diagnostic from NORMAL_CHOLESKY factorization.");
234 ndarray::asEigenMatrix(diagnostic) = _ldlt.vectorD();
237 void updateSolution()
override { ndarray::asEigenMatrix(solution) = _ldlt.solve(rhs); }
239 void updateCovariance()
override {
240 auto cov = ndarray::asEigenMatrix(
covariance);
242 cov = _ldlt.solve(cov);
246 Eigen::LDLT<Eigen::MatrixXd> _ldlt;
249 class SvdSolver :
public LeastSquares::Impl {
251 explicit SvdSolver(
int dimension) : Impl(dimension), _svd(), _tmp(dimension) {}
253 void factor()
override {
254 if (!(state & DESIGN_AND_DATA)) {
255 throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
256 "Cannot initialize DIRECT_SVD solver with normal equations.");
258 _svd.compute(design, Eigen::ComputeThinU | Eigen::ComputeThinV);
259 setRank(_svd.singularValues());
260 LOGL_DEBUG(_log,
"Using direct SVD method; dimension=%d, rank=%d", dimension, rank);
263 void updateRank()
override { setRank(_svd.singularValues()); }
265 void updateDiagnostic()
override {
266 switch (whichDiagnostic) {
268 ndarray::asEigenArray(diagnostic) = _svd.singularValues().array().square();
272 pex::exceptions::LogicError,
273 "Can only compute NORMAL_CHOLESKY diagnostic from DIRECT_SVD factorization.");
275 ndarray::asEigenMatrix(diagnostic) = _svd.singularValues();
280 void updateSolution()
override {
282 ndarray::asEigenMatrix(solution).setZero();
285 _tmp.head(rank) = _svd.matrixU().leftCols(rank).adjoint() *
data;
286 _tmp.head(rank).array() /= _svd.singularValues().head(rank).array();
287 ndarray::asEigenMatrix(solution) = _svd.matrixV().leftCols(rank) * _tmp.head(rank);
290 void updateCovariance()
override {
296 _svd.matrixV().leftCols(rank) *
297 _svd.singularValues().head(rank).array().inverse().square().matrix().asDiagonal() *
298 _svd.matrixV().leftCols(rank).adjoint();
302 Eigen::JacobiSVD<Eigen::MatrixXd> _svd;
303 Eigen::VectorXd _tmp;
309 _impl->threshold = threshold;
310 _impl->state &= ~
Impl::SOLUTION_ARRAY;
311 _impl->state &= ~
Impl::COVARIANCE_ARRAY;
319 return _impl->solution;
324 return _impl->covariance;
331 return ndarray::external(_impl->fisher.data(), ndarray::makeVector(_impl->dimension, _impl->dimension),
332 ndarray::makeVector(_impl->dimension, 1), _impl);
336 if (_impl->whichDiagnostic != factorization) {
337 _impl->state &= ~
Impl::DIAGNOSTIC_ARRAY;
338 _impl->whichDiagnostic = factorization;
341 return _impl->diagnostic;
351 switch (factorization) {
353 _impl = std::make_shared<EigensystemSolver>(dimension);
356 _impl = std::make_shared<CholeskySolver>(dimension);
359 _impl = std::make_shared<SvdSolver>(dimension);
362 _impl->factorization = factorization;
372 Eigen::MatrixXd& LeastSquares::_getDesignMatrix() {
return _impl->design; }
373 Eigen::VectorXd& LeastSquares::_getDataVector() {
return _impl->data; }
375 Eigen::MatrixXd& LeastSquares::_getFisherMatrix() {
return _impl->fisher; }
376 Eigen::VectorXd& LeastSquares::_getRhsVector() {
return _impl->rhs; }
378 void LeastSquares::_factor(
bool haveNormalEquations) {
379 if (haveNormalEquations) {
380 if (_getFisherMatrix().rows() != _impl->dimension) {
381 throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
382 (
boost::format(
"Number of rows of Fisher matrix (%d) does not match"
383 " dimension of LeastSquares solver.") %
384 _getFisherMatrix().rows() % _impl->dimension)
387 if (_getFisherMatrix().cols() != _impl->dimension) {
388 throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
389 (
boost::format(
"Number of columns of Fisher matrix (%d) does not match"
390 " dimension of LeastSquares solver.") %
391 _getFisherMatrix().cols() % _impl->dimension)
394 if (_getRhsVector().size() != _impl->dimension) {
395 throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
396 (
boost::format(
"Number of elements in RHS vector (%d) does not match"
397 " dimension of LeastSquares solver.") %
398 _getRhsVector().size() % _impl->dimension)
403 if (_getDesignMatrix().cols() != _impl->dimension) {
405 pex::exceptions::InvalidParameterError,
406 "Number of columns of design matrix does not match dimension of LeastSquares solver.");
408 if (_getDesignMatrix().rows() != _getDataVector().size()) {
409 throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
410 (
boost::format(
"Number of rows of design matrix (%d) does not match number of "
411 "data points (%d)") %
412 _getDesignMatrix().rows() % _getDataVector().size())
415 if (_getDesignMatrix().cols() > _getDataVector().size()) {
417 pex::exceptions::InvalidParameterError,
418 (
boost::format(
"Number of columns of design matrix (%d) must be smaller than number of "
419 "data points (%d)") %
420 _getDesignMatrix().cols() % _getDataVector().size())