Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame^] | 1 | // Ceres Solver - A fast non-linear least squares minimizer |
| 2 | // Copyright 2015 Google Inc. All rights reserved. |
| 3 | // http://ceres-solver.org/ |
| 4 | // |
| 5 | // Redistribution and use in source and binary forms, with or without |
| 6 | // modification, are permitted provided that the following conditions are met: |
| 7 | // |
| 8 | // * Redistributions of source code must retain the above copyright notice, |
| 9 | // this list of conditions and the following disclaimer. |
| 10 | // * Redistributions in binary form must reproduce the above copyright notice, |
| 11 | // this list of conditions and the following disclaimer in the documentation |
| 12 | // and/or other materials provided with the distribution. |
| 13 | // * Neither the name of Google Inc. nor the names of its contributors may be |
| 14 | // used to endorse or promote products derived from this software without |
| 15 | // specific prior written permission. |
| 16 | // |
| 17 | // THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" |
| 18 | // AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE |
| 19 | // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE |
| 20 | // ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE |
| 21 | // LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR |
| 22 | // CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF |
| 23 | // SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS |
| 24 | // INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN |
| 25 | // CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) |
| 26 | // ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE |
| 27 | // POSSIBILITY OF SUCH DAMAGE. |
| 28 | // |
| 29 | // Author: sameeragarwal@google.com (Sameer Agarwal) |
| 30 | |
| 31 | #include "ceres/trust_region_preprocessor.h" |
| 32 | |
| 33 | #include <numeric> |
| 34 | #include <string> |
| 35 | #include "ceres/callbacks.h" |
| 36 | #include "ceres/context_impl.h" |
| 37 | #include "ceres/evaluator.h" |
| 38 | #include "ceres/linear_solver.h" |
| 39 | #include "ceres/minimizer.h" |
| 40 | #include "ceres/parameter_block.h" |
| 41 | #include "ceres/preconditioner.h" |
| 42 | #include "ceres/preprocessor.h" |
| 43 | #include "ceres/problem_impl.h" |
| 44 | #include "ceres/program.h" |
| 45 | #include "ceres/reorder_program.h" |
| 46 | #include "ceres/suitesparse.h" |
| 47 | #include "ceres/trust_region_strategy.h" |
| 48 | #include "ceres/wall_time.h" |
| 49 | |
| 50 | namespace ceres { |
| 51 | namespace internal { |
| 52 | |
| 53 | using std::vector; |
| 54 | |
| 55 | namespace { |
| 56 | |
| 57 | ParameterBlockOrdering* CreateDefaultLinearSolverOrdering( |
| 58 | const Program& program) { |
| 59 | ParameterBlockOrdering* ordering = new ParameterBlockOrdering; |
| 60 | const vector<ParameterBlock*>& parameter_blocks = |
| 61 | program.parameter_blocks(); |
| 62 | for (int i = 0; i < parameter_blocks.size(); ++i) { |
| 63 | ordering->AddElementToGroup( |
| 64 | const_cast<double*>(parameter_blocks[i]->user_state()), 0); |
| 65 | } |
| 66 | return ordering; |
| 67 | } |
| 68 | |
| 69 | // Check if all the user supplied values in the parameter blocks are |
| 70 | // sane or not, and if the program is feasible or not. |
| 71 | bool IsProgramValid(const Program& program, std::string* error) { |
| 72 | return (program.ParameterBlocksAreFinite(error) && |
| 73 | program.IsFeasible(error)); |
| 74 | } |
| 75 | |
| 76 | void AlternateLinearSolverAndPreconditionerForSchurTypeLinearSolver( |
| 77 | Solver::Options* options) { |
| 78 | if (!IsSchurType(options->linear_solver_type)) { |
| 79 | return; |
| 80 | } |
| 81 | |
| 82 | const LinearSolverType linear_solver_type_given = options->linear_solver_type; |
| 83 | const PreconditionerType preconditioner_type_given = |
| 84 | options->preconditioner_type; |
| 85 | options->linear_solver_type = LinearSolver::LinearSolverForZeroEBlocks( |
| 86 | linear_solver_type_given); |
| 87 | |
| 88 | std::string message; |
| 89 | if (linear_solver_type_given == ITERATIVE_SCHUR) { |
| 90 | options->preconditioner_type = Preconditioner::PreconditionerForZeroEBlocks( |
| 91 | preconditioner_type_given); |
| 92 | |
| 93 | message = |
| 94 | StringPrintf( |
| 95 | "No E blocks. Switching from %s(%s) to %s(%s).", |
| 96 | LinearSolverTypeToString(linear_solver_type_given), |
| 97 | PreconditionerTypeToString(preconditioner_type_given), |
| 98 | LinearSolverTypeToString(options->linear_solver_type), |
| 99 | PreconditionerTypeToString(options->preconditioner_type)); |
| 100 | } else { |
| 101 | message = |
| 102 | StringPrintf( |
| 103 | "No E blocks. Switching from %s to %s.", |
| 104 | LinearSolverTypeToString(linear_solver_type_given), |
| 105 | LinearSolverTypeToString(options->linear_solver_type)); |
| 106 | } |
| 107 | |
| 108 | VLOG_IF(1, options->logging_type != SILENT) << message; |
| 109 | } |
| 110 | |
| 111 | // For Schur type and SPARSE_NORMAL_CHOLESKY linear solvers, reorder |
| 112 | // the program to reduce fill-in and increase cache coherency. |
| 113 | bool ReorderProgram(PreprocessedProblem* pp) { |
| 114 | const Solver::Options& options = pp->options; |
| 115 | if (IsSchurType(options.linear_solver_type)) { |
| 116 | return ReorderProgramForSchurTypeLinearSolver( |
| 117 | options.linear_solver_type, |
| 118 | options.sparse_linear_algebra_library_type, |
| 119 | pp->problem->parameter_map(), |
| 120 | options.linear_solver_ordering.get(), |
| 121 | pp->reduced_program.get(), |
| 122 | &pp->error); |
| 123 | } |
| 124 | |
| 125 | |
| 126 | if (options.linear_solver_type == SPARSE_NORMAL_CHOLESKY && |
| 127 | !options.dynamic_sparsity) { |
| 128 | return ReorderProgramForSparseNormalCholesky( |
| 129 | options.sparse_linear_algebra_library_type, |
| 130 | *options.linear_solver_ordering, |
| 131 | pp->reduced_program.get(), |
| 132 | &pp->error); |
| 133 | } |
| 134 | |
| 135 | return true; |
| 136 | } |
| 137 | |
| 138 | // Configure and create a linear solver object. In doing so, if a |
| 139 | // sparse direct factorization based linear solver is being used, then |
| 140 | // find a fill reducing ordering and reorder the program as needed |
| 141 | // too. |
| 142 | bool SetupLinearSolver(PreprocessedProblem* pp) { |
| 143 | Solver::Options& options = pp->options; |
| 144 | if (options.linear_solver_ordering.get() == NULL) { |
| 145 | // If the user has not supplied a linear solver ordering, then we |
| 146 | // assume that they are giving all the freedom to us in choosing |
| 147 | // the best possible ordering. This intent can be indicated by |
| 148 | // putting all the parameter blocks in the same elimination group. |
| 149 | options.linear_solver_ordering.reset( |
| 150 | CreateDefaultLinearSolverOrdering(*pp->reduced_program)); |
| 151 | } else { |
| 152 | // If the user supplied an ordering, then check if the first |
| 153 | // elimination group is still non-empty after the reduced problem |
| 154 | // has been constructed. |
| 155 | // |
| 156 | // This is important for Schur type linear solvers, where the |
| 157 | // first elimination group is special -- it needs to be an |
| 158 | // independent set. |
| 159 | // |
| 160 | // If the first elimination group is empty, then we cannot use the |
| 161 | // user's requested linear solver (and a preconditioner as the |
| 162 | // case may be) so we must use a different one. |
| 163 | ParameterBlockOrdering* ordering = options.linear_solver_ordering.get(); |
| 164 | const int min_group_id = ordering->MinNonZeroGroup(); |
| 165 | ordering->Remove(pp->removed_parameter_blocks); |
| 166 | if (IsSchurType(options.linear_solver_type) && |
| 167 | min_group_id != ordering->MinNonZeroGroup()) { |
| 168 | AlternateLinearSolverAndPreconditionerForSchurTypeLinearSolver( |
| 169 | &options); |
| 170 | } |
| 171 | } |
| 172 | |
| 173 | // Reorder the program to reduce fill in and improve cache coherency |
| 174 | // of the Jacobian. |
| 175 | if (!ReorderProgram(pp)) { |
| 176 | return false; |
| 177 | } |
| 178 | |
| 179 | // Configure the linear solver. |
| 180 | pp->linear_solver_options = LinearSolver::Options(); |
| 181 | pp->linear_solver_options.min_num_iterations = |
| 182 | options.min_linear_solver_iterations; |
| 183 | pp->linear_solver_options.max_num_iterations = |
| 184 | options.max_linear_solver_iterations; |
| 185 | pp->linear_solver_options.type = options.linear_solver_type; |
| 186 | pp->linear_solver_options.preconditioner_type = options.preconditioner_type; |
| 187 | pp->linear_solver_options.visibility_clustering_type = |
| 188 | options.visibility_clustering_type; |
| 189 | pp->linear_solver_options.sparse_linear_algebra_library_type = |
| 190 | options.sparse_linear_algebra_library_type; |
| 191 | pp->linear_solver_options.dense_linear_algebra_library_type = |
| 192 | options.dense_linear_algebra_library_type; |
| 193 | pp->linear_solver_options.use_explicit_schur_complement = |
| 194 | options.use_explicit_schur_complement; |
| 195 | pp->linear_solver_options.dynamic_sparsity = options.dynamic_sparsity; |
| 196 | pp->linear_solver_options.use_mixed_precision_solves = |
| 197 | options.use_mixed_precision_solves; |
| 198 | pp->linear_solver_options.max_num_refinement_iterations = |
| 199 | options.max_num_refinement_iterations; |
| 200 | pp->linear_solver_options.num_threads = options.num_threads; |
| 201 | pp->linear_solver_options.use_postordering = options.use_postordering; |
| 202 | pp->linear_solver_options.context = pp->problem->context(); |
| 203 | |
| 204 | if (IsSchurType(pp->linear_solver_options.type)) { |
| 205 | OrderingToGroupSizes(options.linear_solver_ordering.get(), |
| 206 | &pp->linear_solver_options.elimination_groups); |
| 207 | |
| 208 | // Schur type solvers expect at least two elimination groups. If |
| 209 | // there is only one elimination group, then it is guaranteed that |
| 210 | // this group only contains e_blocks. Thus we add a dummy |
| 211 | // elimination group with zero blocks in it. |
| 212 | if (pp->linear_solver_options.elimination_groups.size() == 1) { |
| 213 | pp->linear_solver_options.elimination_groups.push_back(0); |
| 214 | } |
| 215 | |
| 216 | if (options.linear_solver_type == SPARSE_SCHUR) { |
| 217 | // When using SPARSE_SCHUR, we ignore the user's postordering |
| 218 | // preferences in certain cases. |
| 219 | // |
| 220 | // 1. SUITE_SPARSE is the sparse linear algebra library requested |
| 221 | // but cholmod_camd is not available. |
| 222 | // 2. CX_SPARSE is the sparse linear algebra library requested. |
| 223 | // |
| 224 | // This ensures that the linear solver does not assume that a |
| 225 | // fill-reducing pre-ordering has been done. |
| 226 | // |
| 227 | // TODO(sameeragarwal): Implement the reordering of parameter |
| 228 | // blocks for CX_SPARSE. |
| 229 | if ((options.sparse_linear_algebra_library_type == SUITE_SPARSE && |
| 230 | !SuiteSparse:: |
| 231 | IsConstrainedApproximateMinimumDegreeOrderingAvailable()) || |
| 232 | (options.sparse_linear_algebra_library_type == CX_SPARSE)) { |
| 233 | pp->linear_solver_options.use_postordering = true; |
| 234 | } |
| 235 | } |
| 236 | } |
| 237 | |
| 238 | pp->linear_solver.reset(LinearSolver::Create(pp->linear_solver_options)); |
| 239 | return (pp->linear_solver.get() != NULL); |
| 240 | } |
| 241 | |
| 242 | // Configure and create the evaluator. |
| 243 | bool SetupEvaluator(PreprocessedProblem* pp) { |
| 244 | const Solver::Options& options = pp->options; |
| 245 | pp->evaluator_options = Evaluator::Options(); |
| 246 | pp->evaluator_options.linear_solver_type = options.linear_solver_type; |
| 247 | pp->evaluator_options.num_eliminate_blocks = 0; |
| 248 | if (IsSchurType(options.linear_solver_type)) { |
| 249 | pp->evaluator_options.num_eliminate_blocks = |
| 250 | options |
| 251 | .linear_solver_ordering |
| 252 | ->group_to_elements().begin() |
| 253 | ->second.size(); |
| 254 | } |
| 255 | |
| 256 | pp->evaluator_options.num_threads = options.num_threads; |
| 257 | pp->evaluator_options.dynamic_sparsity = options.dynamic_sparsity; |
| 258 | pp->evaluator_options.context = pp->problem->context(); |
| 259 | pp->evaluator_options.evaluation_callback = options.evaluation_callback; |
| 260 | pp->evaluator.reset(Evaluator::Create(pp->evaluator_options, |
| 261 | pp->reduced_program.get(), |
| 262 | &pp->error)); |
| 263 | |
| 264 | return (pp->evaluator.get() != NULL); |
| 265 | } |
| 266 | |
| 267 | // If the user requested inner iterations, then find an inner |
| 268 | // iteration ordering as needed and configure and create a |
| 269 | // CoordinateDescentMinimizer object to perform the inner iterations. |
| 270 | bool SetupInnerIterationMinimizer(PreprocessedProblem* pp) { |
| 271 | Solver::Options& options = pp->options; |
| 272 | if (!options.use_inner_iterations) { |
| 273 | return true; |
| 274 | } |
| 275 | |
| 276 | // With just one parameter block, the outer iteration of the trust |
| 277 | // region method and inner iterations are doing exactly the same |
| 278 | // thing, and thus inner iterations are not needed. |
| 279 | if (pp->reduced_program->NumParameterBlocks() == 1) { |
| 280 | LOG(WARNING) << "Reduced problem only contains one parameter block." |
| 281 | << "Disabling inner iterations."; |
| 282 | return true; |
| 283 | } |
| 284 | |
| 285 | if (options.inner_iteration_ordering.get() != NULL) { |
| 286 | // If the user supplied an ordering, then remove the set of |
| 287 | // inactive parameter blocks from it |
| 288 | options.inner_iteration_ordering->Remove(pp->removed_parameter_blocks); |
| 289 | if (options.inner_iteration_ordering->NumElements() == 0) { |
| 290 | LOG(WARNING) << "No remaining elements in the inner iteration ordering."; |
| 291 | return true; |
| 292 | } |
| 293 | |
| 294 | // Validate the reduced ordering. |
| 295 | if (!CoordinateDescentMinimizer::IsOrderingValid( |
| 296 | *pp->reduced_program, |
| 297 | *options.inner_iteration_ordering, |
| 298 | &pp->error)) { |
| 299 | return false; |
| 300 | } |
| 301 | } else { |
| 302 | // The user did not supply an ordering, so create one. |
| 303 | options.inner_iteration_ordering.reset( |
| 304 | CoordinateDescentMinimizer::CreateOrdering(*pp->reduced_program)); |
| 305 | } |
| 306 | |
| 307 | pp->inner_iteration_minimizer.reset( |
| 308 | new CoordinateDescentMinimizer(pp->problem->context())); |
| 309 | return pp->inner_iteration_minimizer->Init(*pp->reduced_program, |
| 310 | pp->problem->parameter_map(), |
| 311 | *options.inner_iteration_ordering, |
| 312 | &pp->error); |
| 313 | } |
| 314 | |
| 315 | // Configure and create a TrustRegionMinimizer object. |
| 316 | void SetupMinimizerOptions(PreprocessedProblem* pp) { |
| 317 | const Solver::Options& options = pp->options; |
| 318 | |
| 319 | SetupCommonMinimizerOptions(pp); |
| 320 | pp->minimizer_options.is_constrained = |
| 321 | pp->reduced_program->IsBoundsConstrained(); |
| 322 | pp->minimizer_options.jacobian.reset(pp->evaluator->CreateJacobian()); |
| 323 | pp->minimizer_options.inner_iteration_minimizer = |
| 324 | pp->inner_iteration_minimizer; |
| 325 | |
| 326 | TrustRegionStrategy::Options strategy_options; |
| 327 | strategy_options.linear_solver = pp->linear_solver.get(); |
| 328 | strategy_options.initial_radius = |
| 329 | options.initial_trust_region_radius; |
| 330 | strategy_options.max_radius = options.max_trust_region_radius; |
| 331 | strategy_options.min_lm_diagonal = options.min_lm_diagonal; |
| 332 | strategy_options.max_lm_diagonal = options.max_lm_diagonal; |
| 333 | strategy_options.trust_region_strategy_type = |
| 334 | options.trust_region_strategy_type; |
| 335 | strategy_options.dogleg_type = options.dogleg_type; |
| 336 | pp->minimizer_options.trust_region_strategy.reset( |
| 337 | TrustRegionStrategy::Create(strategy_options)); |
| 338 | CHECK(pp->minimizer_options.trust_region_strategy != nullptr); |
| 339 | } |
| 340 | |
| 341 | } // namespace |
| 342 | |
| 343 | TrustRegionPreprocessor::~TrustRegionPreprocessor() { |
| 344 | } |
| 345 | |
| 346 | bool TrustRegionPreprocessor::Preprocess(const Solver::Options& options, |
| 347 | ProblemImpl* problem, |
| 348 | PreprocessedProblem* pp) { |
| 349 | CHECK(pp != nullptr); |
| 350 | pp->options = options; |
| 351 | ChangeNumThreadsIfNeeded(&pp->options); |
| 352 | |
| 353 | pp->problem = problem; |
| 354 | Program* program = problem->mutable_program(); |
| 355 | if (!IsProgramValid(*program, &pp->error)) { |
| 356 | return false; |
| 357 | } |
| 358 | |
| 359 | pp->reduced_program.reset( |
| 360 | program->CreateReducedProgram(&pp->removed_parameter_blocks, |
| 361 | &pp->fixed_cost, |
| 362 | &pp->error)); |
| 363 | |
| 364 | if (pp->reduced_program.get() == NULL) { |
| 365 | return false; |
| 366 | } |
| 367 | |
| 368 | if (pp->reduced_program->NumParameterBlocks() == 0) { |
| 369 | // The reduced problem has no parameter or residual blocks. There |
| 370 | // is nothing more to do. |
| 371 | return true; |
| 372 | } |
| 373 | |
| 374 | if (!SetupLinearSolver(pp) || |
| 375 | !SetupEvaluator(pp) || |
| 376 | !SetupInnerIterationMinimizer(pp)) { |
| 377 | return false; |
| 378 | } |
| 379 | |
| 380 | SetupMinimizerOptions(pp); |
| 381 | return true; |
| 382 | } |
| 383 | |
| 384 | } // namespace internal |
| 385 | } // namespace ceres |