photometric_error.h 7.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191
  1. // Ceres Solver - A fast non-linear least squares minimizer
  2. // Copyright 2020 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: nikolaus@nikolaus-demmel.de (Nikolaus Demmel)
  30. //
  31. //
  32. #ifndef CERES_INTERNAL_AUTODIFF_BENCHMARK_PHOTOMETRIC_ERROR_H_
  33. #define CERES_INTERNAL_AUTODIFF_BENCHMARK_PHOTOMETRIC_ERROR_H_
  34. #include <Eigen/Dense>
  35. #include "ceres/cubic_interpolation.h"
  36. namespace ceres {
  37. // Photometric residual that computes the intensity difference for a patch
  38. // between host and target frame. The point is parameterized with inverse
  39. // distance relative to the host frame. The relative pose between host and
  40. // target frame is computed from their respective absolute poses.
  41. //
  42. // The residual is similar to the one defined by Engel et al. [1]. Differences
  43. // include:
  44. //
  45. // 1. Use of a camera model based on spherical projection, namely the enhanced
  46. // unified camera model [2][3]. This is intended to bring some variability to
  47. // the benchmark compared to the SnavelyReprojection that uses a
  48. // polynomial-based distortion model.
  49. //
  50. // 2. To match the camera model, inverse distance parameterization is used for
  51. // points instead of inverse depth [4].
  52. //
  53. // 3. For simplicity, camera intrinsics are assumed constant, and thus host
  54. // frame points are passed as (unprojected) bearing vectors, which avoids the
  55. // need for an 'unproject' function.
  56. //
  57. // 4. Some details of the residual in [1] are omitted for simplicity: The
  58. // brightness transform parameters [a,b], the constant pre-weight w, and the
  59. // per-pixel robust norm.
  60. //
  61. // [1] J. Engel, V. Koltun and D. Cremers, "Direct Sparse Odometry," in IEEE
  62. // Transactions on Pattern Analysis and Machine Intelligence, vol. 40, no. 3,
  63. // pp. 611-625, 1 March 2018.
  64. //
  65. // [2] B. Khomutenko, G. Garcia and P. Martinet, "An Enhanced Unified Camera
  66. // Model," in IEEE Robotics and Automation Letters, vol. 1, no. 1, pp. 137-144,
  67. // Jan. 2016.
  68. //
  69. // [3] V. Usenko, N. Demmel and D. Cremers, "The Double Sphere Camera Model,"
  70. // 2018 International Conference on 3D Vision (3DV), Verona, 2018, pp. 552-560.
  71. //
  72. // [4] H. Matsuki, L. von Stumberg, V. Usenko, J. Stückler and D. Cremers,
  73. // "Omnidirectional DSO: Direct Sparse Odometry With Fisheye Cameras," in IEEE
  74. // Robotics and Automation Letters, vol. 3, no. 4, pp. 3693-3700, Oct. 2018.
  75. template <int PATCH_SIZE_ = 8>
  76. struct PhotometricError {
  77. static constexpr int PATCH_SIZE = PATCH_SIZE_;
  78. static constexpr int POSE_SIZE = 7;
  79. static constexpr int POINT_SIZE = 1;
  80. using Grid = Grid2D<uint8_t, 1>;
  81. using Interpolator = BiCubicInterpolator<Grid>;
  82. using Intrinsics = Eigen::Array<double, 6, 1>;
  83. template <typename T>
  84. using Patch = Eigen::Array<T, PATCH_SIZE, 1>;
  85. template <typename T>
  86. using PatchVectors = Eigen::Matrix<T, 3, PATCH_SIZE>;
  87. PhotometricError(const Patch<double>& intensities_host,
  88. const PatchVectors<double>& bearings_host,
  89. const Interpolator& image_target,
  90. const Intrinsics& intrinsics)
  91. : intensities_host_(intensities_host),
  92. bearings_host_(bearings_host),
  93. image_target_(image_target),
  94. intrinsics_(intrinsics) {}
  95. template <typename T>
  96. inline bool Project(Eigen::Matrix<T, 2, 1>& proj,
  97. const Eigen::Matrix<T, 3, 1>& p) const {
  98. const double& fx = intrinsics_[0];
  99. const double& fy = intrinsics_[1];
  100. const double& cx = intrinsics_[2];
  101. const double& cy = intrinsics_[3];
  102. const double& alpha = intrinsics_[4];
  103. const double& beta = intrinsics_[5];
  104. const T rho2 = beta * (p.x() * p.x() + p.y() * p.y()) + p.z() * p.z();
  105. const T rho = sqrt(rho2);
  106. // Check if 3D point is in domain of projection function.
  107. // See (8) and (17) in [3].
  108. constexpr double NUMERIC_EPSILON = 1e-10;
  109. const double w =
  110. alpha > 0.5 ? (1.0 - alpha) / alpha : alpha / (1.0 - alpha);
  111. if (p.z() <= -w * rho + NUMERIC_EPSILON) {
  112. return false;
  113. }
  114. const T norm = alpha * rho + (1.0 - alpha) * p.z();
  115. const T norm_inv = 1.0 / norm;
  116. const T mx = p.x() * norm_inv;
  117. const T my = p.y() * norm_inv;
  118. proj[0] = fx * mx + cx;
  119. proj[1] = fy * my + cy;
  120. return true;
  121. }
  122. template <typename T>
  123. inline bool operator()(const T* const pose_host_ptr,
  124. const T* const pose_target_ptr,
  125. const T* const idist_ptr,
  126. T* residuals_ptr) const {
  127. Eigen::Map<const Eigen::Quaternion<T>> q_w_h(pose_host_ptr);
  128. Eigen::Map<const Eigen::Matrix<T, 3, 1>> t_w_h(pose_host_ptr + 4);
  129. Eigen::Map<const Eigen::Quaternion<T>> q_w_t(pose_target_ptr);
  130. Eigen::Map<const Eigen::Matrix<T, 3, 1>> t_w_t(pose_target_ptr + 4);
  131. const T& idist = *idist_ptr;
  132. Eigen::Map<Patch<T>> residuals(residuals_ptr);
  133. // Compute relative pose from host to target frame.
  134. const Eigen::Quaternion<T> q_t_h = q_w_t.conjugate() * q_w_h;
  135. const Eigen::Matrix<T, 3, 3> R_t_h = q_t_h.toRotationMatrix();
  136. const Eigen::Matrix<T, 3, 1> t_t_h = q_w_t.conjugate() * (t_w_h - t_w_t);
  137. // Transform points from host to target frame. 3D point in target frame is
  138. // scaled by idist for numerical stability when idist is close to 0
  139. // (projection is invariant to scaling).
  140. PatchVectors<T> p_target_scaled =
  141. (R_t_h * bearings_host_).colwise() + idist * t_t_h;
  142. // Project points and interpolate image.
  143. Patch<T> intensities_target;
  144. for (int i = 0; i < p_target_scaled.cols(); ++i) {
  145. Eigen::Matrix<T, 2, 1> uv;
  146. if (!Project(uv, Eigen::Matrix<T, 3, 1>(p_target_scaled.col(i)))) {
  147. // If any point of the patch is outside the domain of the projection
  148. // function, the residual cannot be evaluated. For the benchmark we want
  149. // to avoid this case and thus throw an exception to indicate
  150. // immediately if it does actually happen after possible future changes.
  151. throw std::runtime_error("Benchmark data leads to invalid projection.");
  152. }
  153. // Mind the order of u and v: Evaluate takes (row, column), but u is
  154. // left-to-right and v top-to-bottom image axis.
  155. image_target_.Evaluate(uv[1], uv[0], &intensities_target[i]);
  156. }
  157. // Residual is intensity difference between host and target frame.
  158. residuals = intensities_target - intensities_host_;
  159. return true;
  160. }
  161. private:
  162. const Patch<double>& intensities_host_;
  163. const PatchVectors<double>& bearings_host_;
  164. const Interpolator& image_target_;
  165. const Intrinsics& intrinsics_;
  166. };
  167. } // namespace ceres
  168. #endif // CERES_INTERNAL_AUTODIFF_BENCHMARK_PHOTOMETRIC_ERROR_H_