problem.h 25 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554
  1. // Ceres Solver - A fast non-linear least squares minimizer
  2. // Copyright 2021 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. // keir@google.com (Keir Mierle)
  31. //
  32. // The Problem object is used to build and hold least squares problems.
  33. #ifndef CERES_PUBLIC_PROBLEM_H_
  34. #define CERES_PUBLIC_PROBLEM_H_
  35. #include <array>
  36. #include <cstddef>
  37. #include <map>
  38. #include <memory>
  39. #include <set>
  40. #include <vector>
  41. #include "ceres/context.h"
  42. #include "ceres/internal/disable_warnings.h"
  43. #include "ceres/internal/export.h"
  44. #include "ceres/internal/port.h"
  45. #include "ceres/types.h"
  46. #include "glog/logging.h"
  47. namespace ceres {
  48. class CostFunction;
  49. class EvaluationCallback;
  50. class LossFunction;
  51. class Manifold;
  52. class Solver;
  53. struct CRSMatrix;
  54. namespace internal {
  55. class Preprocessor;
  56. class ProblemImpl;
  57. class ParameterBlock;
  58. class ResidualBlock;
  59. } // namespace internal
  60. // A ResidualBlockId is an opaque handle clients can use to remove residual
  61. // blocks from a Problem after adding them.
  62. using ResidualBlockId = internal::ResidualBlock*;
  63. // A class to represent non-linear least squares problems. Such
  64. // problems have a cost function that is a sum of error terms (known
  65. // as "residuals"), where each residual is a function of some subset
  66. // of the parameters. The cost function takes the form
  67. //
  68. // N 1
  69. // SUM --- loss( || r_i1, r_i2,..., r_ik ||^2 ),
  70. // i=1 2
  71. //
  72. // where
  73. //
  74. // r_ij is residual number i, component j; the residual is a function of some
  75. // subset of the parameters x1...xk. For example, in a structure from
  76. // motion problem a residual might be the difference between a measured
  77. // point in an image and the reprojected position for the matching
  78. // camera, point pair. The residual would have two components, error in x
  79. // and error in y.
  80. //
  81. // loss(y) is the loss function; for example, squared error or Huber L1
  82. // loss. If loss(y) = y, then the cost function is non-robustified
  83. // least squares.
  84. //
  85. // This class is specifically designed to address the important subset of
  86. // "sparse" least squares problems, where each component of the residual depends
  87. // only on a small number number of parameters, even though the total number of
  88. // residuals and parameters may be very large. This property affords tremendous
  89. // gains in scale, allowing efficient solving of large problems that are
  90. // otherwise inaccessible.
  91. //
  92. // The canonical example of a sparse least squares problem is
  93. // "structure-from-motion" (SFM), where the parameters are points and cameras,
  94. // and residuals are reprojection errors. Typically a single residual will
  95. // depend only on 9 parameters (3 for the point, 6 for the camera).
  96. //
  97. // To create a least squares problem, use the AddResidualBlock() and
  98. // AddParameterBlock() methods, documented below. Here is an example least
  99. // squares problem containing 3 parameter blocks of sizes 3, 4 and 5
  100. // respectively and two residual terms of size 2 and 6:
  101. //
  102. // double x1[] = { 1.0, 2.0, 3.0 };
  103. // double x2[] = { 1.0, 2.0, 3.0, 5.0 };
  104. // double x3[] = { 1.0, 2.0, 3.0, 6.0, 7.0 };
  105. //
  106. // Problem problem;
  107. //
  108. // problem.AddResidualBlock(new MyUnaryCostFunction(...), nullptr, x1);
  109. // problem.AddResidualBlock(new MyBinaryCostFunction(...), nullptr, x2, x3);
  110. //
  111. // Please see cost_function.h for details of the CostFunction object.
  112. class CERES_EXPORT Problem {
  113. public:
  114. struct CERES_EXPORT Options {
  115. // These flags control whether the Problem object owns the CostFunctions,
  116. // LossFunctions, and Manifolds passed into the Problem.
  117. //
  118. // If set to TAKE_OWNERSHIP, then the problem object will delete the
  119. // corresponding object on destruction. The destructor is careful to delete
  120. // the pointers only once, since sharing objects is allowed.
  121. Ownership cost_function_ownership = TAKE_OWNERSHIP;
  122. Ownership loss_function_ownership = TAKE_OWNERSHIP;
  123. Ownership manifold_ownership = TAKE_OWNERSHIP;
  124. // If true, trades memory for faster RemoveResidualBlock() and
  125. // RemoveParameterBlock() operations.
  126. //
  127. // By default, RemoveParameterBlock() and RemoveResidualBlock() take time
  128. // proportional to the size of the entire problem. If you only ever remove
  129. // parameters or residuals from the problem occasionally, this might be
  130. // acceptable. However, if you have memory to spare, enable this option to
  131. // make RemoveParameterBlock() take time proportional to the number of
  132. // residual blocks that depend on it, and RemoveResidualBlock() take (on
  133. // average) constant time.
  134. //
  135. // The increase in memory usage is two-fold: an additional hash set per
  136. // parameter block containing all the residuals that depend on the parameter
  137. // block; and a hash set in the problem containing all residuals.
  138. bool enable_fast_removal = false;
  139. // By default, Ceres performs a variety of safety checks when constructing
  140. // the problem. There is a small but measurable performance penalty to these
  141. // checks, typically around 5% of construction time. If you are sure your
  142. // problem construction is correct, and 5% of the problem construction time
  143. // is truly an overhead you want to avoid, then you can set
  144. // disable_all_safety_checks to true.
  145. //
  146. // WARNING: Do not set this to true, unless you are absolutely sure of what
  147. // you are doing.
  148. bool disable_all_safety_checks = false;
  149. // A Ceres global context to use for solving this problem. This may help to
  150. // reduce computation time as Ceres can reuse expensive objects to create.
  151. // The context object can be nullptr, in which case Ceres may create one.
  152. //
  153. // Ceres does NOT take ownership of the pointer.
  154. Context* context = nullptr;
  155. // Using this callback interface, Ceres can notify you when it is about to
  156. // evaluate the residuals or jacobians. With the callback, you can share
  157. // computation between residual blocks by doing the shared computation in
  158. // EvaluationCallback::PrepareForEvaluation() before Ceres calls
  159. // CostFunction::Evaluate(). It also enables caching results between a pure
  160. // residual evaluation and a residual & jacobian evaluation.
  161. //
  162. // Problem DOES NOT take ownership of the callback.
  163. //
  164. // NOTE: Evaluation callbacks are incompatible with inner iterations. So
  165. // calling Solve with Solver::Options::use_inner_iterations = true on a
  166. // Problem with a non-null evaluation callback is an error.
  167. EvaluationCallback* evaluation_callback = nullptr;
  168. };
  169. // The default constructor is equivalent to the invocation
  170. // Problem(Problem::Options()).
  171. Problem();
  172. explicit Problem(const Options& options);
  173. Problem(Problem&&);
  174. Problem& operator=(Problem&&);
  175. Problem(const Problem&) = delete;
  176. Problem& operator=(const Problem&) = delete;
  177. ~Problem();
  178. // Add a residual block to the overall cost function. The cost function
  179. // carries with its information about the sizes of the parameter blocks it
  180. // expects. The function checks that these match the sizes of the parameter
  181. // blocks listed in parameter_blocks. The program aborts if a mismatch is
  182. // detected. loss_function can be nullptr, in which case the cost of the term
  183. // is just the squared norm of the residuals.
  184. //
  185. // The user has the option of explicitly adding the parameter blocks using
  186. // AddParameterBlock. This causes additional correctness checking; however,
  187. // AddResidualBlock implicitly adds the parameter blocks if they are not
  188. // present, so calling AddParameterBlock explicitly is not required.
  189. //
  190. // The Problem object by default takes ownership of the cost_function and
  191. // loss_function pointers (See Problem::Options to override this behaviour).
  192. // These objects remain live for the life of the Problem object. If the user
  193. // wishes to keep control over the destruction of these objects, then they can
  194. // do this by setting the corresponding enums in the Options struct.
  195. //
  196. // Note: Even though the Problem takes ownership of cost_function and
  197. // loss_function, it does not preclude the user from re-using them in another
  198. // residual block. The destructor takes care to call delete on each
  199. // cost_function or loss_function pointer only once, regardless of how many
  200. // residual blocks refer to them.
  201. //
  202. // Example usage:
  203. //
  204. // double x1[] = {1.0, 2.0, 3.0};
  205. // double x2[] = {1.0, 2.0, 5.0, 6.0};
  206. // double x3[] = {3.0, 6.0, 2.0, 5.0, 1.0};
  207. //
  208. // Problem problem;
  209. //
  210. // problem.AddResidualBlock(new MyUnaryCostFunction(...), nullptr, x1);
  211. // problem.AddResidualBlock(new MyBinaryCostFunction(...), nullptr, x2, x1);
  212. //
  213. // Add a residual block by listing the parameter block pointers directly
  214. // instead of wapping them in a container.
  215. template <typename... Ts>
  216. ResidualBlockId AddResidualBlock(CostFunction* cost_function,
  217. LossFunction* loss_function,
  218. double* x0,
  219. Ts*... xs) {
  220. const std::array<double*, sizeof...(Ts) + 1> parameter_blocks{{x0, xs...}};
  221. return AddResidualBlock(cost_function,
  222. loss_function,
  223. parameter_blocks.data(),
  224. static_cast<int>(parameter_blocks.size()));
  225. }
  226. // Add a residual block by providing a vector of parameter blocks.
  227. ResidualBlockId AddResidualBlock(
  228. CostFunction* cost_function,
  229. LossFunction* loss_function,
  230. const std::vector<double*>& parameter_blocks);
  231. // Add a residual block by providing a pointer to the parameter block array
  232. // and the number of parameter blocks.
  233. ResidualBlockId AddResidualBlock(CostFunction* cost_function,
  234. LossFunction* loss_function,
  235. double* const* const parameter_blocks,
  236. int num_parameter_blocks);
  237. // Add a parameter block with appropriate size to the problem. Repeated calls
  238. // with the same arguments are ignored. Repeated calls with the same double
  239. // pointer but a different size will result in a crash.
  240. void AddParameterBlock(double* values, int size);
  241. // Add a parameter block with appropriate size and Manifold to the
  242. // problem. It is okay for manifold to be nullptr.
  243. //
  244. // Repeated calls with the same arguments are ignored. Repeated calls
  245. // with the same double pointer but a different size results in a crash
  246. // (unless Solver::Options::disable_all_safety_checks is set to true).
  247. //
  248. // Repeated calls with the same double pointer and size but different Manifold
  249. // is equivalent to calling SetManifold(manifold), i.e., any previously
  250. // associated Manifold object will be replaced with the manifold.
  251. void AddParameterBlock(double* values, int size, Manifold* manifold);
  252. // Remove a parameter block from the problem. The Manifold of the parameter
  253. // block, if it exists, will persist until the deletion of the problem
  254. // (similar to cost/loss functions in residual block removal). Any residual
  255. // blocks that depend on the parameter are also removed, as described above
  256. // in RemoveResidualBlock().
  257. //
  258. // If Problem::Options::enable_fast_removal is true, then the removal is fast
  259. // (almost constant time). Otherwise, removing a parameter block will incur a
  260. // scan of the entire Problem object.
  261. //
  262. // WARNING: Removing a residual or parameter block will destroy the implicit
  263. // ordering, rendering the jacobian or residuals returned from the solver
  264. // uninterpretable. If you depend on the evaluated jacobian, do not use
  265. // remove! This may change in a future release.
  266. void RemoveParameterBlock(const double* values);
  267. // Remove a residual block from the problem. Any parameters that the residual
  268. // block depends on are not removed. The cost and loss functions for the
  269. // residual block will not get deleted immediately; won't happen until the
  270. // problem itself is deleted.
  271. //
  272. // WARNING: Removing a residual or parameter block will destroy the implicit
  273. // ordering, rendering the jacobian or residuals returned from the solver
  274. // uninterpretable. If you depend on the evaluated jacobian, do not use
  275. // remove! This may change in a future release.
  276. void RemoveResidualBlock(ResidualBlockId residual_block);
  277. // Hold the indicated parameter block constant during optimization.
  278. void SetParameterBlockConstant(const double* values);
  279. // Allow the indicated parameter block to vary during optimization.
  280. void SetParameterBlockVariable(double* values);
  281. // Returns true if a parameter block is set constant, and false otherwise. A
  282. // parameter block may be set constant in two ways: either by calling
  283. // SetParameterBlockConstant or by associating a Manifold with a zero
  284. // dimensional tangent space with it.
  285. bool IsParameterBlockConstant(const double* values) const;
  286. // Set the Manifold for the parameter block. Calling SetManifold with nullptr
  287. // will clear any previously set Manifold for the parameter block.
  288. //
  289. // Repeated calls will result in any previously associated Manifold object to
  290. // be replaced with the manifold.
  291. //
  292. // The manifold is owned by the Problem by default (See Problem::Options to
  293. // override this behaviour).
  294. //
  295. // It is acceptable to set the same Manifold for multiple parameter blocks.
  296. void SetManifold(double* values, Manifold* manifold);
  297. // Get the Manifold object associated with this parameter block.
  298. //
  299. // If there is no Manifold object associated then nullptr is returned.
  300. const Manifold* GetManifold(const double* values) const;
  301. // Returns true if a Manifold is associated with this parameter block, false
  302. // otherwise.
  303. bool HasManifold(const double* values) const;
  304. // Set the lower/upper bound for the parameter at position "index".
  305. void SetParameterLowerBound(double* values, int index, double lower_bound);
  306. void SetParameterUpperBound(double* values, int index, double upper_bound);
  307. // Get the lower/upper bound for the parameter at position "index". If the
  308. // parameter is not bounded by the user, then its lower bound is
  309. // -std::numeric_limits<double>::max() and upper bound is
  310. // std::numeric_limits<double>::max().
  311. double GetParameterLowerBound(const double* values, int index) const;
  312. double GetParameterUpperBound(const double* values, int index) const;
  313. // Number of parameter blocks in the problem. Always equals
  314. // parameter_blocks().size() and parameter_block_sizes().size().
  315. int NumParameterBlocks() const;
  316. // The size of the parameter vector obtained by summing over the sizes of all
  317. // the parameter blocks.
  318. int NumParameters() const;
  319. // Number of residual blocks in the problem. Always equals
  320. // residual_blocks().size().
  321. int NumResidualBlocks() const;
  322. // The size of the residual vector obtained by summing over the sizes of all
  323. // of the residual blocks.
  324. int NumResiduals() const;
  325. // The size of the parameter block.
  326. int ParameterBlockSize(const double* values) const;
  327. // The dimension of the tangent space of the Manifold for the parameter block.
  328. // If there is no Manifold associated with this parameter block, then
  329. // ParameterBlockTangentSize = ParameterBlockSize.
  330. int ParameterBlockTangentSize(const double* values) const;
  331. // Is the given parameter block present in this problem or not?
  332. bool HasParameterBlock(const double* values) const;
  333. // Fills the passed parameter_blocks vector with pointers to the parameter
  334. // blocks currently in the problem. After this call, parameter_block.size() ==
  335. // NumParameterBlocks.
  336. void GetParameterBlocks(std::vector<double*>* parameter_blocks) const;
  337. // Fills the passed residual_blocks vector with pointers to the residual
  338. // blocks currently in the problem. After this call, residual_blocks.size() ==
  339. // NumResidualBlocks.
  340. void GetResidualBlocks(std::vector<ResidualBlockId>* residual_blocks) const;
  341. // Get all the parameter blocks that depend on the given residual block.
  342. void GetParameterBlocksForResidualBlock(
  343. const ResidualBlockId residual_block,
  344. std::vector<double*>* parameter_blocks) const;
  345. // Get the CostFunction for the given residual block.
  346. const CostFunction* GetCostFunctionForResidualBlock(
  347. const ResidualBlockId residual_block) const;
  348. // Get the LossFunction for the given residual block. Returns nullptr
  349. // if no loss function is associated with this residual block.
  350. const LossFunction* GetLossFunctionForResidualBlock(
  351. const ResidualBlockId residual_block) const;
  352. // Get all the residual blocks that depend on the given parameter block.
  353. //
  354. // If Problem::Options::enable_fast_removal is true, then getting the residual
  355. // blocks is fast and depends only on the number of residual
  356. // blocks. Otherwise, getting the residual blocks for a parameter block will
  357. // incur a scan of the entire Problem object.
  358. void GetResidualBlocksForParameterBlock(
  359. const double* values,
  360. std::vector<ResidualBlockId>* residual_blocks) const;
  361. // Options struct to control Problem::Evaluate.
  362. struct EvaluateOptions {
  363. // The set of parameter blocks for which evaluation should be
  364. // performed. This vector determines the order that parameter blocks occur
  365. // in the gradient vector and in the columns of the jacobian matrix. If
  366. // parameter_blocks is empty, then it is assumed to be equal to vector
  367. // containing ALL the parameter blocks. Generally speaking the parameter
  368. // blocks will occur in the order in which they were added to the
  369. // problem. But, this may change if the user removes any parameter blocks
  370. // from the problem.
  371. //
  372. // NOTE: This vector should contain the same pointers as the ones used to
  373. // add parameter blocks to the Problem. These parameter block should NOT
  374. // point to new memory locations. Bad things will happen otherwise.
  375. std::vector<double*> parameter_blocks;
  376. // The set of residual blocks to evaluate. This vector determines the order
  377. // in which the residuals occur, and how the rows of the jacobian are
  378. // ordered. If residual_blocks is empty, then it is assumed to be equal to
  379. // the vector containing ALL the residual blocks. Generally speaking the
  380. // residual blocks will occur in the order in which they were added to the
  381. // problem. But, this may change if the user removes any residual blocks
  382. // from the problem.
  383. std::vector<ResidualBlockId> residual_blocks;
  384. // Even though the residual blocks in the problem may contain loss
  385. // functions, setting apply_loss_function to false will turn off the
  386. // application of the loss function to the output of the cost function. This
  387. // is of use for example if the user wishes to analyse the solution quality
  388. // by studying the distribution of residuals before and after the solve.
  389. bool apply_loss_function = true;
  390. int num_threads = 1;
  391. };
  392. // Evaluate Problem. Any of the output pointers can be nullptr. Which residual
  393. // blocks and parameter blocks are used is controlled by the EvaluateOptions
  394. // struct above.
  395. //
  396. // Note 1: The evaluation will use the values stored in the memory locations
  397. // pointed to by the parameter block pointers used at the time of the
  398. // construction of the problem. i.e.,
  399. //
  400. // Problem problem;
  401. // double x = 1;
  402. // problem.AddResidualBlock(new MyCostFunction, nullptr, &x);
  403. //
  404. // double cost = 0.0;
  405. // problem.Evaluate(Problem::EvaluateOptions(), &cost,
  406. // nullptr, nullptr, nullptr);
  407. //
  408. // The cost is evaluated at x = 1. If you wish to evaluate the problem at x =
  409. // 2, then
  410. //
  411. // x = 2;
  412. // problem.Evaluate(Problem::EvaluateOptions(), &cost,
  413. // nullptr, nullptr, nullptr);
  414. //
  415. // is the way to do so.
  416. //
  417. // Note 2: If no Manifolds are used, then the size of the gradient vector (and
  418. // the number of columns in the jacobian) is the sum of the sizes of all the
  419. // parameter blocks. If a parameter block has a Manifold, then it contributes
  420. // "TangentSize" entries to the gradient vector (and the number of columns in
  421. // the jacobian).
  422. //
  423. // Note 3: This function cannot be called while the problem is being solved,
  424. // for example it cannot be called from an IterationCallback at the end of an
  425. // iteration during a solve.
  426. //
  427. // Note 4: If an EvaluationCallback is associated with the problem, then its
  428. // PrepareForEvaluation method will be called every time this method is called
  429. // with new_point = true.
  430. bool Evaluate(const EvaluateOptions& options,
  431. double* cost,
  432. std::vector<double>* residuals,
  433. std::vector<double>* gradient,
  434. CRSMatrix* jacobian);
  435. // Evaluates the residual block, storing the scalar cost in *cost, the
  436. // residual components in *residuals, and the jacobians between the parameters
  437. // and residuals in jacobians[i], in row-major order.
  438. //
  439. // If residuals is nullptr, the residuals are not computed.
  440. //
  441. // If jacobians is nullptr, no Jacobians are computed. If jacobians[i] is
  442. // nullptr, then the Jacobian for that parameter block is not computed.
  443. //
  444. // It is not okay to request the Jacobian w.r.t a parameter block that is
  445. // constant.
  446. //
  447. // The return value indicates the success or failure. Even if the function
  448. // returns false, the caller should expect the output memory locations to have
  449. // been modified.
  450. //
  451. // The returned cost and jacobians have had robustification and Manifold
  452. // applied already; for example, the jacobian for a 4-dimensional quaternion
  453. // parameter using the "QuaternionParameterization" is num_residuals by 3
  454. // instead of num_residuals by 4.
  455. //
  456. // apply_loss_function as the name implies allows the user to switch the
  457. // application of the loss function on and off.
  458. //
  459. // If an EvaluationCallback is associated with the problem, then its
  460. // PrepareForEvaluation method will be called every time this method is called
  461. // with new_point = true. This conservatively assumes that the user may have
  462. // changed the parameter values since the previous call to evaluate / solve.
  463. // For improved efficiency, and only if you know that the parameter values
  464. // have not changed between calls, see
  465. // EvaluateResidualBlockAssumingParametersUnchanged().
  466. bool EvaluateResidualBlock(ResidualBlockId residual_block_id,
  467. bool apply_loss_function,
  468. double* cost,
  469. double* residuals,
  470. double** jacobians) const;
  471. // Same as EvaluateResidualBlock except that if an EvaluationCallback is
  472. // associated with the problem, then its PrepareForEvaluation method will be
  473. // called every time this method is called with new_point = false.
  474. //
  475. // This means, if an EvaluationCallback is associated with the problem then it
  476. // is the user's responsibility to call PrepareForEvaluation before calling
  477. // this method if necessary, i.e. iff the parameter values have been changed
  478. // since the last call to evaluate / solve.'
  479. //
  480. // This is because, as the name implies, we assume that the parameter blocks
  481. // did not change since the last time PrepareForEvaluation was called (via
  482. // Solve, Evaluate or EvaluateResidualBlock).
  483. bool EvaluateResidualBlockAssumingParametersUnchanged(
  484. ResidualBlockId residual_block_id,
  485. bool apply_loss_function,
  486. double* cost,
  487. double* residuals,
  488. double** jacobians) const;
  489. // Returns reference to the options with which the Problem was constructed.
  490. const Options& options() const;
  491. // Returns pointer to Problem implementation
  492. internal::ProblemImpl* mutable_impl();
  493. private:
  494. std::unique_ptr<internal::ProblemImpl> impl_;
  495. };
  496. } // namespace ceres
  497. #include "ceres/internal/reenable_warnings.h"
  498. #endif // CERES_PUBLIC_PROBLEM_H_