blob: 8fc1305baa2648a121254e84b4f4cdd3f2892618 [file] [log] [blame]
Keir Mierle8ebb0732012-04-30 23:09:08 -07001// Ceres Solver - A fast non-linear least squares minimizer
2// Copyright 2010, 2011, 2012 Google Inc. All rights reserved.
3// http://code.google.com/p/ceres-solver/
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: keir@google.com (Keir Mierle)
30
31#include "ceres/solver_impl.h"
32
Sameer Agarwal1b7f3b52012-08-09 21:46:19 -070033#include <cstdio>
Keir Mierle8ebb0732012-04-30 23:09:08 -070034#include <iostream> // NOLINT
35#include <numeric>
Sameer Agarwald5b93bf2013-04-26 21:17:49 -070036#include <string>
Sameer Agarwalba8d9672012-10-02 00:48:57 -070037#include "ceres/coordinate_descent_minimizer.h"
Sameer Agarwald5b93bf2013-04-26 21:17:49 -070038#include "ceres/cxsparse.h"
Keir Mierle8ebb0732012-04-30 23:09:08 -070039#include "ceres/evaluator.h"
40#include "ceres/gradient_checking_cost_function.h"
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -070041#include "ceres/iteration_callback.h"
42#include "ceres/levenberg_marquardt_strategy.h"
Sameer Agarwalf4d01642012-11-26 12:55:58 -080043#include "ceres/line_search_minimizer.h"
Sameer Agarwalc8f07902013-04-19 08:01:04 -070044#include "ceres/linear_solver.h"
Keir Mierle8ebb0732012-04-30 23:09:08 -070045#include "ceres/map_util.h"
46#include "ceres/minimizer.h"
Sameer Agarwal2c94eed2012-10-01 16:34:37 -070047#include "ceres/ordered_groups.h"
Keir Mierle8ebb0732012-04-30 23:09:08 -070048#include "ceres/parameter_block.h"
Sameer Agarwalba8d9672012-10-02 00:48:57 -070049#include "ceres/parameter_block_ordering.h"
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -070050#include "ceres/problem.h"
Keir Mierle8ebb0732012-04-30 23:09:08 -070051#include "ceres/problem_impl.h"
52#include "ceres/program.h"
53#include "ceres/residual_block.h"
Keir Mierle8ebb0732012-04-30 23:09:08 -070054#include "ceres/stringprintf.h"
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -070055#include "ceres/suitesparse.h"
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -070056#include "ceres/trust_region_minimizer.h"
Petter Strandmark76533b32012-09-04 22:08:50 -070057#include "ceres/wall_time.h"
Keir Mierle8ebb0732012-04-30 23:09:08 -070058
59namespace ceres {
60namespace internal {
61namespace {
62
Keir Mierle8ebb0732012-04-30 23:09:08 -070063// Callback for updating the user's parameter blocks. Updates are only
64// done if the step is successful.
65class StateUpdatingCallback : public IterationCallback {
66 public:
67 StateUpdatingCallback(Program* program, double* parameters)
68 : program_(program), parameters_(parameters) {}
69
70 CallbackReturnType operator()(const IterationSummary& summary) {
71 if (summary.step_is_successful) {
72 program_->StateVectorToParameterBlocks(parameters_);
73 program_->CopyParameterBlockStateToUserState();
74 }
75 return SOLVER_CONTINUE;
76 }
77
78 private:
79 Program* program_;
80 double* parameters_;
81};
82
Sameer Agarwalf102a682013-02-11 15:08:40 -080083void SetSummaryFinalCost(Solver::Summary* summary) {
84 summary->final_cost = summary->initial_cost;
85 // We need the loop here, instead of just looking at the last
86 // iteration because the minimizer maybe making non-monotonic steps.
87 for (int i = 0; i < summary->iterations.size(); ++i) {
88 const IterationSummary& iteration_summary = summary->iterations[i];
89 summary->final_cost = min(iteration_summary.cost, summary->final_cost);
90 }
91}
92
Keir Mierle8ebb0732012-04-30 23:09:08 -070093// Callback for logging the state of the minimizer to STDERR or STDOUT
94// depending on the user's preferences and logging level.
Sameer Agarwalf4d01642012-11-26 12:55:58 -080095class TrustRegionLoggingCallback : public IterationCallback {
Keir Mierle8ebb0732012-04-30 23:09:08 -070096 public:
Sameer Agarwalf4d01642012-11-26 12:55:58 -080097 explicit TrustRegionLoggingCallback(bool log_to_stdout)
Keir Mierle8ebb0732012-04-30 23:09:08 -070098 : log_to_stdout_(log_to_stdout) {}
99
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800100 ~TrustRegionLoggingCallback() {}
Keir Mierle8ebb0732012-04-30 23:09:08 -0700101
102 CallbackReturnType operator()(const IterationSummary& summary) {
103 const char* kReportRowFormat =
104 "% 4d: f:% 8e d:% 3.2e g:% 3.2e h:% 3.2e "
Sameer Agarwalfa015192012-06-11 14:21:42 -0700105 "rho:% 3.2e mu:% 3.2e li:% 3d it:% 3.2e tt:% 3.2e";
Keir Mierle8ebb0732012-04-30 23:09:08 -0700106 string output = StringPrintf(kReportRowFormat,
107 summary.iteration,
108 summary.cost,
109 summary.cost_change,
110 summary.gradient_max_norm,
111 summary.step_norm,
112 summary.relative_decrease,
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700113 summary.trust_region_radius,
Sameer Agarwalfa015192012-06-11 14:21:42 -0700114 summary.linear_solver_iterations,
115 summary.iteration_time_in_seconds,
116 summary.cumulative_time_in_seconds);
Keir Mierle8ebb0732012-04-30 23:09:08 -0700117 if (log_to_stdout_) {
118 cout << output << endl;
119 } else {
120 VLOG(1) << output;
121 }
122 return SOLVER_CONTINUE;
123 }
124
125 private:
126 const bool log_to_stdout_;
127};
128
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800129// Callback for logging the state of the minimizer to STDERR or STDOUT
130// depending on the user's preferences and logging level.
131class LineSearchLoggingCallback : public IterationCallback {
132 public:
133 explicit LineSearchLoggingCallback(bool log_to_stdout)
134 : log_to_stdout_(log_to_stdout) {}
135
136 ~LineSearchLoggingCallback() {}
137
138 CallbackReturnType operator()(const IterationSummary& summary) {
139 const char* kReportRowFormat =
140 "% 4d: f:% 8e d:% 3.2e g:% 3.2e h:% 3.2e "
141 "s:% 3.2e e:% 3d it:% 3.2e tt:% 3.2e";
142 string output = StringPrintf(kReportRowFormat,
143 summary.iteration,
144 summary.cost,
145 summary.cost_change,
146 summary.gradient_max_norm,
147 summary.step_norm,
148 summary.step_size,
149 summary.line_search_function_evaluations,
150 summary.iteration_time_in_seconds,
151 summary.cumulative_time_in_seconds);
152 if (log_to_stdout_) {
153 cout << output << endl;
154 } else {
155 VLOG(1) << output;
156 }
157 return SOLVER_CONTINUE;
158 }
159
160 private:
161 const bool log_to_stdout_;
162};
163
164
Sameer Agarwal1b7f3b52012-08-09 21:46:19 -0700165// Basic callback to record the execution of the solver to a file for
166// offline analysis.
167class FileLoggingCallback : public IterationCallback {
168 public:
169 explicit FileLoggingCallback(const string& filename)
170 : fptr_(NULL) {
171 fptr_ = fopen(filename.c_str(), "w");
172 CHECK_NOTNULL(fptr_);
173 }
174
175 virtual ~FileLoggingCallback() {
176 if (fptr_ != NULL) {
177 fclose(fptr_);
178 }
179 }
180
181 virtual CallbackReturnType operator()(const IterationSummary& summary) {
182 fprintf(fptr_,
183 "%4d %e %e\n",
184 summary.iteration,
185 summary.cost,
186 summary.cumulative_time_in_seconds);
187 return SOLVER_CONTINUE;
188 }
189 private:
190 FILE* fptr_;
191};
192
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800193// Iterate over each of the groups in order of their priority and fill
194// summary with their sizes.
Sameer Agarwal977be7c2013-01-26 16:01:54 -0800195void SummarizeOrdering(ParameterBlockOrdering* ordering,
196 vector<int>* summary) {
197 CHECK_NOTNULL(summary)->clear();
198 if (ordering == NULL) {
199 return;
200 }
201
202 const map<int, set<double*> >& group_to_elements =
203 ordering->group_to_elements();
204 for (map<int, set<double*> >::const_iterator it = group_to_elements.begin();
205 it != group_to_elements.end();
206 ++it) {
207 summary->push_back(it->second.size());
208 }
209}
210
Keir Mierle8ebb0732012-04-30 23:09:08 -0700211} // namespace
212
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800213void SolverImpl::TrustRegionMinimize(
214 const Solver::Options& options,
215 Program* program,
216 CoordinateDescentMinimizer* inner_iteration_minimizer,
217 Evaluator* evaluator,
218 LinearSolver* linear_solver,
219 double* parameters,
220 Solver::Summary* summary) {
Keir Mierle8ebb0732012-04-30 23:09:08 -0700221 Minimizer::Options minimizer_options(options);
Sameer Agarwal1b7f3b52012-08-09 21:46:19 -0700222
223 // TODO(sameeragarwal): Add support for logging the configuration
224 // and more detailed stats.
225 scoped_ptr<IterationCallback> file_logging_callback;
226 if (!options.solver_log.empty()) {
227 file_logging_callback.reset(new FileLoggingCallback(options.solver_log));
228 minimizer_options.callbacks.insert(minimizer_options.callbacks.begin(),
229 file_logging_callback.get());
230 }
231
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800232 TrustRegionLoggingCallback logging_callback(
233 options.minimizer_progress_to_stdout);
Keir Mierle8ebb0732012-04-30 23:09:08 -0700234 if (options.logging_type != SILENT) {
Keir Mierlef7471832012-06-14 11:31:53 -0700235 minimizer_options.callbacks.insert(minimizer_options.callbacks.begin(),
236 &logging_callback);
Keir Mierle8ebb0732012-04-30 23:09:08 -0700237 }
238
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700239 StateUpdatingCallback updating_callback(program, parameters);
Keir Mierle8ebb0732012-04-30 23:09:08 -0700240 if (options.update_state_every_iteration) {
Keir Mierlef7471832012-06-14 11:31:53 -0700241 // This must get pushed to the front of the callbacks so that it is run
242 // before any of the user callbacks.
243 minimizer_options.callbacks.insert(minimizer_options.callbacks.begin(),
244 &updating_callback);
Keir Mierle8ebb0732012-04-30 23:09:08 -0700245 }
246
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700247 minimizer_options.evaluator = evaluator;
248 scoped_ptr<SparseMatrix> jacobian(evaluator->CreateJacobian());
Sameer Agarwalc58e6dc2013-02-11 16:41:06 -0800249
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700250 minimizer_options.jacobian = jacobian.get();
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700251 minimizer_options.inner_iteration_minimizer = inner_iteration_minimizer;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700252
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700253 TrustRegionStrategy::Options trust_region_strategy_options;
254 trust_region_strategy_options.linear_solver = linear_solver;
255 trust_region_strategy_options.initial_radius =
256 options.initial_trust_region_radius;
257 trust_region_strategy_options.max_radius = options.max_trust_region_radius;
Sameer Agarwaleeedd2e2013-07-07 23:04:31 -0700258 trust_region_strategy_options.min_lm_diagonal = options.min_lm_diagonal;
259 trust_region_strategy_options.max_lm_diagonal = options.max_lm_diagonal;
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700260 trust_region_strategy_options.trust_region_strategy_type =
261 options.trust_region_strategy_type;
Markus Moll51cf7cb2012-08-20 20:10:20 +0200262 trust_region_strategy_options.dogleg_type = options.dogleg_type;
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700263 scoped_ptr<TrustRegionStrategy> strategy(
264 TrustRegionStrategy::Create(trust_region_strategy_options));
265 minimizer_options.trust_region_strategy = strategy.get();
266
267 TrustRegionMinimizer minimizer;
Petter Strandmark76533b32012-09-04 22:08:50 -0700268 double minimizer_start_time = WallTimeInSeconds();
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700269 minimizer.Minimize(minimizer_options, parameters, summary);
Petter Strandmark76533b32012-09-04 22:08:50 -0700270 summary->minimizer_time_in_seconds =
271 WallTimeInSeconds() - minimizer_start_time;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700272}
273
Sameer Agarwal8140f0f2013-03-12 09:45:08 -0700274#ifndef CERES_NO_LINE_SEARCH_MINIMIZER
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800275void SolverImpl::LineSearchMinimize(
276 const Solver::Options& options,
277 Program* program,
278 Evaluator* evaluator,
279 double* parameters,
280 Solver::Summary* summary) {
281 Minimizer::Options minimizer_options(options);
282
283 // TODO(sameeragarwal): Add support for logging the configuration
284 // and more detailed stats.
285 scoped_ptr<IterationCallback> file_logging_callback;
286 if (!options.solver_log.empty()) {
287 file_logging_callback.reset(new FileLoggingCallback(options.solver_log));
288 minimizer_options.callbacks.insert(minimizer_options.callbacks.begin(),
289 file_logging_callback.get());
290 }
291
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800292 LineSearchLoggingCallback logging_callback(
293 options.minimizer_progress_to_stdout);
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800294 if (options.logging_type != SILENT) {
295 minimizer_options.callbacks.insert(minimizer_options.callbacks.begin(),
296 &logging_callback);
297 }
298
299 StateUpdatingCallback updating_callback(program, parameters);
300 if (options.update_state_every_iteration) {
301 // This must get pushed to the front of the callbacks so that it is run
302 // before any of the user callbacks.
303 minimizer_options.callbacks.insert(minimizer_options.callbacks.begin(),
304 &updating_callback);
305 }
306
307 minimizer_options.evaluator = evaluator;
308
309 LineSearchMinimizer minimizer;
310 double minimizer_start_time = WallTimeInSeconds();
311 minimizer.Minimize(minimizer_options, parameters, summary);
312 summary->minimizer_time_in_seconds =
313 WallTimeInSeconds() - minimizer_start_time;
314}
Sameer Agarwal8140f0f2013-03-12 09:45:08 -0700315#endif // CERES_NO_LINE_SEARCH_MINIMIZER
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800316
317void SolverImpl::Solve(const Solver::Options& options,
318 ProblemImpl* problem_impl,
Keir Mierle8ebb0732012-04-30 23:09:08 -0700319 Solver::Summary* summary) {
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800320 if (options.minimizer_type == TRUST_REGION) {
321 TrustRegionSolve(options, problem_impl, summary);
322 } else {
Sameer Agarwal8140f0f2013-03-12 09:45:08 -0700323#ifndef CERES_NO_LINE_SEARCH_MINIMIZER
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800324 LineSearchSolve(options, problem_impl, summary);
Sameer Agarwal8140f0f2013-03-12 09:45:08 -0700325#else
326 LOG(FATAL) << "Ceres Solver was compiled with -DLINE_SEARCH_MINIMIZER=OFF";
327#endif
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800328 }
329}
330
331void SolverImpl::TrustRegionSolve(const Solver::Options& original_options,
332 ProblemImpl* original_problem_impl,
333 Solver::Summary* summary) {
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800334 EventLogger event_logger("TrustRegionSolve");
Petter Strandmark76533b32012-09-04 22:08:50 -0700335 double solver_start_time = WallTimeInSeconds();
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700336
Sameer Agarwal4997cbc2012-07-02 12:44:34 -0700337 Program* original_program = original_problem_impl->mutable_program();
338 ProblemImpl* problem_impl = original_problem_impl;
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700339
Sameer Agarwal4997cbc2012-07-02 12:44:34 -0700340 // Reset the summary object to its default values.
341 *CHECK_NOTNULL(summary) = Solver::Summary();
342
Sameer Agarwald010de52013-02-15 14:26:56 -0800343 summary->minimizer_type = TRUST_REGION;
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700344 summary->num_parameter_blocks = problem_impl->NumParameterBlocks();
345 summary->num_parameters = problem_impl->NumParameters();
Keir Mierleba944212013-02-25 12:46:44 -0800346 summary->num_effective_parameters =
347 original_program->NumEffectiveParameters();
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700348 summary->num_residual_blocks = problem_impl->NumResidualBlocks();
349 summary->num_residuals = problem_impl->NumResiduals();
350
351 // Empty programs are usually a user error.
352 if (summary->num_parameter_blocks == 0) {
353 summary->error = "Problem contains no parameter blocks.";
354 LOG(ERROR) << summary->error;
355 return;
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700356 }
357
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700358 if (summary->num_residual_blocks == 0) {
359 summary->error = "Problem contains no residual blocks.";
360 LOG(ERROR) << summary->error;
361 return;
362 }
363
Sameer Agarwal977be7c2013-01-26 16:01:54 -0800364 SummarizeOrdering(original_options.linear_solver_ordering,
365 &(summary->linear_solver_ordering_given));
366
367 SummarizeOrdering(original_options.inner_iteration_ordering,
368 &(summary->inner_iteration_ordering_given));
369
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700370 Solver::Options options(original_options);
371 options.linear_solver_ordering = NULL;
372 options.inner_iteration_ordering = NULL;
373
Keir Mierle8ebb0732012-04-30 23:09:08 -0700374#ifndef CERES_USE_OPENMP
375 if (options.num_threads > 1) {
376 LOG(WARNING)
377 << "OpenMP support is not compiled into this binary; "
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700378 << "only options.num_threads=1 is supported. Switching "
Keir Mierle8ebb0732012-04-30 23:09:08 -0700379 << "to single threaded mode.";
380 options.num_threads = 1;
381 }
382 if (options.num_linear_solver_threads > 1) {
383 LOG(WARNING)
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700384 << "OpenMP support is not compiled into this binary; "
385 << "only options.num_linear_solver_threads=1 is supported. Switching "
Keir Mierle8ebb0732012-04-30 23:09:08 -0700386 << "to single threaded mode.";
387 options.num_linear_solver_threads = 1;
388 }
389#endif
390
Keir Mierle8ebb0732012-04-30 23:09:08 -0700391 summary->num_threads_given = original_options.num_threads;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700392 summary->num_threads_used = options.num_threads;
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700393
Sameer Agarwalc4a32912013-06-13 22:00:48 -0700394 if (options.trust_region_minimizer_iterations_to_dump.size() > 0 &&
395 options.trust_region_problem_dump_format_type != CONSOLE &&
396 options.trust_region_problem_dump_directory.empty()) {
397 summary->error =
398 "Solver::Options::trust_region_problem_dump_directory is empty.";
399 LOG(ERROR) << summary->error;
400 return;
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700401 }
402
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800403 event_logger.AddEvent("Init");
Sameer Agarwalf102a682013-02-11 15:08:40 -0800404
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700405 original_program->SetParameterBlockStatePtrsToUserStatePtrs();
Sameer Agarwalc58e6dc2013-02-11 16:41:06 -0800406 event_logger.AddEvent("SetParameterBlockPtrs");
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700407
408 // If the user requests gradient checking, construct a new
409 // ProblemImpl by wrapping the CostFunctions of problem_impl inside
410 // GradientCheckingCostFunction and replacing problem_impl with
411 // gradient_checking_problem_impl.
412 scoped_ptr<ProblemImpl> gradient_checking_problem_impl;
413 if (options.check_gradients) {
414 VLOG(1) << "Checking Gradients";
415 gradient_checking_problem_impl.reset(
416 CreateGradientCheckingProblemImpl(
417 problem_impl,
418 options.numeric_derivative_relative_step_size,
419 options.gradient_check_relative_precision));
420
421 // From here on, problem_impl will point to the gradient checking
422 // version.
423 problem_impl = gradient_checking_problem_impl.get();
424 }
Keir Mierle8ebb0732012-04-30 23:09:08 -0700425
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700426 if (original_options.linear_solver_ordering != NULL) {
Sameer Agarwal65625f72012-09-17 12:06:57 -0700427 if (!IsOrderingValid(original_options, problem_impl, &summary->error)) {
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700428 LOG(ERROR) << summary->error;
Sameer Agarwal65625f72012-09-17 12:06:57 -0700429 return;
430 }
Sameer Agarwalc58e6dc2013-02-11 16:41:06 -0800431 event_logger.AddEvent("CheckOrdering");
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700432 options.linear_solver_ordering =
433 new ParameterBlockOrdering(*original_options.linear_solver_ordering);
Sameer Agarwalc58e6dc2013-02-11 16:41:06 -0800434 event_logger.AddEvent("CopyOrdering");
Sameer Agarwal65625f72012-09-17 12:06:57 -0700435 } else {
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700436 options.linear_solver_ordering = new ParameterBlockOrdering;
Sameer Agarwal65625f72012-09-17 12:06:57 -0700437 const ProblemImpl::ParameterMap& parameter_map =
438 problem_impl->parameter_map();
439 for (ProblemImpl::ParameterMap::const_iterator it = parameter_map.begin();
440 it != parameter_map.end();
441 ++it) {
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700442 options.linear_solver_ordering->AddElementToGroup(it->first, 0);
Sameer Agarwal65625f72012-09-17 12:06:57 -0700443 }
Sameer Agarwalc58e6dc2013-02-11 16:41:06 -0800444 event_logger.AddEvent("ConstructOrdering");
Sameer Agarwal65625f72012-09-17 12:06:57 -0700445 }
446
Keir Mierle8ebb0732012-04-30 23:09:08 -0700447 // Create the three objects needed to minimize: the transformed program, the
448 // evaluator, and the linear solver.
Keir Mierle31c1e782012-08-17 16:16:32 -0700449 scoped_ptr<Program> reduced_program(CreateReducedProgram(&options,
450 problem_impl,
451 &summary->fixed_cost,
452 &summary->error));
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800453
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800454 event_logger.AddEvent("CreateReducedProgram");
Keir Mierle8ebb0732012-04-30 23:09:08 -0700455 if (reduced_program == NULL) {
456 return;
457 }
458
Sameer Agarwal977be7c2013-01-26 16:01:54 -0800459 SummarizeOrdering(options.linear_solver_ordering,
460 &(summary->linear_solver_ordering_used));
461
Keir Mierle8ebb0732012-04-30 23:09:08 -0700462 summary->num_parameter_blocks_reduced = reduced_program->NumParameterBlocks();
463 summary->num_parameters_reduced = reduced_program->NumParameters();
Keir Mierleba944212013-02-25 12:46:44 -0800464 summary->num_effective_parameters_reduced =
465 reduced_program->NumEffectiveParameters();
Keir Mierle8ebb0732012-04-30 23:09:08 -0700466 summary->num_residual_blocks_reduced = reduced_program->NumResidualBlocks();
467 summary->num_residuals_reduced = reduced_program->NumResiduals();
468
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700469 if (summary->num_parameter_blocks_reduced == 0) {
470 summary->preprocessor_time_in_seconds =
471 WallTimeInSeconds() - solver_start_time;
472
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800473 double post_process_start_time = WallTimeInSeconds();
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700474 LOG(INFO) << "Terminating: FUNCTION_TOLERANCE reached. "
475 << "No non-constant parameter blocks found.";
476
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800477 summary->initial_cost = summary->fixed_cost;
478 summary->final_cost = summary->fixed_cost;
479
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700480 // FUNCTION_TOLERANCE is the right convergence here, as we know
481 // that the objective function is constant and cannot be changed
482 // any further.
483 summary->termination_type = FUNCTION_TOLERANCE;
484
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700485 // Ensure the program state is set to the user parameters on the way out.
486 original_program->SetParameterBlockStatePtrsToUserStatePtrs();
487
488 summary->postprocessor_time_in_seconds =
489 WallTimeInSeconds() - post_process_start_time;
490 return;
491 }
492
Keir Mierle8ebb0732012-04-30 23:09:08 -0700493 scoped_ptr<LinearSolver>
494 linear_solver(CreateLinearSolver(&options, &summary->error));
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800495 event_logger.AddEvent("CreateLinearSolver");
Keir Mierle8ebb0732012-04-30 23:09:08 -0700496 if (linear_solver == NULL) {
497 return;
498 }
499
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700500 summary->linear_solver_type_given = original_options.linear_solver_type;
501 summary->linear_solver_type_used = options.linear_solver_type;
502
503 summary->preconditioner_type = options.preconditioner_type;
504
505 summary->num_linear_solver_threads_given =
506 original_options.num_linear_solver_threads;
507 summary->num_linear_solver_threads_used = options.num_linear_solver_threads;
508
509 summary->sparse_linear_algebra_library =
510 options.sparse_linear_algebra_library;
511
512 summary->trust_region_strategy_type = options.trust_region_strategy_type;
513 summary->dogleg_type = options.dogleg_type;
514
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700515 scoped_ptr<Evaluator> evaluator(CreateEvaluator(options,
516 problem_impl->parameter_map(),
517 reduced_program.get(),
518 &summary->error));
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800519
520 event_logger.AddEvent("CreateEvaluator");
521
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700522 if (evaluator == NULL) {
Keir Mierle8ebb0732012-04-30 23:09:08 -0700523 return;
524 }
525
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700526 scoped_ptr<CoordinateDescentMinimizer> inner_iteration_minimizer;
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700527 if (options.use_inner_iterations) {
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700528 if (reduced_program->parameter_blocks().size() < 2) {
529 LOG(WARNING) << "Reduced problem only contains one parameter block."
530 << "Disabling inner iterations.";
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700531 } else {
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700532 inner_iteration_minimizer.reset(
533 CreateInnerIterationMinimizer(original_options,
534 *reduced_program,
535 problem_impl->parameter_map(),
Sameer Agarwal977be7c2013-01-26 16:01:54 -0800536 summary));
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700537 if (inner_iteration_minimizer == NULL) {
538 LOG(ERROR) << summary->error;
539 return;
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700540 }
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700541 }
Keir Mierle8ebb0732012-04-30 23:09:08 -0700542 }
Sameer Agarwal9f9488b2013-05-23 09:40:40 -0700543 event_logger.AddEvent("CreateInnerIterationMinimizer");
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800544
Keir Mierle8ebb0732012-04-30 23:09:08 -0700545 // The optimizer works on contiguous parameter vectors; allocate some.
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700546 Vector parameters(reduced_program->NumParameters());
Keir Mierle8ebb0732012-04-30 23:09:08 -0700547
548 // Collect the discontiguous parameters into a contiguous state vector.
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700549 reduced_program->ParameterBlocksToStateVector(parameters.data());
Keir Mierle8ebb0732012-04-30 23:09:08 -0700550
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700551 Vector original_parameters = parameters;
552
Petter Strandmark76533b32012-09-04 22:08:50 -0700553 double minimizer_start_time = WallTimeInSeconds();
Sameer Agarwalfa015192012-06-11 14:21:42 -0700554 summary->preprocessor_time_in_seconds =
555 minimizer_start_time - solver_start_time;
556
Keir Mierle8ebb0732012-04-30 23:09:08 -0700557 // Run the optimization.
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800558 TrustRegionMinimize(options,
559 reduced_program.get(),
560 inner_iteration_minimizer.get(),
561 evaluator.get(),
562 linear_solver.get(),
563 parameters.data(),
564 summary);
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800565 event_logger.AddEvent("Minimize");
Keir Mierle8ebb0732012-04-30 23:09:08 -0700566
Sameer Agarwalf102a682013-02-11 15:08:40 -0800567 SetSummaryFinalCost(summary);
568
Keir Mierle8ebb0732012-04-30 23:09:08 -0700569 // If the user aborted mid-optimization or the optimization
570 // terminated because of a numerical failure, then return without
571 // updating user state.
572 if (summary->termination_type == USER_ABORT ||
573 summary->termination_type == NUMERICAL_FAILURE) {
574 return;
575 }
576
Petter Strandmark76533b32012-09-04 22:08:50 -0700577 double post_process_start_time = WallTimeInSeconds();
Sameer Agarwal4997cbc2012-07-02 12:44:34 -0700578
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800579 // Push the contiguous optimized parameters back to the user's
580 // parameters.
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700581 reduced_program->StateVectorToParameterBlocks(parameters.data());
Keir Mierle8ebb0732012-04-30 23:09:08 -0700582 reduced_program->CopyParameterBlockStateToUserState();
583
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800584 // Ensure the program state is set to the user parameters on the way
585 // out.
Sameer Agarwal4997cbc2012-07-02 12:44:34 -0700586 original_program->SetParameterBlockStatePtrsToUserStatePtrs();
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700587
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800588 const map<string, double>& linear_solver_time_statistics =
589 linear_solver->TimeStatistics();
590 summary->linear_solver_time_in_seconds =
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800591 FindWithDefault(linear_solver_time_statistics,
592 "LinearSolver::Solve",
593 0.0);
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800594
595 const map<string, double>& evaluator_time_statistics =
596 evaluator->TimeStatistics();
597
598 summary->residual_evaluation_time_in_seconds =
599 FindWithDefault(evaluator_time_statistics, "Evaluator::Residual", 0.0);
600 summary->jacobian_evaluation_time_in_seconds =
601 FindWithDefault(evaluator_time_statistics, "Evaluator::Jacobian", 0.0);
602
Keir Mierle8ebb0732012-04-30 23:09:08 -0700603 // Stick a fork in it, we're done.
Petter Strandmark76533b32012-09-04 22:08:50 -0700604 summary->postprocessor_time_in_seconds =
605 WallTimeInSeconds() - post_process_start_time;
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800606 event_logger.AddEvent("PostProcess");
Keir Mierle8ebb0732012-04-30 23:09:08 -0700607}
608
Sameer Agarwal8140f0f2013-03-12 09:45:08 -0700609
610#ifndef CERES_NO_LINE_SEARCH_MINIMIZER
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800611void SolverImpl::LineSearchSolve(const Solver::Options& original_options,
612 ProblemImpl* original_problem_impl,
613 Solver::Summary* summary) {
614 double solver_start_time = WallTimeInSeconds();
615
616 Program* original_program = original_problem_impl->mutable_program();
617 ProblemImpl* problem_impl = original_problem_impl;
618
619 // Reset the summary object to its default values.
620 *CHECK_NOTNULL(summary) = Solver::Summary();
621
Sameer Agarwald010de52013-02-15 14:26:56 -0800622 summary->minimizer_type = LINE_SEARCH;
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800623 summary->line_search_direction_type =
624 original_options.line_search_direction_type;
Sameer Agarwald010de52013-02-15 14:26:56 -0800625 summary->max_lbfgs_rank = original_options.max_lbfgs_rank;
626 summary->line_search_type = original_options.line_search_type;
Sameer Agarwal4f010b22013-07-01 08:01:01 -0700627 summary->line_search_interpolation_type =
628 original_options.line_search_interpolation_type;
629 summary->nonlinear_conjugate_gradient_type =
630 original_options.nonlinear_conjugate_gradient_type;
631
Sameer Agarwal1c70ae92013-06-30 12:50:43 -0700632 summary->num_parameter_blocks = original_program->NumParameterBlocks();
633 summary->num_parameters = original_program->NumParameters();
634 summary->num_residual_blocks = original_program->NumResidualBlocks();
635 summary->num_residuals = original_program->NumResiduals();
Sameer Agarwal7a8f7972013-07-03 09:03:55 -0700636 summary->num_effective_parameters =
637 original_program->NumEffectiveParameters();
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800638
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100639 // Validate values for configuration parameters supplied by user.
640 if ((original_options.line_search_direction_type == ceres::BFGS ||
641 original_options.line_search_direction_type == ceres::LBFGS) &&
642 original_options.line_search_type != ceres::WOLFE) {
643 summary->error =
644 string("Invalid configuration: require line_search_type == "
645 "ceres::WOLFE when using (L)BFGS to ensure that underlying "
646 "assumptions are guaranteed to be satisfied.");
647 LOG(ERROR) << summary->error;
648 return;
649 }
650 if (original_options.max_lbfgs_rank == 0) {
651 summary->error =
652 string("Invalid configuration: require max_lbfgs_rank > 0");
653 LOG(ERROR) << summary->error;
654 return;
655 }
656 if (original_options.min_line_search_step_size <= 0.0) {
657 summary->error = "Invalid configuration: min_line_search_step_size <= 0.0.";
658 LOG(ERROR) << summary->error;
659 return;
660 }
661 if (original_options.line_search_sufficient_function_decrease <= 0.0) {
662 summary->error =
663 string("Invalid configuration: require ") +
664 string("line_search_sufficient_function_decrease <= 0.0.");
665 LOG(ERROR) << summary->error;
666 return;
667 }
668 if (original_options.max_line_search_step_contraction <= 0.0 ||
669 original_options.max_line_search_step_contraction >= 1.0) {
670 summary->error = string("Invalid configuration: require ") +
671 string("0.0 < max_line_search_step_contraction < 1.0.");
672 LOG(ERROR) << summary->error;
673 return;
674 }
675 if (original_options.min_line_search_step_contraction <=
676 original_options.max_line_search_step_contraction ||
677 original_options.min_line_search_step_contraction > 1.0) {
678 summary->error = string("Invalid configuration: require ") +
679 string("max_line_search_step_contraction < ") +
680 string("min_line_search_step_contraction <= 1.0.");
681 LOG(ERROR) << summary->error;
682 return;
683 }
684 // Warn user if they have requested BISECTION interpolation, but constraints
685 // on max/min step size change during line search prevent bisection scaling
686 // from occurring. Warn only, as this is likely a user mistake, but one which
687 // does not prevent us from continuing.
688 LOG_IF(WARNING,
689 (original_options.line_search_interpolation_type == ceres::BISECTION &&
690 (original_options.max_line_search_step_contraction > 0.5 ||
691 original_options.min_line_search_step_contraction < 0.5)))
692 << "Line search interpolation type is BISECTION, but specified "
693 << "max_line_search_step_contraction: "
694 << original_options.max_line_search_step_contraction << ", and "
695 << "min_line_search_step_contraction: "
696 << original_options.min_line_search_step_contraction
697 << ", prevent bisection (0.5) scaling, continuing with solve regardless.";
698 if (original_options.max_num_line_search_step_size_iterations <= 0) {
699 summary->error = string("Invalid configuration: require ") +
700 string("max_num_line_search_step_size_iterations > 0.");
701 LOG(ERROR) << summary->error;
702 return;
703 }
704 if (original_options.line_search_sufficient_curvature_decrease <=
705 original_options.line_search_sufficient_function_decrease ||
706 original_options.line_search_sufficient_curvature_decrease > 1.0) {
707 summary->error = string("Invalid configuration: require ") +
708 string("line_search_sufficient_function_decrease < ") +
709 string("line_search_sufficient_curvature_decrease < 1.0.");
710 LOG(ERROR) << summary->error;
711 return;
712 }
713 if (original_options.max_line_search_step_expansion <= 1.0) {
714 summary->error = string("Invalid configuration: require ") +
715 string("max_line_search_step_expansion > 1.0.");
716 LOG(ERROR) << summary->error;
717 return;
718 }
719
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800720 // Empty programs are usually a user error.
721 if (summary->num_parameter_blocks == 0) {
722 summary->error = "Problem contains no parameter blocks.";
723 LOG(ERROR) << summary->error;
724 return;
725 }
726
727 if (summary->num_residual_blocks == 0) {
728 summary->error = "Problem contains no residual blocks.";
729 LOG(ERROR) << summary->error;
730 return;
731 }
732
733 Solver::Options options(original_options);
734
735 // This ensures that we get a Block Jacobian Evaluator along with
736 // none of the Schur nonsense. This file will have to be extensively
737 // refactored to deal with the various bits of cleanups related to
738 // line search.
739 options.linear_solver_type = CGNR;
740
741 options.linear_solver_ordering = NULL;
742 options.inner_iteration_ordering = NULL;
743
744#ifndef CERES_USE_OPENMP
745 if (options.num_threads > 1) {
746 LOG(WARNING)
747 << "OpenMP support is not compiled into this binary; "
748 << "only options.num_threads=1 is supported. Switching "
749 << "to single threaded mode.";
750 options.num_threads = 1;
751 }
Sameer Agarwal8140f0f2013-03-12 09:45:08 -0700752#endif // CERES_USE_OPENMP
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800753
754 summary->num_threads_given = original_options.num_threads;
755 summary->num_threads_used = options.num_threads;
756
757 if (original_options.linear_solver_ordering != NULL) {
758 if (!IsOrderingValid(original_options, problem_impl, &summary->error)) {
759 LOG(ERROR) << summary->error;
760 return;
761 }
762 options.linear_solver_ordering =
763 new ParameterBlockOrdering(*original_options.linear_solver_ordering);
764 } else {
765 options.linear_solver_ordering = new ParameterBlockOrdering;
766 const ProblemImpl::ParameterMap& parameter_map =
767 problem_impl->parameter_map();
768 for (ProblemImpl::ParameterMap::const_iterator it = parameter_map.begin();
769 it != parameter_map.end();
770 ++it) {
771 options.linear_solver_ordering->AddElementToGroup(it->first, 0);
772 }
773 }
774
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800775 original_program->SetParameterBlockStatePtrsToUserStatePtrs();
776
777 // If the user requests gradient checking, construct a new
778 // ProblemImpl by wrapping the CostFunctions of problem_impl inside
779 // GradientCheckingCostFunction and replacing problem_impl with
780 // gradient_checking_problem_impl.
781 scoped_ptr<ProblemImpl> gradient_checking_problem_impl;
782 if (options.check_gradients) {
783 VLOG(1) << "Checking Gradients";
784 gradient_checking_problem_impl.reset(
785 CreateGradientCheckingProblemImpl(
786 problem_impl,
787 options.numeric_derivative_relative_step_size,
788 options.gradient_check_relative_precision));
789
790 // From here on, problem_impl will point to the gradient checking
791 // version.
792 problem_impl = gradient_checking_problem_impl.get();
793 }
794
795 // Create the three objects needed to minimize: the transformed program, the
796 // evaluator, and the linear solver.
797 scoped_ptr<Program> reduced_program(CreateReducedProgram(&options,
798 problem_impl,
799 &summary->fixed_cost,
800 &summary->error));
801 if (reduced_program == NULL) {
802 return;
803 }
804
805 summary->num_parameter_blocks_reduced = reduced_program->NumParameterBlocks();
806 summary->num_parameters_reduced = reduced_program->NumParameters();
807 summary->num_residual_blocks_reduced = reduced_program->NumResidualBlocks();
Sameer Agarwal1c70ae92013-06-30 12:50:43 -0700808 summary->num_effective_parameters_reduced =
809 reduced_program->NumEffectiveParameters();
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800810 summary->num_residuals_reduced = reduced_program->NumResiduals();
811
812 if (summary->num_parameter_blocks_reduced == 0) {
813 summary->preprocessor_time_in_seconds =
814 WallTimeInSeconds() - solver_start_time;
815
816 LOG(INFO) << "Terminating: FUNCTION_TOLERANCE reached. "
817 << "No non-constant parameter blocks found.";
818
819 // FUNCTION_TOLERANCE is the right convergence here, as we know
820 // that the objective function is constant and cannot be changed
821 // any further.
822 summary->termination_type = FUNCTION_TOLERANCE;
823
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800824 const double post_process_start_time = WallTimeInSeconds();
Sameer Agarwalf102a682013-02-11 15:08:40 -0800825
826 SetSummaryFinalCost(summary);
827
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800828 // Ensure the program state is set to the user parameters on the way out.
829 original_program->SetParameterBlockStatePtrsToUserStatePtrs();
Sameer Agarwal956ed7e2013-02-19 07:06:15 -0800830 summary->postprocessor_time_in_seconds =
831 WallTimeInSeconds() - post_process_start_time;
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800832 return;
833 }
834
835 scoped_ptr<Evaluator> evaluator(CreateEvaluator(options,
836 problem_impl->parameter_map(),
837 reduced_program.get(),
838 &summary->error));
839 if (evaluator == NULL) {
840 return;
841 }
842
843 // The optimizer works on contiguous parameter vectors; allocate some.
844 Vector parameters(reduced_program->NumParameters());
845
846 // Collect the discontiguous parameters into a contiguous state vector.
847 reduced_program->ParameterBlocksToStateVector(parameters.data());
848
849 Vector original_parameters = parameters;
850
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800851 const double minimizer_start_time = WallTimeInSeconds();
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800852 summary->preprocessor_time_in_seconds =
853 minimizer_start_time - solver_start_time;
854
855 // Run the optimization.
856 LineSearchMinimize(options,
857 reduced_program.get(),
858 evaluator.get(),
859 parameters.data(),
860 summary);
861
862 // If the user aborted mid-optimization or the optimization
863 // terminated because of a numerical failure, then return without
864 // updating user state.
865 if (summary->termination_type == USER_ABORT ||
866 summary->termination_type == NUMERICAL_FAILURE) {
867 return;
868 }
869
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800870 const double post_process_start_time = WallTimeInSeconds();
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800871
872 // Push the contiguous optimized parameters back to the user's parameters.
873 reduced_program->StateVectorToParameterBlocks(parameters.data());
874 reduced_program->CopyParameterBlockStateToUserState();
875
Sameer Agarwalf102a682013-02-11 15:08:40 -0800876 SetSummaryFinalCost(summary);
877
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800878 // Ensure the program state is set to the user parameters on the way out.
879 original_program->SetParameterBlockStatePtrsToUserStatePtrs();
880
Sameer Agarwald010de52013-02-15 14:26:56 -0800881 const map<string, double>& evaluator_time_statistics =
882 evaluator->TimeStatistics();
883
884 summary->residual_evaluation_time_in_seconds =
885 FindWithDefault(evaluator_time_statistics, "Evaluator::Residual", 0.0);
886 summary->jacobian_evaluation_time_in_seconds =
887 FindWithDefault(evaluator_time_statistics, "Evaluator::Jacobian", 0.0);
888
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800889 // Stick a fork in it, we're done.
890 summary->postprocessor_time_in_seconds =
891 WallTimeInSeconds() - post_process_start_time;
892}
Sameer Agarwal700d50d2013-03-12 16:12:42 -0700893#endif // CERES_NO_LINE_SEARCH_MINIMIZER
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800894
Sameer Agarwal65625f72012-09-17 12:06:57 -0700895bool SolverImpl::IsOrderingValid(const Solver::Options& options,
896 const ProblemImpl* problem_impl,
897 string* error) {
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700898 if (options.linear_solver_ordering->NumElements() !=
Sameer Agarwal65625f72012-09-17 12:06:57 -0700899 problem_impl->NumParameterBlocks()) {
900 *error = "Number of parameter blocks in user supplied ordering "
901 "does not match the number of parameter blocks in the problem";
902 return false;
903 }
904
905 const Program& program = problem_impl->program();
906 const vector<ParameterBlock*>& parameter_blocks = program.parameter_blocks();
Sameer Agarwal65625f72012-09-17 12:06:57 -0700907 for (vector<ParameterBlock*>::const_iterator it = parameter_blocks.begin();
908 it != parameter_blocks.end();
909 ++it) {
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700910 if (!options.linear_solver_ordering
911 ->IsMember(const_cast<double*>((*it)->user_state()))) {
Sameer Agarwal65625f72012-09-17 12:06:57 -0700912 *error = "Problem contains a parameter block that is not in "
913 "the user specified ordering.";
914 return false;
915 }
916 }
917
918 if (IsSchurType(options.linear_solver_type) &&
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700919 options.linear_solver_ordering->NumGroups() > 1) {
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700920 const vector<ResidualBlock*>& residual_blocks = program.residual_blocks();
Sameer Agarwal65625f72012-09-17 12:06:57 -0700921 const set<double*>& e_blocks =
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700922 options.linear_solver_ordering->group_to_elements().begin()->second;
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700923 if (!IsParameterBlockSetIndependent(e_blocks, residual_blocks)) {
924 *error = "The user requested the use of a Schur type solver. "
925 "But the first elimination group in the ordering is not an "
926 "independent set.";
927 return false;
Sameer Agarwal65625f72012-09-17 12:06:57 -0700928 }
929 }
930 return true;
931}
932
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800933bool SolverImpl::IsParameterBlockSetIndependent(
934 const set<double*>& parameter_block_ptrs,
935 const vector<ResidualBlock*>& residual_blocks) {
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700936 // Loop over each residual block and ensure that no two parameter
937 // blocks in the same residual block are part of
938 // parameter_block_ptrs as that would violate the assumption that it
939 // is an independent set in the Hessian matrix.
940 for (vector<ResidualBlock*>::const_iterator it = residual_blocks.begin();
941 it != residual_blocks.end();
942 ++it) {
943 ParameterBlock* const* parameter_blocks = (*it)->parameter_blocks();
944 const int num_parameter_blocks = (*it)->NumParameterBlocks();
945 int count = 0;
946 for (int i = 0; i < num_parameter_blocks; ++i) {
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700947 count += parameter_block_ptrs.count(
948 parameter_blocks[i]->mutable_user_state());
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700949 }
950 if (count > 1) {
951 return false;
952 }
953 }
954 return true;
955}
956
957
Keir Mierle8ebb0732012-04-30 23:09:08 -0700958// Strips varying parameters and residuals, maintaining order, and updating
959// num_eliminate_blocks.
960bool SolverImpl::RemoveFixedBlocksFromProgram(Program* program,
Sameer Agarwal2c94eed2012-10-01 16:34:37 -0700961 ParameterBlockOrdering* ordering,
Markus Moll8de78db2012-07-14 15:49:11 +0200962 double* fixed_cost,
Keir Mierle8ebb0732012-04-30 23:09:08 -0700963 string* error) {
Keir Mierle8ebb0732012-04-30 23:09:08 -0700964 vector<ParameterBlock*>* parameter_blocks =
965 program->mutable_parameter_blocks();
966
Markus Moll8de78db2012-07-14 15:49:11 +0200967 scoped_array<double> residual_block_evaluate_scratch;
968 if (fixed_cost != NULL) {
Keir Mierlec161a9d2012-07-18 14:01:48 -0700969 residual_block_evaluate_scratch.reset(
Markus Moll8de78db2012-07-14 15:49:11 +0200970 new double[program->MaxScratchDoublesNeededForEvaluate()]);
971 *fixed_cost = 0.0;
972 }
973
Keir Mierle8ebb0732012-04-30 23:09:08 -0700974 // Mark all the parameters as unused. Abuse the index member of the parameter
975 // blocks for the marking.
976 for (int i = 0; i < parameter_blocks->size(); ++i) {
977 (*parameter_blocks)[i]->set_index(-1);
978 }
979
980 // Filter out residual that have all-constant parameters, and mark all the
981 // parameter blocks that appear in residuals.
982 {
983 vector<ResidualBlock*>* residual_blocks =
984 program->mutable_residual_blocks();
985 int j = 0;
986 for (int i = 0; i < residual_blocks->size(); ++i) {
987 ResidualBlock* residual_block = (*residual_blocks)[i];
988 int num_parameter_blocks = residual_block->NumParameterBlocks();
989
990 // Determine if the residual block is fixed, and also mark varying
991 // parameters that appear in the residual block.
992 bool all_constant = true;
993 for (int k = 0; k < num_parameter_blocks; k++) {
994 ParameterBlock* parameter_block = residual_block->parameter_blocks()[k];
995 if (!parameter_block->IsConstant()) {
996 all_constant = false;
997 parameter_block->set_index(1);
998 }
999 }
1000
1001 if (!all_constant) {
1002 (*residual_blocks)[j++] = (*residual_blocks)[i];
Markus Moll8de78db2012-07-14 15:49:11 +02001003 } else if (fixed_cost != NULL) {
1004 // The residual is constant and will be removed, so its cost is
1005 // added to the variable fixed_cost.
1006 double cost = 0.0;
Sameer Agarwal039ff072013-02-26 09:15:39 -08001007 if (!residual_block->Evaluate(true,
1008 &cost,
1009 NULL,
1010 NULL,
1011 residual_block_evaluate_scratch.get())) {
Markus Moll8de78db2012-07-14 15:49:11 +02001012 *error = StringPrintf("Evaluation of the residual %d failed during "
1013 "removal of fixed residual blocks.", i);
1014 return false;
1015 }
1016 *fixed_cost += cost;
Keir Mierle8ebb0732012-04-30 23:09:08 -07001017 }
1018 }
1019 residual_blocks->resize(j);
1020 }
1021
1022 // Filter out unused or fixed parameter blocks, and update
Sameer Agarwal65625f72012-09-17 12:06:57 -07001023 // the ordering.
Keir Mierle8ebb0732012-04-30 23:09:08 -07001024 {
1025 vector<ParameterBlock*>* parameter_blocks =
1026 program->mutable_parameter_blocks();
1027 int j = 0;
1028 for (int i = 0; i < parameter_blocks->size(); ++i) {
1029 ParameterBlock* parameter_block = (*parameter_blocks)[i];
1030 if (parameter_block->index() == 1) {
1031 (*parameter_blocks)[j++] = parameter_block;
Sameer Agarwal65625f72012-09-17 12:06:57 -07001032 } else {
Sameer Agarwal2c94eed2012-10-01 16:34:37 -07001033 ordering->Remove(parameter_block->mutable_user_state());
Keir Mierle8ebb0732012-04-30 23:09:08 -07001034 }
1035 }
1036 parameter_blocks->resize(j);
1037 }
1038
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001039 if (!(((program->NumResidualBlocks() == 0) &&
Keir Mierle8ebb0732012-04-30 23:09:08 -07001040 (program->NumParameterBlocks() == 0)) ||
1041 ((program->NumResidualBlocks() != 0) &&
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001042 (program->NumParameterBlocks() != 0)))) {
1043 *error = "Congratulations, you found a bug in Ceres. Please report it.";
1044 return false;
1045 }
1046
Keir Mierle8ebb0732012-04-30 23:09:08 -07001047 return true;
1048}
1049
1050Program* SolverImpl::CreateReducedProgram(Solver::Options* options,
1051 ProblemImpl* problem_impl,
Markus Moll8de78db2012-07-14 15:49:11 +02001052 double* fixed_cost,
Keir Mierle8ebb0732012-04-30 23:09:08 -07001053 string* error) {
Sameer Agarwal68b32a92012-10-06 23:10:51 -07001054 CHECK_NOTNULL(options->linear_solver_ordering);
Keir Mierle8ebb0732012-04-30 23:09:08 -07001055 Program* original_program = problem_impl->mutable_program();
1056 scoped_ptr<Program> transformed_program(new Program(*original_program));
Sameer Agarwal42a84b82013-02-01 12:22:53 -08001057
Sameer Agarwal68b32a92012-10-06 23:10:51 -07001058 ParameterBlockOrdering* linear_solver_ordering =
1059 options->linear_solver_ordering;
Sameer Agarwal65625f72012-09-17 12:06:57 -07001060 const int min_group_id =
Sameer Agarwal68b32a92012-10-06 23:10:51 -07001061 linear_solver_ordering->group_to_elements().begin()->first;
Keir Mierle8ebb0732012-04-30 23:09:08 -07001062
1063 if (!RemoveFixedBlocksFromProgram(transformed_program.get(),
Sameer Agarwal68b32a92012-10-06 23:10:51 -07001064 linear_solver_ordering,
Markus Moll8de78db2012-07-14 15:49:11 +02001065 fixed_cost,
Keir Mierle8ebb0732012-04-30 23:09:08 -07001066 error)) {
1067 return NULL;
1068 }
1069
1070 if (transformed_program->NumParameterBlocks() == 0) {
1071 LOG(WARNING) << "No varying parameter blocks to optimize; "
1072 << "bailing early.";
1073 return transformed_program.release();
1074 }
1075
Sameer Agarwal65625f72012-09-17 12:06:57 -07001076 if (IsSchurType(options->linear_solver_type) &&
Sameer Agarwal68b32a92012-10-06 23:10:51 -07001077 linear_solver_ordering->GroupSize(min_group_id) == 0) {
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001078 // If the user requested the use of a Schur type solver, and
1079 // supplied a non-NULL linear_solver_ordering object with more than
1080 // one elimination group, then it can happen that after all the
1081 // parameter blocks which are fixed or unused have been removed from
1082 // the program and the ordering, there are no more parameter blocks
1083 // in the first elimination group.
1084 //
1085 // In such a case, the use of a Schur type solver is not possible,
1086 // as they assume there is at least one e_block. Thus, we
1087 // automatically switch to the closest solver to the one indicated
1088 // by the user.
1089 AlternateLinearSolverForSchurTypeLinearSolver(options);
Sameer Agarwal65625f72012-09-17 12:06:57 -07001090 }
Keir Mierle8ebb0732012-04-30 23:09:08 -07001091
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001092 if (IsSchurType(options->linear_solver_type)) {
Sameer Agarwalac626962013-05-06 07:04:26 -07001093 if (!ReorderProgramForSchurTypeLinearSolver(
1094 options->linear_solver_type,
1095 options->sparse_linear_algebra_library,
1096 problem_impl->parameter_map(),
1097 linear_solver_ordering,
1098 transformed_program.get(),
1099 error)) {
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001100 return NULL;
1101 }
1102 return transformed_program.release();
1103 }
Sameer Agarwal42a84b82013-02-01 12:22:53 -08001104
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001105 if (options->linear_solver_type == SPARSE_NORMAL_CHOLESKY) {
1106 if (!ReorderProgramForSparseNormalCholesky(
1107 options->sparse_linear_algebra_library,
1108 linear_solver_ordering,
1109 transformed_program.get(),
1110 error)) {
1111 return NULL;
1112 }
1113
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001114 return transformed_program.release();
1115 }
1116
Keir Mierle8ebb0732012-04-30 23:09:08 -07001117 transformed_program->SetParameterOffsetsAndIndex();
1118 return transformed_program.release();
1119}
1120
1121LinearSolver* SolverImpl::CreateLinearSolver(Solver::Options* options,
1122 string* error) {
Sameer Agarwal65625f72012-09-17 12:06:57 -07001123 CHECK_NOTNULL(options);
Sameer Agarwal68b32a92012-10-06 23:10:51 -07001124 CHECK_NOTNULL(options->linear_solver_ordering);
Sameer Agarwal65625f72012-09-17 12:06:57 -07001125 CHECK_NOTNULL(error);
1126
Sameer Agarwal5ecd2512012-06-17 16:34:03 -07001127 if (options->trust_region_strategy_type == DOGLEG) {
1128 if (options->linear_solver_type == ITERATIVE_SCHUR ||
1129 options->linear_solver_type == CGNR) {
1130 *error = "DOGLEG only supports exact factorization based linear "
1131 "solvers. If you want to use an iterative solver please "
1132 "use LEVENBERG_MARQUARDT as the trust_region_strategy_type";
1133 return NULL;
1134 }
1135 }
1136
Keir Mierle8ebb0732012-04-30 23:09:08 -07001137#ifdef CERES_NO_SUITESPARSE
Sameer Agarwalb0518732012-05-29 00:27:57 -07001138 if (options->linear_solver_type == SPARSE_NORMAL_CHOLESKY &&
1139 options->sparse_linear_algebra_library == SUITE_SPARSE) {
1140 *error = "Can't use SPARSE_NORMAL_CHOLESKY with SUITESPARSE because "
1141 "SuiteSparse was not enabled when Ceres was built.";
Keir Mierle8ebb0732012-04-30 23:09:08 -07001142 return NULL;
1143 }
Sameer Agarwal65625f72012-09-17 12:06:57 -07001144
Petter Strandmarkdd1a2762012-09-19 12:55:05 +02001145 if (options->preconditioner_type == CLUSTER_JACOBI) {
Sameer Agarwal65625f72012-09-17 12:06:57 -07001146 *error = "CLUSTER_JACOBI preconditioner not suppored. Please build Ceres "
1147 "with SuiteSparse support.";
1148 return NULL;
1149 }
1150
Petter Strandmarkdd1a2762012-09-19 12:55:05 +02001151 if (options->preconditioner_type == CLUSTER_TRIDIAGONAL) {
Sameer Agarwal65625f72012-09-17 12:06:57 -07001152 *error = "CLUSTER_TRIDIAGONAL preconditioner not suppored. Please build "
1153 "Ceres with SuiteSparse support.";
1154 return NULL;
1155 }
Sameer Agarwalb0518732012-05-29 00:27:57 -07001156#endif
1157
1158#ifdef CERES_NO_CXSPARSE
1159 if (options->linear_solver_type == SPARSE_NORMAL_CHOLESKY &&
1160 options->sparse_linear_algebra_library == CX_SPARSE) {
1161 *error = "Can't use SPARSE_NORMAL_CHOLESKY with CXSPARSE because "
1162 "CXSparse was not enabled when Ceres was built.";
1163 return NULL;
1164 }
1165#endif
1166
Sameer Agarwal65625f72012-09-17 12:06:57 -07001167#if defined(CERES_NO_SUITESPARSE) && defined(CERES_NO_CXSPARSE)
Keir Mierle27dd0d32012-10-15 13:52:36 -07001168 if (options->linear_solver_type == SPARSE_SCHUR) {
Sameer Agarwal65625f72012-09-17 12:06:57 -07001169 *error = "Can't use SPARSE_SCHUR because neither SuiteSparse nor"
1170 "CXSparse was enabled when Ceres was compiled.";
1171 return NULL;
1172 }
1173#endif
Keir Mierle8ebb0732012-04-30 23:09:08 -07001174
Sameer Agarwaleeedd2e2013-07-07 23:04:31 -07001175 if (options->max_linear_solver_iterations <= 0) {
1176 *error = "Solver::Options::max_linear_solver_iterations is not positive.";
Keir Mierle8ebb0732012-04-30 23:09:08 -07001177 return NULL;
1178 }
Sameer Agarwaleeedd2e2013-07-07 23:04:31 -07001179 if (options->min_linear_solver_iterations <= 0) {
1180 *error = "Solver::Options::min_linear_solver_iterations is not positive.";
Keir Mierle8ebb0732012-04-30 23:09:08 -07001181 return NULL;
1182 }
Sameer Agarwaleeedd2e2013-07-07 23:04:31 -07001183 if (options->min_linear_solver_iterations >
1184 options->max_linear_solver_iterations) {
1185 *error = "Solver::Options::min_linear_solver_iterations > "
1186 "Solver::Options::max_linear_solver_iterations.";
Keir Mierle8ebb0732012-04-30 23:09:08 -07001187 return NULL;
1188 }
1189
1190 LinearSolver::Options linear_solver_options;
Keir Mierle8ebb0732012-04-30 23:09:08 -07001191 linear_solver_options.min_num_iterations =
Sameer Agarwaleeedd2e2013-07-07 23:04:31 -07001192 options->min_linear_solver_iterations;
Keir Mierle8ebb0732012-04-30 23:09:08 -07001193 linear_solver_options.max_num_iterations =
Sameer Agarwaleeedd2e2013-07-07 23:04:31 -07001194 options->max_linear_solver_iterations;
Keir Mierle8ebb0732012-04-30 23:09:08 -07001195 linear_solver_options.type = options->linear_solver_type;
1196 linear_solver_options.preconditioner_type = options->preconditioner_type;
Sameer Agarwalb0518732012-05-29 00:27:57 -07001197 linear_solver_options.sparse_linear_algebra_library =
1198 options->sparse_linear_algebra_library;
Sameer Agarwal9189f4e2013-04-19 17:09:49 -07001199 linear_solver_options.use_postordering = options->use_postordering;
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001200
1201 // Ignore user's postordering preferences and force it to be true if
1202 // cholmod_camd is not available. This ensures that the linear
1203 // solver does not assume that a fill-reducing pre-ordering has been
1204 // done.
1205#if !defined(CERES_NO_SUITESPARSE) && defined(CERES_NO_CAMD)
1206 if (IsSchurType(linear_solver_options.type) &&
1207 linear_solver_options.sparse_linear_algebra_library == SUITE_SPARSE) {
1208 linear_solver_options.use_postordering = true;
1209 }
1210#endif
1211
Keir Mierle8ebb0732012-04-30 23:09:08 -07001212 linear_solver_options.num_threads = options->num_linear_solver_threads;
Keir Mierle8ebb0732012-04-30 23:09:08 -07001213 options->num_linear_solver_threads = linear_solver_options.num_threads;
1214
Sameer Agarwal65625f72012-09-17 12:06:57 -07001215 const map<int, set<double*> >& groups =
Sameer Agarwal68b32a92012-10-06 23:10:51 -07001216 options->linear_solver_ordering->group_to_elements();
Sameer Agarwal65625f72012-09-17 12:06:57 -07001217 for (map<int, set<double*> >::const_iterator it = groups.begin();
1218 it != groups.end();
1219 ++it) {
1220 linear_solver_options.elimination_groups.push_back(it->second.size());
1221 }
1222 // Schur type solvers, expect at least two elimination groups. If
Sameer Agarwal1dc544a2012-09-18 18:09:03 -07001223 // there is only one elimination group, then CreateReducedProgram
Sameer Agarwal65625f72012-09-17 12:06:57 -07001224 // guarantees that this group only contains e_blocks. Thus we add a
1225 // dummy elimination group with zero blocks in it.
1226 if (IsSchurType(linear_solver_options.type) &&
1227 linear_solver_options.elimination_groups.size() == 1) {
1228 linear_solver_options.elimination_groups.push_back(0);
1229 }
1230
Keir Mierle8ebb0732012-04-30 23:09:08 -07001231 return LinearSolver::Create(linear_solver_options);
1232}
1233
Keir Mierle8ebb0732012-04-30 23:09:08 -07001234
1235// Find the minimum index of any parameter block to the given residual.
1236// Parameter blocks that have indices greater than num_eliminate_blocks are
1237// considered to have an index equal to num_eliminate_blocks.
Sergey Sharybinb53c9662013-02-25 01:14:26 +06001238static int MinParameterBlock(const ResidualBlock* residual_block,
1239 int num_eliminate_blocks) {
Keir Mierle8ebb0732012-04-30 23:09:08 -07001240 int min_parameter_block_position = num_eliminate_blocks;
1241 for (int i = 0; i < residual_block->NumParameterBlocks(); ++i) {
1242 ParameterBlock* parameter_block = residual_block->parameter_blocks()[i];
Keir Mierle32de18d2012-05-13 16:45:05 -07001243 if (!parameter_block->IsConstant()) {
1244 CHECK_NE(parameter_block->index(), -1)
1245 << "Did you forget to call Program::SetParameterOffsetsAndIndex()? "
1246 << "This is a Ceres bug; please contact the developers!";
1247 min_parameter_block_position = std::min(parameter_block->index(),
1248 min_parameter_block_position);
1249 }
Keir Mierle8ebb0732012-04-30 23:09:08 -07001250 }
1251 return min_parameter_block_position;
1252}
1253
1254// Reorder the residuals for program, if necessary, so that the residuals
1255// involving each E block occur together. This is a necessary condition for the
1256// Schur eliminator, which works on these "row blocks" in the jacobian.
Sameer Agarwal509f68c2013-02-20 01:39:03 -08001257bool SolverImpl::LexicographicallyOrderResidualBlocks(
1258 const int num_eliminate_blocks,
1259 Program* program,
1260 string* error) {
Sameer Agarwal9123e2f2012-09-18 21:49:06 -07001261 CHECK_GE(num_eliminate_blocks, 1)
1262 << "Congratulations, you found a Ceres bug! Please report this error "
1263 << "to the developers.";
Keir Mierle8ebb0732012-04-30 23:09:08 -07001264
1265 // Create a histogram of the number of residuals for each E block. There is an
1266 // extra bucket at the end to catch all non-eliminated F blocks.
Sameer Agarwal65625f72012-09-17 12:06:57 -07001267 vector<int> residual_blocks_per_e_block(num_eliminate_blocks + 1);
Keir Mierle8ebb0732012-04-30 23:09:08 -07001268 vector<ResidualBlock*>* residual_blocks = program->mutable_residual_blocks();
1269 vector<int> min_position_per_residual(residual_blocks->size());
1270 for (int i = 0; i < residual_blocks->size(); ++i) {
1271 ResidualBlock* residual_block = (*residual_blocks)[i];
Sameer Agarwal65625f72012-09-17 12:06:57 -07001272 int position = MinParameterBlock(residual_block, num_eliminate_blocks);
Keir Mierle8ebb0732012-04-30 23:09:08 -07001273 min_position_per_residual[i] = position;
Sameer Agarwal65625f72012-09-17 12:06:57 -07001274 DCHECK_LE(position, num_eliminate_blocks);
Keir Mierle8ebb0732012-04-30 23:09:08 -07001275 residual_blocks_per_e_block[position]++;
1276 }
1277
1278 // Run a cumulative sum on the histogram, to obtain offsets to the start of
1279 // each histogram bucket (where each bucket is for the residuals for that
1280 // E-block).
Sameer Agarwal65625f72012-09-17 12:06:57 -07001281 vector<int> offsets(num_eliminate_blocks + 1);
Keir Mierle8ebb0732012-04-30 23:09:08 -07001282 std::partial_sum(residual_blocks_per_e_block.begin(),
1283 residual_blocks_per_e_block.end(),
1284 offsets.begin());
1285 CHECK_EQ(offsets.back(), residual_blocks->size())
1286 << "Congratulations, you found a Ceres bug! Please report this error "
1287 << "to the developers.";
1288
1289 CHECK(find(residual_blocks_per_e_block.begin(),
1290 residual_blocks_per_e_block.end() - 1, 0) !=
1291 residual_blocks_per_e_block.end())
1292 << "Congratulations, you found a Ceres bug! Please report this error "
1293 << "to the developers.";
1294
1295 // Fill in each bucket with the residual blocks for its corresponding E block.
1296 // Each bucket is individually filled from the back of the bucket to the front
1297 // of the bucket. The filling order among the buckets is dictated by the
1298 // residual blocks. This loop uses the offsets as counters; subtracting one
1299 // from each offset as a residual block is placed in the bucket. When the
1300 // filling is finished, the offset pointerts should have shifted down one
1301 // entry (this is verified below).
1302 vector<ResidualBlock*> reordered_residual_blocks(
1303 (*residual_blocks).size(), static_cast<ResidualBlock*>(NULL));
1304 for (int i = 0; i < residual_blocks->size(); ++i) {
1305 int bucket = min_position_per_residual[i];
1306
1307 // Decrement the cursor, which should now point at the next empty position.
1308 offsets[bucket]--;
1309
1310 // Sanity.
1311 CHECK(reordered_residual_blocks[offsets[bucket]] == NULL)
1312 << "Congratulations, you found a Ceres bug! Please report this error "
1313 << "to the developers.";
1314
1315 reordered_residual_blocks[offsets[bucket]] = (*residual_blocks)[i];
1316 }
1317
1318 // Sanity check #1: The difference in bucket offsets should match the
1319 // histogram sizes.
Sameer Agarwal65625f72012-09-17 12:06:57 -07001320 for (int i = 0; i < num_eliminate_blocks; ++i) {
Keir Mierle8ebb0732012-04-30 23:09:08 -07001321 CHECK_EQ(residual_blocks_per_e_block[i], offsets[i + 1] - offsets[i])
1322 << "Congratulations, you found a Ceres bug! Please report this error "
1323 << "to the developers.";
1324 }
1325 // Sanity check #2: No NULL's left behind.
1326 for (int i = 0; i < reordered_residual_blocks.size(); ++i) {
1327 CHECK(reordered_residual_blocks[i] != NULL)
1328 << "Congratulations, you found a Ceres bug! Please report this error "
1329 << "to the developers.";
1330 }
1331
1332 // Now that the residuals are collected by E block, swap them in place.
1333 swap(*program->mutable_residual_blocks(), reordered_residual_blocks);
1334 return true;
1335}
1336
Sameer Agarwal509f68c2013-02-20 01:39:03 -08001337Evaluator* SolverImpl::CreateEvaluator(
1338 const Solver::Options& options,
1339 const ProblemImpl::ParameterMap& parameter_map,
1340 Program* program,
1341 string* error) {
Keir Mierle8ebb0732012-04-30 23:09:08 -07001342 Evaluator::Options evaluator_options;
1343 evaluator_options.linear_solver_type = options.linear_solver_type;
Sameer Agarwal65625f72012-09-17 12:06:57 -07001344 evaluator_options.num_eliminate_blocks =
Sameer Agarwal68b32a92012-10-06 23:10:51 -07001345 (options.linear_solver_ordering->NumGroups() > 0 &&
Sameer Agarwal65625f72012-09-17 12:06:57 -07001346 IsSchurType(options.linear_solver_type))
Sameer Agarwal68b32a92012-10-06 23:10:51 -07001347 ? (options.linear_solver_ordering
1348 ->group_to_elements().begin()
1349 ->second.size())
Sameer Agarwal9123e2f2012-09-18 21:49:06 -07001350 : 0;
Keir Mierle8ebb0732012-04-30 23:09:08 -07001351 evaluator_options.num_threads = options.num_threads;
1352 return Evaluator::Create(evaluator_options, program, error);
1353}
1354
Sameer Agarwal67a107b2012-10-08 14:35:41 -07001355CoordinateDescentMinimizer* SolverImpl::CreateInnerIterationMinimizer(
1356 const Solver::Options& options,
1357 const Program& program,
1358 const ProblemImpl::ParameterMap& parameter_map,
Sameer Agarwal977be7c2013-01-26 16:01:54 -08001359 Solver::Summary* summary) {
Sameer Agarwal9f9488b2013-05-23 09:40:40 -07001360 summary->inner_iterations_given = true;
1361
Sameer Agarwal67a107b2012-10-08 14:35:41 -07001362 scoped_ptr<CoordinateDescentMinimizer> inner_iteration_minimizer(
1363 new CoordinateDescentMinimizer);
1364 scoped_ptr<ParameterBlockOrdering> inner_iteration_ordering;
1365 ParameterBlockOrdering* ordering_ptr = NULL;
1366
1367 if (options.inner_iteration_ordering == NULL) {
1368 // Find a recursive decomposition of the Hessian matrix as a set
1369 // of independent sets of decreasing size and invert it. This
1370 // seems to work better in practice, i.e., Cameras before
1371 // points.
1372 inner_iteration_ordering.reset(new ParameterBlockOrdering);
1373 ComputeRecursiveIndependentSetOrdering(program,
1374 inner_iteration_ordering.get());
1375 inner_iteration_ordering->Reverse();
1376 ordering_ptr = inner_iteration_ordering.get();
1377 } else {
1378 const map<int, set<double*> >& group_to_elements =
1379 options.inner_iteration_ordering->group_to_elements();
1380
1381 // Iterate over each group and verify that it is an independent
1382 // set.
1383 map<int, set<double*> >::const_iterator it = group_to_elements.begin();
Sameer Agarwal509f68c2013-02-20 01:39:03 -08001384 for ( ; it != group_to_elements.end(); ++it) {
Sameer Agarwal67a107b2012-10-08 14:35:41 -07001385 if (!IsParameterBlockSetIndependent(it->second,
1386 program.residual_blocks())) {
Sameer Agarwal977be7c2013-01-26 16:01:54 -08001387 summary->error =
Sameer Agarwal67a107b2012-10-08 14:35:41 -07001388 StringPrintf("The user-provided "
1389 "parameter_blocks_for_inner_iterations does not "
1390 "form an independent set. Group Id: %d", it->first);
1391 return NULL;
1392 }
1393 }
1394 ordering_ptr = options.inner_iteration_ordering;
1395 }
1396
1397 if (!inner_iteration_minimizer->Init(program,
1398 parameter_map,
1399 *ordering_ptr,
Sameer Agarwal977be7c2013-01-26 16:01:54 -08001400 &summary->error)) {
Sameer Agarwal67a107b2012-10-08 14:35:41 -07001401 return NULL;
1402 }
1403
Sameer Agarwal9f9488b2013-05-23 09:40:40 -07001404 summary->inner_iterations_used = true;
1405 summary->inner_iteration_time_in_seconds = 0.0;
Sameer Agarwal977be7c2013-01-26 16:01:54 -08001406 SummarizeOrdering(ordering_ptr, &(summary->inner_iteration_ordering_used));
Sameer Agarwal67a107b2012-10-08 14:35:41 -07001407 return inner_iteration_minimizer.release();
1408}
1409
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001410void SolverImpl::AlternateLinearSolverForSchurTypeLinearSolver(
1411 Solver::Options* options) {
1412 if (!IsSchurType(options->linear_solver_type)) {
1413 return;
1414 }
1415
1416 string msg = "No e_blocks remaining. Switching from ";
1417 if (options->linear_solver_type == SPARSE_SCHUR) {
1418 options->linear_solver_type = SPARSE_NORMAL_CHOLESKY;
1419 msg += "SPARSE_SCHUR to SPARSE_NORMAL_CHOLESKY.";
1420 } else if (options->linear_solver_type == DENSE_SCHUR) {
1421 // TODO(sameeragarwal): This is probably not a great choice.
1422 // Ideally, we should have a DENSE_NORMAL_CHOLESKY, that can
1423 // take a BlockSparseMatrix as input.
1424 options->linear_solver_type = DENSE_QR;
1425 msg += "DENSE_SCHUR to DENSE_QR.";
1426 } else if (options->linear_solver_type == ITERATIVE_SCHUR) {
1427 options->linear_solver_type = CGNR;
1428 if (options->preconditioner_type != IDENTITY) {
1429 msg += StringPrintf("ITERATIVE_SCHUR with %s preconditioner "
1430 "to CGNR with JACOBI preconditioner.",
1431 PreconditionerTypeToString(
1432 options->preconditioner_type));
1433 // CGNR currently only supports the JACOBI preconditioner.
1434 options->preconditioner_type = JACOBI;
1435 } else {
Sameer Agarwalcbdeb792013-04-22 10:18:18 -07001436 msg += "ITERATIVE_SCHUR with IDENTITY preconditioner"
1437 "to CGNR with IDENTITY preconditioner.";
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001438 }
1439 }
1440 LOG(WARNING) << msg;
1441}
1442
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001443bool SolverImpl::ApplyUserOrdering(
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001444 const ProblemImpl::ParameterMap& parameter_map,
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001445 const ParameterBlockOrdering* parameter_block_ordering,
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001446 Program* program,
1447 string* error) {
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001448 const int num_parameter_blocks = program->NumParameterBlocks();
1449 if (parameter_block_ordering->NumElements() != num_parameter_blocks) {
1450 *error = StringPrintf("User specified ordering does not have the same "
1451 "number of parameters as the problem. The problem"
1452 "has %d blocks while the ordering has %d blocks.",
1453 num_parameter_blocks,
1454 parameter_block_ordering->NumElements());
1455 return false;
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001456 }
1457
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001458 vector<ParameterBlock*>* parameter_blocks =
1459 program->mutable_parameter_blocks();
1460 parameter_blocks->clear();
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001461
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001462 const map<int, set<double*> >& groups =
1463 parameter_block_ordering->group_to_elements();
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001464
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001465 for (map<int, set<double*> >::const_iterator group_it = groups.begin();
1466 group_it != groups.end();
1467 ++group_it) {
1468 const set<double*>& group = group_it->second;
1469 for (set<double*>::const_iterator parameter_block_ptr_it = group.begin();
1470 parameter_block_ptr_it != group.end();
1471 ++parameter_block_ptr_it) {
1472 ProblemImpl::ParameterMap::const_iterator parameter_block_it =
1473 parameter_map.find(*parameter_block_ptr_it);
1474 if (parameter_block_it == parameter_map.end()) {
1475 *error = StringPrintf("User specified ordering contains a pointer "
1476 "to a double that is not a parameter block in "
1477 "the problem. The invalid double is in group: %d",
1478 group_it->first);
1479 return false;
1480 }
1481 parameter_blocks->push_back(parameter_block_it->second);
1482 }
1483 }
1484 return true;
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001485}
1486
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001487
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001488TripletSparseMatrix* SolverImpl::CreateJacobianBlockSparsityTranspose(
1489 const Program* program) {
1490
Sameer Agarwalcbdeb792013-04-22 10:18:18 -07001491 // Matrix to store the block sparsity structure of the Jacobian.
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001492 TripletSparseMatrix* tsm =
1493 new TripletSparseMatrix(program->NumParameterBlocks(),
1494 program->NumResidualBlocks(),
1495 10 * program->NumResidualBlocks());
1496 int num_nonzeros = 0;
1497 int* rows = tsm->mutable_rows();
1498 int* cols = tsm->mutable_cols();
1499 double* values = tsm->mutable_values();
1500
1501 const vector<ResidualBlock*>& residual_blocks = program->residual_blocks();
1502 for (int c = 0; c < residual_blocks.size(); ++c) {
1503 const ResidualBlock* residual_block = residual_blocks[c];
1504 const int num_parameter_blocks = residual_block->NumParameterBlocks();
1505 ParameterBlock* const* parameter_blocks =
1506 residual_block->parameter_blocks();
1507
1508 for (int j = 0; j < num_parameter_blocks; ++j) {
1509 if (parameter_blocks[j]->IsConstant()) {
1510 continue;
1511 }
1512
1513 // Re-size the matrix if needed.
1514 if (num_nonzeros >= tsm->max_num_nonzeros()) {
Sameer Agarwal5f433c82013-06-13 06:57:58 -07001515 tsm->set_num_nonzeros(num_nonzeros);
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001516 tsm->Reserve(2 * num_nonzeros);
1517 rows = tsm->mutable_rows();
1518 cols = tsm->mutable_cols();
1519 values = tsm->mutable_values();
1520 }
1521 CHECK_LT(num_nonzeros, tsm->max_num_nonzeros());
1522
1523 const int r = parameter_blocks[j]->index();
1524 rows[num_nonzeros] = r;
1525 cols[num_nonzeros] = c;
1526 values[num_nonzeros] = 1.0;
1527 ++num_nonzeros;
1528 }
1529 }
1530
1531 tsm->set_num_nonzeros(num_nonzeros);
1532 return tsm;
1533}
1534
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001535bool SolverImpl::ReorderProgramForSchurTypeLinearSolver(
1536 const LinearSolverType linear_solver_type,
1537 const SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type,
1538 const ProblemImpl::ParameterMap& parameter_map,
1539 ParameterBlockOrdering* parameter_block_ordering,
1540 Program* program,
1541 string* error) {
1542 if (parameter_block_ordering->NumGroups() == 1) {
1543 // If the user supplied an parameter_block_ordering with just one
1544 // group, it is equivalent to the user supplying NULL as an
1545 // parameter_block_ordering. Ceres is completely free to choose the
1546 // parameter block ordering as it sees fit. For Schur type solvers,
1547 // this means that the user wishes for Ceres to identify the
1548 // e_blocks, which we do by computing a maximal independent set.
1549 vector<ParameterBlock*> schur_ordering;
Sameer Agarwala427c872013-06-24 17:50:56 -07001550 const int num_eliminate_blocks =
1551 ComputeStableSchurOrdering(*program, &schur_ordering);
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001552
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001553 CHECK_EQ(schur_ordering.size(), program->NumParameterBlocks())
1554 << "Congratulations, you found a Ceres bug! Please report this error "
1555 << "to the developers.";
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001556
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001557 // Update the parameter_block_ordering object.
1558 for (int i = 0; i < schur_ordering.size(); ++i) {
1559 double* parameter_block = schur_ordering[i]->mutable_user_state();
1560 const int group_id = (i < num_eliminate_blocks) ? 0 : 1;
1561 parameter_block_ordering->AddElementToGroup(parameter_block, group_id);
1562 }
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001563
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001564 // We could call ApplyUserOrdering but this is cheaper and
1565 // simpler.
1566 swap(*program->mutable_parameter_blocks(), schur_ordering);
1567 } else {
1568 // The user provided an ordering with more than one elimination
1569 // group. Trust the user and apply the ordering.
1570 if (!ApplyUserOrdering(parameter_map,
1571 parameter_block_ordering,
1572 program,
1573 error)) {
1574 return false;
1575 }
1576 }
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001577
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001578 // Pre-order the columns corresponding to the schur complement if
1579 // possible.
1580#if !defined(CERES_NO_SUITESPARSE) && !defined(CERES_NO_CAMD)
1581 if (linear_solver_type == SPARSE_SCHUR &&
1582 sparse_linear_algebra_library_type == SUITE_SPARSE) {
1583 vector<int> constraints;
1584 vector<ParameterBlock*>& parameter_blocks =
1585 *(program->mutable_parameter_blocks());
1586
1587 for (int i = 0; i < parameter_blocks.size(); ++i) {
1588 constraints.push_back(
1589 parameter_block_ordering->GroupId(
1590 parameter_blocks[i]->mutable_user_state()));
1591 }
1592
1593 // Set the offsets and index for CreateJacobianSparsityTranspose.
1594 program->SetParameterOffsetsAndIndex();
1595 // Compute a block sparse presentation of J'.
1596 scoped_ptr<TripletSparseMatrix> tsm_block_jacobian_transpose(
1597 SolverImpl::CreateJacobianBlockSparsityTranspose(program));
1598
1599 SuiteSparse ss;
1600 cholmod_sparse* block_jacobian_transpose =
1601 ss.CreateSparseMatrix(tsm_block_jacobian_transpose.get());
1602
1603 vector<int> ordering(parameter_blocks.size(), 0);
1604 ss.ConstrainedApproximateMinimumDegreeOrdering(block_jacobian_transpose,
1605 &constraints[0],
1606 &ordering[0]);
1607 ss.Free(block_jacobian_transpose);
1608
1609 const vector<ParameterBlock*> parameter_blocks_copy(parameter_blocks);
1610 for (int i = 0; i < program->NumParameterBlocks(); ++i) {
1611 parameter_blocks[i] = parameter_blocks_copy[ordering[i]];
1612 }
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001613 }
1614#endif
1615
1616 program->SetParameterOffsetsAndIndex();
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001617 // Schur type solvers also require that their residual blocks be
1618 // lexicographically ordered.
1619 const int num_eliminate_blocks =
1620 parameter_block_ordering->group_to_elements().begin()->second.size();
1621 return LexicographicallyOrderResidualBlocks(num_eliminate_blocks,
1622 program,
1623 error);
1624}
1625
1626bool SolverImpl::ReorderProgramForSparseNormalCholesky(
1627 const SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type,
1628 const ParameterBlockOrdering* parameter_block_ordering,
1629 Program* program,
1630 string* error) {
1631 // Set the offsets and index for CreateJacobianSparsityTranspose.
1632 program->SetParameterOffsetsAndIndex();
1633 // Compute a block sparse presentation of J'.
1634 scoped_ptr<TripletSparseMatrix> tsm_block_jacobian_transpose(
1635 SolverImpl::CreateJacobianBlockSparsityTranspose(program));
1636
1637 vector<int> ordering(program->NumParameterBlocks(), 0);
1638 vector<ParameterBlock*>& parameter_blocks =
1639 *(program->mutable_parameter_blocks());
1640
1641 if (sparse_linear_algebra_library_type == SUITE_SPARSE) {
1642#ifdef CERES_NO_SUITESPARSE
1643 *error = "Can't use SPARSE_NORMAL_CHOLESKY with SUITE_SPARSE because "
1644 "SuiteSparse was not enabled when Ceres was built.";
1645 return false;
1646#else
1647 SuiteSparse ss;
1648 cholmod_sparse* block_jacobian_transpose =
1649 ss.CreateSparseMatrix(tsm_block_jacobian_transpose.get());
1650
1651# ifdef CERES_NO_CAMD
1652 // No cholmod_camd, so ignore user's parameter_block_ordering and
1653 // use plain old AMD.
1654 ss.ApproximateMinimumDegreeOrdering(block_jacobian_transpose, &ordering[0]);
1655# else
1656 if (parameter_block_ordering->NumGroups() > 1) {
1657 // If the user specified more than one elimination groups use them
1658 // to constrain the ordering.
1659 vector<int> constraints;
1660 for (int i = 0; i < parameter_blocks.size(); ++i) {
1661 constraints.push_back(
1662 parameter_block_ordering->GroupId(
1663 parameter_blocks[i]->mutable_user_state()));
1664 }
1665 ss.ConstrainedApproximateMinimumDegreeOrdering(
1666 block_jacobian_transpose,
1667 &constraints[0],
1668 &ordering[0]);
1669 } else {
1670 ss.ApproximateMinimumDegreeOrdering(block_jacobian_transpose,
1671 &ordering[0]);
1672 }
1673# endif // CERES_NO_CAMD
1674
1675 ss.Free(block_jacobian_transpose);
1676#endif // CERES_NO_SUITESPARSE
1677
1678 } else if (sparse_linear_algebra_library_type == CX_SPARSE) {
1679#ifndef CERES_NO_CXSPARSE
1680
1681 // CXSparse works with J'J instead of J'. So compute the block
1682 // sparsity for J'J and compute an approximate minimum degree
1683 // ordering.
1684 CXSparse cxsparse;
1685 cs_di* block_jacobian_transpose;
1686 block_jacobian_transpose =
1687 cxsparse.CreateSparseMatrix(tsm_block_jacobian_transpose.get());
1688 cs_di* block_jacobian = cxsparse.TransposeMatrix(block_jacobian_transpose);
1689 cs_di* block_hessian =
1690 cxsparse.MatrixMatrixMultiply(block_jacobian_transpose, block_jacobian);
1691 cxsparse.Free(block_jacobian);
1692 cxsparse.Free(block_jacobian_transpose);
1693
1694 cxsparse.ApproximateMinimumDegreeOrdering(block_hessian, &ordering[0]);
1695 cxsparse.Free(block_hessian);
1696#else // CERES_NO_CXSPARSE
1697 *error = "Can't use SPARSE_NORMAL_CHOLESKY with CX_SPARSE because "
1698 "CXSparse was not enabled when Ceres was built.";
1699 return false;
1700#endif // CERES_NO_CXSPARSE
1701 } else {
1702 *error = "Unknown sparse linear algebra library.";
1703 return false;
1704 }
1705
1706 // Apply ordering.
1707 const vector<ParameterBlock*> parameter_blocks_copy(parameter_blocks);
1708 for (int i = 0; i < program->NumParameterBlocks(); ++i) {
1709 parameter_blocks[i] = parameter_blocks_copy[ordering[i]];
1710 }
1711
1712 program->SetParameterOffsetsAndIndex();
1713 return true;
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001714}
1715
Keir Mierle8ebb0732012-04-30 23:09:08 -07001716} // namespace internal
1717} // namespace ceres