Change ownership of pointers in Minimizer::Options.
This is a intermediate change to clean things up
in preparation for a broader refactoring of the SolverImpl.
Essentially we are replacing raw pointers in Minimizer::Options
with shared_ptr objects. For now this only makes things a bit
more complicated looking inside solver_impl.cc, but going
forward this will lead to considerable simplifications in
tracking ownership of various pointers.
Change-Id: I21db8fc6763c29b0d15e834d7c968a0f514042a0
diff --git a/internal/ceres/coordinate_descent_minimizer.cc b/internal/ceres/coordinate_descent_minimizer.cc
index 3b0553e..1d55458 100644
--- a/internal/ceres/coordinate_descent_minimizer.cc
+++ b/internal/ceres/coordinate_descent_minimizer.cc
@@ -210,23 +210,16 @@
summary->final_cost = 0.0;
string error;
- scoped_ptr<Evaluator> evaluator(
- Evaluator::Create(evaluator_options_, program, &error));
- CHECK_NOTNULL(evaluator.get());
-
- scoped_ptr<SparseMatrix> jacobian(evaluator->CreateJacobian());
- CHECK_NOTNULL(jacobian.get());
+ Minimizer::Options minimizer_options;
+ minimizer_options.evaluator.reset(
+ CHECK_NOTNULL(Evaluator::Create(evaluator_options_, program, &error)));
+ minimizer_options.jacobian.reset(
+ CHECK_NOTNULL(minimizer_options.evaluator->CreateJacobian()));
TrustRegionStrategy::Options trs_options;
trs_options.linear_solver = linear_solver;
-
- scoped_ptr<TrustRegionStrategy>trust_region_strategy(
+ minimizer_options.trust_region_strategy.reset(
CHECK_NOTNULL(TrustRegionStrategy::Create(trs_options)));
-
- Minimizer::Options minimizer_options;
- minimizer_options.evaluator = evaluator.get();
- minimizer_options.jacobian = jacobian.get();
- minimizer_options.trust_region_strategy = trust_region_strategy.get();
minimizer_options.is_silent = true;
TrustRegionMinimizer minimizer;
diff --git a/internal/ceres/coordinate_descent_minimizer.h b/internal/ceres/coordinate_descent_minimizer.h
index e324b38..c1f8ffc 100644
--- a/internal/ceres/coordinate_descent_minimizer.h
+++ b/internal/ceres/coordinate_descent_minimizer.h
@@ -43,6 +43,7 @@
namespace internal {
class Program;
+class LinearSolver;
// Given a Program, and a ParameterBlockOrdering which partitions
// (non-exhaustively) the Hessian matrix into independent sets,
diff --git a/internal/ceres/line_search_minimizer.cc b/internal/ceres/line_search_minimizer.cc
index ae77a73..f494eda 100644
--- a/internal/ceres/line_search_minimizer.cc
+++ b/internal/ceres/line_search_minimizer.cc
@@ -103,7 +103,7 @@
double start_time = WallTimeInSeconds();
double iteration_start_time = start_time;
- Evaluator* evaluator = CHECK_NOTNULL(options.evaluator);
+ Evaluator* evaluator = CHECK_NOTNULL(options.evaluator.get());
const int num_parameters = evaluator->NumParameters();
const int num_effective_parameters = evaluator->NumEffectiveParameters();
diff --git a/internal/ceres/minimizer.cc b/internal/ceres/minimizer.cc
index 6c3b68d..558921b 100644
--- a/internal/ceres/minimizer.cc
+++ b/internal/ceres/minimizer.cc
@@ -28,13 +28,29 @@
//
// Author: sameeragarwal@google.com (Sameer Agarwal)
+#include "ceres/line_search_minimizer.h"
#include "ceres/minimizer.h"
+#include "ceres/trust_region_minimizer.h"
#include "ceres/types.h"
#include "glog/logging.h"
namespace ceres {
namespace internal {
+Minimizer* Minimizer::Create(MinimizerType minimizer_type) {
+ if (minimizer_type == TRUST_REGION) {
+ return new TrustRegionMinimizer;
+ }
+
+ if (minimizer_type == LINE_SEARCH) {
+ return new LineSearchMinimizer;
+ }
+
+ LOG(FATAL) << "Unknown minimizer_type: " << minimizer_type;
+ return NULL;
+}
+
+
Minimizer::~Minimizer() {}
bool Minimizer::RunCallbacks(const Minimizer::Options& options,
diff --git a/internal/ceres/minimizer.h b/internal/ceres/minimizer.h
index f1da3f7..dabf07e 100644
--- a/internal/ceres/minimizer.h
+++ b/internal/ceres/minimizer.h
@@ -41,9 +41,10 @@
namespace internal {
class Evaluator;
-class LinearSolver;
class SparseMatrix;
class TrustRegionStrategy;
+class CoordinateDescentMinimizer;
+class LinearSolver;
// Interface for non-linear least squares solvers.
class Minimizer {
@@ -107,14 +108,10 @@
options.line_search_sufficient_curvature_decrease;
max_line_search_step_expansion =
options.max_line_search_step_expansion;
- is_silent = (options.logging_type == SILENT);
- evaluator = NULL;
- trust_region_strategy = NULL;
- jacobian = NULL;
- callbacks = options.callbacks;
- inner_iteration_minimizer = NULL;
inner_iteration_tolerance = options.inner_iteration_tolerance;
+ is_silent = (options.logging_type == SILENT);
is_constrained = false;
+ callbacks = options.callbacks;
}
int max_num_iterations;
@@ -154,10 +151,14 @@
int max_num_line_search_direction_restarts;
double line_search_sufficient_curvature_decrease;
double max_line_search_step_expansion;
+ double inner_iteration_tolerance;
// If true, then all logging is disabled.
bool is_silent;
+ // Use a bounds constrained optimization algorithm.
+ bool is_constrained;
+
// List of callbacks that are executed by the Minimizer at the end
// of each iteration.
//
@@ -165,27 +166,23 @@
vector<IterationCallback*> callbacks;
// Object responsible for evaluating the cost, residuals and
- // Jacobian matrix. The Options struct does not own this pointer.
- Evaluator* evaluator;
+ // Jacobian matrix.
+ shared_ptr<Evaluator> evaluator;
// Object responsible for actually computing the trust region
- // step, and sizing the trust region radius. The Options struct
- // does not own this pointer.
- TrustRegionStrategy* trust_region_strategy;
+ // step, and sizing the trust region radius.
+ shared_ptr<TrustRegionStrategy> trust_region_strategy;
// Object holding the Jacobian matrix. It is assumed that the
// sparsity structure of the matrix has already been initialized
// and will remain constant for the life time of the
- // optimization. The Options struct does not own this pointer.
- SparseMatrix* jacobian;
+ // optimization.
+ shared_ptr<SparseMatrix> jacobian;
- Minimizer* inner_iteration_minimizer;
- double inner_iteration_tolerance;
-
- // Use a bounds constrained optimization algorithm.
- bool is_constrained;
+ shared_ptr<CoordinateDescentMinimizer> inner_iteration_minimizer;
};
+ static Minimizer* Create(MinimizerType minimizer_type);
static bool RunCallbacks(const Options& options,
const IterationSummary& iteration_summary,
Solver::Summary* summary);
diff --git a/internal/ceres/solver_impl.cc b/internal/ceres/solver_impl.cc
index a1cf4ca..41dbcde 100644
--- a/internal/ceres/solver_impl.cc
+++ b/internal/ceres/solver_impl.cc
@@ -67,8 +67,8 @@
void SolverImpl::TrustRegionMinimize(
const Solver::Options& options,
Program* program,
- CoordinateDescentMinimizer* inner_iteration_minimizer,
- Evaluator* evaluator,
+ shared_ptr<CoordinateDescentMinimizer> inner_iteration_minimizer,
+ shared_ptr<Evaluator> evaluator,
LinearSolver* linear_solver,
Solver::Summary* summary) {
Minimizer::Options minimizer_options(options);
@@ -98,9 +98,7 @@
}
minimizer_options.evaluator = evaluator;
- scoped_ptr<SparseMatrix> jacobian(evaluator->CreateJacobian());
-
- minimizer_options.jacobian = jacobian.get();
+ minimizer_options.jacobian.reset(evaluator->CreateJacobian());
minimizer_options.inner_iteration_minimizer = inner_iteration_minimizer;
TrustRegionStrategy::Options trust_region_strategy_options;
@@ -113,9 +111,7 @@
trust_region_strategy_options.trust_region_strategy_type =
options.trust_region_strategy_type;
trust_region_strategy_options.dogleg_type = options.dogleg_type;
- scoped_ptr<TrustRegionStrategy> strategy(
- TrustRegionStrategy::Create(trust_region_strategy_options));
- minimizer_options.trust_region_strategy = strategy.get();
+ minimizer_options.trust_region_strategy.reset(TrustRegionStrategy::Create(trust_region_strategy_options));
TrustRegionMinimizer minimizer;
double minimizer_start_time = WallTimeInSeconds();
@@ -137,7 +133,7 @@
void SolverImpl::LineSearchMinimize(
const Solver::Options& options,
Program* program,
- Evaluator* evaluator,
+ shared_ptr<Evaluator> evaluator,
Solver::Summary* summary) {
Minimizer::Options minimizer_options(options);
@@ -364,18 +360,18 @@
summary->trust_region_strategy_type = options.trust_region_strategy_type;
summary->dogleg_type = options.dogleg_type;
- scoped_ptr<Evaluator> evaluator(CreateEvaluator(options,
+ shared_ptr<Evaluator> evaluator(CreateEvaluator(options,
problem_impl->parameter_map(),
reduced_program.get(),
&summary->message));
event_logger.AddEvent("CreateEvaluator");
- if (evaluator == NULL) {
+ if (evaluator.get() == NULL) {
return;
}
- scoped_ptr<CoordinateDescentMinimizer> inner_iteration_minimizer;
+ shared_ptr<CoordinateDescentMinimizer> inner_iteration_minimizer;
if (options.use_inner_iterations) {
if (reduced_program->parameter_blocks().size() < 2) {
LOG(WARNING) << "Reduced problem only contains one parameter block."
@@ -386,7 +382,7 @@
*reduced_program,
problem_impl->parameter_map(),
summary));
- if (inner_iteration_minimizer == NULL) {
+ if (inner_iteration_minimizer.get() == NULL) {
LOG(ERROR) << summary->message;
return;
}
@@ -401,8 +397,8 @@
// Run the optimization.
TrustRegionMinimize(options,
reduced_program.get(),
- inner_iteration_minimizer.get(),
- evaluator.get(),
+ inner_iteration_minimizer,
+ evaluator,
linear_solver.get(),
summary);
event_logger.AddEvent("Minimize");
@@ -561,11 +557,11 @@
return;
}
- scoped_ptr<Evaluator> evaluator(CreateEvaluator(options,
+ shared_ptr<Evaluator> evaluator(CreateEvaluator(options,
problem_impl->parameter_map(),
reduced_program.get(),
&summary->message));
- if (evaluator == NULL) {
+ if (evaluator.get() == NULL) {
return;
}
@@ -574,7 +570,7 @@
minimizer_start_time - solver_start_time;
// Run the optimization.
- LineSearchMinimize(options, reduced_program.get(), evaluator.get(), summary);
+ LineSearchMinimize(options, reduced_program.get(), evaluator, summary);
const double post_process_start_time = WallTimeInSeconds();
diff --git a/internal/ceres/solver_impl.h b/internal/ceres/solver_impl.h
index c42c32a..5cefdb1 100644
--- a/internal/ceres/solver_impl.h
+++ b/internal/ceres/solver_impl.h
@@ -64,8 +64,8 @@
static void TrustRegionMinimize(
const Solver::Options &options,
Program* program,
- CoordinateDescentMinimizer* inner_iteration_minimizer,
- Evaluator* evaluator,
+ shared_ptr<CoordinateDescentMinimizer> inner_iteration_minimizer,
+ shared_ptr<Evaluator> evaluator,
LinearSolver* linear_solver,
Solver::Summary* summary);
@@ -76,7 +76,7 @@
// Run the LineSearchMinimizer for the given evaluator and configuration.
static void LineSearchMinimize(const Solver::Options &options,
Program* program,
- Evaluator* evaluator,
+ shared_ptr<Evaluator> evaluator,
Solver::Summary* summary);
// Create the transformed Program, which has all the fixed blocks
diff --git a/internal/ceres/trust_region_minimizer.cc b/internal/ceres/trust_region_minimizer.cc
index 4be5619..52e3b23 100644
--- a/internal/ceres/trust_region_minimizer.cc
+++ b/internal/ceres/trust_region_minimizer.cc
@@ -1,5 +1,5 @@
// Ceres Solver - A fast non-linear least squares minimizer
-// Copyright 2012 Google Inc. All rights reserved.
+// Copyright 2014 Google Inc. All rights reserved.
// http://code.google.com/p/ceres-solver/
//
// Redistribution and use in source and binary forms, with or without
@@ -40,6 +40,7 @@
#include "Eigen/Core"
#include "ceres/array_utils.h"
+#include "ceres/coordinate_descent_minimizer.h"
#include "ceres/evaluator.h"
#include "ceres/file.h"
#include "ceres/internal/eigen.h"
@@ -128,9 +129,9 @@
double iteration_start_time = start_time;
Init(options);
- Evaluator* evaluator = CHECK_NOTNULL(options_.evaluator);
- SparseMatrix* jacobian = CHECK_NOTNULL(options_.jacobian);
- TrustRegionStrategy* strategy = CHECK_NOTNULL(options_.trust_region_strategy);
+ Evaluator* evaluator = CHECK_NOTNULL(options_.evaluator.get());
+ SparseMatrix* jacobian = CHECK_NOTNULL(options_.jacobian.get());
+ TrustRegionStrategy* strategy = CHECK_NOTNULL(options_.trust_region_strategy.get());
const bool is_not_silent = !options.is_silent;
@@ -253,7 +254,8 @@
double candidate_cost = cost;
double accumulated_candidate_model_cost_change = 0.0;
int num_consecutive_invalid_steps = 0;
- bool inner_iterations_are_enabled = options.inner_iteration_minimizer != NULL;
+ bool inner_iterations_are_enabled =
+ options.inner_iteration_minimizer.get() != NULL;
while (true) {
bool inner_iterations_were_useful = false;
if (!RunCallbacks(options, iteration_summary, summary)) {
diff --git a/internal/ceres/trust_region_minimizer_test.cc b/internal/ceres/trust_region_minimizer_test.cc
index ef49206..94be541 100644
--- a/internal/ceres/trust_region_minimizer_test.cc
+++ b/internal/ceres/trust_region_minimizer_test.cc
@@ -223,15 +223,12 @@
parameters[2] = (col3 ? parameters[2] : 0.0);
parameters[3] = (col4 ? parameters[3] : 0.0);
- PowellEvaluator2<col1, col2, col3, col4> powell_evaluator;
- scoped_ptr<SparseMatrix> jacobian(powell_evaluator.CreateJacobian());
-
Minimizer::Options minimizer_options(solver_options);
minimizer_options.gradient_tolerance = 1e-26;
minimizer_options.function_tolerance = 1e-26;
minimizer_options.parameter_tolerance = 1e-26;
- minimizer_options.evaluator = &powell_evaluator;
- minimizer_options.jacobian = jacobian.get();
+ minimizer_options.evaluator.reset(new PowellEvaluator2<col1, col2, col3, col4>);
+ minimizer_options.jacobian.reset(minimizer_options.evaluator->CreateJacobian());
TrustRegionStrategy::Options trust_region_strategy_options;
trust_region_strategy_options.trust_region_strategy_type = strategy_type;
@@ -240,9 +237,8 @@
trust_region_strategy_options.max_radius = 1e20;
trust_region_strategy_options.min_lm_diagonal = 1e-6;
trust_region_strategy_options.max_lm_diagonal = 1e32;
- scoped_ptr<TrustRegionStrategy> strategy(
+ minimizer_options.trust_region_strategy.reset(
TrustRegionStrategy::Create(trust_region_strategy_options));
- minimizer_options.trust_region_strategy = strategy.get();
TrustRegionMinimizer minimizer;
Solver::Summary summary;