blob: 0943edbba855a69233f741bba00d235c8e9cf168 [file] [log] [blame]
Austin Schuh70cc9552019-01-21 19:46:48 -08001// 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#include "ceres/trust_region_preprocessor.h"
32
33#include <numeric>
34#include <string>
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080035
Austin Schuh70cc9552019-01-21 19:46:48 -080036#include "ceres/callbacks.h"
37#include "ceres/context_impl.h"
38#include "ceres/evaluator.h"
39#include "ceres/linear_solver.h"
40#include "ceres/minimizer.h"
41#include "ceres/parameter_block.h"
42#include "ceres/preconditioner.h"
43#include "ceres/preprocessor.h"
44#include "ceres/problem_impl.h"
45#include "ceres/program.h"
46#include "ceres/reorder_program.h"
47#include "ceres/suitesparse.h"
48#include "ceres/trust_region_strategy.h"
49#include "ceres/wall_time.h"
50
51namespace ceres {
52namespace internal {
53
54using std::vector;
55
56namespace {
57
58ParameterBlockOrdering* CreateDefaultLinearSolverOrdering(
59 const Program& program) {
60 ParameterBlockOrdering* ordering = new ParameterBlockOrdering;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080061 const vector<ParameterBlock*>& parameter_blocks = program.parameter_blocks();
Austin Schuh70cc9552019-01-21 19:46:48 -080062 for (int i = 0; i < parameter_blocks.size(); ++i) {
63 ordering->AddElementToGroup(
64 const_cast<double*>(parameter_blocks[i]->user_state()), 0);
65 }
66 return ordering;
67}
68
69// Check if all the user supplied values in the parameter blocks are
70// sane or not, and if the program is feasible or not.
71bool IsProgramValid(const Program& program, std::string* error) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080072 return (program.ParameterBlocksAreFinite(error) && program.IsFeasible(error));
Austin Schuh70cc9552019-01-21 19:46:48 -080073}
74
75void AlternateLinearSolverAndPreconditionerForSchurTypeLinearSolver(
76 Solver::Options* options) {
77 if (!IsSchurType(options->linear_solver_type)) {
78 return;
79 }
80
81 const LinearSolverType linear_solver_type_given = options->linear_solver_type;
82 const PreconditionerType preconditioner_type_given =
83 options->preconditioner_type;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080084 options->linear_solver_type =
85 LinearSolver::LinearSolverForZeroEBlocks(linear_solver_type_given);
Austin Schuh70cc9552019-01-21 19:46:48 -080086
87 std::string message;
88 if (linear_solver_type_given == ITERATIVE_SCHUR) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080089 options->preconditioner_type =
90 Preconditioner::PreconditionerForZeroEBlocks(preconditioner_type_given);
Austin Schuh70cc9552019-01-21 19:46:48 -080091
92 message =
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080093 StringPrintf("No E blocks. Switching from %s(%s) to %s(%s).",
94 LinearSolverTypeToString(linear_solver_type_given),
95 PreconditionerTypeToString(preconditioner_type_given),
96 LinearSolverTypeToString(options->linear_solver_type),
97 PreconditionerTypeToString(options->preconditioner_type));
Austin Schuh70cc9552019-01-21 19:46:48 -080098 } else {
99 message =
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800100 StringPrintf("No E blocks. Switching from %s to %s.",
101 LinearSolverTypeToString(linear_solver_type_given),
102 LinearSolverTypeToString(options->linear_solver_type));
Austin Schuh70cc9552019-01-21 19:46:48 -0800103 }
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800104 if (options->logging_type != SILENT) {
105 VLOG(1) << message;
106 }
Austin Schuh70cc9552019-01-21 19:46:48 -0800107}
108
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800109// Reorder the program to reduce fill-in and increase cache coherency.
Austin Schuh70cc9552019-01-21 19:46:48 -0800110bool ReorderProgram(PreprocessedProblem* pp) {
111 const Solver::Options& options = pp->options;
112 if (IsSchurType(options.linear_solver_type)) {
113 return ReorderProgramForSchurTypeLinearSolver(
114 options.linear_solver_type,
115 options.sparse_linear_algebra_library_type,
116 pp->problem->parameter_map(),
117 options.linear_solver_ordering.get(),
118 pp->reduced_program.get(),
119 &pp->error);
120 }
121
Austin Schuh70cc9552019-01-21 19:46:48 -0800122 if (options.linear_solver_type == SPARSE_NORMAL_CHOLESKY &&
123 !options.dynamic_sparsity) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800124 return ReorderProgramForSparseCholesky(
Austin Schuh70cc9552019-01-21 19:46:48 -0800125 options.sparse_linear_algebra_library_type,
126 *options.linear_solver_ordering,
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800127 0, /* use all the rows of the jacobian */
128 pp->reduced_program.get(),
129 &pp->error);
130 }
131
132 if (options.linear_solver_type == CGNR &&
133 options.preconditioner_type == SUBSET) {
134 pp->linear_solver_options.subset_preconditioner_start_row_block =
135 ReorderResidualBlocksByPartition(
136 options.residual_blocks_for_subset_preconditioner,
137 pp->reduced_program.get());
138
139 return ReorderProgramForSparseCholesky(
140 options.sparse_linear_algebra_library_type,
141 *options.linear_solver_ordering,
142 pp->linear_solver_options.subset_preconditioner_start_row_block,
Austin Schuh70cc9552019-01-21 19:46:48 -0800143 pp->reduced_program.get(),
144 &pp->error);
145 }
146
147 return true;
148}
149
150// Configure and create a linear solver object. In doing so, if a
151// sparse direct factorization based linear solver is being used, then
152// find a fill reducing ordering and reorder the program as needed
153// too.
154bool SetupLinearSolver(PreprocessedProblem* pp) {
155 Solver::Options& options = pp->options;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800156 pp->linear_solver_options = LinearSolver::Options();
157
158 if (!options.linear_solver_ordering) {
Austin Schuh70cc9552019-01-21 19:46:48 -0800159 // If the user has not supplied a linear solver ordering, then we
160 // assume that they are giving all the freedom to us in choosing
161 // the best possible ordering. This intent can be indicated by
162 // putting all the parameter blocks in the same elimination group.
163 options.linear_solver_ordering.reset(
164 CreateDefaultLinearSolverOrdering(*pp->reduced_program));
165 } else {
166 // If the user supplied an ordering, then check if the first
167 // elimination group is still non-empty after the reduced problem
168 // has been constructed.
169 //
170 // This is important for Schur type linear solvers, where the
171 // first elimination group is special -- it needs to be an
172 // independent set.
173 //
174 // If the first elimination group is empty, then we cannot use the
175 // user's requested linear solver (and a preconditioner as the
176 // case may be) so we must use a different one.
177 ParameterBlockOrdering* ordering = options.linear_solver_ordering.get();
178 const int min_group_id = ordering->MinNonZeroGroup();
179 ordering->Remove(pp->removed_parameter_blocks);
180 if (IsSchurType(options.linear_solver_type) &&
181 min_group_id != ordering->MinNonZeroGroup()) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800182 AlternateLinearSolverAndPreconditionerForSchurTypeLinearSolver(&options);
Austin Schuh70cc9552019-01-21 19:46:48 -0800183 }
184 }
185
186 // Reorder the program to reduce fill in and improve cache coherency
187 // of the Jacobian.
188 if (!ReorderProgram(pp)) {
189 return false;
190 }
191
192 // Configure the linear solver.
Austin Schuh70cc9552019-01-21 19:46:48 -0800193 pp->linear_solver_options.min_num_iterations =
194 options.min_linear_solver_iterations;
195 pp->linear_solver_options.max_num_iterations =
196 options.max_linear_solver_iterations;
197 pp->linear_solver_options.type = options.linear_solver_type;
198 pp->linear_solver_options.preconditioner_type = options.preconditioner_type;
199 pp->linear_solver_options.visibility_clustering_type =
200 options.visibility_clustering_type;
201 pp->linear_solver_options.sparse_linear_algebra_library_type =
202 options.sparse_linear_algebra_library_type;
203 pp->linear_solver_options.dense_linear_algebra_library_type =
204 options.dense_linear_algebra_library_type;
205 pp->linear_solver_options.use_explicit_schur_complement =
206 options.use_explicit_schur_complement;
207 pp->linear_solver_options.dynamic_sparsity = options.dynamic_sparsity;
208 pp->linear_solver_options.use_mixed_precision_solves =
209 options.use_mixed_precision_solves;
210 pp->linear_solver_options.max_num_refinement_iterations =
211 options.max_num_refinement_iterations;
212 pp->linear_solver_options.num_threads = options.num_threads;
213 pp->linear_solver_options.use_postordering = options.use_postordering;
214 pp->linear_solver_options.context = pp->problem->context();
215
216 if (IsSchurType(pp->linear_solver_options.type)) {
217 OrderingToGroupSizes(options.linear_solver_ordering.get(),
218 &pp->linear_solver_options.elimination_groups);
219
220 // Schur type solvers expect at least two elimination groups. If
221 // there is only one elimination group, then it is guaranteed that
222 // this group only contains e_blocks. Thus we add a dummy
223 // elimination group with zero blocks in it.
224 if (pp->linear_solver_options.elimination_groups.size() == 1) {
225 pp->linear_solver_options.elimination_groups.push_back(0);
226 }
227
228 if (options.linear_solver_type == SPARSE_SCHUR) {
229 // When using SPARSE_SCHUR, we ignore the user's postordering
230 // preferences in certain cases.
231 //
232 // 1. SUITE_SPARSE is the sparse linear algebra library requested
233 // but cholmod_camd is not available.
234 // 2. CX_SPARSE is the sparse linear algebra library requested.
235 //
236 // This ensures that the linear solver does not assume that a
237 // fill-reducing pre-ordering has been done.
238 //
239 // TODO(sameeragarwal): Implement the reordering of parameter
240 // blocks for CX_SPARSE.
241 if ((options.sparse_linear_algebra_library_type == SUITE_SPARSE &&
242 !SuiteSparse::
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800243 IsConstrainedApproximateMinimumDegreeOrderingAvailable()) ||
Austin Schuh70cc9552019-01-21 19:46:48 -0800244 (options.sparse_linear_algebra_library_type == CX_SPARSE)) {
245 pp->linear_solver_options.use_postordering = true;
246 }
247 }
248 }
249
250 pp->linear_solver.reset(LinearSolver::Create(pp->linear_solver_options));
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800251 return (pp->linear_solver != nullptr);
Austin Schuh70cc9552019-01-21 19:46:48 -0800252}
253
254// Configure and create the evaluator.
255bool SetupEvaluator(PreprocessedProblem* pp) {
256 const Solver::Options& options = pp->options;
257 pp->evaluator_options = Evaluator::Options();
258 pp->evaluator_options.linear_solver_type = options.linear_solver_type;
259 pp->evaluator_options.num_eliminate_blocks = 0;
260 if (IsSchurType(options.linear_solver_type)) {
261 pp->evaluator_options.num_eliminate_blocks =
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800262 options.linear_solver_ordering->group_to_elements()
263 .begin()
264 ->second.size();
Austin Schuh70cc9552019-01-21 19:46:48 -0800265 }
266
267 pp->evaluator_options.num_threads = options.num_threads;
268 pp->evaluator_options.dynamic_sparsity = options.dynamic_sparsity;
269 pp->evaluator_options.context = pp->problem->context();
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800270 pp->evaluator_options.evaluation_callback =
271 pp->reduced_program->mutable_evaluation_callback();
272 pp->evaluator.reset(Evaluator::Create(
273 pp->evaluator_options, pp->reduced_program.get(), &pp->error));
Austin Schuh70cc9552019-01-21 19:46:48 -0800274
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800275 return (pp->evaluator != nullptr);
Austin Schuh70cc9552019-01-21 19:46:48 -0800276}
277
278// If the user requested inner iterations, then find an inner
279// iteration ordering as needed and configure and create a
280// CoordinateDescentMinimizer object to perform the inner iterations.
281bool SetupInnerIterationMinimizer(PreprocessedProblem* pp) {
282 Solver::Options& options = pp->options;
283 if (!options.use_inner_iterations) {
284 return true;
285 }
286
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800287 if (pp->reduced_program->mutable_evaluation_callback()) {
288 pp->error = "Inner iterations cannot be used with EvaluationCallbacks";
289 return false;
290 }
291
Austin Schuh70cc9552019-01-21 19:46:48 -0800292 // With just one parameter block, the outer iteration of the trust
293 // region method and inner iterations are doing exactly the same
294 // thing, and thus inner iterations are not needed.
295 if (pp->reduced_program->NumParameterBlocks() == 1) {
296 LOG(WARNING) << "Reduced problem only contains one parameter block."
297 << "Disabling inner iterations.";
298 return true;
299 }
300
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800301 if (options.inner_iteration_ordering != nullptr) {
Austin Schuh70cc9552019-01-21 19:46:48 -0800302 // If the user supplied an ordering, then remove the set of
303 // inactive parameter blocks from it
304 options.inner_iteration_ordering->Remove(pp->removed_parameter_blocks);
305 if (options.inner_iteration_ordering->NumElements() == 0) {
306 LOG(WARNING) << "No remaining elements in the inner iteration ordering.";
307 return true;
308 }
309
310 // Validate the reduced ordering.
311 if (!CoordinateDescentMinimizer::IsOrderingValid(
312 *pp->reduced_program,
313 *options.inner_iteration_ordering,
314 &pp->error)) {
315 return false;
316 }
317 } else {
318 // The user did not supply an ordering, so create one.
319 options.inner_iteration_ordering.reset(
320 CoordinateDescentMinimizer::CreateOrdering(*pp->reduced_program));
321 }
322
323 pp->inner_iteration_minimizer.reset(
324 new CoordinateDescentMinimizer(pp->problem->context()));
325 return pp->inner_iteration_minimizer->Init(*pp->reduced_program,
326 pp->problem->parameter_map(),
327 *options.inner_iteration_ordering,
328 &pp->error);
329}
330
331// Configure and create a TrustRegionMinimizer object.
332void SetupMinimizerOptions(PreprocessedProblem* pp) {
333 const Solver::Options& options = pp->options;
334
335 SetupCommonMinimizerOptions(pp);
336 pp->minimizer_options.is_constrained =
337 pp->reduced_program->IsBoundsConstrained();
338 pp->minimizer_options.jacobian.reset(pp->evaluator->CreateJacobian());
339 pp->minimizer_options.inner_iteration_minimizer =
340 pp->inner_iteration_minimizer;
341
342 TrustRegionStrategy::Options strategy_options;
343 strategy_options.linear_solver = pp->linear_solver.get();
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800344 strategy_options.initial_radius = options.initial_trust_region_radius;
Austin Schuh70cc9552019-01-21 19:46:48 -0800345 strategy_options.max_radius = options.max_trust_region_radius;
346 strategy_options.min_lm_diagonal = options.min_lm_diagonal;
347 strategy_options.max_lm_diagonal = options.max_lm_diagonal;
348 strategy_options.trust_region_strategy_type =
349 options.trust_region_strategy_type;
350 strategy_options.dogleg_type = options.dogleg_type;
351 pp->minimizer_options.trust_region_strategy.reset(
352 TrustRegionStrategy::Create(strategy_options));
353 CHECK(pp->minimizer_options.trust_region_strategy != nullptr);
354}
355
356} // namespace
357
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800358TrustRegionPreprocessor::~TrustRegionPreprocessor() {}
Austin Schuh70cc9552019-01-21 19:46:48 -0800359
360bool TrustRegionPreprocessor::Preprocess(const Solver::Options& options,
361 ProblemImpl* problem,
362 PreprocessedProblem* pp) {
363 CHECK(pp != nullptr);
364 pp->options = options;
365 ChangeNumThreadsIfNeeded(&pp->options);
366
367 pp->problem = problem;
368 Program* program = problem->mutable_program();
369 if (!IsProgramValid(*program, &pp->error)) {
370 return false;
371 }
372
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800373 pp->reduced_program.reset(program->CreateReducedProgram(
374 &pp->removed_parameter_blocks, &pp->fixed_cost, &pp->error));
Austin Schuh70cc9552019-01-21 19:46:48 -0800375
376 if (pp->reduced_program.get() == NULL) {
377 return false;
378 }
379
380 if (pp->reduced_program->NumParameterBlocks() == 0) {
381 // The reduced problem has no parameter or residual blocks. There
382 // is nothing more to do.
383 return true;
384 }
385
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800386 if (!SetupLinearSolver(pp) || !SetupEvaluator(pp) ||
Austin Schuh70cc9552019-01-21 19:46:48 -0800387 !SetupInnerIterationMinimizer(pp)) {
388 return false;
389 }
390
391 SetupMinimizerOptions(pp);
392 return true;
393}
394
395} // namespace internal
396} // namespace ceres