Welcome to mirror list, hosted at ThFree Co, Russian Federation.

dynamic_sparse_normal_cholesky_solver.cc « ceres « internal « ceres « extern - git.blender.org/blender.git - Unnamed repository; edit this file 'description' to name the repository.
summaryrefslogtreecommitdiff
blob: d31c4228f055fc0c33f67d8ba64d8e5fb02ece45 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
// Ceres Solver - A fast non-linear least squares minimizer
// Copyright 2017 Google Inc. All rights reserved.
// http://ceres-solver.org/
//
// Redistribution and use in source and binary forms, with or without
// modification, are permitted provided that the following conditions are met:
//
// * Redistributions of source code must retain the above copyright notice,
//   this list of conditions and the following disclaimer.
// * Redistributions in binary form must reproduce the above copyright notice,
//   this list of conditions and the following disclaimer in the documentation
//   and/or other materials provided with the distribution.
// * Neither the name of Google Inc. nor the names of its contributors may be
//   used to endorse or promote products derived from this software without
//   specific prior written permission.
//
// THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
// AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
// ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
// LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
// CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
// SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
// INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
// CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
// ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
// POSSIBILITY OF SUCH DAMAGE.
//
// Author: sameeragarwal@google.com (Sameer Agarwal)

#include "ceres/dynamic_sparse_normal_cholesky_solver.h"

#include <algorithm>
#include <cstring>
#include <ctime>
#include <memory>
#include <sstream>

#include "Eigen/SparseCore"
#include "ceres/compressed_row_sparse_matrix.h"
#include "ceres/cxsparse.h"
#include "ceres/internal/eigen.h"
#include "ceres/linear_solver.h"
#include "ceres/suitesparse.h"
#include "ceres/triplet_sparse_matrix.h"
#include "ceres/types.h"
#include "ceres/wall_time.h"

#ifdef CERES_USE_EIGEN_SPARSE
#include "Eigen/SparseCholesky"
#endif

namespace ceres {
namespace internal {

DynamicSparseNormalCholeskySolver::DynamicSparseNormalCholeskySolver(
    const LinearSolver::Options& options)
    : options_(options) {}

LinearSolver::Summary DynamicSparseNormalCholeskySolver::SolveImpl(
    CompressedRowSparseMatrix* A,
    const double* b,
    const LinearSolver::PerSolveOptions& per_solve_options,
    double* x) {
  const int num_cols = A->num_cols();
  VectorRef(x, num_cols).setZero();
  A->LeftMultiply(b, x);

  if (per_solve_options.D != nullptr) {
    // Temporarily append a diagonal block to the A matrix, but undo
    // it before returning the matrix to the user.
    std::unique_ptr<CompressedRowSparseMatrix> regularizer;
    if (!A->col_blocks().empty()) {
      regularizer.reset(CompressedRowSparseMatrix::CreateBlockDiagonalMatrix(
          per_solve_options.D, A->col_blocks()));
    } else {
      regularizer.reset(
          new CompressedRowSparseMatrix(per_solve_options.D, num_cols));
    }
    A->AppendRows(*regularizer);
  }

  LinearSolver::Summary summary;
  switch (options_.sparse_linear_algebra_library_type) {
    case SUITE_SPARSE:
      summary = SolveImplUsingSuiteSparse(A, x);
      break;
    case CX_SPARSE:
      summary = SolveImplUsingCXSparse(A, x);
      break;
    case EIGEN_SPARSE:
      summary = SolveImplUsingEigen(A, x);
      break;
    default:
      LOG(FATAL) << "Unsupported sparse linear algebra library for "
                 << "dynamic sparsity: "
                 << SparseLinearAlgebraLibraryTypeToString(
                        options_.sparse_linear_algebra_library_type);
  }

  if (per_solve_options.D != nullptr) {
    A->DeleteRows(num_cols);
  }

  return summary;
}

LinearSolver::Summary DynamicSparseNormalCholeskySolver::SolveImplUsingEigen(
    CompressedRowSparseMatrix* A, double* rhs_and_solution) {
#ifndef CERES_USE_EIGEN_SPARSE

  LinearSolver::Summary summary;
  summary.num_iterations = 0;
  summary.termination_type = LINEAR_SOLVER_FATAL_ERROR;
  summary.message =
      "SPARSE_NORMAL_CHOLESKY cannot be used with EIGEN_SPARSE "
      "because Ceres was not built with support for "
      "Eigen's SimplicialLDLT decomposition. "
      "This requires enabling building with -DEIGENSPARSE=ON.";
  return summary;

#else

  EventLogger event_logger("DynamicSparseNormalCholeskySolver::Eigen::Solve");

  Eigen::MappedSparseMatrix<double, Eigen::RowMajor> a(A->num_rows(),
                                                       A->num_cols(),
                                                       A->num_nonzeros(),
                                                       A->mutable_rows(),
                                                       A->mutable_cols(),
                                                       A->mutable_values());

  Eigen::SparseMatrix<double> lhs = a.transpose() * a;
  Eigen::SimplicialLDLT<Eigen::SparseMatrix<double>> solver;

  LinearSolver::Summary summary;
  summary.num_iterations = 1;
  summary.termination_type = LINEAR_SOLVER_SUCCESS;
  summary.message = "Success.";

  solver.analyzePattern(lhs);
  if (VLOG_IS_ON(2)) {
    std::stringstream ss;
    solver.dumpMemory(ss);
    VLOG(2) << "Symbolic Analysis\n" << ss.str();
  }

  event_logger.AddEvent("Analyze");
  if (solver.info() != Eigen::Success) {
    summary.termination_type = LINEAR_SOLVER_FATAL_ERROR;
    summary.message = "Eigen failure. Unable to find symbolic factorization.";
    return summary;
  }

  solver.factorize(lhs);
  event_logger.AddEvent("Factorize");
  if (solver.info() != Eigen::Success) {
    summary.termination_type = LINEAR_SOLVER_FAILURE;
    summary.message = "Eigen failure. Unable to find numeric factorization.";
    return summary;
  }

  const Vector rhs = VectorRef(rhs_and_solution, lhs.cols());
  VectorRef(rhs_and_solution, lhs.cols()) = solver.solve(rhs);
  event_logger.AddEvent("Solve");
  if (solver.info() != Eigen::Success) {
    summary.termination_type = LINEAR_SOLVER_FAILURE;
    summary.message = "Eigen failure. Unable to do triangular solve.";
    return summary;
  }

  return summary;
#endif  // CERES_USE_EIGEN_SPARSE
}

LinearSolver::Summary DynamicSparseNormalCholeskySolver::SolveImplUsingCXSparse(
    CompressedRowSparseMatrix* A, double* rhs_and_solution) {
#ifdef CERES_NO_CXSPARSE

  LinearSolver::Summary summary;
  summary.num_iterations = 0;
  summary.termination_type = LINEAR_SOLVER_FATAL_ERROR;
  summary.message =
      "SPARSE_NORMAL_CHOLESKY cannot be used with CX_SPARSE "
      "because Ceres was not built with support for CXSparse. "
      "This requires enabling building with -DCXSPARSE=ON.";

  return summary;

#else
  EventLogger event_logger(
      "DynamicSparseNormalCholeskySolver::CXSparse::Solve");

  LinearSolver::Summary summary;
  summary.num_iterations = 1;
  summary.termination_type = LINEAR_SOLVER_SUCCESS;
  summary.message = "Success.";

  CXSparse cxsparse;

  // Wrap the augmented Jacobian in a compressed sparse column matrix.
  cs_di a_transpose = cxsparse.CreateSparseMatrixTransposeView(A);

  // Compute the normal equations. J'J delta = J'f and solve them
  // using a sparse Cholesky factorization. Notice that when compared
  // to SuiteSparse we have to explicitly compute the transpose of Jt,
  // and then the normal equations before they can be
  // factorized. CHOLMOD/SuiteSparse on the other hand can just work
  // off of Jt to compute the Cholesky factorization of the normal
  // equations.
  cs_di* a = cxsparse.TransposeMatrix(&a_transpose);
  cs_di* lhs = cxsparse.MatrixMatrixMultiply(&a_transpose, a);
  cxsparse.Free(a);
  event_logger.AddEvent("NormalEquations");

  if (!cxsparse.SolveCholesky(lhs, rhs_and_solution)) {
    summary.termination_type = LINEAR_SOLVER_FAILURE;
    summary.message = "CXSparse::SolveCholesky failed";
  }
  event_logger.AddEvent("Solve");

  cxsparse.Free(lhs);
  event_logger.AddEvent("TearDown");
  return summary;
#endif
}

LinearSolver::Summary
DynamicSparseNormalCholeskySolver::SolveImplUsingSuiteSparse(
    CompressedRowSparseMatrix* A, double* rhs_and_solution) {
#ifdef CERES_NO_SUITESPARSE

  LinearSolver::Summary summary;
  summary.num_iterations = 0;
  summary.termination_type = LINEAR_SOLVER_FATAL_ERROR;
  summary.message =
      "SPARSE_NORMAL_CHOLESKY cannot be used with SUITE_SPARSE "
      "because Ceres was not built with support for SuiteSparse. "
      "This requires enabling building with -DSUITESPARSE=ON.";
  return summary;

#else

  EventLogger event_logger(
      "DynamicSparseNormalCholeskySolver::SuiteSparse::Solve");
  LinearSolver::Summary summary;
  summary.termination_type = LINEAR_SOLVER_SUCCESS;
  summary.num_iterations = 1;
  summary.message = "Success.";

  SuiteSparse ss;
  const int num_cols = A->num_cols();
  cholmod_sparse lhs = ss.CreateSparseMatrixTransposeView(A);
  event_logger.AddEvent("Setup");
  cholmod_factor* factor = ss.AnalyzeCholesky(&lhs, &summary.message);
  event_logger.AddEvent("Analysis");

  if (factor == nullptr) {
    summary.termination_type = LINEAR_SOLVER_FATAL_ERROR;
    return summary;
  }

  summary.termination_type = ss.Cholesky(&lhs, factor, &summary.message);
  if (summary.termination_type == LINEAR_SOLVER_SUCCESS) {
    cholmod_dense cholmod_rhs =
        ss.CreateDenseVectorView(rhs_and_solution, num_cols);
    cholmod_dense* solution = ss.Solve(factor, &cholmod_rhs, &summary.message);
    event_logger.AddEvent("Solve");
    if (solution != nullptr) {
      memcpy(
          rhs_and_solution, solution->x, num_cols * sizeof(*rhs_and_solution));
      ss.Free(solution);
    } else {
      summary.termination_type = LINEAR_SOLVER_FAILURE;
    }
  }

  ss.Free(factor);
  event_logger.AddEvent("Teardown");
  return summary;

#endif
}

}  // namespace internal
}  // namespace ceres