Imported Upstream version ceres 1.13.0
[platform/upstream/ceres-solver.git] / internal / ceres / dense_qr_solver.h
1 // Ceres Solver - A fast non-linear least squares minimizer
2 // Copyright 2015 Google Inc. All rights reserved.
3 // http://ceres-solver.org/
4 //
5 // Redistribution and use in source and binary forms, with or without
6 // modification, are permitted provided that the following conditions are met:
7 //
8 // * Redistributions of source code must retain the above copyright notice,
9 //   this list of conditions and the following disclaimer.
10 // * Redistributions in binary form must reproduce the above copyright notice,
11 //   this list of conditions and the following disclaimer in the documentation
12 //   and/or other materials provided with the distribution.
13 // * Neither the name of Google Inc. nor the names of its contributors may be
14 //   used to endorse or promote products derived from this software without
15 //   specific prior written permission.
16 //
17 // THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
18 // AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
19 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
20 // ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
21 // LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
22 // CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
23 // SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
24 // INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
25 // CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
26 // ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
27 // POSSIBILITY OF SUCH DAMAGE.
28 //
29 // Author: sameeragarwal@google.com (Sameer Agarwal)
30 //
31 // Solve dense rectangular systems Ax = b using the QR factorization.
32 #ifndef CERES_INTERNAL_DENSE_QR_SOLVER_H_
33 #define CERES_INTERNAL_DENSE_QR_SOLVER_H_
34
35 #include "ceres/linear_solver.h"
36 #include "ceres/internal/eigen.h"
37 #include "ceres/internal/macros.h"
38
39 namespace ceres {
40 namespace internal {
41
42 class DenseSparseMatrix;
43
44 // This class implements the LinearSolver interface for solving
45 // rectangular/unsymmetric (well constrained) linear systems of the
46 // form
47 //
48 //   Ax = b
49 //
50 // Since there does not usually exist a solution that satisfies these
51 // equations, the solver instead solves the linear least squares
52 // problem
53 //
54 //   min_x |Ax - b|^2
55 //
56 // The solution strategy is based on computing the QR decomposition of
57 // A, i.e.
58 //
59 //   A = QR
60 //
61 // Where Q is an orthonormal matrix and R is an upper triangular
62 // matrix. Then
63 //
64 //     Ax = b
65 //    QRx = b
66 //  Q'QRx = Q'b
67 //     Rx = Q'b
68 //      x = R^{-1} Q'b
69 //
70 // If the PerSolveOptions struct has a non-null array D, then the
71 // augmented/regularized linear system
72 //
73 //   [    A    ]x = [b]
74 //   [ diag(D) ]    [0]
75 //
76 // is solved.
77 //
78 // This class uses the dense QR factorization routines from the Eigen
79 // library. This solver always returns a solution, it is the user's
80 // responsibility to judge if the solution is good enough for their
81 // purposes.
82 class DenseQRSolver: public DenseSparseMatrixSolver {
83  public:
84   explicit DenseQRSolver(const LinearSolver::Options& options);
85
86  private:
87   virtual LinearSolver::Summary SolveImpl(
88       DenseSparseMatrix* A,
89       const double* b,
90       const LinearSolver::PerSolveOptions& per_solve_options,
91       double* x);
92
93   LinearSolver::Summary SolveUsingEigen(
94       DenseSparseMatrix* A,
95       const double* b,
96       const LinearSolver::PerSolveOptions& per_solve_options,
97       double* x);
98
99   LinearSolver::Summary SolveUsingLAPACK(
100       DenseSparseMatrix* A,
101       const double* b,
102       const LinearSolver::PerSolveOptions& per_solve_options,
103       double* x);
104
105   const LinearSolver::Options options_;
106   ColMajorMatrix lhs_;
107   Vector rhs_;
108   Vector work_;
109   CERES_DISALLOW_COPY_AND_ASSIGN(DenseQRSolver);
110 };
111
112 }  // namespace internal
113 }  // namespace ceres
114
115 #endif  // CERES_INTERNAL_DENSE_QR_SOLVER_H_