blob: 8b8395cd8eb340ce47932e9d86438475fdbd02da [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;
258 trust_region_strategy_options.lm_min_diagonal = options.lm_min_diagonal;
259 trust_region_strategy_options.lm_max_diagonal = options.lm_max_diagonal;
260 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 Agarwalf4d01642012-11-26 12:55:58 -0800627 summary->num_parameter_blocks = problem_impl->NumParameterBlocks();
628 summary->num_parameters = problem_impl->NumParameters();
629 summary->num_residual_blocks = problem_impl->NumResidualBlocks();
630 summary->num_residuals = problem_impl->NumResiduals();
631
632 // Empty programs are usually a user error.
633 if (summary->num_parameter_blocks == 0) {
634 summary->error = "Problem contains no parameter blocks.";
635 LOG(ERROR) << summary->error;
636 return;
637 }
638
639 if (summary->num_residual_blocks == 0) {
640 summary->error = "Problem contains no residual blocks.";
641 LOG(ERROR) << summary->error;
642 return;
643 }
644
645 Solver::Options options(original_options);
646
647 // This ensures that we get a Block Jacobian Evaluator along with
648 // none of the Schur nonsense. This file will have to be extensively
649 // refactored to deal with the various bits of cleanups related to
650 // line search.
651 options.linear_solver_type = CGNR;
652
653 options.linear_solver_ordering = NULL;
654 options.inner_iteration_ordering = NULL;
655
656#ifndef CERES_USE_OPENMP
657 if (options.num_threads > 1) {
658 LOG(WARNING)
659 << "OpenMP support is not compiled into this binary; "
660 << "only options.num_threads=1 is supported. Switching "
661 << "to single threaded mode.";
662 options.num_threads = 1;
663 }
Sameer Agarwal8140f0f2013-03-12 09:45:08 -0700664#endif // CERES_USE_OPENMP
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800665
666 summary->num_threads_given = original_options.num_threads;
667 summary->num_threads_used = options.num_threads;
668
669 if (original_options.linear_solver_ordering != NULL) {
670 if (!IsOrderingValid(original_options, problem_impl, &summary->error)) {
671 LOG(ERROR) << summary->error;
672 return;
673 }
674 options.linear_solver_ordering =
675 new ParameterBlockOrdering(*original_options.linear_solver_ordering);
676 } else {
677 options.linear_solver_ordering = new ParameterBlockOrdering;
678 const ProblemImpl::ParameterMap& parameter_map =
679 problem_impl->parameter_map();
680 for (ProblemImpl::ParameterMap::const_iterator it = parameter_map.begin();
681 it != parameter_map.end();
682 ++it) {
683 options.linear_solver_ordering->AddElementToGroup(it->first, 0);
684 }
685 }
686
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800687 original_program->SetParameterBlockStatePtrsToUserStatePtrs();
688
689 // If the user requests gradient checking, construct a new
690 // ProblemImpl by wrapping the CostFunctions of problem_impl inside
691 // GradientCheckingCostFunction and replacing problem_impl with
692 // gradient_checking_problem_impl.
693 scoped_ptr<ProblemImpl> gradient_checking_problem_impl;
694 if (options.check_gradients) {
695 VLOG(1) << "Checking Gradients";
696 gradient_checking_problem_impl.reset(
697 CreateGradientCheckingProblemImpl(
698 problem_impl,
699 options.numeric_derivative_relative_step_size,
700 options.gradient_check_relative_precision));
701
702 // From here on, problem_impl will point to the gradient checking
703 // version.
704 problem_impl = gradient_checking_problem_impl.get();
705 }
706
707 // Create the three objects needed to minimize: the transformed program, the
708 // evaluator, and the linear solver.
709 scoped_ptr<Program> reduced_program(CreateReducedProgram(&options,
710 problem_impl,
711 &summary->fixed_cost,
712 &summary->error));
713 if (reduced_program == NULL) {
714 return;
715 }
716
717 summary->num_parameter_blocks_reduced = reduced_program->NumParameterBlocks();
718 summary->num_parameters_reduced = reduced_program->NumParameters();
719 summary->num_residual_blocks_reduced = reduced_program->NumResidualBlocks();
720 summary->num_residuals_reduced = reduced_program->NumResiduals();
721
722 if (summary->num_parameter_blocks_reduced == 0) {
723 summary->preprocessor_time_in_seconds =
724 WallTimeInSeconds() - solver_start_time;
725
726 LOG(INFO) << "Terminating: FUNCTION_TOLERANCE reached. "
727 << "No non-constant parameter blocks found.";
728
729 // FUNCTION_TOLERANCE is the right convergence here, as we know
730 // that the objective function is constant and cannot be changed
731 // any further.
732 summary->termination_type = FUNCTION_TOLERANCE;
733
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800734 const double post_process_start_time = WallTimeInSeconds();
Sameer Agarwalf102a682013-02-11 15:08:40 -0800735
736 SetSummaryFinalCost(summary);
737
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800738 // Ensure the program state is set to the user parameters on the way out.
739 original_program->SetParameterBlockStatePtrsToUserStatePtrs();
Sameer Agarwal956ed7e2013-02-19 07:06:15 -0800740 summary->postprocessor_time_in_seconds =
741 WallTimeInSeconds() - post_process_start_time;
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800742 return;
743 }
744
745 scoped_ptr<Evaluator> evaluator(CreateEvaluator(options,
746 problem_impl->parameter_map(),
747 reduced_program.get(),
748 &summary->error));
749 if (evaluator == NULL) {
750 return;
751 }
752
753 // The optimizer works on contiguous parameter vectors; allocate some.
754 Vector parameters(reduced_program->NumParameters());
755
756 // Collect the discontiguous parameters into a contiguous state vector.
757 reduced_program->ParameterBlocksToStateVector(parameters.data());
758
759 Vector original_parameters = parameters;
760
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800761 const double minimizer_start_time = WallTimeInSeconds();
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800762 summary->preprocessor_time_in_seconds =
763 minimizer_start_time - solver_start_time;
764
765 // Run the optimization.
766 LineSearchMinimize(options,
767 reduced_program.get(),
768 evaluator.get(),
769 parameters.data(),
770 summary);
771
772 // If the user aborted mid-optimization or the optimization
773 // terminated because of a numerical failure, then return without
774 // updating user state.
775 if (summary->termination_type == USER_ABORT ||
776 summary->termination_type == NUMERICAL_FAILURE) {
777 return;
778 }
779
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800780 const double post_process_start_time = WallTimeInSeconds();
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800781
782 // Push the contiguous optimized parameters back to the user's parameters.
783 reduced_program->StateVectorToParameterBlocks(parameters.data());
784 reduced_program->CopyParameterBlockStateToUserState();
785
Sameer Agarwalf102a682013-02-11 15:08:40 -0800786 SetSummaryFinalCost(summary);
787
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800788 // Ensure the program state is set to the user parameters on the way out.
789 original_program->SetParameterBlockStatePtrsToUserStatePtrs();
790
Sameer Agarwald010de52013-02-15 14:26:56 -0800791 const map<string, double>& evaluator_time_statistics =
792 evaluator->TimeStatistics();
793
794 summary->residual_evaluation_time_in_seconds =
795 FindWithDefault(evaluator_time_statistics, "Evaluator::Residual", 0.0);
796 summary->jacobian_evaluation_time_in_seconds =
797 FindWithDefault(evaluator_time_statistics, "Evaluator::Jacobian", 0.0);
798
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800799 // Stick a fork in it, we're done.
800 summary->postprocessor_time_in_seconds =
801 WallTimeInSeconds() - post_process_start_time;
802}
Sameer Agarwal700d50d2013-03-12 16:12:42 -0700803#endif // CERES_NO_LINE_SEARCH_MINIMIZER
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800804
Sameer Agarwal65625f72012-09-17 12:06:57 -0700805bool SolverImpl::IsOrderingValid(const Solver::Options& options,
806 const ProblemImpl* problem_impl,
807 string* error) {
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700808 if (options.linear_solver_ordering->NumElements() !=
Sameer Agarwal65625f72012-09-17 12:06:57 -0700809 problem_impl->NumParameterBlocks()) {
810 *error = "Number of parameter blocks in user supplied ordering "
811 "does not match the number of parameter blocks in the problem";
812 return false;
813 }
814
815 const Program& program = problem_impl->program();
816 const vector<ParameterBlock*>& parameter_blocks = program.parameter_blocks();
Sameer Agarwal65625f72012-09-17 12:06:57 -0700817 for (vector<ParameterBlock*>::const_iterator it = parameter_blocks.begin();
818 it != parameter_blocks.end();
819 ++it) {
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700820 if (!options.linear_solver_ordering
821 ->IsMember(const_cast<double*>((*it)->user_state()))) {
Sameer Agarwal65625f72012-09-17 12:06:57 -0700822 *error = "Problem contains a parameter block that is not in "
823 "the user specified ordering.";
824 return false;
825 }
826 }
827
828 if (IsSchurType(options.linear_solver_type) &&
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700829 options.linear_solver_ordering->NumGroups() > 1) {
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700830 const vector<ResidualBlock*>& residual_blocks = program.residual_blocks();
Sameer Agarwal65625f72012-09-17 12:06:57 -0700831 const set<double*>& e_blocks =
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700832 options.linear_solver_ordering->group_to_elements().begin()->second;
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700833 if (!IsParameterBlockSetIndependent(e_blocks, residual_blocks)) {
834 *error = "The user requested the use of a Schur type solver. "
835 "But the first elimination group in the ordering is not an "
836 "independent set.";
837 return false;
Sameer Agarwal65625f72012-09-17 12:06:57 -0700838 }
839 }
840 return true;
841}
842
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800843bool SolverImpl::IsParameterBlockSetIndependent(
844 const set<double*>& parameter_block_ptrs,
845 const vector<ResidualBlock*>& residual_blocks) {
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700846 // Loop over each residual block and ensure that no two parameter
847 // blocks in the same residual block are part of
848 // parameter_block_ptrs as that would violate the assumption that it
849 // is an independent set in the Hessian matrix.
850 for (vector<ResidualBlock*>::const_iterator it = residual_blocks.begin();
851 it != residual_blocks.end();
852 ++it) {
853 ParameterBlock* const* parameter_blocks = (*it)->parameter_blocks();
854 const int num_parameter_blocks = (*it)->NumParameterBlocks();
855 int count = 0;
856 for (int i = 0; i < num_parameter_blocks; ++i) {
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700857 count += parameter_block_ptrs.count(
858 parameter_blocks[i]->mutable_user_state());
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700859 }
860 if (count > 1) {
861 return false;
862 }
863 }
864 return true;
865}
866
867
Keir Mierle8ebb0732012-04-30 23:09:08 -0700868// Strips varying parameters and residuals, maintaining order, and updating
869// num_eliminate_blocks.
870bool SolverImpl::RemoveFixedBlocksFromProgram(Program* program,
Sameer Agarwal2c94eed2012-10-01 16:34:37 -0700871 ParameterBlockOrdering* ordering,
Markus Moll8de78db2012-07-14 15:49:11 +0200872 double* fixed_cost,
Keir Mierle8ebb0732012-04-30 23:09:08 -0700873 string* error) {
Keir Mierle8ebb0732012-04-30 23:09:08 -0700874 vector<ParameterBlock*>* parameter_blocks =
875 program->mutable_parameter_blocks();
876
Markus Moll8de78db2012-07-14 15:49:11 +0200877 scoped_array<double> residual_block_evaluate_scratch;
878 if (fixed_cost != NULL) {
Keir Mierlec161a9d2012-07-18 14:01:48 -0700879 residual_block_evaluate_scratch.reset(
Markus Moll8de78db2012-07-14 15:49:11 +0200880 new double[program->MaxScratchDoublesNeededForEvaluate()]);
881 *fixed_cost = 0.0;
882 }
883
Keir Mierle8ebb0732012-04-30 23:09:08 -0700884 // Mark all the parameters as unused. Abuse the index member of the parameter
885 // blocks for the marking.
886 for (int i = 0; i < parameter_blocks->size(); ++i) {
887 (*parameter_blocks)[i]->set_index(-1);
888 }
889
890 // Filter out residual that have all-constant parameters, and mark all the
891 // parameter blocks that appear in residuals.
892 {
893 vector<ResidualBlock*>* residual_blocks =
894 program->mutable_residual_blocks();
895 int j = 0;
896 for (int i = 0; i < residual_blocks->size(); ++i) {
897 ResidualBlock* residual_block = (*residual_blocks)[i];
898 int num_parameter_blocks = residual_block->NumParameterBlocks();
899
900 // Determine if the residual block is fixed, and also mark varying
901 // parameters that appear in the residual block.
902 bool all_constant = true;
903 for (int k = 0; k < num_parameter_blocks; k++) {
904 ParameterBlock* parameter_block = residual_block->parameter_blocks()[k];
905 if (!parameter_block->IsConstant()) {
906 all_constant = false;
907 parameter_block->set_index(1);
908 }
909 }
910
911 if (!all_constant) {
912 (*residual_blocks)[j++] = (*residual_blocks)[i];
Markus Moll8de78db2012-07-14 15:49:11 +0200913 } else if (fixed_cost != NULL) {
914 // The residual is constant and will be removed, so its cost is
915 // added to the variable fixed_cost.
916 double cost = 0.0;
Sameer Agarwal039ff072013-02-26 09:15:39 -0800917 if (!residual_block->Evaluate(true,
918 &cost,
919 NULL,
920 NULL,
921 residual_block_evaluate_scratch.get())) {
Markus Moll8de78db2012-07-14 15:49:11 +0200922 *error = StringPrintf("Evaluation of the residual %d failed during "
923 "removal of fixed residual blocks.", i);
924 return false;
925 }
926 *fixed_cost += cost;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700927 }
928 }
929 residual_blocks->resize(j);
930 }
931
932 // Filter out unused or fixed parameter blocks, and update
Sameer Agarwal65625f72012-09-17 12:06:57 -0700933 // the ordering.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700934 {
935 vector<ParameterBlock*>* parameter_blocks =
936 program->mutable_parameter_blocks();
937 int j = 0;
938 for (int i = 0; i < parameter_blocks->size(); ++i) {
939 ParameterBlock* parameter_block = (*parameter_blocks)[i];
940 if (parameter_block->index() == 1) {
941 (*parameter_blocks)[j++] = parameter_block;
Sameer Agarwal65625f72012-09-17 12:06:57 -0700942 } else {
Sameer Agarwal2c94eed2012-10-01 16:34:37 -0700943 ordering->Remove(parameter_block->mutable_user_state());
Keir Mierle8ebb0732012-04-30 23:09:08 -0700944 }
945 }
946 parameter_blocks->resize(j);
947 }
948
Sameer Agarwalc8f07902013-04-19 08:01:04 -0700949 if (!(((program->NumResidualBlocks() == 0) &&
Keir Mierle8ebb0732012-04-30 23:09:08 -0700950 (program->NumParameterBlocks() == 0)) ||
951 ((program->NumResidualBlocks() != 0) &&
Sameer Agarwalc8f07902013-04-19 08:01:04 -0700952 (program->NumParameterBlocks() != 0)))) {
953 *error = "Congratulations, you found a bug in Ceres. Please report it.";
954 return false;
955 }
956
Keir Mierle8ebb0732012-04-30 23:09:08 -0700957 return true;
958}
959
960Program* SolverImpl::CreateReducedProgram(Solver::Options* options,
961 ProblemImpl* problem_impl,
Markus Moll8de78db2012-07-14 15:49:11 +0200962 double* fixed_cost,
Keir Mierle8ebb0732012-04-30 23:09:08 -0700963 string* error) {
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700964 CHECK_NOTNULL(options->linear_solver_ordering);
Keir Mierle8ebb0732012-04-30 23:09:08 -0700965 Program* original_program = problem_impl->mutable_program();
966 scoped_ptr<Program> transformed_program(new Program(*original_program));
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800967
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700968 ParameterBlockOrdering* linear_solver_ordering =
969 options->linear_solver_ordering;
Sameer Agarwal65625f72012-09-17 12:06:57 -0700970 const int min_group_id =
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700971 linear_solver_ordering->group_to_elements().begin()->first;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700972
973 if (!RemoveFixedBlocksFromProgram(transformed_program.get(),
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700974 linear_solver_ordering,
Markus Moll8de78db2012-07-14 15:49:11 +0200975 fixed_cost,
Keir Mierle8ebb0732012-04-30 23:09:08 -0700976 error)) {
977 return NULL;
978 }
979
980 if (transformed_program->NumParameterBlocks() == 0) {
981 LOG(WARNING) << "No varying parameter blocks to optimize; "
982 << "bailing early.";
983 return transformed_program.release();
984 }
985
Sameer Agarwal65625f72012-09-17 12:06:57 -0700986 if (IsSchurType(options->linear_solver_type) &&
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700987 linear_solver_ordering->GroupSize(min_group_id) == 0) {
Sameer Agarwalc8f07902013-04-19 08:01:04 -0700988 // If the user requested the use of a Schur type solver, and
989 // supplied a non-NULL linear_solver_ordering object with more than
990 // one elimination group, then it can happen that after all the
991 // parameter blocks which are fixed or unused have been removed from
992 // the program and the ordering, there are no more parameter blocks
993 // in the first elimination group.
994 //
995 // In such a case, the use of a Schur type solver is not possible,
996 // as they assume there is at least one e_block. Thus, we
997 // automatically switch to the closest solver to the one indicated
998 // by the user.
999 AlternateLinearSolverForSchurTypeLinearSolver(options);
Sameer Agarwal65625f72012-09-17 12:06:57 -07001000 }
Keir Mierle8ebb0732012-04-30 23:09:08 -07001001
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001002 if (IsSchurType(options->linear_solver_type)) {
Sameer Agarwalac626962013-05-06 07:04:26 -07001003 if (!ReorderProgramForSchurTypeLinearSolver(
1004 options->linear_solver_type,
1005 options->sparse_linear_algebra_library,
1006 problem_impl->parameter_map(),
1007 linear_solver_ordering,
1008 transformed_program.get(),
1009 error)) {
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001010 return NULL;
1011 }
1012 return transformed_program.release();
1013 }
Sameer Agarwal42a84b82013-02-01 12:22:53 -08001014
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001015 if (options->linear_solver_type == SPARSE_NORMAL_CHOLESKY) {
1016 if (!ReorderProgramForSparseNormalCholesky(
1017 options->sparse_linear_algebra_library,
1018 linear_solver_ordering,
1019 transformed_program.get(),
1020 error)) {
1021 return NULL;
1022 }
1023
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001024 return transformed_program.release();
1025 }
1026
Keir Mierle8ebb0732012-04-30 23:09:08 -07001027 transformed_program->SetParameterOffsetsAndIndex();
1028 return transformed_program.release();
1029}
1030
1031LinearSolver* SolverImpl::CreateLinearSolver(Solver::Options* options,
1032 string* error) {
Sameer Agarwal65625f72012-09-17 12:06:57 -07001033 CHECK_NOTNULL(options);
Sameer Agarwal68b32a92012-10-06 23:10:51 -07001034 CHECK_NOTNULL(options->linear_solver_ordering);
Sameer Agarwal65625f72012-09-17 12:06:57 -07001035 CHECK_NOTNULL(error);
1036
Sameer Agarwal5ecd2512012-06-17 16:34:03 -07001037 if (options->trust_region_strategy_type == DOGLEG) {
1038 if (options->linear_solver_type == ITERATIVE_SCHUR ||
1039 options->linear_solver_type == CGNR) {
1040 *error = "DOGLEG only supports exact factorization based linear "
1041 "solvers. If you want to use an iterative solver please "
1042 "use LEVENBERG_MARQUARDT as the trust_region_strategy_type";
1043 return NULL;
1044 }
1045 }
1046
Keir Mierle8ebb0732012-04-30 23:09:08 -07001047#ifdef CERES_NO_SUITESPARSE
Sameer Agarwalb0518732012-05-29 00:27:57 -07001048 if (options->linear_solver_type == SPARSE_NORMAL_CHOLESKY &&
1049 options->sparse_linear_algebra_library == SUITE_SPARSE) {
1050 *error = "Can't use SPARSE_NORMAL_CHOLESKY with SUITESPARSE because "
1051 "SuiteSparse was not enabled when Ceres was built.";
Keir Mierle8ebb0732012-04-30 23:09:08 -07001052 return NULL;
1053 }
Sameer Agarwal65625f72012-09-17 12:06:57 -07001054
Petter Strandmarkdd1a2762012-09-19 12:55:05 +02001055 if (options->preconditioner_type == CLUSTER_JACOBI) {
Sameer Agarwal65625f72012-09-17 12:06:57 -07001056 *error = "CLUSTER_JACOBI preconditioner not suppored. Please build Ceres "
1057 "with SuiteSparse support.";
1058 return NULL;
1059 }
1060
Petter Strandmarkdd1a2762012-09-19 12:55:05 +02001061 if (options->preconditioner_type == CLUSTER_TRIDIAGONAL) {
Sameer Agarwal65625f72012-09-17 12:06:57 -07001062 *error = "CLUSTER_TRIDIAGONAL preconditioner not suppored. Please build "
1063 "Ceres with SuiteSparse support.";
1064 return NULL;
1065 }
Sameer Agarwalb0518732012-05-29 00:27:57 -07001066#endif
1067
1068#ifdef CERES_NO_CXSPARSE
1069 if (options->linear_solver_type == SPARSE_NORMAL_CHOLESKY &&
1070 options->sparse_linear_algebra_library == CX_SPARSE) {
1071 *error = "Can't use SPARSE_NORMAL_CHOLESKY with CXSPARSE because "
1072 "CXSparse was not enabled when Ceres was built.";
1073 return NULL;
1074 }
1075#endif
1076
Sameer Agarwal65625f72012-09-17 12:06:57 -07001077#if defined(CERES_NO_SUITESPARSE) && defined(CERES_NO_CXSPARSE)
Keir Mierle27dd0d32012-10-15 13:52:36 -07001078 if (options->linear_solver_type == SPARSE_SCHUR) {
Sameer Agarwal65625f72012-09-17 12:06:57 -07001079 *error = "Can't use SPARSE_SCHUR because neither SuiteSparse nor"
1080 "CXSparse was enabled when Ceres was compiled.";
1081 return NULL;
1082 }
1083#endif
Keir Mierle8ebb0732012-04-30 23:09:08 -07001084
1085 if (options->linear_solver_max_num_iterations <= 0) {
1086 *error = "Solver::Options::linear_solver_max_num_iterations is 0.";
1087 return NULL;
1088 }
1089 if (options->linear_solver_min_num_iterations <= 0) {
1090 *error = "Solver::Options::linear_solver_min_num_iterations is 0.";
1091 return NULL;
1092 }
1093 if (options->linear_solver_min_num_iterations >
1094 options->linear_solver_max_num_iterations) {
1095 *error = "Solver::Options::linear_solver_min_num_iterations > "
1096 "Solver::Options::linear_solver_max_num_iterations.";
1097 return NULL;
1098 }
1099
1100 LinearSolver::Options linear_solver_options;
Keir Mierle8ebb0732012-04-30 23:09:08 -07001101 linear_solver_options.min_num_iterations =
1102 options->linear_solver_min_num_iterations;
1103 linear_solver_options.max_num_iterations =
1104 options->linear_solver_max_num_iterations;
1105 linear_solver_options.type = options->linear_solver_type;
1106 linear_solver_options.preconditioner_type = options->preconditioner_type;
Sameer Agarwalb0518732012-05-29 00:27:57 -07001107 linear_solver_options.sparse_linear_algebra_library =
1108 options->sparse_linear_algebra_library;
Sameer Agarwal9189f4e2013-04-19 17:09:49 -07001109 linear_solver_options.use_postordering = options->use_postordering;
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001110
1111 // Ignore user's postordering preferences and force it to be true if
1112 // cholmod_camd is not available. This ensures that the linear
1113 // solver does not assume that a fill-reducing pre-ordering has been
1114 // done.
1115#if !defined(CERES_NO_SUITESPARSE) && defined(CERES_NO_CAMD)
1116 if (IsSchurType(linear_solver_options.type) &&
1117 linear_solver_options.sparse_linear_algebra_library == SUITE_SPARSE) {
1118 linear_solver_options.use_postordering = true;
1119 }
1120#endif
1121
Keir Mierle8ebb0732012-04-30 23:09:08 -07001122 linear_solver_options.num_threads = options->num_linear_solver_threads;
Keir Mierle8ebb0732012-04-30 23:09:08 -07001123 options->num_linear_solver_threads = linear_solver_options.num_threads;
1124
Sameer Agarwal65625f72012-09-17 12:06:57 -07001125 const map<int, set<double*> >& groups =
Sameer Agarwal68b32a92012-10-06 23:10:51 -07001126 options->linear_solver_ordering->group_to_elements();
Sameer Agarwal65625f72012-09-17 12:06:57 -07001127 for (map<int, set<double*> >::const_iterator it = groups.begin();
1128 it != groups.end();
1129 ++it) {
1130 linear_solver_options.elimination_groups.push_back(it->second.size());
1131 }
1132 // Schur type solvers, expect at least two elimination groups. If
Sameer Agarwal1dc544a2012-09-18 18:09:03 -07001133 // there is only one elimination group, then CreateReducedProgram
Sameer Agarwal65625f72012-09-17 12:06:57 -07001134 // guarantees that this group only contains e_blocks. Thus we add a
1135 // dummy elimination group with zero blocks in it.
1136 if (IsSchurType(linear_solver_options.type) &&
1137 linear_solver_options.elimination_groups.size() == 1) {
1138 linear_solver_options.elimination_groups.push_back(0);
1139 }
1140
Keir Mierle8ebb0732012-04-30 23:09:08 -07001141 return LinearSolver::Create(linear_solver_options);
1142}
1143
Keir Mierle8ebb0732012-04-30 23:09:08 -07001144
1145// Find the minimum index of any parameter block to the given residual.
1146// Parameter blocks that have indices greater than num_eliminate_blocks are
1147// considered to have an index equal to num_eliminate_blocks.
Sergey Sharybinb53c9662013-02-25 01:14:26 +06001148static int MinParameterBlock(const ResidualBlock* residual_block,
1149 int num_eliminate_blocks) {
Keir Mierle8ebb0732012-04-30 23:09:08 -07001150 int min_parameter_block_position = num_eliminate_blocks;
1151 for (int i = 0; i < residual_block->NumParameterBlocks(); ++i) {
1152 ParameterBlock* parameter_block = residual_block->parameter_blocks()[i];
Keir Mierle32de18d2012-05-13 16:45:05 -07001153 if (!parameter_block->IsConstant()) {
1154 CHECK_NE(parameter_block->index(), -1)
1155 << "Did you forget to call Program::SetParameterOffsetsAndIndex()? "
1156 << "This is a Ceres bug; please contact the developers!";
1157 min_parameter_block_position = std::min(parameter_block->index(),
1158 min_parameter_block_position);
1159 }
Keir Mierle8ebb0732012-04-30 23:09:08 -07001160 }
1161 return min_parameter_block_position;
1162}
1163
1164// Reorder the residuals for program, if necessary, so that the residuals
1165// involving each E block occur together. This is a necessary condition for the
1166// Schur eliminator, which works on these "row blocks" in the jacobian.
Sameer Agarwal509f68c2013-02-20 01:39:03 -08001167bool SolverImpl::LexicographicallyOrderResidualBlocks(
1168 const int num_eliminate_blocks,
1169 Program* program,
1170 string* error) {
Sameer Agarwal9123e2f2012-09-18 21:49:06 -07001171 CHECK_GE(num_eliminate_blocks, 1)
1172 << "Congratulations, you found a Ceres bug! Please report this error "
1173 << "to the developers.";
Keir Mierle8ebb0732012-04-30 23:09:08 -07001174
1175 // Create a histogram of the number of residuals for each E block. There is an
1176 // extra bucket at the end to catch all non-eliminated F blocks.
Sameer Agarwal65625f72012-09-17 12:06:57 -07001177 vector<int> residual_blocks_per_e_block(num_eliminate_blocks + 1);
Keir Mierle8ebb0732012-04-30 23:09:08 -07001178 vector<ResidualBlock*>* residual_blocks = program->mutable_residual_blocks();
1179 vector<int> min_position_per_residual(residual_blocks->size());
1180 for (int i = 0; i < residual_blocks->size(); ++i) {
1181 ResidualBlock* residual_block = (*residual_blocks)[i];
Sameer Agarwal65625f72012-09-17 12:06:57 -07001182 int position = MinParameterBlock(residual_block, num_eliminate_blocks);
Keir Mierle8ebb0732012-04-30 23:09:08 -07001183 min_position_per_residual[i] = position;
Sameer Agarwal65625f72012-09-17 12:06:57 -07001184 DCHECK_LE(position, num_eliminate_blocks);
Keir Mierle8ebb0732012-04-30 23:09:08 -07001185 residual_blocks_per_e_block[position]++;
1186 }
1187
1188 // Run a cumulative sum on the histogram, to obtain offsets to the start of
1189 // each histogram bucket (where each bucket is for the residuals for that
1190 // E-block).
Sameer Agarwal65625f72012-09-17 12:06:57 -07001191 vector<int> offsets(num_eliminate_blocks + 1);
Keir Mierle8ebb0732012-04-30 23:09:08 -07001192 std::partial_sum(residual_blocks_per_e_block.begin(),
1193 residual_blocks_per_e_block.end(),
1194 offsets.begin());
1195 CHECK_EQ(offsets.back(), residual_blocks->size())
1196 << "Congratulations, you found a Ceres bug! Please report this error "
1197 << "to the developers.";
1198
1199 CHECK(find(residual_blocks_per_e_block.begin(),
1200 residual_blocks_per_e_block.end() - 1, 0) !=
1201 residual_blocks_per_e_block.end())
1202 << "Congratulations, you found a Ceres bug! Please report this error "
1203 << "to the developers.";
1204
1205 // Fill in each bucket with the residual blocks for its corresponding E block.
1206 // Each bucket is individually filled from the back of the bucket to the front
1207 // of the bucket. The filling order among the buckets is dictated by the
1208 // residual blocks. This loop uses the offsets as counters; subtracting one
1209 // from each offset as a residual block is placed in the bucket. When the
1210 // filling is finished, the offset pointerts should have shifted down one
1211 // entry (this is verified below).
1212 vector<ResidualBlock*> reordered_residual_blocks(
1213 (*residual_blocks).size(), static_cast<ResidualBlock*>(NULL));
1214 for (int i = 0; i < residual_blocks->size(); ++i) {
1215 int bucket = min_position_per_residual[i];
1216
1217 // Decrement the cursor, which should now point at the next empty position.
1218 offsets[bucket]--;
1219
1220 // Sanity.
1221 CHECK(reordered_residual_blocks[offsets[bucket]] == NULL)
1222 << "Congratulations, you found a Ceres bug! Please report this error "
1223 << "to the developers.";
1224
1225 reordered_residual_blocks[offsets[bucket]] = (*residual_blocks)[i];
1226 }
1227
1228 // Sanity check #1: The difference in bucket offsets should match the
1229 // histogram sizes.
Sameer Agarwal65625f72012-09-17 12:06:57 -07001230 for (int i = 0; i < num_eliminate_blocks; ++i) {
Keir Mierle8ebb0732012-04-30 23:09:08 -07001231 CHECK_EQ(residual_blocks_per_e_block[i], offsets[i + 1] - offsets[i])
1232 << "Congratulations, you found a Ceres bug! Please report this error "
1233 << "to the developers.";
1234 }
1235 // Sanity check #2: No NULL's left behind.
1236 for (int i = 0; i < reordered_residual_blocks.size(); ++i) {
1237 CHECK(reordered_residual_blocks[i] != NULL)
1238 << "Congratulations, you found a Ceres bug! Please report this error "
1239 << "to the developers.";
1240 }
1241
1242 // Now that the residuals are collected by E block, swap them in place.
1243 swap(*program->mutable_residual_blocks(), reordered_residual_blocks);
1244 return true;
1245}
1246
Sameer Agarwal509f68c2013-02-20 01:39:03 -08001247Evaluator* SolverImpl::CreateEvaluator(
1248 const Solver::Options& options,
1249 const ProblemImpl::ParameterMap& parameter_map,
1250 Program* program,
1251 string* error) {
Keir Mierle8ebb0732012-04-30 23:09:08 -07001252 Evaluator::Options evaluator_options;
1253 evaluator_options.linear_solver_type = options.linear_solver_type;
Sameer Agarwal65625f72012-09-17 12:06:57 -07001254 evaluator_options.num_eliminate_blocks =
Sameer Agarwal68b32a92012-10-06 23:10:51 -07001255 (options.linear_solver_ordering->NumGroups() > 0 &&
Sameer Agarwal65625f72012-09-17 12:06:57 -07001256 IsSchurType(options.linear_solver_type))
Sameer Agarwal68b32a92012-10-06 23:10:51 -07001257 ? (options.linear_solver_ordering
1258 ->group_to_elements().begin()
1259 ->second.size())
Sameer Agarwal9123e2f2012-09-18 21:49:06 -07001260 : 0;
Keir Mierle8ebb0732012-04-30 23:09:08 -07001261 evaluator_options.num_threads = options.num_threads;
1262 return Evaluator::Create(evaluator_options, program, error);
1263}
1264
Sameer Agarwal67a107b2012-10-08 14:35:41 -07001265CoordinateDescentMinimizer* SolverImpl::CreateInnerIterationMinimizer(
1266 const Solver::Options& options,
1267 const Program& program,
1268 const ProblemImpl::ParameterMap& parameter_map,
Sameer Agarwal977be7c2013-01-26 16:01:54 -08001269 Solver::Summary* summary) {
Sameer Agarwal9f9488b2013-05-23 09:40:40 -07001270 summary->inner_iterations_given = true;
1271
Sameer Agarwal67a107b2012-10-08 14:35:41 -07001272 scoped_ptr<CoordinateDescentMinimizer> inner_iteration_minimizer(
1273 new CoordinateDescentMinimizer);
1274 scoped_ptr<ParameterBlockOrdering> inner_iteration_ordering;
1275 ParameterBlockOrdering* ordering_ptr = NULL;
1276
1277 if (options.inner_iteration_ordering == NULL) {
1278 // Find a recursive decomposition of the Hessian matrix as a set
1279 // of independent sets of decreasing size and invert it. This
1280 // seems to work better in practice, i.e., Cameras before
1281 // points.
1282 inner_iteration_ordering.reset(new ParameterBlockOrdering);
1283 ComputeRecursiveIndependentSetOrdering(program,
1284 inner_iteration_ordering.get());
1285 inner_iteration_ordering->Reverse();
1286 ordering_ptr = inner_iteration_ordering.get();
1287 } else {
1288 const map<int, set<double*> >& group_to_elements =
1289 options.inner_iteration_ordering->group_to_elements();
1290
1291 // Iterate over each group and verify that it is an independent
1292 // set.
1293 map<int, set<double*> >::const_iterator it = group_to_elements.begin();
Sameer Agarwal509f68c2013-02-20 01:39:03 -08001294 for ( ; it != group_to_elements.end(); ++it) {
Sameer Agarwal67a107b2012-10-08 14:35:41 -07001295 if (!IsParameterBlockSetIndependent(it->second,
1296 program.residual_blocks())) {
Sameer Agarwal977be7c2013-01-26 16:01:54 -08001297 summary->error =
Sameer Agarwal67a107b2012-10-08 14:35:41 -07001298 StringPrintf("The user-provided "
1299 "parameter_blocks_for_inner_iterations does not "
1300 "form an independent set. Group Id: %d", it->first);
1301 return NULL;
1302 }
1303 }
1304 ordering_ptr = options.inner_iteration_ordering;
1305 }
1306
1307 if (!inner_iteration_minimizer->Init(program,
1308 parameter_map,
1309 *ordering_ptr,
Sameer Agarwal977be7c2013-01-26 16:01:54 -08001310 &summary->error)) {
Sameer Agarwal67a107b2012-10-08 14:35:41 -07001311 return NULL;
1312 }
1313
Sameer Agarwal9f9488b2013-05-23 09:40:40 -07001314 summary->inner_iterations_used = true;
1315 summary->inner_iteration_time_in_seconds = 0.0;
Sameer Agarwal977be7c2013-01-26 16:01:54 -08001316 SummarizeOrdering(ordering_ptr, &(summary->inner_iteration_ordering_used));
Sameer Agarwal67a107b2012-10-08 14:35:41 -07001317 return inner_iteration_minimizer.release();
1318}
1319
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001320void SolverImpl::AlternateLinearSolverForSchurTypeLinearSolver(
1321 Solver::Options* options) {
1322 if (!IsSchurType(options->linear_solver_type)) {
1323 return;
1324 }
1325
1326 string msg = "No e_blocks remaining. Switching from ";
1327 if (options->linear_solver_type == SPARSE_SCHUR) {
1328 options->linear_solver_type = SPARSE_NORMAL_CHOLESKY;
1329 msg += "SPARSE_SCHUR to SPARSE_NORMAL_CHOLESKY.";
1330 } else if (options->linear_solver_type == DENSE_SCHUR) {
1331 // TODO(sameeragarwal): This is probably not a great choice.
1332 // Ideally, we should have a DENSE_NORMAL_CHOLESKY, that can
1333 // take a BlockSparseMatrix as input.
1334 options->linear_solver_type = DENSE_QR;
1335 msg += "DENSE_SCHUR to DENSE_QR.";
1336 } else if (options->linear_solver_type == ITERATIVE_SCHUR) {
1337 options->linear_solver_type = CGNR;
1338 if (options->preconditioner_type != IDENTITY) {
1339 msg += StringPrintf("ITERATIVE_SCHUR with %s preconditioner "
1340 "to CGNR with JACOBI preconditioner.",
1341 PreconditionerTypeToString(
1342 options->preconditioner_type));
1343 // CGNR currently only supports the JACOBI preconditioner.
1344 options->preconditioner_type = JACOBI;
1345 } else {
Sameer Agarwalcbdeb792013-04-22 10:18:18 -07001346 msg += "ITERATIVE_SCHUR with IDENTITY preconditioner"
1347 "to CGNR with IDENTITY preconditioner.";
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001348 }
1349 }
1350 LOG(WARNING) << msg;
1351}
1352
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001353bool SolverImpl::ApplyUserOrdering(
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001354 const ProblemImpl::ParameterMap& parameter_map,
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001355 const ParameterBlockOrdering* parameter_block_ordering,
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001356 Program* program,
1357 string* error) {
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001358 const int num_parameter_blocks = program->NumParameterBlocks();
1359 if (parameter_block_ordering->NumElements() != num_parameter_blocks) {
1360 *error = StringPrintf("User specified ordering does not have the same "
1361 "number of parameters as the problem. The problem"
1362 "has %d blocks while the ordering has %d blocks.",
1363 num_parameter_blocks,
1364 parameter_block_ordering->NumElements());
1365 return false;
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001366 }
1367
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001368 vector<ParameterBlock*>* parameter_blocks =
1369 program->mutable_parameter_blocks();
1370 parameter_blocks->clear();
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001371
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001372 const map<int, set<double*> >& groups =
1373 parameter_block_ordering->group_to_elements();
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001374
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001375 for (map<int, set<double*> >::const_iterator group_it = groups.begin();
1376 group_it != groups.end();
1377 ++group_it) {
1378 const set<double*>& group = group_it->second;
1379 for (set<double*>::const_iterator parameter_block_ptr_it = group.begin();
1380 parameter_block_ptr_it != group.end();
1381 ++parameter_block_ptr_it) {
1382 ProblemImpl::ParameterMap::const_iterator parameter_block_it =
1383 parameter_map.find(*parameter_block_ptr_it);
1384 if (parameter_block_it == parameter_map.end()) {
1385 *error = StringPrintf("User specified ordering contains a pointer "
1386 "to a double that is not a parameter block in "
1387 "the problem. The invalid double is in group: %d",
1388 group_it->first);
1389 return false;
1390 }
1391 parameter_blocks->push_back(parameter_block_it->second);
1392 }
1393 }
1394 return true;
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001395}
1396
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001397
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001398TripletSparseMatrix* SolverImpl::CreateJacobianBlockSparsityTranspose(
1399 const Program* program) {
1400
Sameer Agarwalcbdeb792013-04-22 10:18:18 -07001401 // Matrix to store the block sparsity structure of the Jacobian.
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001402 TripletSparseMatrix* tsm =
1403 new TripletSparseMatrix(program->NumParameterBlocks(),
1404 program->NumResidualBlocks(),
1405 10 * program->NumResidualBlocks());
1406 int num_nonzeros = 0;
1407 int* rows = tsm->mutable_rows();
1408 int* cols = tsm->mutable_cols();
1409 double* values = tsm->mutable_values();
1410
1411 const vector<ResidualBlock*>& residual_blocks = program->residual_blocks();
1412 for (int c = 0; c < residual_blocks.size(); ++c) {
1413 const ResidualBlock* residual_block = residual_blocks[c];
1414 const int num_parameter_blocks = residual_block->NumParameterBlocks();
1415 ParameterBlock* const* parameter_blocks =
1416 residual_block->parameter_blocks();
1417
1418 for (int j = 0; j < num_parameter_blocks; ++j) {
1419 if (parameter_blocks[j]->IsConstant()) {
1420 continue;
1421 }
1422
1423 // Re-size the matrix if needed.
1424 if (num_nonzeros >= tsm->max_num_nonzeros()) {
Sameer Agarwal5f433c82013-06-13 06:57:58 -07001425 tsm->set_num_nonzeros(num_nonzeros);
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001426 tsm->Reserve(2 * num_nonzeros);
1427 rows = tsm->mutable_rows();
1428 cols = tsm->mutable_cols();
1429 values = tsm->mutable_values();
1430 }
1431 CHECK_LT(num_nonzeros, tsm->max_num_nonzeros());
1432
1433 const int r = parameter_blocks[j]->index();
1434 rows[num_nonzeros] = r;
1435 cols[num_nonzeros] = c;
1436 values[num_nonzeros] = 1.0;
1437 ++num_nonzeros;
1438 }
1439 }
1440
1441 tsm->set_num_nonzeros(num_nonzeros);
1442 return tsm;
1443}
1444
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001445bool SolverImpl::ReorderProgramForSchurTypeLinearSolver(
1446 const LinearSolverType linear_solver_type,
1447 const SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type,
1448 const ProblemImpl::ParameterMap& parameter_map,
1449 ParameterBlockOrdering* parameter_block_ordering,
1450 Program* program,
1451 string* error) {
1452 if (parameter_block_ordering->NumGroups() == 1) {
1453 // If the user supplied an parameter_block_ordering with just one
1454 // group, it is equivalent to the user supplying NULL as an
1455 // parameter_block_ordering. Ceres is completely free to choose the
1456 // parameter block ordering as it sees fit. For Schur type solvers,
1457 // this means that the user wishes for Ceres to identify the
1458 // e_blocks, which we do by computing a maximal independent set.
1459 vector<ParameterBlock*> schur_ordering;
Sameer Agarwala427c872013-06-24 17:50:56 -07001460 const int num_eliminate_blocks =
1461 ComputeStableSchurOrdering(*program, &schur_ordering);
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001462
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001463 CHECK_EQ(schur_ordering.size(), program->NumParameterBlocks())
1464 << "Congratulations, you found a Ceres bug! Please report this error "
1465 << "to the developers.";
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001466
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001467 // Update the parameter_block_ordering object.
1468 for (int i = 0; i < schur_ordering.size(); ++i) {
1469 double* parameter_block = schur_ordering[i]->mutable_user_state();
1470 const int group_id = (i < num_eliminate_blocks) ? 0 : 1;
1471 parameter_block_ordering->AddElementToGroup(parameter_block, group_id);
1472 }
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001473
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001474 // We could call ApplyUserOrdering but this is cheaper and
1475 // simpler.
1476 swap(*program->mutable_parameter_blocks(), schur_ordering);
1477 } else {
1478 // The user provided an ordering with more than one elimination
1479 // group. Trust the user and apply the ordering.
1480 if (!ApplyUserOrdering(parameter_map,
1481 parameter_block_ordering,
1482 program,
1483 error)) {
1484 return false;
1485 }
1486 }
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001487
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001488 // Pre-order the columns corresponding to the schur complement if
1489 // possible.
1490#if !defined(CERES_NO_SUITESPARSE) && !defined(CERES_NO_CAMD)
1491 if (linear_solver_type == SPARSE_SCHUR &&
1492 sparse_linear_algebra_library_type == SUITE_SPARSE) {
1493 vector<int> constraints;
1494 vector<ParameterBlock*>& parameter_blocks =
1495 *(program->mutable_parameter_blocks());
1496
1497 for (int i = 0; i < parameter_blocks.size(); ++i) {
1498 constraints.push_back(
1499 parameter_block_ordering->GroupId(
1500 parameter_blocks[i]->mutable_user_state()));
1501 }
1502
1503 // Set the offsets and index for CreateJacobianSparsityTranspose.
1504 program->SetParameterOffsetsAndIndex();
1505 // Compute a block sparse presentation of J'.
1506 scoped_ptr<TripletSparseMatrix> tsm_block_jacobian_transpose(
1507 SolverImpl::CreateJacobianBlockSparsityTranspose(program));
1508
1509 SuiteSparse ss;
1510 cholmod_sparse* block_jacobian_transpose =
1511 ss.CreateSparseMatrix(tsm_block_jacobian_transpose.get());
1512
1513 vector<int> ordering(parameter_blocks.size(), 0);
1514 ss.ConstrainedApproximateMinimumDegreeOrdering(block_jacobian_transpose,
1515 &constraints[0],
1516 &ordering[0]);
1517 ss.Free(block_jacobian_transpose);
1518
1519 const vector<ParameterBlock*> parameter_blocks_copy(parameter_blocks);
1520 for (int i = 0; i < program->NumParameterBlocks(); ++i) {
1521 parameter_blocks[i] = parameter_blocks_copy[ordering[i]];
1522 }
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001523 }
1524#endif
1525
1526 program->SetParameterOffsetsAndIndex();
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001527 // Schur type solvers also require that their residual blocks be
1528 // lexicographically ordered.
1529 const int num_eliminate_blocks =
1530 parameter_block_ordering->group_to_elements().begin()->second.size();
1531 return LexicographicallyOrderResidualBlocks(num_eliminate_blocks,
1532 program,
1533 error);
1534}
1535
1536bool SolverImpl::ReorderProgramForSparseNormalCholesky(
1537 const SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type,
1538 const ParameterBlockOrdering* parameter_block_ordering,
1539 Program* program,
1540 string* error) {
1541 // Set the offsets and index for CreateJacobianSparsityTranspose.
1542 program->SetParameterOffsetsAndIndex();
1543 // Compute a block sparse presentation of J'.
1544 scoped_ptr<TripletSparseMatrix> tsm_block_jacobian_transpose(
1545 SolverImpl::CreateJacobianBlockSparsityTranspose(program));
1546
1547 vector<int> ordering(program->NumParameterBlocks(), 0);
1548 vector<ParameterBlock*>& parameter_blocks =
1549 *(program->mutable_parameter_blocks());
1550
1551 if (sparse_linear_algebra_library_type == SUITE_SPARSE) {
1552#ifdef CERES_NO_SUITESPARSE
1553 *error = "Can't use SPARSE_NORMAL_CHOLESKY with SUITE_SPARSE because "
1554 "SuiteSparse was not enabled when Ceres was built.";
1555 return false;
1556#else
1557 SuiteSparse ss;
1558 cholmod_sparse* block_jacobian_transpose =
1559 ss.CreateSparseMatrix(tsm_block_jacobian_transpose.get());
1560
1561# ifdef CERES_NO_CAMD
1562 // No cholmod_camd, so ignore user's parameter_block_ordering and
1563 // use plain old AMD.
1564 ss.ApproximateMinimumDegreeOrdering(block_jacobian_transpose, &ordering[0]);
1565# else
1566 if (parameter_block_ordering->NumGroups() > 1) {
1567 // If the user specified more than one elimination groups use them
1568 // to constrain the ordering.
1569 vector<int> constraints;
1570 for (int i = 0; i < parameter_blocks.size(); ++i) {
1571 constraints.push_back(
1572 parameter_block_ordering->GroupId(
1573 parameter_blocks[i]->mutable_user_state()));
1574 }
1575 ss.ConstrainedApproximateMinimumDegreeOrdering(
1576 block_jacobian_transpose,
1577 &constraints[0],
1578 &ordering[0]);
1579 } else {
1580 ss.ApproximateMinimumDegreeOrdering(block_jacobian_transpose,
1581 &ordering[0]);
1582 }
1583# endif // CERES_NO_CAMD
1584
1585 ss.Free(block_jacobian_transpose);
1586#endif // CERES_NO_SUITESPARSE
1587
1588 } else if (sparse_linear_algebra_library_type == CX_SPARSE) {
1589#ifndef CERES_NO_CXSPARSE
1590
1591 // CXSparse works with J'J instead of J'. So compute the block
1592 // sparsity for J'J and compute an approximate minimum degree
1593 // ordering.
1594 CXSparse cxsparse;
1595 cs_di* block_jacobian_transpose;
1596 block_jacobian_transpose =
1597 cxsparse.CreateSparseMatrix(tsm_block_jacobian_transpose.get());
1598 cs_di* block_jacobian = cxsparse.TransposeMatrix(block_jacobian_transpose);
1599 cs_di* block_hessian =
1600 cxsparse.MatrixMatrixMultiply(block_jacobian_transpose, block_jacobian);
1601 cxsparse.Free(block_jacobian);
1602 cxsparse.Free(block_jacobian_transpose);
1603
1604 cxsparse.ApproximateMinimumDegreeOrdering(block_hessian, &ordering[0]);
1605 cxsparse.Free(block_hessian);
1606#else // CERES_NO_CXSPARSE
1607 *error = "Can't use SPARSE_NORMAL_CHOLESKY with CX_SPARSE because "
1608 "CXSparse was not enabled when Ceres was built.";
1609 return false;
1610#endif // CERES_NO_CXSPARSE
1611 } else {
1612 *error = "Unknown sparse linear algebra library.";
1613 return false;
1614 }
1615
1616 // Apply ordering.
1617 const vector<ParameterBlock*> parameter_blocks_copy(parameter_blocks);
1618 for (int i = 0; i < program->NumParameterBlocks(); ++i) {
1619 parameter_blocks[i] = parameter_blocks_copy[ordering[i]];
1620 }
1621
1622 program->SetParameterOffsetsAndIndex();
1623 return true;
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001624}
1625
Keir Mierle8ebb0732012-04-30 23:09:08 -07001626} // namespace internal
1627} // namespace ceres