1 // Ceres Solver - A fast non-linear least squares minimizer
2 // Copyright 2014 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: sameeragarwal@google.com (Sameer Agarwal)
30 
31 #ifndef CERES_INTERNAL_REORDER_PROGRAM_H_
32 #define CERES_INTERNAL_REORDER_PROGRAM_H_
33 
34 #include <string>
35 #include "ceres/internal/port.h"
36 #include "ceres/parameter_block_ordering.h"
37 #include "ceres/problem_impl.h"
38 #include "ceres/types.h"
39 
40 namespace ceres {
41 namespace internal {
42 
43 class Program;
44 
45 // Reorder the parameter blocks in program using the ordering
46 bool ApplyOrdering(const ProblemImpl::ParameterMap& parameter_map,
47                    const ParameterBlockOrdering& ordering,
48                    Program* program,
49                    string* error);
50 
51 // Reorder the residuals for program, if necessary, so that the residuals
52 // involving each E block occur together. This is a necessary condition for the
53 // Schur eliminator, which works on these "row blocks" in the jacobian.
54 bool LexicographicallyOrderResidualBlocks(int num_eliminate_blocks,
55                                           Program* program,
56                                           string* error);
57 
58 // Schur type solvers require that all parameter blocks eliminated
59 // by the Schur eliminator occur before others and the residuals be
60 // sorted in lexicographic order of their parameter blocks.
61 //
62 // If the parameter_block_ordering only contains one elimination
63 // group then a maximal independent set is computed and used as the
64 // first elimination group, otherwise the user's ordering is used.
65 //
66 // If the linear solver type is SPARSE_SCHUR and support for
67 // constrained fill-reducing ordering is available in the sparse
68 // linear algebra library (SuiteSparse version >= 4.2.0) then
69 // columns of the schur complement matrix are ordered to reduce the
70 // fill-in the Cholesky factorization.
71 //
72 // Upon return, ordering contains the parameter block ordering that
73 // was used to order the program.
74 bool ReorderProgramForSchurTypeLinearSolver(
75     LinearSolverType linear_solver_type,
76     SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type,
77     const ProblemImpl::ParameterMap& parameter_map,
78     ParameterBlockOrdering* parameter_block_ordering,
79     Program* program,
80     string* error);
81 
82 // Sparse cholesky factorization routines when doing the sparse
83 // cholesky factorization of the Jacobian matrix, reorders its
84 // columns to reduce the fill-in. Compute this permutation and
85 // re-order the parameter blocks.
86 //
87 // When using SuiteSparse, if the parameter_block_ordering contains
88 // more than one elimination group and support for constrained
89 // fill-reducing ordering is available in the sparse linear algebra
90 // library (SuiteSparse version >= 4.2.0) then the fill reducing
91 // ordering will take it into account, otherwise it will be ignored.
92 bool ReorderProgramForSparseNormalCholesky(
93     SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type,
94     const ParameterBlockOrdering& parameter_block_ordering,
95     Program* program,
96     string* error);
97 
98 }  // namespace internal
99 }  // namespace ceres
100 
101 #endif  // CERES_INTERNAL_REORDER_PROGRAM_
102