LSST Applications  21.0.0+04719a4bac,21.0.0-1-ga51b5d4+f5e6047307,21.0.0-11-g2b59f77+a9c1acf22d,21.0.0-11-ga42c5b2+86977b0b17,21.0.0-12-gf4ce030+76814010d2,21.0.0-13-g1721dae+760e7a6536,21.0.0-13-g3a573fe+768d78a30a,21.0.0-15-g5a7caf0+f21cbc5713,21.0.0-16-g0fb55c1+b60e2d390c,21.0.0-19-g4cded4ca+71a93a33c0,21.0.0-2-g103fe59+bb20972958,21.0.0-2-g45278ab+04719a4bac,21.0.0-2-g5242d73+3ad5d60fb1,21.0.0-2-g7f82c8f+8babb168e8,21.0.0-2-g8f08a60+06509c8b61,21.0.0-2-g8faa9b5+616205b9df,21.0.0-2-ga326454+8babb168e8,21.0.0-2-gde069b7+5e4aea9c2f,21.0.0-2-gecfae73+1d3a86e577,21.0.0-2-gfc62afb+3ad5d60fb1,21.0.0-25-g1d57be3cd+e73869a214,21.0.0-3-g357aad2+ed88757d29,21.0.0-3-g4a4ce7f+3ad5d60fb1,21.0.0-3-g4be5c26+3ad5d60fb1,21.0.0-3-g65f322c+e0b24896a3,21.0.0-3-g7d9da8d+616205b9df,21.0.0-3-ge02ed75+a9c1acf22d,21.0.0-4-g591bb35+a9c1acf22d,21.0.0-4-g65b4814+b60e2d390c,21.0.0-4-gccdca77+0de219a2bc,21.0.0-4-ge8a399c+6c55c39e83,21.0.0-5-gd00fb1e+05fce91b99,21.0.0-6-gc675373+3ad5d60fb1,21.0.0-64-g1122c245+4fb2b8f86e,21.0.0-7-g04766d7+cd19d05db2,21.0.0-7-gdf92d54+04719a4bac,21.0.0-8-g5674e7b+d1bd76f71f,master-gac4afde19b+a9c1acf22d,w.2021.13
LSST Data Management Base Package
LeastSquares.h
Go to the documentation of this file.
1 // -*- LSST-C++ -*-
2 
3 /*
4  * LSST Data Management System
5  * Copyright 2008-2013 LSST Corporation.
6  *
7  * This product includes software developed by the
8  * LSST Project (http://www.lsst.org/).
9  *
10  * This program is free software: you can redistribute it and/or modify
11  * it under the terms of the GNU General Public License as published by
12  * the Free Software Foundation, either version 3 of the License, or
13  * (at your option) any later version.
14  *
15  * This program is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18  * GNU General Public License for more details.
19  *
20  * You should have received a copy of the LSST License Statement and
21  * the GNU General Public License along with this program. If not,
22  * see <http://www.lsstcorp.org/LegalNotices/>.
23  */
24 
25 #ifndef LSST_AFW_MATH_LeastSquares_h_INCLUDED
26 #define LSST_AFW_MATH_LeastSquares_h_INCLUDED
27 
28 #include "lsst/base.h"
29 #include "ndarray/eigen.h"
30 
31 namespace lsst {
32 namespace afw {
33 namespace math {
34 
68 public:
69  class Impl;
70 
88  DIRECT_SVD
96  };
97 
99  template <typename T1, typename T2, int C1, int C2>
100  static LeastSquares fromDesignMatrix(ndarray::Array<T1, 2, C1> const& design,
101  ndarray::Array<T2, 1, C2> const& data,
102  Factorization factorization = NORMAL_EIGENSYSTEM) {
103  LeastSquares r(factorization, design.template getSize<1>());
104  r.setDesignMatrix(design, data);
105  return r;
106  }
107 
109  template <typename D1, typename D2>
110  static LeastSquares fromDesignMatrix(Eigen::MatrixBase<D1> const& design,
111  Eigen::MatrixBase<D2> const& data,
112  Factorization factorization = NORMAL_EIGENSYSTEM) {
113  LeastSquares r(factorization, design.cols());
114  r.setDesignMatrix(design, data);
115  return r;
116  }
117 
119  template <typename T1, typename T2, int C1, int C2>
120  void setDesignMatrix(ndarray::Array<T1, 2, C1> const& design, ndarray::Array<T2, 1, C2> const& data) {
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  }
129 
131  template <typename D1, typename D2>
132  void setDesignMatrix(Eigen::MatrixBase<D1> const& design, Eigen::MatrixBase<D2> const& data) {
133  _getDesignMatrix() = design.template cast<double>();
134  _getDataVector() = data.template cast<double>();
135  _factor(false);
136  }
137 
139  template <typename T1, int C1>
140  void setDesignMatrix(ndarray::Array<T1, 2, C1> const& design) {
141  _getDesignMatrix() = ndarray::asEigenMatrix(design).template cast<double>();
142  _factor(false);
143  }
144 
146  template <typename D1, typename D2>
147  void setDesignMatrix(Eigen::MatrixBase<D1> const& design) {
148  _getDesignMatrix() = design.template cast<double>();
149  _factor(false);
150  }
151 
153  template <typename T1, typename T2, int C1, int C2>
154  static LeastSquares fromNormalEquations(ndarray::Array<T1, 2, C1> const& fisher,
155  ndarray::Array<T2, 1, C2> const& rhs,
156  Factorization factorization = NORMAL_EIGENSYSTEM) {
157  LeastSquares r(factorization, fisher.template getSize<0>());
158  r.setNormalEquations(fisher, rhs);
159  return r;
160  }
161 
163  template <typename D1, typename D2>
164  static LeastSquares fromNormalEquations(Eigen::MatrixBase<D1> const& fisher,
165  Eigen::MatrixBase<D2> const& rhs,
166  Factorization factorization = NORMAL_EIGENSYSTEM) {
167  LeastSquares r(factorization, fisher.rows());
168  r.setNormalEquations(fisher, rhs);
169  return r;
170  }
171 
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>();
177  } else {
178  _getFisherMatrix() = ndarray::asEigenMatrix(fisher).transpose().template cast<double>();
179  }
180  _getRhsVector() = ndarray::asEigenMatrix(rhs).template cast<double>();
181  _factor(true);
182  }
183 
185  template <typename D1, typename D2>
186  void setNormalEquations(Eigen::MatrixBase<D1> const& fisher, Eigen::MatrixBase<D2> const& rhs) {
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  }
195 
214  void setThreshold(double threshold);
215 
217  double getThreshold() const;
218 
231  ndarray::Array<double const, 1, 1> getSolution();
232 
245  ndarray::Array<double const, 2, 2> getCovariance();
246 
259  ndarray::Array<double const, 2, 2> getFisherMatrix();
260 
286  ndarray::Array<double const, 1, 1> getDiagnostic(Factorization factorization);
287 
289  int getDimension() const;
290 
297  int getRank() const;
298 
301 
308  LeastSquares(Factorization factorization, int dimension);
309 
314 
315  // Need to define dtor in source file so it can see Impl declaration.
317 
318 private:
319  // We want a column-major design matrix so the self-adjoint product is cache-friendly, so we always
320  // copy a (possibly row-major) design matrix into a col-major one. This is an unnecessarily and
321  // cache-unfriendly operation when solver is DIRECT_SVD, but right now it doesn't seem to be worth
322  // special-casing the design for that case. In other cases it's a cache-unfriendly op that avoids
323  // an even worse one, and it can always be avoided by using a column-major design matrix.
324  Eigen::MatrixXd& _getDesignMatrix();
325  Eigen::VectorXd& _getDataVector();
326 
327  // Storage order matters less here, so we just go with what Eigen is most accustomed to.
328  // Because the Fisher matrix is symmetric, we can just transpose it before copying to avoid doing
329  // any expensive copies between different storage orders.
330  Eigen::MatrixXd& _getFisherMatrix();
331  Eigen::VectorXd& _getRhsVector();
332 
333  // Do the actual high-level linear algrebra factorization; the appropriate input matrices and
334  // vectors must already be set before this is called. The solution and other quantities
335  // will not be computed until they are first requested.
336  void _factor(bool haveNormalEquations);
337 
338  std::shared_ptr<Impl> _impl;
339 };
340 } // namespace math
341 } // namespace afw
342 } // namespace lsst
343 
344 #endif // !LSST_AFW_MATH_LeastSquares_h_INCLUDED
char * data
Definition: BaseRecord.cc:62
Basic LSST definitions.
Solver for linear least-squares problems.
Definition: LeastSquares.h:67
LeastSquares & operator=(LeastSquares &&)
void setDesignMatrix(ndarray::Array< T1, 2, C1 > const &design)
Reset the design matrix given as an ndarray; dimension and data are not changed.
Definition: LeastSquares.h:140
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.
Definition: LeastSquares.h:100
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.
Definition: LeastSquares.h:71
@ 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
LeastSquares(Factorization factorization, int dimension)
Construct a least-squares object for the given factorization and dimensionality.
int getDimension() const
Return the number of parameters.
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.
Definition: LeastSquares.h:132
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.
Definition: LeastSquares.h:164
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...
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.
Definition: LeastSquares.h:186
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.
Definition: LeastSquares.h:120
double getThreshold() const
Get the threshold used to determine when to truncate Eigenvalues.
LeastSquares(LeastSquares &&)
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.
Definition: LeastSquares.h:174
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.
Definition: LeastSquares.h:110
void setThreshold(double threshold)
Set the threshold used to determine when to truncate Eigenvalues.
void setDesignMatrix(Eigen::MatrixBase< D1 > const &design)
Reset the design matrix given as an Eigen object; dimension and data are not changed.
Definition: LeastSquares.h:147
LeastSquares(LeastSquares const &)
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.
Definition: LeastSquares.h:154
class[[deprecated("Removed with no replacement (but see lsst::afw::image::TransmissionCurve). Will be " "removed after v22.")]] FilterProperty final
Describe the properties of a Filter (e.g.
Definition: Filter.h:53
A base class for image defects.