LSST Applications  21.0.0-172-gfb10e10a+18fedfabac,22.0.0+297cba6710,22.0.0+80564b0ff1,22.0.0+8d77f4f51a,22.0.0+a28f4c53b1,22.0.0+dcf3732eb2,22.0.1-1-g7d6de66+2a20fdde0d,22.0.1-1-g8e32f31+297cba6710,22.0.1-1-geca5380+7fa3b7d9b6,22.0.1-12-g44dc1dc+2a20fdde0d,22.0.1-15-g6a90155+515f58c32b,22.0.1-16-g9282f48+790f5f2caa,22.0.1-2-g92698f7+dcf3732eb2,22.0.1-2-ga9b0f51+7fa3b7d9b6,22.0.1-2-gd1925c9+bf4f0e694f,22.0.1-24-g1ad7a390+a9625a72a8,22.0.1-25-g5bf6245+3ad8ecd50b,22.0.1-25-gb120d7b+8b5510f75f,22.0.1-27-g97737f7+2a20fdde0d,22.0.1-32-gf62ce7b1+aa4237961e,22.0.1-4-g0b3f228+2a20fdde0d,22.0.1-4-g243d05b+871c1b8305,22.0.1-4-g3a563be+32dcf1063f,22.0.1-4-g44f2e3d+9e4ab0f4fa,22.0.1-42-gca6935d93+ba5e5ca3eb,22.0.1-5-g15c806e+85460ae5f3,22.0.1-5-g58711c4+611d128589,22.0.1-5-g75bb458+99c117b92f,22.0.1-6-g1c63a23+7fa3b7d9b6,22.0.1-6-g50866e6+84ff5a128b,22.0.1-6-g8d3140d+720564cf76,22.0.1-6-gd805d02+cc5644f571,22.0.1-8-ge5750ce+85460ae5f3,master-g6e05de7fdc+babf819c66,master-g99da0e417a+8d77f4f51a,w.2021.48
LSST Data Management Base Package
Classes | Public Types | Public Member Functions | Static Public Member Functions | List of all members
lsst::afw::math::LeastSquares Class Referencefinal

Solver for linear least-squares problems. More...

#include <LeastSquares.h>

Classes

class  Impl
 

Public Types

enum  Factorization { NORMAL_EIGENSYSTEM , NORMAL_CHOLESKY , DIRECT_SVD }
 Private implementation; forward-declared publicly so we can inherit from it in .cc. More...
 

Public Member Functions

template<typename T1 , typename T2 , int C1, int C2>
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. More...
 
template<typename D1 , typename D2 >
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. More...
 
template<typename T1 , int C1>
void setDesignMatrix (ndarray::Array< T1, 2, C1 > const &design)
 Reset the design matrix given as an ndarray; dimension and data are not changed. More...
 
template<typename D1 , typename D2 >
void setDesignMatrix (Eigen::MatrixBase< D1 > const &design)
 Reset the design matrix given as an Eigen object; dimension and data are not changed. More...
 
template<typename T1 , typename T2 , int C1, int C2>
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. More...
 
template<typename D1 , typename D2 >
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. More...
 
void setThreshold (double threshold)
 Set the threshold used to determine when to truncate Eigenvalues. More...
 
double getThreshold () const
 Get the threshold used to determine when to truncate Eigenvalues. More...
 
ndarray::Array< double const, 1, 1 > getSolution ()
 Return the vector solution to the least squares problem. More...
 
ndarray::Array< double const, 2, 2 > getCovariance ()
 Return the covariance matrix of the least squares problem. More...
 
ndarray::Array< double const, 2, 2 > getFisherMatrix ()
 Return the Fisher matrix (inverse of the covariance) of the parameters. More...
 
ndarray::Array< double const, 1, 1 > getDiagnostic (Factorization factorization)
 Return a factorization-dependent vector that can be used to characterize the stability of the solution. More...
 
int getDimension () const
 Return the number of parameters. More...
 
int getRank () const
 Return the rank of the problem (number of nonzero Eigenvalues). More...
 
Factorization getFactorization () const
 Retun the type of factorization used by the solver. More...
 
 LeastSquares (Factorization factorization, int dimension)
 Construct a least-squares object for the given factorization and dimensionality. More...
 
 LeastSquares (LeastSquares const &)
 
 LeastSquares (LeastSquares &&)
 
LeastSquaresoperator= (LeastSquares const &)
 
LeastSquaresoperator= (LeastSquares &&)
 
 ~LeastSquares ()
 

Static Public Member Functions

template<typename T1 , typename T2 , int C1, int C2>
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. More...
 
template<typename D1 , typename D2 >
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. More...
 
template<typename T1 , typename T2 , int C1, int C2>
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. More...
 
template<typename D1 , typename D2 >
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. More...
 

Detailed Description

Solver for linear least-squares problems.

Linear least-squares problems are defined as finding the vector \(x\) that minimizes \(\left|A x - b\right|_2\), with the number of rows of \(A\) generally greater than the number of columns. We call \(A\) the design matrix, \(b\) the data vector, and \(x\) the solution vector. When the rank of \(A\) is equal to the number of columns, we can obtain using the solution using the normal equations:

\[ A^T A x = A^T b \]

If \(A\) is not full-rank, the problem is underconstrained, and we usually wish to solve the minimum-norm least-squares problem, which also minimizes \(|x|_2\). This can be done by computing a pseudo-inverse of \(A\) using an SVD or complete orthogonal factorization of \(A\), or by performing an Eigen decomposition of \(A^T A\).

This class can be constructed from the design matrix and data vector, or from the two terms in the normal equations (below, we call the matrix \(A^TA\) the Fisher matrix, and the vector \(A^T b\) simply the "right-hand side" (RHS) vector). If initialized with the design matrix and data vector, we can still use the normal equations to solve it. The solution via the normal equations is more susceptible to round-off error, but it is also faster, and if the normal equation terms can be computed directly it can be significantly less expensive in terms of memory. The Fisher matrix is a symmetric matrix, and it should be exactly symmetric when provided as input, because which triangle will be used is an implementation detail that is subject to change and may depend on the storage order.

The solver always operates in double precision, and returns all results in double precision. However, it can be initialized from single precision inputs. It isn't usually a good idea to construct the normal equations in single precision, however, even when the data are single precision.

Definition at line 67 of file LeastSquares.h.

Member Enumeration Documentation

◆ Factorization

Private implementation; forward-declared publicly so we can inherit from it in .cc.

Enumerator
NORMAL_EIGENSYSTEM 

Use the normal equations with a symmetric Eigensystem decomposition.

This method is fully robust and computes the minimum-norm solution when the problem does not have full rank. It is affected slightly more by round-off error than the SVD method, but because it can handle singular matrices this usually isn't a problem.

NORMAL_CHOLESKY 

Use the normal equations with a Cholesky decomposition.

While this method uses a robust LDL^T decomposition that does not require square roots, it is not appropriate for problems that do not have full rank, and cannot be used to determine whether a problem has full rank. It is the fastest decomposition.

DIRECT_SVD 

Use a thin singular value decomposition of the design matrix.

This method is the most robust and computes the minimum-norm solution the problem does not have full rank. However, it is also the slowest and is not available when the solver is initialized with the normal equations.

Definition at line 71 of file LeastSquares.h.

71  {
88  DIRECT_SVD
96  };
@ NORMAL_EIGENSYSTEM
Use the normal equations with a symmetric Eigensystem decomposition.
Definition: LeastSquares.h:72
@ NORMAL_CHOLESKY
Use the normal equations with a Cholesky decomposition.
Definition: LeastSquares.h:80
@ DIRECT_SVD
Use a thin singular value decomposition of the design matrix.
Definition: LeastSquares.h:88

Constructor & Destructor Documentation

◆ LeastSquares() [1/3]

lsst::afw::math::LeastSquares::LeastSquares ( Factorization  factorization,
int  dimension 
)

Construct a least-squares object for the given factorization and dimensionality.

One of the set* member functions must be called before any other operations can be performed on a LeastSquares object initialized this way.

Definition at line 350 of file LeastSquares.cc.

350  {
351  switch (factorization) {
352  case NORMAL_EIGENSYSTEM:
353  _impl = std::make_shared<EigensystemSolver>(dimension);
354  break;
355  case NORMAL_CHOLESKY:
356  _impl = std::make_shared<CholeskySolver>(dimension);
357  break;
358  case DIRECT_SVD:
359  _impl = std::make_shared<SvdSolver>(dimension);
360  break;
361  }
362  _impl->factorization = factorization;
363 }

◆ LeastSquares() [2/3]

lsst::afw::math::LeastSquares::LeastSquares ( LeastSquares const &  )
default

◆ LeastSquares() [3/3]

lsst::afw::math::LeastSquares::LeastSquares ( LeastSquares &&  )
default

◆ ~LeastSquares()

lsst::afw::math::LeastSquares::~LeastSquares ( )
default

Member Function Documentation

◆ fromDesignMatrix() [1/2]

template<typename D1 , typename D2 >
static LeastSquares lsst::afw::math::LeastSquares::fromDesignMatrix ( Eigen::MatrixBase< D1 > const &  design,
Eigen::MatrixBase< D2 > const &  data,
Factorization  factorization = NORMAL_EIGENSYSTEM 
)
inlinestatic

Initialize from the design matrix and data vector given as an Eigen objects.

Definition at line 110 of file LeastSquares.h.

112  {
113  LeastSquares r(factorization, design.cols());
114  r.setDesignMatrix(design, data);
115  return r;
116  }
char * data
Definition: BaseRecord.cc:61
LeastSquares(Factorization factorization, int dimension)
Construct a least-squares object for the given factorization and dimensionality.

◆ fromDesignMatrix() [2/2]

template<typename T1 , typename T2 , int C1, int C2>
static LeastSquares lsst::afw::math::LeastSquares::fromDesignMatrix ( ndarray::Array< T1, 2, C1 > const &  design,
ndarray::Array< T2, 1, C2 > const &  data,
Factorization  factorization = NORMAL_EIGENSYSTEM 
)
inlinestatic

Initialize from the design matrix and data vector given as ndarrays.

Definition at line 100 of file LeastSquares.h.

102  {
103  LeastSquares r(factorization, design.template getSize<1>());
104  r.setDesignMatrix(design, data);
105  return r;
106  }

◆ fromNormalEquations() [1/2]

template<typename D1 , typename D2 >
static LeastSquares lsst::afw::math::LeastSquares::fromNormalEquations ( Eigen::MatrixBase< D1 > const &  fisher,
Eigen::MatrixBase< D2 > const &  rhs,
Factorization  factorization = NORMAL_EIGENSYSTEM 
)
inlinestatic

Initialize from the terms in the normal equations, given as Eigen objects.

Definition at line 164 of file LeastSquares.h.

166  {
167  LeastSquares r(factorization, fisher.rows());
168  r.setNormalEquations(fisher, rhs);
169  return r;
170  }

◆ fromNormalEquations() [2/2]

template<typename T1 , typename T2 , int C1, int C2>
static LeastSquares lsst::afw::math::LeastSquares::fromNormalEquations ( ndarray::Array< T1, 2, C1 > const &  fisher,
ndarray::Array< T2, 1, C2 > const &  rhs,
Factorization  factorization = NORMAL_EIGENSYSTEM 
)
inlinestatic

Initialize from the terms in the normal equations, given as ndarrays.

Definition at line 154 of file LeastSquares.h.

156  {
157  LeastSquares r(factorization, fisher.template getSize<0>());
158  r.setNormalEquations(fisher, rhs);
159  return r;
160  }

◆ getCovariance()

ndarray::Array< double const, 2, 2 > lsst::afw::math::LeastSquares::getCovariance ( )

Return the covariance matrix of the least squares problem.

The returned array is owned by the LeastSquares object and may be modified in-place by future calls to LeastSquares member functions, so it's best to promptly copy the result elsewhere.

If you want an Eigen object instead, just use ndarray::asEigenMatrix(getCovariance()).

The covariance is cached the first time this member function is called, and will be reused unless the matrices are reset or the threshold is changed.

Definition at line 322 of file LeastSquares.cc.

322  {
324  return _impl->covariance;
325 }
ndarray::Array< double, 2, 2 > covariance
Definition: LeastSquares.cc:68

◆ getDiagnostic()

ndarray::Array< double const, 1, 1 > lsst::afw::math::LeastSquares::getDiagnostic ( Factorization  factorization)

Return a factorization-dependent vector that can be used to characterize the stability of the solution.

The returned array's size is always equal to getDimension(). When getRank() is less than getDimension(), some elements of the array were considered effectively zero (see setThreshold).

For the NORMAL_EIGENSYSTEM method, this is the vector of Eigenvalues of the Fisher matrix, including those rejected as being below the threshold, sorted in descending order (all values are nonnegative). This is the square of the singular values of the design matrix, and is only available when the factorization is either NORMAL_EIGENSYSTEM or DIRECT_SVD.

For the DIRECT_SVD method, this is the vector of singular values of the design matrix, including those rejected as being below the threshold, sorted in descending order (all values are nonnegative). This is the square root of the Eigenvalues of the Fisher matrix, and is only available when the factorization is either NORMAL_EIGENSYSTEM or DIRECT_SVD.

For the NORMAL_CHOLESKY method, this is \(D\) in the pivoted Cholesky factorization \(P L D L^T P^T\) of the Fisher matrix. This does not provide a reliable way to test the stability of the problem, but it does provide a way to compute the determinant of the Fisher matrix. It is only available when the factorization is NORMAL_CHOLESKY.

Definition at line 335 of file LeastSquares.cc.

335  {
336  if (_impl->whichDiagnostic != factorization) {
337  _impl->state &= ~Impl::DIAGNOSTIC_ARRAY;
338  _impl->whichDiagnostic = factorization;
339  }
341  return _impl->diagnostic;
342 }
ndarray::Array< double, 1, 1 > diagnostic
Definition: LeastSquares.cc:69

◆ getDimension()

int lsst::afw::math::LeastSquares::getDimension ( ) const

Return the number of parameters.

Definition at line 344 of file LeastSquares.cc.

344 { return _impl->dimension; }

◆ getFactorization()

LeastSquares::Factorization lsst::afw::math::LeastSquares::getFactorization ( ) const

Retun the type of factorization used by the solver.

Definition at line 348 of file LeastSquares.cc.

348 { return _impl->factorization; }

◆ getFisherMatrix()

ndarray::Array< double const, 2, 2 > lsst::afw::math::LeastSquares::getFisherMatrix ( )

Return the Fisher matrix (inverse of the covariance) of the parameters.

Note that the Fisher matrix is exactly the same as the matrix on the lhs of the normal equations.

The returned array is owned by the LeastSquares object and may be modified in-place by future calls to LeastSquares member functions, so it's best to promptly copy the result elsewhere.

If you want an Eigen object instead, just use ndarray::asEigenMatrix(getFisherMatrix()).

Definition at line 327 of file LeastSquares.cc.

327  {
329  // Wrap the Eigen::MatrixXd in an ndarray::Array, using _impl as the reference-counted owner.
330  // Doesn't matter if we swap strides, because it's symmetric.
331  return ndarray::external(_impl->fisher.data(), ndarray::makeVector(_impl->dimension, _impl->dimension),
332  ndarray::makeVector(_impl->dimension, 1), _impl);
333 }

◆ getRank()

int lsst::afw::math::LeastSquares::getRank ( ) const

Return the rank of the problem (number of nonzero Eigenvalues).

The returned value is always the same as getDimension() when the factorization is NORMAL_CHOLESKY (which may be incorrect, because a Cholesky decomposition is not rank-revealing).

Definition at line 346 of file LeastSquares.cc.

346 { return _impl->rank; }

◆ getSolution()

ndarray::Array< double const, 1, 1 > lsst::afw::math::LeastSquares::getSolution ( )

Return the vector solution to the least squares problem.

The returned array is owned by the LeastSquares object and may be modified in-place by future calls to LeastSquares member functions, so it's best to promptly copy the result elsewhere.

If you want an Eigen object instead, just use ndarray::asEigenMatrix(getSolution()).

The solution is cached the first time this member function is called, and will be reused unless the matrices are reset or the threshold is changed.

Definition at line 317 of file LeastSquares.cc.

317  {
319  return _impl->solution;
320 }
ndarray::Array< double, 1, 1 > solution
Definition: LeastSquares.cc:67

◆ getThreshold()

double lsst::afw::math::LeastSquares::getThreshold ( ) const

Get the threshold used to determine when to truncate Eigenvalues.

Definition at line 315 of file LeastSquares.cc.

315 { return _impl->threshold; }

◆ operator=() [1/2]

LeastSquares & lsst::afw::math::LeastSquares::operator= ( LeastSquares &&  )
default

◆ operator=() [2/2]

LeastSquares & lsst::afw::math::LeastSquares::operator= ( LeastSquares const &  )
default

◆ setDesignMatrix() [1/4]

template<typename D1 , typename D2 >
void lsst::afw::math::LeastSquares::setDesignMatrix ( Eigen::MatrixBase< D1 > const &  design)
inline

Reset the design matrix given as an Eigen object; dimension and data are not changed.

Definition at line 147 of file LeastSquares.h.

147  {
148  _getDesignMatrix() = design.template cast<double>();
149  _factor(false);
150  }

◆ setDesignMatrix() [2/4]

template<typename D1 , typename D2 >
void lsst::afw::math::LeastSquares::setDesignMatrix ( Eigen::MatrixBase< D1 > const &  design,
Eigen::MatrixBase< D2 > const &  data 
)
inline

Reset the design matrix and data vector given as Eigen objects; dimension must not change.

Definition at line 132 of file LeastSquares.h.

132  {
133  _getDesignMatrix() = design.template cast<double>();
134  _getDataVector() = data.template cast<double>();
135  _factor(false);
136  }

◆ setDesignMatrix() [3/4]

template<typename T1 , int C1>
void lsst::afw::math::LeastSquares::setDesignMatrix ( ndarray::Array< T1, 2, C1 > const &  design)
inline

Reset the design matrix given as an ndarray; dimension and data are not changed.

Definition at line 140 of file LeastSquares.h.

140  {
141  _getDesignMatrix() = ndarray::asEigenMatrix(design).template cast<double>();
142  _factor(false);
143  }

◆ setDesignMatrix() [4/4]

template<typename T1 , typename T2 , int C1, int C2>
void lsst::afw::math::LeastSquares::setDesignMatrix ( ndarray::Array< T1, 2, C1 > const &  design,
ndarray::Array< T2, 1, C2 > const &  data 
)
inline

Reset the design matrix and data vector given as ndarrays; dimension must not change.

Definition at line 120 of file LeastSquares.h.

120  {
121  // n.b. "template cast<T>" below is not a special kind of cast; it's just
122  // the weird C++ syntax required for calling a templated member function
123  // called "cast" in this context; see
124  // http://eigen.tuxfamily.org/dox-devel/TopicTemplateKeyword.html
125  _getDesignMatrix() = ndarray::asEigenMatrix(design).template cast<double>();
126  _getDataVector() = ndarray::asEigenMatrix(data).template cast<double>();
127  _factor(false);
128  }

◆ setNormalEquations() [1/2]

template<typename D1 , typename D2 >
void lsst::afw::math::LeastSquares::setNormalEquations ( Eigen::MatrixBase< D1 > const &  fisher,
Eigen::MatrixBase< D2 > const &  rhs 
)
inline

Reset the terms in the normal equations given as Eigen objects; dimension must not change.

Definition at line 186 of file LeastSquares.h.

186  {
187  if (bool(Eigen::MatrixBase<D1>::IsRowMajor) == bool(Eigen::MatrixXd::IsRowMajor)) {
188  _getFisherMatrix() = fisher.template cast<double>();
189  } else {
190  _getFisherMatrix() = fisher.transpose().template cast<double>();
191  }
192  _getRhsVector() = rhs.template cast<double>();
193  _factor(true);
194  }

◆ setNormalEquations() [2/2]

template<typename T1 , typename T2 , int C1, int C2>
void lsst::afw::math::LeastSquares::setNormalEquations ( ndarray::Array< T1, 2, C1 > const &  fisher,
ndarray::Array< T2, 1, C2 > const &  rhs 
)
inline

Reset the terms in the normal equations given as ndarrays; dimension must not change.

Definition at line 174 of file LeastSquares.h.

174  {
175  if ((C1 > 0) == bool(Eigen::MatrixXd::IsRowMajor)) {
176  _getFisherMatrix() = ndarray::asEigenMatrix(fisher).template cast<double>();
177  } else {
178  _getFisherMatrix() = ndarray::asEigenMatrix(fisher).transpose().template cast<double>();
179  }
180  _getRhsVector() = ndarray::asEigenMatrix(rhs).template cast<double>();
181  _factor(true);
182  }

◆ setThreshold()

void lsst::afw::math::LeastSquares::setThreshold ( double  threshold)

Set the threshold used to determine when to truncate Eigenvalues.

The rank of the matrix is determined by comparing the product of this threshold and the first (largest) element of the array returned by getDiagnostic() to all other elements of that array. Elements smaller than this product are ignored and reduce the rank.

In the DIRECT_SVD case, the diagnostic array contains the singular values of the design matrix, while in the NORMAL_EIGENSYSTEM case the diagnostic array holds the Eigenvalues of the Fisher matrix, and the latter are the square of the former. The default threshold is the same (std::numeric_limits<double>::epsilon()) in both cases, reflecting the fact that using the normal equations squares the condition number of the problem.

The NORMAL_CHOLESKY method does not use the threshold and assumes the problem is full-rank.

Definition at line 308 of file LeastSquares.cc.

308  {
309  _impl->threshold = threshold;
310  _impl->state &= ~Impl::SOLUTION_ARRAY;
311  _impl->state &= ~Impl::COVARIANCE_ARRAY;
312  _impl->updateRank();
313 }

The documentation for this class was generated from the following files: