1 // 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 (at) google.com (Keir Mierle) 30 31 #ifndef CERES_INTERNAL_SOLVER_IMPL_H_ 32 #define CERES_INTERNAL_SOLVER_IMPL_H_ 33 34 #include <set> 35 #include <string> 36 #include <vector> 37 #include "ceres/internal/port.h" 38 #include "ceres/ordered_groups.h" 39 #include "ceres/problem_impl.h" 40 #include "ceres/solver.h" 41 42 namespace ceres { 43 namespace internal { 44 45 class CoordinateDescentMinimizer; 46 class Evaluator; 47 class LinearSolver; 48 class Program; 49 class TripletSparseMatrix; 50 51 class SolverImpl { 52 public: 53 // Mirrors the interface in solver.h, but exposes implementation 54 // details for testing internally. 55 static void Solve(const Solver::Options& options, 56 ProblemImpl* problem_impl, 57 Solver::Summary* summary); 58 59 static void TrustRegionSolve(const Solver::Options& options, 60 ProblemImpl* problem_impl, 61 Solver::Summary* summary); 62 63 // Run the TrustRegionMinimizer for the given evaluator and configuration. 64 static void TrustRegionMinimize( 65 const Solver::Options &options, 66 Program* program, 67 CoordinateDescentMinimizer* inner_iteration_minimizer, 68 Evaluator* evaluator, 69 LinearSolver* linear_solver, 70 double* parameters, 71 Solver::Summary* summary); 72 73 #ifndef CERES_NO_LINE_SEARCH_MINIMIZER 74 static void LineSearchSolve(const Solver::Options& options, 75 ProblemImpl* problem_impl, 76 Solver::Summary* summary); 77 78 // Run the LineSearchMinimizer for the given evaluator and configuration. 79 static void LineSearchMinimize(const Solver::Options &options, 80 Program* program, 81 Evaluator* evaluator, 82 double* parameters, 83 Solver::Summary* summary); 84 #endif // CERES_NO_LINE_SEARCH_MINIMIZER 85 86 // Create the transformed Program, which has all the fixed blocks 87 // and residuals eliminated, and in the case of automatic schur 88 // ordering, has the E blocks first in the resulting program, with 89 // options.num_eliminate_blocks set appropriately. 90 // 91 // If fixed_cost is not NULL, the residual blocks that are removed 92 // are evaluated and the sum of their cost is returned in fixed_cost. 93 static Program* CreateReducedProgram(Solver::Options* options, 94 ProblemImpl* problem_impl, 95 double* fixed_cost, 96 string* error); 97 98 // Create the appropriate linear solver, taking into account any 99 // config changes decided by CreateTransformedProgram(). The 100 // selected linear solver, which may be different from what the user 101 // selected; consider the case that the remaining elimininated 102 // blocks is zero after removing fixed blocks. 103 static LinearSolver* CreateLinearSolver(Solver::Options* options, 104 string* error); 105 106 // Reorder the residuals for program, if necessary, so that the 107 // residuals involving e block (i.e., the first num_eliminate_block 108 // parameter blocks) occur together. This is a necessary condition 109 // for the Schur eliminator. 110 static bool LexicographicallyOrderResidualBlocks( 111 const int num_eliminate_blocks, 112 Program* program, 113 string* error); 114 115 // Create the appropriate evaluator for the transformed program. 116 static Evaluator* CreateEvaluator( 117 const Solver::Options& options, 118 const ProblemImpl::ParameterMap& parameter_map, 119 Program* program, 120 string* error); 121 122 // Remove the fixed or unused parameter blocks and residuals 123 // depending only on fixed parameters from the problem. Also updates 124 // num_eliminate_blocks, since removed parameters changes the point 125 // at which the eliminated blocks is valid. If fixed_cost is not 126 // NULL, the residual blocks that are removed are evaluated and the 127 // sum of their cost is returned in fixed_cost. 128 static bool RemoveFixedBlocksFromProgram(Program* program, 129 ParameterBlockOrdering* ordering, 130 double* fixed_cost, 131 string* error); 132 133 static bool IsOrderingValid(const Solver::Options& options, 134 const ProblemImpl* problem_impl, 135 string* error); 136 137 static bool IsParameterBlockSetIndependent( 138 const set<double*>& parameter_block_ptrs, 139 const vector<ResidualBlock*>& residual_blocks); 140 141 static CoordinateDescentMinimizer* CreateInnerIterationMinimizer( 142 const Solver::Options& options, 143 const Program& program, 144 const ProblemImpl::ParameterMap& parameter_map, 145 Solver::Summary* summary); 146 147 // If the linear solver is of Schur type, then replace it with the 148 // closest equivalent linear solver. This is done when the user 149 // requested a Schur type solver but the problem structure makes it 150 // impossible to use one. 151 // 152 // If the linear solver is not of Schur type, the function is a 153 // no-op. 154 static void AlternateLinearSolverForSchurTypeLinearSolver( 155 Solver::Options* options); 156 157 // Create a TripletSparseMatrix which contains the zero-one 158 // structure corresponding to the block sparsity of the transpose of 159 // the Jacobian matrix. 160 // 161 // Caller owns the result. 162 static TripletSparseMatrix* CreateJacobianBlockSparsityTranspose( 163 const Program* program); 164 165 // Reorder the parameter blocks in program using the ordering 166 static bool ApplyUserOrdering( 167 const ProblemImpl::ParameterMap& parameter_map, 168 const ParameterBlockOrdering* parameter_block_ordering, 169 Program* program, 170 string* error); 171 172 // Sparse cholesky factorization routines when doing the sparse 173 // cholesky factorization of the Jacobian matrix, reorders its 174 // columns to reduce the fill-in. Compute this permutation and 175 // re-order the parameter blocks. 176 // 177 // If the parameter_block_ordering contains more than one 178 // elimination group and support for constrained fill-reducing 179 // ordering is available in the sparse linear algebra library 180 // (SuiteSparse version >= 4.2.0) then the fill reducing 181 // ordering will take it into account, otherwise it will be ignored. 182 static bool ReorderProgramForSparseNormalCholesky( 183 const SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type, 184 const ParameterBlockOrdering* parameter_block_ordering, 185 Program* program, 186 string* error); 187 188 // Schur type solvers require that all parameter blocks eliminated 189 // by the Schur eliminator occur before others and the residuals be 190 // sorted in lexicographic order of their parameter blocks. 191 // 192 // If the parameter_block_ordering only contains one elimination 193 // group then a maximal independent set is computed and used as the 194 // first elimination group, otherwise the user's ordering is used. 195 // 196 // If the linear solver type is SPARSE_SCHUR and support for 197 // constrained fill-reducing ordering is available in the sparse 198 // linear algebra library (SuiteSparse version >= 4.2.0) then 199 // columns of the schur complement matrix are ordered to reduce the 200 // fill-in the Cholesky factorization. 201 // 202 // Upon return, ordering contains the parameter block ordering that 203 // was used to order the program. 204 static bool ReorderProgramForSchurTypeLinearSolver( 205 const LinearSolverType linear_solver_type, 206 const SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type, 207 const ProblemImpl::ParameterMap& parameter_map, 208 ParameterBlockOrdering* parameter_block_ordering, 209 Program* program, 210 string* error); 211 212 // array contains a list of (possibly repeating) non-negative 213 // integers. Let us assume that we have constructed another array 214 // `p` by sorting and uniqueing the entries of array. 215 // CompactifyArray replaces each entry in "array" with its position 216 // in `p`. 217 static void CompactifyArray(vector<int>* array); 218 }; 219 220 } // namespace internal 221 } // namespace ceres 222 223 #endif // CERES_INTERNAL_SOLVER_IMPL_H_ 224