blob: 56993c874240b2da6ee4fbf1e664aa02a98f395e [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
394 if (options.lsqp_iterations_to_dump.size() > 0) {
395 LOG(WARNING) << "Dumping linear least squares problems to disk is"
396 " currently broken. Ignoring Solver::Options::lsqp_iterations_to_dump";
397 }
398
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800399 event_logger.AddEvent("Init");
Sameer Agarwalf102a682013-02-11 15:08:40 -0800400
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700401 original_program->SetParameterBlockStatePtrsToUserStatePtrs();
Sameer Agarwalc58e6dc2013-02-11 16:41:06 -0800402 event_logger.AddEvent("SetParameterBlockPtrs");
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700403
404 // If the user requests gradient checking, construct a new
405 // ProblemImpl by wrapping the CostFunctions of problem_impl inside
406 // GradientCheckingCostFunction and replacing problem_impl with
407 // gradient_checking_problem_impl.
408 scoped_ptr<ProblemImpl> gradient_checking_problem_impl;
409 if (options.check_gradients) {
410 VLOG(1) << "Checking Gradients";
411 gradient_checking_problem_impl.reset(
412 CreateGradientCheckingProblemImpl(
413 problem_impl,
414 options.numeric_derivative_relative_step_size,
415 options.gradient_check_relative_precision));
416
417 // From here on, problem_impl will point to the gradient checking
418 // version.
419 problem_impl = gradient_checking_problem_impl.get();
420 }
Keir Mierle8ebb0732012-04-30 23:09:08 -0700421
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700422 if (original_options.linear_solver_ordering != NULL) {
Sameer Agarwal65625f72012-09-17 12:06:57 -0700423 if (!IsOrderingValid(original_options, problem_impl, &summary->error)) {
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700424 LOG(ERROR) << summary->error;
Sameer Agarwal65625f72012-09-17 12:06:57 -0700425 return;
426 }
Sameer Agarwalc58e6dc2013-02-11 16:41:06 -0800427 event_logger.AddEvent("CheckOrdering");
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700428 options.linear_solver_ordering =
429 new ParameterBlockOrdering(*original_options.linear_solver_ordering);
Sameer Agarwalc58e6dc2013-02-11 16:41:06 -0800430 event_logger.AddEvent("CopyOrdering");
Sameer Agarwal65625f72012-09-17 12:06:57 -0700431 } else {
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700432 options.linear_solver_ordering = new ParameterBlockOrdering;
Sameer Agarwal65625f72012-09-17 12:06:57 -0700433 const ProblemImpl::ParameterMap& parameter_map =
434 problem_impl->parameter_map();
435 for (ProblemImpl::ParameterMap::const_iterator it = parameter_map.begin();
436 it != parameter_map.end();
437 ++it) {
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700438 options.linear_solver_ordering->AddElementToGroup(it->first, 0);
Sameer Agarwal65625f72012-09-17 12:06:57 -0700439 }
Sameer Agarwalc58e6dc2013-02-11 16:41:06 -0800440 event_logger.AddEvent("ConstructOrdering");
Sameer Agarwal65625f72012-09-17 12:06:57 -0700441 }
442
Keir Mierle8ebb0732012-04-30 23:09:08 -0700443 // Create the three objects needed to minimize: the transformed program, the
444 // evaluator, and the linear solver.
Keir Mierle31c1e782012-08-17 16:16:32 -0700445 scoped_ptr<Program> reduced_program(CreateReducedProgram(&options,
446 problem_impl,
447 &summary->fixed_cost,
448 &summary->error));
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800449
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800450 event_logger.AddEvent("CreateReducedProgram");
Keir Mierle8ebb0732012-04-30 23:09:08 -0700451 if (reduced_program == NULL) {
452 return;
453 }
454
Sameer Agarwal977be7c2013-01-26 16:01:54 -0800455 SummarizeOrdering(options.linear_solver_ordering,
456 &(summary->linear_solver_ordering_used));
457
Keir Mierle8ebb0732012-04-30 23:09:08 -0700458 summary->num_parameter_blocks_reduced = reduced_program->NumParameterBlocks();
459 summary->num_parameters_reduced = reduced_program->NumParameters();
Keir Mierleba944212013-02-25 12:46:44 -0800460 summary->num_effective_parameters_reduced =
461 reduced_program->NumEffectiveParameters();
Keir Mierle8ebb0732012-04-30 23:09:08 -0700462 summary->num_residual_blocks_reduced = reduced_program->NumResidualBlocks();
463 summary->num_residuals_reduced = reduced_program->NumResiduals();
464
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700465 if (summary->num_parameter_blocks_reduced == 0) {
466 summary->preprocessor_time_in_seconds =
467 WallTimeInSeconds() - solver_start_time;
468
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800469 double post_process_start_time = WallTimeInSeconds();
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700470 LOG(INFO) << "Terminating: FUNCTION_TOLERANCE reached. "
471 << "No non-constant parameter blocks found.";
472
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800473 summary->initial_cost = summary->fixed_cost;
474 summary->final_cost = summary->fixed_cost;
475
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700476 // FUNCTION_TOLERANCE is the right convergence here, as we know
477 // that the objective function is constant and cannot be changed
478 // any further.
479 summary->termination_type = FUNCTION_TOLERANCE;
480
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700481 // Ensure the program state is set to the user parameters on the way out.
482 original_program->SetParameterBlockStatePtrsToUserStatePtrs();
483
484 summary->postprocessor_time_in_seconds =
485 WallTimeInSeconds() - post_process_start_time;
486 return;
487 }
488
Keir Mierle8ebb0732012-04-30 23:09:08 -0700489 scoped_ptr<LinearSolver>
490 linear_solver(CreateLinearSolver(&options, &summary->error));
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800491 event_logger.AddEvent("CreateLinearSolver");
Keir Mierle8ebb0732012-04-30 23:09:08 -0700492 if (linear_solver == NULL) {
493 return;
494 }
495
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700496 summary->linear_solver_type_given = original_options.linear_solver_type;
497 summary->linear_solver_type_used = options.linear_solver_type;
498
499 summary->preconditioner_type = options.preconditioner_type;
500
501 summary->num_linear_solver_threads_given =
502 original_options.num_linear_solver_threads;
503 summary->num_linear_solver_threads_used = options.num_linear_solver_threads;
504
505 summary->sparse_linear_algebra_library =
506 options.sparse_linear_algebra_library;
507
508 summary->trust_region_strategy_type = options.trust_region_strategy_type;
509 summary->dogleg_type = options.dogleg_type;
510
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700511 scoped_ptr<Evaluator> evaluator(CreateEvaluator(options,
512 problem_impl->parameter_map(),
513 reduced_program.get(),
514 &summary->error));
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800515
516 event_logger.AddEvent("CreateEvaluator");
517
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700518 if (evaluator == NULL) {
Keir Mierle8ebb0732012-04-30 23:09:08 -0700519 return;
520 }
521
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700522 scoped_ptr<CoordinateDescentMinimizer> inner_iteration_minimizer;
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700523 if (options.use_inner_iterations) {
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700524 if (reduced_program->parameter_blocks().size() < 2) {
525 LOG(WARNING) << "Reduced problem only contains one parameter block."
526 << "Disabling inner iterations.";
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700527 } else {
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700528 inner_iteration_minimizer.reset(
529 CreateInnerIterationMinimizer(original_options,
530 *reduced_program,
531 problem_impl->parameter_map(),
Sameer Agarwal977be7c2013-01-26 16:01:54 -0800532 summary));
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700533 if (inner_iteration_minimizer == NULL) {
534 LOG(ERROR) << summary->error;
535 return;
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700536 }
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700537 }
Keir Mierle8ebb0732012-04-30 23:09:08 -0700538 }
539
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800540 event_logger.AddEvent("CreateIIM");
541
Keir Mierle8ebb0732012-04-30 23:09:08 -0700542 // The optimizer works on contiguous parameter vectors; allocate some.
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700543 Vector parameters(reduced_program->NumParameters());
Keir Mierle8ebb0732012-04-30 23:09:08 -0700544
545 // Collect the discontiguous parameters into a contiguous state vector.
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700546 reduced_program->ParameterBlocksToStateVector(parameters.data());
Keir Mierle8ebb0732012-04-30 23:09:08 -0700547
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700548 Vector original_parameters = parameters;
549
Petter Strandmark76533b32012-09-04 22:08:50 -0700550 double minimizer_start_time = WallTimeInSeconds();
Sameer Agarwalfa015192012-06-11 14:21:42 -0700551 summary->preprocessor_time_in_seconds =
552 minimizer_start_time - solver_start_time;
553
Keir Mierle8ebb0732012-04-30 23:09:08 -0700554 // Run the optimization.
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800555 TrustRegionMinimize(options,
556 reduced_program.get(),
557 inner_iteration_minimizer.get(),
558 evaluator.get(),
559 linear_solver.get(),
560 parameters.data(),
561 summary);
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800562 event_logger.AddEvent("Minimize");
Keir Mierle8ebb0732012-04-30 23:09:08 -0700563
Sameer Agarwalf102a682013-02-11 15:08:40 -0800564 SetSummaryFinalCost(summary);
565
Keir Mierle8ebb0732012-04-30 23:09:08 -0700566 // If the user aborted mid-optimization or the optimization
567 // terminated because of a numerical failure, then return without
568 // updating user state.
569 if (summary->termination_type == USER_ABORT ||
570 summary->termination_type == NUMERICAL_FAILURE) {
571 return;
572 }
573
Petter Strandmark76533b32012-09-04 22:08:50 -0700574 double post_process_start_time = WallTimeInSeconds();
Sameer Agarwal4997cbc2012-07-02 12:44:34 -0700575
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800576 // Push the contiguous optimized parameters back to the user's
577 // parameters.
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700578 reduced_program->StateVectorToParameterBlocks(parameters.data());
Keir Mierle8ebb0732012-04-30 23:09:08 -0700579 reduced_program->CopyParameterBlockStateToUserState();
580
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800581 // Ensure the program state is set to the user parameters on the way
582 // out.
Sameer Agarwal4997cbc2012-07-02 12:44:34 -0700583 original_program->SetParameterBlockStatePtrsToUserStatePtrs();
Sameer Agarwal67a107b2012-10-08 14:35:41 -0700584
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800585 const map<string, double>& linear_solver_time_statistics =
586 linear_solver->TimeStatistics();
587 summary->linear_solver_time_in_seconds =
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800588 FindWithDefault(linear_solver_time_statistics,
589 "LinearSolver::Solve",
590 0.0);
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800591
592 const map<string, double>& evaluator_time_statistics =
593 evaluator->TimeStatistics();
594
595 summary->residual_evaluation_time_in_seconds =
596 FindWithDefault(evaluator_time_statistics, "Evaluator::Residual", 0.0);
597 summary->jacobian_evaluation_time_in_seconds =
598 FindWithDefault(evaluator_time_statistics, "Evaluator::Jacobian", 0.0);
599
Keir Mierle8ebb0732012-04-30 23:09:08 -0700600 // Stick a fork in it, we're done.
Petter Strandmark76533b32012-09-04 22:08:50 -0700601 summary->postprocessor_time_in_seconds =
602 WallTimeInSeconds() - post_process_start_time;
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800603 event_logger.AddEvent("PostProcess");
Keir Mierle8ebb0732012-04-30 23:09:08 -0700604}
605
Sameer Agarwal8140f0f2013-03-12 09:45:08 -0700606
607#ifndef CERES_NO_LINE_SEARCH_MINIMIZER
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800608void SolverImpl::LineSearchSolve(const Solver::Options& original_options,
609 ProblemImpl* original_problem_impl,
610 Solver::Summary* summary) {
611 double solver_start_time = WallTimeInSeconds();
612
613 Program* original_program = original_problem_impl->mutable_program();
614 ProblemImpl* problem_impl = original_problem_impl;
615
616 // Reset the summary object to its default values.
617 *CHECK_NOTNULL(summary) = Solver::Summary();
618
Sameer Agarwald010de52013-02-15 14:26:56 -0800619 summary->minimizer_type = LINE_SEARCH;
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800620 summary->line_search_direction_type =
621 original_options.line_search_direction_type;
Sameer Agarwald010de52013-02-15 14:26:56 -0800622 summary->max_lbfgs_rank = original_options.max_lbfgs_rank;
623 summary->line_search_type = original_options.line_search_type;
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800624 summary->num_parameter_blocks = problem_impl->NumParameterBlocks();
625 summary->num_parameters = problem_impl->NumParameters();
626 summary->num_residual_blocks = problem_impl->NumResidualBlocks();
627 summary->num_residuals = problem_impl->NumResiduals();
628
629 // Empty programs are usually a user error.
630 if (summary->num_parameter_blocks == 0) {
631 summary->error = "Problem contains no parameter blocks.";
632 LOG(ERROR) << summary->error;
633 return;
634 }
635
636 if (summary->num_residual_blocks == 0) {
637 summary->error = "Problem contains no residual blocks.";
638 LOG(ERROR) << summary->error;
639 return;
640 }
641
642 Solver::Options options(original_options);
643
644 // This ensures that we get a Block Jacobian Evaluator along with
645 // none of the Schur nonsense. This file will have to be extensively
646 // refactored to deal with the various bits of cleanups related to
647 // line search.
648 options.linear_solver_type = CGNR;
649
650 options.linear_solver_ordering = NULL;
651 options.inner_iteration_ordering = NULL;
652
653#ifndef CERES_USE_OPENMP
654 if (options.num_threads > 1) {
655 LOG(WARNING)
656 << "OpenMP support is not compiled into this binary; "
657 << "only options.num_threads=1 is supported. Switching "
658 << "to single threaded mode.";
659 options.num_threads = 1;
660 }
Sameer Agarwal8140f0f2013-03-12 09:45:08 -0700661#endif // CERES_USE_OPENMP
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800662
663 summary->num_threads_given = original_options.num_threads;
664 summary->num_threads_used = options.num_threads;
665
666 if (original_options.linear_solver_ordering != NULL) {
667 if (!IsOrderingValid(original_options, problem_impl, &summary->error)) {
668 LOG(ERROR) << summary->error;
669 return;
670 }
671 options.linear_solver_ordering =
672 new ParameterBlockOrdering(*original_options.linear_solver_ordering);
673 } else {
674 options.linear_solver_ordering = new ParameterBlockOrdering;
675 const ProblemImpl::ParameterMap& parameter_map =
676 problem_impl->parameter_map();
677 for (ProblemImpl::ParameterMap::const_iterator it = parameter_map.begin();
678 it != parameter_map.end();
679 ++it) {
680 options.linear_solver_ordering->AddElementToGroup(it->first, 0);
681 }
682 }
683
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800684 original_program->SetParameterBlockStatePtrsToUserStatePtrs();
685
686 // If the user requests gradient checking, construct a new
687 // ProblemImpl by wrapping the CostFunctions of problem_impl inside
688 // GradientCheckingCostFunction and replacing problem_impl with
689 // gradient_checking_problem_impl.
690 scoped_ptr<ProblemImpl> gradient_checking_problem_impl;
691 if (options.check_gradients) {
692 VLOG(1) << "Checking Gradients";
693 gradient_checking_problem_impl.reset(
694 CreateGradientCheckingProblemImpl(
695 problem_impl,
696 options.numeric_derivative_relative_step_size,
697 options.gradient_check_relative_precision));
698
699 // From here on, problem_impl will point to the gradient checking
700 // version.
701 problem_impl = gradient_checking_problem_impl.get();
702 }
703
704 // Create the three objects needed to minimize: the transformed program, the
705 // evaluator, and the linear solver.
706 scoped_ptr<Program> reduced_program(CreateReducedProgram(&options,
707 problem_impl,
708 &summary->fixed_cost,
709 &summary->error));
710 if (reduced_program == NULL) {
711 return;
712 }
713
714 summary->num_parameter_blocks_reduced = reduced_program->NumParameterBlocks();
715 summary->num_parameters_reduced = reduced_program->NumParameters();
716 summary->num_residual_blocks_reduced = reduced_program->NumResidualBlocks();
717 summary->num_residuals_reduced = reduced_program->NumResiduals();
718
719 if (summary->num_parameter_blocks_reduced == 0) {
720 summary->preprocessor_time_in_seconds =
721 WallTimeInSeconds() - solver_start_time;
722
723 LOG(INFO) << "Terminating: FUNCTION_TOLERANCE reached. "
724 << "No non-constant parameter blocks found.";
725
726 // FUNCTION_TOLERANCE is the right convergence here, as we know
727 // that the objective function is constant and cannot be changed
728 // any further.
729 summary->termination_type = FUNCTION_TOLERANCE;
730
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800731 const double post_process_start_time = WallTimeInSeconds();
Sameer Agarwalf102a682013-02-11 15:08:40 -0800732
733 SetSummaryFinalCost(summary);
734
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800735 // Ensure the program state is set to the user parameters on the way out.
736 original_program->SetParameterBlockStatePtrsToUserStatePtrs();
Sameer Agarwal956ed7e2013-02-19 07:06:15 -0800737 summary->postprocessor_time_in_seconds =
738 WallTimeInSeconds() - post_process_start_time;
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800739 return;
740 }
741
742 scoped_ptr<Evaluator> evaluator(CreateEvaluator(options,
743 problem_impl->parameter_map(),
744 reduced_program.get(),
745 &summary->error));
746 if (evaluator == NULL) {
747 return;
748 }
749
750 // The optimizer works on contiguous parameter vectors; allocate some.
751 Vector parameters(reduced_program->NumParameters());
752
753 // Collect the discontiguous parameters into a contiguous state vector.
754 reduced_program->ParameterBlocksToStateVector(parameters.data());
755
756 Vector original_parameters = parameters;
757
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800758 const double minimizer_start_time = WallTimeInSeconds();
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800759 summary->preprocessor_time_in_seconds =
760 minimizer_start_time - solver_start_time;
761
762 // Run the optimization.
763 LineSearchMinimize(options,
764 reduced_program.get(),
765 evaluator.get(),
766 parameters.data(),
767 summary);
768
769 // If the user aborted mid-optimization or the optimization
770 // terminated because of a numerical failure, then return without
771 // updating user state.
772 if (summary->termination_type == USER_ABORT ||
773 summary->termination_type == NUMERICAL_FAILURE) {
774 return;
775 }
776
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800777 const double post_process_start_time = WallTimeInSeconds();
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800778
779 // Push the contiguous optimized parameters back to the user's parameters.
780 reduced_program->StateVectorToParameterBlocks(parameters.data());
781 reduced_program->CopyParameterBlockStateToUserState();
782
Sameer Agarwalf102a682013-02-11 15:08:40 -0800783 SetSummaryFinalCost(summary);
784
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800785 // Ensure the program state is set to the user parameters on the way out.
786 original_program->SetParameterBlockStatePtrsToUserStatePtrs();
787
Sameer Agarwald010de52013-02-15 14:26:56 -0800788 const map<string, double>& evaluator_time_statistics =
789 evaluator->TimeStatistics();
790
791 summary->residual_evaluation_time_in_seconds =
792 FindWithDefault(evaluator_time_statistics, "Evaluator::Residual", 0.0);
793 summary->jacobian_evaluation_time_in_seconds =
794 FindWithDefault(evaluator_time_statistics, "Evaluator::Jacobian", 0.0);
795
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800796 // Stick a fork in it, we're done.
797 summary->postprocessor_time_in_seconds =
798 WallTimeInSeconds() - post_process_start_time;
799}
Sameer Agarwal700d50d2013-03-12 16:12:42 -0700800#endif // CERES_NO_LINE_SEARCH_MINIMIZER
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800801
Sameer Agarwal65625f72012-09-17 12:06:57 -0700802bool SolverImpl::IsOrderingValid(const Solver::Options& options,
803 const ProblemImpl* problem_impl,
804 string* error) {
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700805 if (options.linear_solver_ordering->NumElements() !=
Sameer Agarwal65625f72012-09-17 12:06:57 -0700806 problem_impl->NumParameterBlocks()) {
807 *error = "Number of parameter blocks in user supplied ordering "
808 "does not match the number of parameter blocks in the problem";
809 return false;
810 }
811
812 const Program& program = problem_impl->program();
813 const vector<ParameterBlock*>& parameter_blocks = program.parameter_blocks();
Sameer Agarwal65625f72012-09-17 12:06:57 -0700814 for (vector<ParameterBlock*>::const_iterator it = parameter_blocks.begin();
815 it != parameter_blocks.end();
816 ++it) {
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700817 if (!options.linear_solver_ordering
818 ->IsMember(const_cast<double*>((*it)->user_state()))) {
Sameer Agarwal65625f72012-09-17 12:06:57 -0700819 *error = "Problem contains a parameter block that is not in "
820 "the user specified ordering.";
821 return false;
822 }
823 }
824
825 if (IsSchurType(options.linear_solver_type) &&
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700826 options.linear_solver_ordering->NumGroups() > 1) {
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700827 const vector<ResidualBlock*>& residual_blocks = program.residual_blocks();
Sameer Agarwal65625f72012-09-17 12:06:57 -0700828 const set<double*>& e_blocks =
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700829 options.linear_solver_ordering->group_to_elements().begin()->second;
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700830 if (!IsParameterBlockSetIndependent(e_blocks, residual_blocks)) {
831 *error = "The user requested the use of a Schur type solver. "
832 "But the first elimination group in the ordering is not an "
833 "independent set.";
834 return false;
Sameer Agarwal65625f72012-09-17 12:06:57 -0700835 }
836 }
837 return true;
838}
839
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800840bool SolverImpl::IsParameterBlockSetIndependent(
841 const set<double*>& parameter_block_ptrs,
842 const vector<ResidualBlock*>& residual_blocks) {
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700843 // Loop over each residual block and ensure that no two parameter
844 // blocks in the same residual block are part of
845 // parameter_block_ptrs as that would violate the assumption that it
846 // is an independent set in the Hessian matrix.
847 for (vector<ResidualBlock*>::const_iterator it = residual_blocks.begin();
848 it != residual_blocks.end();
849 ++it) {
850 ParameterBlock* const* parameter_blocks = (*it)->parameter_blocks();
851 const int num_parameter_blocks = (*it)->NumParameterBlocks();
852 int count = 0;
853 for (int i = 0; i < num_parameter_blocks; ++i) {
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700854 count += parameter_block_ptrs.count(
855 parameter_blocks[i]->mutable_user_state());
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700856 }
857 if (count > 1) {
858 return false;
859 }
860 }
861 return true;
862}
863
864
Keir Mierle8ebb0732012-04-30 23:09:08 -0700865// Strips varying parameters and residuals, maintaining order, and updating
866// num_eliminate_blocks.
867bool SolverImpl::RemoveFixedBlocksFromProgram(Program* program,
Sameer Agarwal2c94eed2012-10-01 16:34:37 -0700868 ParameterBlockOrdering* ordering,
Markus Moll8de78db2012-07-14 15:49:11 +0200869 double* fixed_cost,
Keir Mierle8ebb0732012-04-30 23:09:08 -0700870 string* error) {
Keir Mierle8ebb0732012-04-30 23:09:08 -0700871 vector<ParameterBlock*>* parameter_blocks =
872 program->mutable_parameter_blocks();
873
Markus Moll8de78db2012-07-14 15:49:11 +0200874 scoped_array<double> residual_block_evaluate_scratch;
875 if (fixed_cost != NULL) {
Keir Mierlec161a9d2012-07-18 14:01:48 -0700876 residual_block_evaluate_scratch.reset(
Markus Moll8de78db2012-07-14 15:49:11 +0200877 new double[program->MaxScratchDoublesNeededForEvaluate()]);
878 *fixed_cost = 0.0;
879 }
880
Keir Mierle8ebb0732012-04-30 23:09:08 -0700881 // Mark all the parameters as unused. Abuse the index member of the parameter
882 // blocks for the marking.
883 for (int i = 0; i < parameter_blocks->size(); ++i) {
884 (*parameter_blocks)[i]->set_index(-1);
885 }
886
887 // Filter out residual that have all-constant parameters, and mark all the
888 // parameter blocks that appear in residuals.
889 {
890 vector<ResidualBlock*>* residual_blocks =
891 program->mutable_residual_blocks();
892 int j = 0;
893 for (int i = 0; i < residual_blocks->size(); ++i) {
894 ResidualBlock* residual_block = (*residual_blocks)[i];
895 int num_parameter_blocks = residual_block->NumParameterBlocks();
896
897 // Determine if the residual block is fixed, and also mark varying
898 // parameters that appear in the residual block.
899 bool all_constant = true;
900 for (int k = 0; k < num_parameter_blocks; k++) {
901 ParameterBlock* parameter_block = residual_block->parameter_blocks()[k];
902 if (!parameter_block->IsConstant()) {
903 all_constant = false;
904 parameter_block->set_index(1);
905 }
906 }
907
908 if (!all_constant) {
909 (*residual_blocks)[j++] = (*residual_blocks)[i];
Markus Moll8de78db2012-07-14 15:49:11 +0200910 } else if (fixed_cost != NULL) {
911 // The residual is constant and will be removed, so its cost is
912 // added to the variable fixed_cost.
913 double cost = 0.0;
Sameer Agarwal039ff072013-02-26 09:15:39 -0800914 if (!residual_block->Evaluate(true,
915 &cost,
916 NULL,
917 NULL,
918 residual_block_evaluate_scratch.get())) {
Markus Moll8de78db2012-07-14 15:49:11 +0200919 *error = StringPrintf("Evaluation of the residual %d failed during "
920 "removal of fixed residual blocks.", i);
921 return false;
922 }
923 *fixed_cost += cost;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700924 }
925 }
926 residual_blocks->resize(j);
927 }
928
929 // Filter out unused or fixed parameter blocks, and update
Sameer Agarwal65625f72012-09-17 12:06:57 -0700930 // the ordering.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700931 {
932 vector<ParameterBlock*>* parameter_blocks =
933 program->mutable_parameter_blocks();
934 int j = 0;
935 for (int i = 0; i < parameter_blocks->size(); ++i) {
936 ParameterBlock* parameter_block = (*parameter_blocks)[i];
937 if (parameter_block->index() == 1) {
938 (*parameter_blocks)[j++] = parameter_block;
Sameer Agarwal65625f72012-09-17 12:06:57 -0700939 } else {
Sameer Agarwal2c94eed2012-10-01 16:34:37 -0700940 ordering->Remove(parameter_block->mutable_user_state());
Keir Mierle8ebb0732012-04-30 23:09:08 -0700941 }
942 }
943 parameter_blocks->resize(j);
944 }
945
Sameer Agarwalc8f07902013-04-19 08:01:04 -0700946 if (!(((program->NumResidualBlocks() == 0) &&
Keir Mierle8ebb0732012-04-30 23:09:08 -0700947 (program->NumParameterBlocks() == 0)) ||
948 ((program->NumResidualBlocks() != 0) &&
Sameer Agarwalc8f07902013-04-19 08:01:04 -0700949 (program->NumParameterBlocks() != 0)))) {
950 *error = "Congratulations, you found a bug in Ceres. Please report it.";
951 return false;
952 }
953
Keir Mierle8ebb0732012-04-30 23:09:08 -0700954 return true;
955}
956
957Program* SolverImpl::CreateReducedProgram(Solver::Options* options,
958 ProblemImpl* problem_impl,
Markus Moll8de78db2012-07-14 15:49:11 +0200959 double* fixed_cost,
Keir Mierle8ebb0732012-04-30 23:09:08 -0700960 string* error) {
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700961 CHECK_NOTNULL(options->linear_solver_ordering);
Keir Mierle8ebb0732012-04-30 23:09:08 -0700962 Program* original_program = problem_impl->mutable_program();
963 scoped_ptr<Program> transformed_program(new Program(*original_program));
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800964
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700965 ParameterBlockOrdering* linear_solver_ordering =
966 options->linear_solver_ordering;
Sameer Agarwal65625f72012-09-17 12:06:57 -0700967 const int min_group_id =
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700968 linear_solver_ordering->group_to_elements().begin()->first;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700969
970 if (!RemoveFixedBlocksFromProgram(transformed_program.get(),
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700971 linear_solver_ordering,
Markus Moll8de78db2012-07-14 15:49:11 +0200972 fixed_cost,
Keir Mierle8ebb0732012-04-30 23:09:08 -0700973 error)) {
974 return NULL;
975 }
976
977 if (transformed_program->NumParameterBlocks() == 0) {
978 LOG(WARNING) << "No varying parameter blocks to optimize; "
979 << "bailing early.";
980 return transformed_program.release();
981 }
982
Sameer Agarwal65625f72012-09-17 12:06:57 -0700983 if (IsSchurType(options->linear_solver_type) &&
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700984 linear_solver_ordering->GroupSize(min_group_id) == 0) {
Sameer Agarwalc8f07902013-04-19 08:01:04 -0700985 // If the user requested the use of a Schur type solver, and
986 // supplied a non-NULL linear_solver_ordering object with more than
987 // one elimination group, then it can happen that after all the
988 // parameter blocks which are fixed or unused have been removed from
989 // the program and the ordering, there are no more parameter blocks
990 // in the first elimination group.
991 //
992 // In such a case, the use of a Schur type solver is not possible,
993 // as they assume there is at least one e_block. Thus, we
994 // automatically switch to the closest solver to the one indicated
995 // by the user.
996 AlternateLinearSolverForSchurTypeLinearSolver(options);
Sameer Agarwal65625f72012-09-17 12:06:57 -0700997 }
Keir Mierle8ebb0732012-04-30 23:09:08 -0700998
Sameer Agarwalc8f07902013-04-19 08:01:04 -0700999 if (IsSchurType(options->linear_solver_type)) {
Sameer Agarwalac626962013-05-06 07:04:26 -07001000 if (!ReorderProgramForSchurTypeLinearSolver(
1001 options->linear_solver_type,
1002 options->sparse_linear_algebra_library,
1003 problem_impl->parameter_map(),
1004 linear_solver_ordering,
1005 transformed_program.get(),
1006 error)) {
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001007 return NULL;
1008 }
1009 return transformed_program.release();
1010 }
Sameer Agarwal42a84b82013-02-01 12:22:53 -08001011
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001012 if (options->linear_solver_type == SPARSE_NORMAL_CHOLESKY) {
1013 if (!ReorderProgramForSparseNormalCholesky(
1014 options->sparse_linear_algebra_library,
1015 linear_solver_ordering,
1016 transformed_program.get(),
1017 error)) {
1018 return NULL;
1019 }
1020
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001021 return transformed_program.release();
1022 }
1023
Keir Mierle8ebb0732012-04-30 23:09:08 -07001024 transformed_program->SetParameterOffsetsAndIndex();
1025 return transformed_program.release();
1026}
1027
1028LinearSolver* SolverImpl::CreateLinearSolver(Solver::Options* options,
1029 string* error) {
Sameer Agarwal65625f72012-09-17 12:06:57 -07001030 CHECK_NOTNULL(options);
Sameer Agarwal68b32a92012-10-06 23:10:51 -07001031 CHECK_NOTNULL(options->linear_solver_ordering);
Sameer Agarwal65625f72012-09-17 12:06:57 -07001032 CHECK_NOTNULL(error);
1033
Sameer Agarwal5ecd2512012-06-17 16:34:03 -07001034 if (options->trust_region_strategy_type == DOGLEG) {
1035 if (options->linear_solver_type == ITERATIVE_SCHUR ||
1036 options->linear_solver_type == CGNR) {
1037 *error = "DOGLEG only supports exact factorization based linear "
1038 "solvers. If you want to use an iterative solver please "
1039 "use LEVENBERG_MARQUARDT as the trust_region_strategy_type";
1040 return NULL;
1041 }
1042 }
1043
Keir Mierle8ebb0732012-04-30 23:09:08 -07001044#ifdef CERES_NO_SUITESPARSE
Sameer Agarwalb0518732012-05-29 00:27:57 -07001045 if (options->linear_solver_type == SPARSE_NORMAL_CHOLESKY &&
1046 options->sparse_linear_algebra_library == SUITE_SPARSE) {
1047 *error = "Can't use SPARSE_NORMAL_CHOLESKY with SUITESPARSE because "
1048 "SuiteSparse was not enabled when Ceres was built.";
Keir Mierle8ebb0732012-04-30 23:09:08 -07001049 return NULL;
1050 }
Sameer Agarwal65625f72012-09-17 12:06:57 -07001051
Petter Strandmarkdd1a2762012-09-19 12:55:05 +02001052 if (options->preconditioner_type == CLUSTER_JACOBI) {
Sameer Agarwal65625f72012-09-17 12:06:57 -07001053 *error = "CLUSTER_JACOBI preconditioner not suppored. Please build Ceres "
1054 "with SuiteSparse support.";
1055 return NULL;
1056 }
1057
Petter Strandmarkdd1a2762012-09-19 12:55:05 +02001058 if (options->preconditioner_type == CLUSTER_TRIDIAGONAL) {
Sameer Agarwal65625f72012-09-17 12:06:57 -07001059 *error = "CLUSTER_TRIDIAGONAL preconditioner not suppored. Please build "
1060 "Ceres with SuiteSparse support.";
1061 return NULL;
1062 }
Sameer Agarwalb0518732012-05-29 00:27:57 -07001063#endif
1064
1065#ifdef CERES_NO_CXSPARSE
1066 if (options->linear_solver_type == SPARSE_NORMAL_CHOLESKY &&
1067 options->sparse_linear_algebra_library == CX_SPARSE) {
1068 *error = "Can't use SPARSE_NORMAL_CHOLESKY with CXSPARSE because "
1069 "CXSparse was not enabled when Ceres was built.";
1070 return NULL;
1071 }
1072#endif
1073
Sameer Agarwal65625f72012-09-17 12:06:57 -07001074#if defined(CERES_NO_SUITESPARSE) && defined(CERES_NO_CXSPARSE)
Keir Mierle27dd0d32012-10-15 13:52:36 -07001075 if (options->linear_solver_type == SPARSE_SCHUR) {
Sameer Agarwal65625f72012-09-17 12:06:57 -07001076 *error = "Can't use SPARSE_SCHUR because neither SuiteSparse nor"
1077 "CXSparse was enabled when Ceres was compiled.";
1078 return NULL;
1079 }
1080#endif
Keir Mierle8ebb0732012-04-30 23:09:08 -07001081
1082 if (options->linear_solver_max_num_iterations <= 0) {
1083 *error = "Solver::Options::linear_solver_max_num_iterations is 0.";
1084 return NULL;
1085 }
1086 if (options->linear_solver_min_num_iterations <= 0) {
1087 *error = "Solver::Options::linear_solver_min_num_iterations is 0.";
1088 return NULL;
1089 }
1090 if (options->linear_solver_min_num_iterations >
1091 options->linear_solver_max_num_iterations) {
1092 *error = "Solver::Options::linear_solver_min_num_iterations > "
1093 "Solver::Options::linear_solver_max_num_iterations.";
1094 return NULL;
1095 }
1096
1097 LinearSolver::Options linear_solver_options;
Keir Mierle8ebb0732012-04-30 23:09:08 -07001098 linear_solver_options.min_num_iterations =
1099 options->linear_solver_min_num_iterations;
1100 linear_solver_options.max_num_iterations =
1101 options->linear_solver_max_num_iterations;
1102 linear_solver_options.type = options->linear_solver_type;
1103 linear_solver_options.preconditioner_type = options->preconditioner_type;
Sameer Agarwalb0518732012-05-29 00:27:57 -07001104 linear_solver_options.sparse_linear_algebra_library =
1105 options->sparse_linear_algebra_library;
Sameer Agarwal9189f4e2013-04-19 17:09:49 -07001106 linear_solver_options.use_postordering = options->use_postordering;
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001107
1108 // Ignore user's postordering preferences and force it to be true if
1109 // cholmod_camd is not available. This ensures that the linear
1110 // solver does not assume that a fill-reducing pre-ordering has been
1111 // done.
1112#if !defined(CERES_NO_SUITESPARSE) && defined(CERES_NO_CAMD)
1113 if (IsSchurType(linear_solver_options.type) &&
1114 linear_solver_options.sparse_linear_algebra_library == SUITE_SPARSE) {
1115 linear_solver_options.use_postordering = true;
1116 }
1117#endif
1118
Keir Mierle8ebb0732012-04-30 23:09:08 -07001119 linear_solver_options.num_threads = options->num_linear_solver_threads;
Keir Mierle8ebb0732012-04-30 23:09:08 -07001120 options->num_linear_solver_threads = linear_solver_options.num_threads;
1121
Sameer Agarwal65625f72012-09-17 12:06:57 -07001122 const map<int, set<double*> >& groups =
Sameer Agarwal68b32a92012-10-06 23:10:51 -07001123 options->linear_solver_ordering->group_to_elements();
Sameer Agarwal65625f72012-09-17 12:06:57 -07001124 for (map<int, set<double*> >::const_iterator it = groups.begin();
1125 it != groups.end();
1126 ++it) {
1127 linear_solver_options.elimination_groups.push_back(it->second.size());
1128 }
1129 // Schur type solvers, expect at least two elimination groups. If
Sameer Agarwal1dc544a2012-09-18 18:09:03 -07001130 // there is only one elimination group, then CreateReducedProgram
Sameer Agarwal65625f72012-09-17 12:06:57 -07001131 // guarantees that this group only contains e_blocks. Thus we add a
1132 // dummy elimination group with zero blocks in it.
1133 if (IsSchurType(linear_solver_options.type) &&
1134 linear_solver_options.elimination_groups.size() == 1) {
1135 linear_solver_options.elimination_groups.push_back(0);
1136 }
1137
Keir Mierle8ebb0732012-04-30 23:09:08 -07001138 return LinearSolver::Create(linear_solver_options);
1139}
1140
Keir Mierle8ebb0732012-04-30 23:09:08 -07001141
1142// Find the minimum index of any parameter block to the given residual.
1143// Parameter blocks that have indices greater than num_eliminate_blocks are
1144// considered to have an index equal to num_eliminate_blocks.
Sergey Sharybinb53c9662013-02-25 01:14:26 +06001145static int MinParameterBlock(const ResidualBlock* residual_block,
1146 int num_eliminate_blocks) {
Keir Mierle8ebb0732012-04-30 23:09:08 -07001147 int min_parameter_block_position = num_eliminate_blocks;
1148 for (int i = 0; i < residual_block->NumParameterBlocks(); ++i) {
1149 ParameterBlock* parameter_block = residual_block->parameter_blocks()[i];
Keir Mierle32de18d2012-05-13 16:45:05 -07001150 if (!parameter_block->IsConstant()) {
1151 CHECK_NE(parameter_block->index(), -1)
1152 << "Did you forget to call Program::SetParameterOffsetsAndIndex()? "
1153 << "This is a Ceres bug; please contact the developers!";
1154 min_parameter_block_position = std::min(parameter_block->index(),
1155 min_parameter_block_position);
1156 }
Keir Mierle8ebb0732012-04-30 23:09:08 -07001157 }
1158 return min_parameter_block_position;
1159}
1160
1161// Reorder the residuals for program, if necessary, so that the residuals
1162// involving each E block occur together. This is a necessary condition for the
1163// Schur eliminator, which works on these "row blocks" in the jacobian.
Sameer Agarwal509f68c2013-02-20 01:39:03 -08001164bool SolverImpl::LexicographicallyOrderResidualBlocks(
1165 const int num_eliminate_blocks,
1166 Program* program,
1167 string* error) {
Sameer Agarwal9123e2f2012-09-18 21:49:06 -07001168 CHECK_GE(num_eliminate_blocks, 1)
1169 << "Congratulations, you found a Ceres bug! Please report this error "
1170 << "to the developers.";
Keir Mierle8ebb0732012-04-30 23:09:08 -07001171
1172 // Create a histogram of the number of residuals for each E block. There is an
1173 // extra bucket at the end to catch all non-eliminated F blocks.
Sameer Agarwal65625f72012-09-17 12:06:57 -07001174 vector<int> residual_blocks_per_e_block(num_eliminate_blocks + 1);
Keir Mierle8ebb0732012-04-30 23:09:08 -07001175 vector<ResidualBlock*>* residual_blocks = program->mutable_residual_blocks();
1176 vector<int> min_position_per_residual(residual_blocks->size());
1177 for (int i = 0; i < residual_blocks->size(); ++i) {
1178 ResidualBlock* residual_block = (*residual_blocks)[i];
Sameer Agarwal65625f72012-09-17 12:06:57 -07001179 int position = MinParameterBlock(residual_block, num_eliminate_blocks);
Keir Mierle8ebb0732012-04-30 23:09:08 -07001180 min_position_per_residual[i] = position;
Sameer Agarwal65625f72012-09-17 12:06:57 -07001181 DCHECK_LE(position, num_eliminate_blocks);
Keir Mierle8ebb0732012-04-30 23:09:08 -07001182 residual_blocks_per_e_block[position]++;
1183 }
1184
1185 // Run a cumulative sum on the histogram, to obtain offsets to the start of
1186 // each histogram bucket (where each bucket is for the residuals for that
1187 // E-block).
Sameer Agarwal65625f72012-09-17 12:06:57 -07001188 vector<int> offsets(num_eliminate_blocks + 1);
Keir Mierle8ebb0732012-04-30 23:09:08 -07001189 std::partial_sum(residual_blocks_per_e_block.begin(),
1190 residual_blocks_per_e_block.end(),
1191 offsets.begin());
1192 CHECK_EQ(offsets.back(), residual_blocks->size())
1193 << "Congratulations, you found a Ceres bug! Please report this error "
1194 << "to the developers.";
1195
1196 CHECK(find(residual_blocks_per_e_block.begin(),
1197 residual_blocks_per_e_block.end() - 1, 0) !=
1198 residual_blocks_per_e_block.end())
1199 << "Congratulations, you found a Ceres bug! Please report this error "
1200 << "to the developers.";
1201
1202 // Fill in each bucket with the residual blocks for its corresponding E block.
1203 // Each bucket is individually filled from the back of the bucket to the front
1204 // of the bucket. The filling order among the buckets is dictated by the
1205 // residual blocks. This loop uses the offsets as counters; subtracting one
1206 // from each offset as a residual block is placed in the bucket. When the
1207 // filling is finished, the offset pointerts should have shifted down one
1208 // entry (this is verified below).
1209 vector<ResidualBlock*> reordered_residual_blocks(
1210 (*residual_blocks).size(), static_cast<ResidualBlock*>(NULL));
1211 for (int i = 0; i < residual_blocks->size(); ++i) {
1212 int bucket = min_position_per_residual[i];
1213
1214 // Decrement the cursor, which should now point at the next empty position.
1215 offsets[bucket]--;
1216
1217 // Sanity.
1218 CHECK(reordered_residual_blocks[offsets[bucket]] == NULL)
1219 << "Congratulations, you found a Ceres bug! Please report this error "
1220 << "to the developers.";
1221
1222 reordered_residual_blocks[offsets[bucket]] = (*residual_blocks)[i];
1223 }
1224
1225 // Sanity check #1: The difference in bucket offsets should match the
1226 // histogram sizes.
Sameer Agarwal65625f72012-09-17 12:06:57 -07001227 for (int i = 0; i < num_eliminate_blocks; ++i) {
Keir Mierle8ebb0732012-04-30 23:09:08 -07001228 CHECK_EQ(residual_blocks_per_e_block[i], offsets[i + 1] - offsets[i])
1229 << "Congratulations, you found a Ceres bug! Please report this error "
1230 << "to the developers.";
1231 }
1232 // Sanity check #2: No NULL's left behind.
1233 for (int i = 0; i < reordered_residual_blocks.size(); ++i) {
1234 CHECK(reordered_residual_blocks[i] != NULL)
1235 << "Congratulations, you found a Ceres bug! Please report this error "
1236 << "to the developers.";
1237 }
1238
1239 // Now that the residuals are collected by E block, swap them in place.
1240 swap(*program->mutable_residual_blocks(), reordered_residual_blocks);
1241 return true;
1242}
1243
Sameer Agarwal509f68c2013-02-20 01:39:03 -08001244Evaluator* SolverImpl::CreateEvaluator(
1245 const Solver::Options& options,
1246 const ProblemImpl::ParameterMap& parameter_map,
1247 Program* program,
1248 string* error) {
Keir Mierle8ebb0732012-04-30 23:09:08 -07001249 Evaluator::Options evaluator_options;
1250 evaluator_options.linear_solver_type = options.linear_solver_type;
Sameer Agarwal65625f72012-09-17 12:06:57 -07001251 evaluator_options.num_eliminate_blocks =
Sameer Agarwal68b32a92012-10-06 23:10:51 -07001252 (options.linear_solver_ordering->NumGroups() > 0 &&
Sameer Agarwal65625f72012-09-17 12:06:57 -07001253 IsSchurType(options.linear_solver_type))
Sameer Agarwal68b32a92012-10-06 23:10:51 -07001254 ? (options.linear_solver_ordering
1255 ->group_to_elements().begin()
1256 ->second.size())
Sameer Agarwal9123e2f2012-09-18 21:49:06 -07001257 : 0;
Keir Mierle8ebb0732012-04-30 23:09:08 -07001258 evaluator_options.num_threads = options.num_threads;
1259 return Evaluator::Create(evaluator_options, program, error);
1260}
1261
Sameer Agarwal67a107b2012-10-08 14:35:41 -07001262CoordinateDescentMinimizer* SolverImpl::CreateInnerIterationMinimizer(
1263 const Solver::Options& options,
1264 const Program& program,
1265 const ProblemImpl::ParameterMap& parameter_map,
Sameer Agarwal977be7c2013-01-26 16:01:54 -08001266 Solver::Summary* summary) {
Sameer Agarwal67a107b2012-10-08 14:35:41 -07001267 scoped_ptr<CoordinateDescentMinimizer> inner_iteration_minimizer(
1268 new CoordinateDescentMinimizer);
1269 scoped_ptr<ParameterBlockOrdering> inner_iteration_ordering;
1270 ParameterBlockOrdering* ordering_ptr = NULL;
1271
1272 if (options.inner_iteration_ordering == NULL) {
1273 // Find a recursive decomposition of the Hessian matrix as a set
1274 // of independent sets of decreasing size and invert it. This
1275 // seems to work better in practice, i.e., Cameras before
1276 // points.
1277 inner_iteration_ordering.reset(new ParameterBlockOrdering);
1278 ComputeRecursiveIndependentSetOrdering(program,
1279 inner_iteration_ordering.get());
1280 inner_iteration_ordering->Reverse();
1281 ordering_ptr = inner_iteration_ordering.get();
1282 } else {
1283 const map<int, set<double*> >& group_to_elements =
1284 options.inner_iteration_ordering->group_to_elements();
1285
1286 // Iterate over each group and verify that it is an independent
1287 // set.
1288 map<int, set<double*> >::const_iterator it = group_to_elements.begin();
Sameer Agarwal509f68c2013-02-20 01:39:03 -08001289 for ( ; it != group_to_elements.end(); ++it) {
Sameer Agarwal67a107b2012-10-08 14:35:41 -07001290 if (!IsParameterBlockSetIndependent(it->second,
1291 program.residual_blocks())) {
Sameer Agarwal977be7c2013-01-26 16:01:54 -08001292 summary->error =
Sameer Agarwal67a107b2012-10-08 14:35:41 -07001293 StringPrintf("The user-provided "
1294 "parameter_blocks_for_inner_iterations does not "
1295 "form an independent set. Group Id: %d", it->first);
1296 return NULL;
1297 }
1298 }
1299 ordering_ptr = options.inner_iteration_ordering;
1300 }
1301
1302 if (!inner_iteration_minimizer->Init(program,
1303 parameter_map,
1304 *ordering_ptr,
Sameer Agarwal977be7c2013-01-26 16:01:54 -08001305 &summary->error)) {
Sameer Agarwal67a107b2012-10-08 14:35:41 -07001306 return NULL;
1307 }
1308
Sameer Agarwal977be7c2013-01-26 16:01:54 -08001309 summary->inner_iterations = true;
1310 SummarizeOrdering(ordering_ptr, &(summary->inner_iteration_ordering_used));
1311
Sameer Agarwal67a107b2012-10-08 14:35:41 -07001312 return inner_iteration_minimizer.release();
1313}
1314
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001315void SolverImpl::AlternateLinearSolverForSchurTypeLinearSolver(
1316 Solver::Options* options) {
1317 if (!IsSchurType(options->linear_solver_type)) {
1318 return;
1319 }
1320
1321 string msg = "No e_blocks remaining. Switching from ";
1322 if (options->linear_solver_type == SPARSE_SCHUR) {
1323 options->linear_solver_type = SPARSE_NORMAL_CHOLESKY;
1324 msg += "SPARSE_SCHUR to SPARSE_NORMAL_CHOLESKY.";
1325 } else if (options->linear_solver_type == DENSE_SCHUR) {
1326 // TODO(sameeragarwal): This is probably not a great choice.
1327 // Ideally, we should have a DENSE_NORMAL_CHOLESKY, that can
1328 // take a BlockSparseMatrix as input.
1329 options->linear_solver_type = DENSE_QR;
1330 msg += "DENSE_SCHUR to DENSE_QR.";
1331 } else if (options->linear_solver_type == ITERATIVE_SCHUR) {
1332 options->linear_solver_type = CGNR;
1333 if (options->preconditioner_type != IDENTITY) {
1334 msg += StringPrintf("ITERATIVE_SCHUR with %s preconditioner "
1335 "to CGNR with JACOBI preconditioner.",
1336 PreconditionerTypeToString(
1337 options->preconditioner_type));
1338 // CGNR currently only supports the JACOBI preconditioner.
1339 options->preconditioner_type = JACOBI;
1340 } else {
Sameer Agarwalcbdeb792013-04-22 10:18:18 -07001341 msg += "ITERATIVE_SCHUR with IDENTITY preconditioner"
1342 "to CGNR with IDENTITY preconditioner.";
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001343 }
1344 }
1345 LOG(WARNING) << msg;
1346}
1347
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001348bool SolverImpl::ApplyUserOrdering(
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001349 const ProblemImpl::ParameterMap& parameter_map,
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001350 const ParameterBlockOrdering* parameter_block_ordering,
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001351 Program* program,
1352 string* error) {
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001353 const int num_parameter_blocks = program->NumParameterBlocks();
1354 if (parameter_block_ordering->NumElements() != num_parameter_blocks) {
1355 *error = StringPrintf("User specified ordering does not have the same "
1356 "number of parameters as the problem. The problem"
1357 "has %d blocks while the ordering has %d blocks.",
1358 num_parameter_blocks,
1359 parameter_block_ordering->NumElements());
1360 return false;
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001361 }
1362
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001363 vector<ParameterBlock*>* parameter_blocks =
1364 program->mutable_parameter_blocks();
1365 parameter_blocks->clear();
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001366
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001367 const map<int, set<double*> >& groups =
1368 parameter_block_ordering->group_to_elements();
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001369
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001370 for (map<int, set<double*> >::const_iterator group_it = groups.begin();
1371 group_it != groups.end();
1372 ++group_it) {
1373 const set<double*>& group = group_it->second;
1374 for (set<double*>::const_iterator parameter_block_ptr_it = group.begin();
1375 parameter_block_ptr_it != group.end();
1376 ++parameter_block_ptr_it) {
1377 ProblemImpl::ParameterMap::const_iterator parameter_block_it =
1378 parameter_map.find(*parameter_block_ptr_it);
1379 if (parameter_block_it == parameter_map.end()) {
1380 *error = StringPrintf("User specified ordering contains a pointer "
1381 "to a double that is not a parameter block in "
1382 "the problem. The invalid double is in group: %d",
1383 group_it->first);
1384 return false;
1385 }
1386 parameter_blocks->push_back(parameter_block_it->second);
1387 }
1388 }
1389 return true;
Sameer Agarwalc8f07902013-04-19 08:01:04 -07001390}
1391
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001392
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001393TripletSparseMatrix* SolverImpl::CreateJacobianBlockSparsityTranspose(
1394 const Program* program) {
1395
Sameer Agarwalcbdeb792013-04-22 10:18:18 -07001396 // Matrix to store the block sparsity structure of the Jacobian.
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001397 TripletSparseMatrix* tsm =
1398 new TripletSparseMatrix(program->NumParameterBlocks(),
1399 program->NumResidualBlocks(),
1400 10 * program->NumResidualBlocks());
1401 int num_nonzeros = 0;
1402 int* rows = tsm->mutable_rows();
1403 int* cols = tsm->mutable_cols();
1404 double* values = tsm->mutable_values();
1405
1406 const vector<ResidualBlock*>& residual_blocks = program->residual_blocks();
1407 for (int c = 0; c < residual_blocks.size(); ++c) {
1408 const ResidualBlock* residual_block = residual_blocks[c];
1409 const int num_parameter_blocks = residual_block->NumParameterBlocks();
1410 ParameterBlock* const* parameter_blocks =
1411 residual_block->parameter_blocks();
1412
1413 for (int j = 0; j < num_parameter_blocks; ++j) {
1414 if (parameter_blocks[j]->IsConstant()) {
1415 continue;
1416 }
1417
1418 // Re-size the matrix if needed.
1419 if (num_nonzeros >= tsm->max_num_nonzeros()) {
1420 tsm->Reserve(2 * num_nonzeros);
1421 rows = tsm->mutable_rows();
1422 cols = tsm->mutable_cols();
1423 values = tsm->mutable_values();
1424 }
1425 CHECK_LT(num_nonzeros, tsm->max_num_nonzeros());
1426
1427 const int r = parameter_blocks[j]->index();
1428 rows[num_nonzeros] = r;
1429 cols[num_nonzeros] = c;
1430 values[num_nonzeros] = 1.0;
1431 ++num_nonzeros;
1432 }
1433 }
1434
1435 tsm->set_num_nonzeros(num_nonzeros);
1436 return tsm;
1437}
1438
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001439bool SolverImpl::ReorderProgramForSchurTypeLinearSolver(
1440 const LinearSolverType linear_solver_type,
1441 const SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type,
1442 const ProblemImpl::ParameterMap& parameter_map,
1443 ParameterBlockOrdering* parameter_block_ordering,
1444 Program* program,
1445 string* error) {
1446 if (parameter_block_ordering->NumGroups() == 1) {
1447 // If the user supplied an parameter_block_ordering with just one
1448 // group, it is equivalent to the user supplying NULL as an
1449 // parameter_block_ordering. Ceres is completely free to choose the
1450 // parameter block ordering as it sees fit. For Schur type solvers,
1451 // this means that the user wishes for Ceres to identify the
1452 // e_blocks, which we do by computing a maximal independent set.
1453 vector<ParameterBlock*> schur_ordering;
Sameer Agarwal36c73c22013-05-17 22:52:21 -07001454 const int num_eliminate_blocks = ComputeStableSchurOrdering(*program,
1455 &schur_ordering);
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001456
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001457 CHECK_EQ(schur_ordering.size(), program->NumParameterBlocks())
1458 << "Congratulations, you found a Ceres bug! Please report this error "
1459 << "to the developers.";
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001460
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001461 // Update the parameter_block_ordering object.
1462 for (int i = 0; i < schur_ordering.size(); ++i) {
1463 double* parameter_block = schur_ordering[i]->mutable_user_state();
1464 const int group_id = (i < num_eliminate_blocks) ? 0 : 1;
1465 parameter_block_ordering->AddElementToGroup(parameter_block, group_id);
1466 }
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001467
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001468 // We could call ApplyUserOrdering but this is cheaper and
1469 // simpler.
1470 swap(*program->mutable_parameter_blocks(), schur_ordering);
1471 } else {
1472 // The user provided an ordering with more than one elimination
1473 // group. Trust the user and apply the ordering.
1474 if (!ApplyUserOrdering(parameter_map,
1475 parameter_block_ordering,
1476 program,
1477 error)) {
1478 return false;
1479 }
1480 }
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001481
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001482 // Pre-order the columns corresponding to the schur complement if
1483 // possible.
1484#if !defined(CERES_NO_SUITESPARSE) && !defined(CERES_NO_CAMD)
1485 if (linear_solver_type == SPARSE_SCHUR &&
1486 sparse_linear_algebra_library_type == SUITE_SPARSE) {
1487 vector<int> constraints;
1488 vector<ParameterBlock*>& parameter_blocks =
1489 *(program->mutable_parameter_blocks());
1490
1491 for (int i = 0; i < parameter_blocks.size(); ++i) {
1492 constraints.push_back(
1493 parameter_block_ordering->GroupId(
1494 parameter_blocks[i]->mutable_user_state()));
1495 }
1496
1497 // Set the offsets and index for CreateJacobianSparsityTranspose.
1498 program->SetParameterOffsetsAndIndex();
1499 // Compute a block sparse presentation of J'.
1500 scoped_ptr<TripletSparseMatrix> tsm_block_jacobian_transpose(
1501 SolverImpl::CreateJacobianBlockSparsityTranspose(program));
1502
1503 SuiteSparse ss;
1504 cholmod_sparse* block_jacobian_transpose =
1505 ss.CreateSparseMatrix(tsm_block_jacobian_transpose.get());
1506
1507 vector<int> ordering(parameter_blocks.size(), 0);
1508 ss.ConstrainedApproximateMinimumDegreeOrdering(block_jacobian_transpose,
1509 &constraints[0],
1510 &ordering[0]);
1511 ss.Free(block_jacobian_transpose);
1512
1513 const vector<ParameterBlock*> parameter_blocks_copy(parameter_blocks);
1514 for (int i = 0; i < program->NumParameterBlocks(); ++i) {
1515 parameter_blocks[i] = parameter_blocks_copy[ordering[i]];
1516 }
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001517 }
1518#endif
1519
1520 program->SetParameterOffsetsAndIndex();
Sameer Agarwald5b93bf2013-04-26 21:17:49 -07001521 // Schur type solvers also require that their residual blocks be
1522 // lexicographically ordered.
1523 const int num_eliminate_blocks =
1524 parameter_block_ordering->group_to_elements().begin()->second.size();
1525 return LexicographicallyOrderResidualBlocks(num_eliminate_blocks,
1526 program,
1527 error);
1528}
1529
1530bool SolverImpl::ReorderProgramForSparseNormalCholesky(
1531 const SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type,
1532 const ParameterBlockOrdering* parameter_block_ordering,
1533 Program* program,
1534 string* error) {
1535 // Set the offsets and index for CreateJacobianSparsityTranspose.
1536 program->SetParameterOffsetsAndIndex();
1537 // Compute a block sparse presentation of J'.
1538 scoped_ptr<TripletSparseMatrix> tsm_block_jacobian_transpose(
1539 SolverImpl::CreateJacobianBlockSparsityTranspose(program));
1540
1541 vector<int> ordering(program->NumParameterBlocks(), 0);
1542 vector<ParameterBlock*>& parameter_blocks =
1543 *(program->mutable_parameter_blocks());
1544
1545 if (sparse_linear_algebra_library_type == SUITE_SPARSE) {
1546#ifdef CERES_NO_SUITESPARSE
1547 *error = "Can't use SPARSE_NORMAL_CHOLESKY with SUITE_SPARSE because "
1548 "SuiteSparse was not enabled when Ceres was built.";
1549 return false;
1550#else
1551 SuiteSparse ss;
1552 cholmod_sparse* block_jacobian_transpose =
1553 ss.CreateSparseMatrix(tsm_block_jacobian_transpose.get());
1554
1555# ifdef CERES_NO_CAMD
1556 // No cholmod_camd, so ignore user's parameter_block_ordering and
1557 // use plain old AMD.
1558 ss.ApproximateMinimumDegreeOrdering(block_jacobian_transpose, &ordering[0]);
1559# else
1560 if (parameter_block_ordering->NumGroups() > 1) {
1561 // If the user specified more than one elimination groups use them
1562 // to constrain the ordering.
1563 vector<int> constraints;
1564 for (int i = 0; i < parameter_blocks.size(); ++i) {
1565 constraints.push_back(
1566 parameter_block_ordering->GroupId(
1567 parameter_blocks[i]->mutable_user_state()));
1568 }
1569 ss.ConstrainedApproximateMinimumDegreeOrdering(
1570 block_jacobian_transpose,
1571 &constraints[0],
1572 &ordering[0]);
1573 } else {
1574 ss.ApproximateMinimumDegreeOrdering(block_jacobian_transpose,
1575 &ordering[0]);
1576 }
1577# endif // CERES_NO_CAMD
1578
1579 ss.Free(block_jacobian_transpose);
1580#endif // CERES_NO_SUITESPARSE
1581
1582 } else if (sparse_linear_algebra_library_type == CX_SPARSE) {
1583#ifndef CERES_NO_CXSPARSE
1584
1585 // CXSparse works with J'J instead of J'. So compute the block
1586 // sparsity for J'J and compute an approximate minimum degree
1587 // ordering.
1588 CXSparse cxsparse;
1589 cs_di* block_jacobian_transpose;
1590 block_jacobian_transpose =
1591 cxsparse.CreateSparseMatrix(tsm_block_jacobian_transpose.get());
1592 cs_di* block_jacobian = cxsparse.TransposeMatrix(block_jacobian_transpose);
1593 cs_di* block_hessian =
1594 cxsparse.MatrixMatrixMultiply(block_jacobian_transpose, block_jacobian);
1595 cxsparse.Free(block_jacobian);
1596 cxsparse.Free(block_jacobian_transpose);
1597
1598 cxsparse.ApproximateMinimumDegreeOrdering(block_hessian, &ordering[0]);
1599 cxsparse.Free(block_hessian);
1600#else // CERES_NO_CXSPARSE
1601 *error = "Can't use SPARSE_NORMAL_CHOLESKY with CX_SPARSE because "
1602 "CXSparse was not enabled when Ceres was built.";
1603 return false;
1604#endif // CERES_NO_CXSPARSE
1605 } else {
1606 *error = "Unknown sparse linear algebra library.";
1607 return false;
1608 }
1609
1610 // Apply ordering.
1611 const vector<ParameterBlock*> parameter_blocks_copy(parameter_blocks);
1612 for (int i = 0; i < program->NumParameterBlocks(); ++i) {
1613 parameter_blocks[i] = parameter_blocks_copy[ordering[i]];
1614 }
1615
1616 program->SetParameterOffsetsAndIndex();
1617 return true;
Sameer Agarwalf7ed22e2013-04-19 14:24:48 -07001618}
1619
Keir Mierle8ebb0732012-04-30 23:09:08 -07001620} // namespace internal
1621} // namespace ceres