program.h 8.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201
  1. // Ceres Solver - A fast non-linear least squares minimizer
  2. // Copyright 2023 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: keir@google.com (Keir Mierle)
  30. #ifndef CERES_INTERNAL_PROGRAM_H_
  31. #define CERES_INTERNAL_PROGRAM_H_
  32. #include <memory>
  33. #include <set>
  34. #include <string>
  35. #include <vector>
  36. #include "ceres/evaluation_callback.h"
  37. #include "ceres/internal/disable_warnings.h"
  38. #include "ceres/internal/export.h"
  39. namespace ceres::internal {
  40. class ParameterBlock;
  41. class ProblemImpl;
  42. class ResidualBlock;
  43. class TripletSparseMatrix;
  44. class ContextImpl;
  45. // A nonlinear least squares optimization problem. This is different from the
  46. // similarly-named "Problem" object, which offers a mutation interface for
  47. // adding and modifying parameters and residuals. The Program contains the core
  48. // part of the Problem, which is the parameters and the residuals, stored in a
  49. // particular ordering. The ordering is critical, since it defines the mapping
  50. // between (residual, parameter) pairs and a position in the jacobian of the
  51. // objective function. Various parts of Ceres transform one Program into
  52. // another; for example, the first stage of solving involves stripping all
  53. // constant parameters and residuals. This is in contrast with Problem, which is
  54. // not built for transformation.
  55. class CERES_NO_EXPORT Program {
  56. public:
  57. // The ordered parameter and residual blocks for the program.
  58. const std::vector<ParameterBlock*>& parameter_blocks() const;
  59. const std::vector<ResidualBlock*>& residual_blocks() const;
  60. std::vector<ParameterBlock*>* mutable_parameter_blocks();
  61. std::vector<ResidualBlock*>* mutable_residual_blocks();
  62. EvaluationCallback* mutable_evaluation_callback();
  63. // Serialize to/from the program and update states.
  64. //
  65. // NOTE: Setting the state of a parameter block can trigger the
  66. // computation of the Jacobian of its manifold. If this computation fails for
  67. // some reason, then this method returns false and the state of the parameter
  68. // blocks cannot be trusted.
  69. bool StateVectorToParameterBlocks(const double* state);
  70. void ParameterBlocksToStateVector(double* state) const;
  71. // Copy internal state to the user's parameters.
  72. void CopyParameterBlockStateToUserState();
  73. // Set the parameter block pointers to the user pointers. Since this
  74. // runs parameter block set state internally, which may call manifold, this
  75. // can fail. False is returned on failure.
  76. bool SetParameterBlockStatePtrsToUserStatePtrs();
  77. // Update a state vector for the program given a delta.
  78. bool Plus(const double* state,
  79. const double* delta,
  80. double* state_plus_delta,
  81. ContextImpl* context,
  82. int num_threads) const;
  83. // Set the parameter indices and offsets. This permits mapping backward
  84. // from a ParameterBlock* to an index in the parameter_blocks() vector. For
  85. // any parameter block p, after calling SetParameterOffsetsAndIndex(), it
  86. // is true that
  87. //
  88. // parameter_blocks()[p->index()] == p
  89. //
  90. // If a parameter appears in a residual but not in the parameter block, then
  91. // it will have an index of -1.
  92. //
  93. // This also updates p->state_offset() and p->delta_offset(), which are the
  94. // position of the parameter in the state and delta vector respectively.
  95. void SetParameterOffsetsAndIndex();
  96. // Check if the internal state of the program (the indexing and the
  97. // offsets) are correct.
  98. bool IsValid() const;
  99. bool ParameterBlocksAreFinite(std::string* message) const;
  100. // Returns true if the program has any non-constant parameter blocks
  101. // which have non-trivial bounds constraints.
  102. bool IsBoundsConstrained() const;
  103. // Returns false, if the program has any constant parameter blocks
  104. // which are not feasible, or any variable parameter blocks which
  105. // have a lower bound greater than or equal to the upper bound.
  106. bool IsFeasible(std::string* message) const;
  107. // Loop over each residual block and ensure that no two parameter
  108. // blocks in the same residual block are part of
  109. // parameter_blocks as that would violate the assumption that it
  110. // is an independent set in the Hessian matrix.
  111. bool IsParameterBlockSetIndependent(
  112. const std::set<double*>& independent_set) const;
  113. // Create a TripletSparseMatrix which contains the zero-one
  114. // structure corresponding to the block sparsity of the transpose of
  115. // the Jacobian matrix.
  116. //
  117. // start_residual_block which allows the user to ignore the first
  118. // start_residual_block residuals.
  119. std::unique_ptr<TripletSparseMatrix> CreateJacobianBlockSparsityTranspose(
  120. int start_residual_block = 0) const;
  121. // Create a copy of this program and removes constant parameter
  122. // blocks and residual blocks with no varying parameter blocks while
  123. // preserving their relative order.
  124. //
  125. // removed_parameter_blocks on exit will contain the list of
  126. // parameter blocks that were removed.
  127. //
  128. // fixed_cost will be equal to the sum of the costs of the residual
  129. // blocks that were removed.
  130. //
  131. // If there was a problem, then the function will return a nullptr
  132. // pointer and error will contain a human readable description of
  133. // the problem.
  134. std::unique_ptr<Program> CreateReducedProgram(
  135. std::vector<double*>* removed_parameter_blocks,
  136. double* fixed_cost,
  137. std::string* error) const;
  138. // See problem.h for what these do.
  139. int NumParameterBlocks() const;
  140. int NumParameters() const;
  141. int NumEffectiveParameters() const;
  142. int NumResidualBlocks() const;
  143. int NumResiduals() const;
  144. int MaxScratchDoublesNeededForEvaluate() const;
  145. int MaxDerivativesPerResidualBlock() const;
  146. int MaxParametersPerResidualBlock() const;
  147. int MaxResidualsPerResidualBlock() const;
  148. // A human-readable dump of the parameter blocks for debugging.
  149. // TODO(keir): If necessary, also dump the residual blocks.
  150. std::string ToString() const;
  151. private:
  152. // Remove constant parameter blocks and residual blocks with no
  153. // varying parameter blocks while preserving their relative order.
  154. //
  155. // removed_parameter_blocks on exit will contain the list of
  156. // parameter blocks that were removed.
  157. //
  158. // fixed_cost will be equal to the sum of the costs of the residual
  159. // blocks that were removed.
  160. //
  161. // If there was a problem, then the function will return false and
  162. // error will contain a human readable description of the problem.
  163. bool RemoveFixedBlocks(std::vector<double*>* removed_parameter_blocks,
  164. double* fixed_cost,
  165. std::string* message);
  166. // The Program does not own the ParameterBlock or ResidualBlock objects.
  167. std::vector<ParameterBlock*> parameter_blocks_;
  168. std::vector<ResidualBlock*> residual_blocks_;
  169. EvaluationCallback* evaluation_callback_ = nullptr;
  170. friend class ProblemImpl;
  171. };
  172. } // namespace ceres::internal
  173. #include "ceres/internal/reenable_warnings.h"
  174. #endif // CERES_INTERNAL_PROGRAM_H_