blob: 368a6edd27ec98e037b7d50dc62d94b0baeaef68 [file] [log] [blame]
Austin Schuh70cc9552019-01-21 19:46:48 -08001// Ceres Solver - A fast non-linear least squares minimizer
Austin Schuh3de38b02024-06-25 18:25:10 -07002// Copyright 2023 Google Inc. All rights reserved.
Austin Schuh70cc9552019-01-21 19:46:48 -08003// 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#ifndef CERES_INTERNAL_REORDER_PROGRAM_H_
32#define CERES_INTERNAL_REORDER_PROGRAM_H_
33
34#include <string>
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080035
Austin Schuh3de38b02024-06-25 18:25:10 -070036#include "ceres/internal/disable_warnings.h"
37#include "ceres/internal/export.h"
38#include "ceres/linear_solver.h"
Austin Schuh70cc9552019-01-21 19:46:48 -080039#include "ceres/parameter_block_ordering.h"
40#include "ceres/problem_impl.h"
41#include "ceres/types.h"
42
Austin Schuh3de38b02024-06-25 18:25:10 -070043namespace ceres::internal {
Austin Schuh70cc9552019-01-21 19:46:48 -080044
45class Program;
46
47// Reorder the parameter blocks in program using the ordering
Austin Schuh3de38b02024-06-25 18:25:10 -070048CERES_NO_EXPORT bool ApplyOrdering(
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080049 const ProblemImpl::ParameterMap& parameter_map,
50 const ParameterBlockOrdering& ordering,
51 Program* program,
52 std::string* error);
Austin Schuh70cc9552019-01-21 19:46:48 -080053
54// Reorder the residuals for program, if necessary, so that the residuals
55// involving each E block occur together. This is a necessary condition for the
56// Schur eliminator, which works on these "row blocks" in the jacobian.
Austin Schuh3de38b02024-06-25 18:25:10 -070057CERES_NO_EXPORT bool LexicographicallyOrderResidualBlocks(
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080058 int size_of_first_elimination_group, Program* program, std::string* error);
Austin Schuh70cc9552019-01-21 19:46:48 -080059
60// Schur type solvers require that all parameter blocks eliminated
61// by the Schur eliminator occur before others and the residuals be
62// sorted in lexicographic order of their parameter blocks.
63//
64// If the parameter_block_ordering only contains one elimination
65// group then a maximal independent set is computed and used as the
66// first elimination group, otherwise the user's ordering is used.
67//
68// If the linear solver type is SPARSE_SCHUR and support for
69// constrained fill-reducing ordering is available in the sparse
70// linear algebra library (SuiteSparse version >= 4.2.0) then
71// columns of the schur complement matrix are ordered to reduce the
72// fill-in the Cholesky factorization.
73//
74// Upon return, ordering contains the parameter block ordering that
75// was used to order the program.
Austin Schuh3de38b02024-06-25 18:25:10 -070076CERES_NO_EXPORT bool ReorderProgramForSchurTypeLinearSolver(
Austin Schuh70cc9552019-01-21 19:46:48 -080077 LinearSolverType linear_solver_type,
78 SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type,
Austin Schuh3de38b02024-06-25 18:25:10 -070079 LinearSolverOrderingType linear_solver_ordering_type,
Austin Schuh70cc9552019-01-21 19:46:48 -080080 const ProblemImpl::ParameterMap& parameter_map,
81 ParameterBlockOrdering* parameter_block_ordering,
82 Program* program,
83 std::string* error);
84
85// Sparse cholesky factorization routines when doing the sparse
86// cholesky factorization of the Jacobian matrix, reorders its
87// columns to reduce the fill-in. Compute this permutation and
88// re-order the parameter blocks.
89//
90// When using SuiteSparse, if the parameter_block_ordering contains
91// more than one elimination group and support for constrained
92// fill-reducing ordering is available in the sparse linear algebra
93// library (SuiteSparse version >= 4.2.0) then the fill reducing
94// ordering will take it into account, otherwise it will be ignored.
Austin Schuh3de38b02024-06-25 18:25:10 -070095CERES_NO_EXPORT bool ReorderProgramForSparseCholesky(
Austin Schuh70cc9552019-01-21 19:46:48 -080096 SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type,
Austin Schuh3de38b02024-06-25 18:25:10 -070097 LinearSolverOrderingType linear_solver_ordering_type,
Austin Schuh70cc9552019-01-21 19:46:48 -080098 const ParameterBlockOrdering& parameter_block_ordering,
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080099 int start_row_block,
Austin Schuh70cc9552019-01-21 19:46:48 -0800100 Program* program,
101 std::string* error);
102
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800103// Reorder the residual blocks in the program so that all the residual
104// blocks in bottom_residual_blocks are at the bottom. The return
105// value is the number of residual blocks in the program in "top" part
106// of the Program, i.e., the ones not included in
107// bottom_residual_blocks.
108//
109// This number can be different from program->NumResidualBlocks() -
110// bottom_residual_blocks.size() because we allow
111// bottom_residual_blocks to contain residual blocks not present in
112// the Program.
Austin Schuh3de38b02024-06-25 18:25:10 -0700113CERES_NO_EXPORT int ReorderResidualBlocksByPartition(
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800114 const std::unordered_set<ResidualBlockId>& bottom_residual_blocks,
115 Program* program);
116
Austin Schuh3de38b02024-06-25 18:25:10 -0700117// The return value of this function indicates whether the columns of
118// the Jacobian can be reordered using a fill reducing ordering.
119CERES_NO_EXPORT bool AreJacobianColumnsOrdered(
120 LinearSolverType linear_solver_type,
121 PreconditionerType preconditioner_type,
122 SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type,
123 LinearSolverOrderingType linear_solver_ordering_type);
Austin Schuh70cc9552019-01-21 19:46:48 -0800124
Austin Schuh3de38b02024-06-25 18:25:10 -0700125} // namespace ceres::internal
126
127#include "ceres/internal/reenable_warnings.h"
128
129#endif // CERES_INTERNAL_REORDER_PROGRAM_H_