You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

133 lines
4.7 KiB

  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud@inria.fr>
  5. //
  6. // This Source Code Form is subject to the terms of the Mozilla
  7. // Public License v. 2.0. If a copy of the MPL was not distributed
  8. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  9. #include "main.h"
  10. #include <Eigen/Geometry>
  11. #include <Eigen/LU>
  12. #include <Eigen/SVD>
  13. /* this test covers the following files:
  14. Geometry/OrthoMethods.h
  15. */
  16. template<typename Scalar> void orthomethods_3()
  17. {
  18. typedef typename NumTraits<Scalar>::Real RealScalar;
  19. typedef Matrix<Scalar,3,3> Matrix3;
  20. typedef Matrix<Scalar,3,1> Vector3;
  21. typedef Matrix<Scalar,4,1> Vector4;
  22. Vector3 v0 = Vector3::Random(),
  23. v1 = Vector3::Random(),
  24. v2 = Vector3::Random();
  25. // cross product
  26. VERIFY_IS_MUCH_SMALLER_THAN(v1.cross(v2).dot(v1), Scalar(1));
  27. VERIFY_IS_MUCH_SMALLER_THAN(v1.dot(v1.cross(v2)), Scalar(1));
  28. VERIFY_IS_MUCH_SMALLER_THAN(v1.cross(v2).dot(v2), Scalar(1));
  29. VERIFY_IS_MUCH_SMALLER_THAN(v2.dot(v1.cross(v2)), Scalar(1));
  30. VERIFY_IS_MUCH_SMALLER_THAN(v1.cross(Vector3::Random()).dot(v1), Scalar(1));
  31. Matrix3 mat3;
  32. mat3 << v0.normalized(),
  33. (v0.cross(v1)).normalized(),
  34. (v0.cross(v1).cross(v0)).normalized();
  35. VERIFY(mat3.isUnitary());
  36. mat3.setRandom();
  37. VERIFY_IS_APPROX(v0.cross(mat3*v1), -(mat3*v1).cross(v0));
  38. VERIFY_IS_APPROX(v0.cross(mat3.lazyProduct(v1)), -(mat3.lazyProduct(v1)).cross(v0));
  39. // colwise/rowwise cross product
  40. mat3.setRandom();
  41. Vector3 vec3 = Vector3::Random();
  42. Matrix3 mcross;
  43. int i = internal::random<int>(0,2);
  44. mcross = mat3.colwise().cross(vec3);
  45. VERIFY_IS_APPROX(mcross.col(i), mat3.col(i).cross(vec3));
  46. VERIFY_IS_MUCH_SMALLER_THAN((mat3.adjoint() * mat3.colwise().cross(vec3)).diagonal().cwiseAbs().sum(), Scalar(1));
  47. VERIFY_IS_MUCH_SMALLER_THAN((mat3.adjoint() * mat3.colwise().cross(Vector3::Random())).diagonal().cwiseAbs().sum(), Scalar(1));
  48. VERIFY_IS_MUCH_SMALLER_THAN((vec3.adjoint() * mat3.colwise().cross(vec3)).cwiseAbs().sum(), Scalar(1));
  49. VERIFY_IS_MUCH_SMALLER_THAN((vec3.adjoint() * Matrix3::Random().colwise().cross(vec3)).cwiseAbs().sum(), Scalar(1));
  50. mcross = mat3.rowwise().cross(vec3);
  51. VERIFY_IS_APPROX(mcross.row(i), mat3.row(i).cross(vec3));
  52. // cross3
  53. Vector4 v40 = Vector4::Random(),
  54. v41 = Vector4::Random(),
  55. v42 = Vector4::Random();
  56. v40.w() = v41.w() = v42.w() = 0;
  57. v42.template head<3>() = v40.template head<3>().cross(v41.template head<3>());
  58. VERIFY_IS_APPROX(v40.cross3(v41), v42);
  59. VERIFY_IS_MUCH_SMALLER_THAN(v40.cross3(Vector4::Random()).dot(v40), Scalar(1));
  60. // check mixed product
  61. typedef Matrix<RealScalar, 3, 1> RealVector3;
  62. RealVector3 rv1 = RealVector3::Random();
  63. VERIFY_IS_APPROX(v1.cross(rv1.template cast<Scalar>()), v1.cross(rv1));
  64. VERIFY_IS_APPROX(rv1.template cast<Scalar>().cross(v1), rv1.cross(v1));
  65. }
  66. template<typename Scalar, int Size> void orthomethods(int size=Size)
  67. {
  68. typedef typename NumTraits<Scalar>::Real RealScalar;
  69. typedef Matrix<Scalar,Size,1> VectorType;
  70. typedef Matrix<Scalar,3,Size> Matrix3N;
  71. typedef Matrix<Scalar,Size,3> MatrixN3;
  72. typedef Matrix<Scalar,3,1> Vector3;
  73. VectorType v0 = VectorType::Random(size);
  74. // unitOrthogonal
  75. VERIFY_IS_MUCH_SMALLER_THAN(v0.unitOrthogonal().dot(v0), Scalar(1));
  76. VERIFY_IS_APPROX(v0.unitOrthogonal().norm(), RealScalar(1));
  77. if (size>=3)
  78. {
  79. v0.template head<2>().setZero();
  80. v0.tail(size-2).setRandom();
  81. VERIFY_IS_MUCH_SMALLER_THAN(v0.unitOrthogonal().dot(v0), Scalar(1));
  82. VERIFY_IS_APPROX(v0.unitOrthogonal().norm(), RealScalar(1));
  83. }
  84. // colwise/rowwise cross product
  85. Vector3 vec3 = Vector3::Random();
  86. int i = internal::random<int>(0,size-1);
  87. Matrix3N mat3N(3,size), mcross3N(3,size);
  88. mat3N.setRandom();
  89. mcross3N = mat3N.colwise().cross(vec3);
  90. VERIFY_IS_APPROX(mcross3N.col(i), mat3N.col(i).cross(vec3));
  91. MatrixN3 matN3(size,3), mcrossN3(size,3);
  92. matN3.setRandom();
  93. mcrossN3 = matN3.rowwise().cross(vec3);
  94. VERIFY_IS_APPROX(mcrossN3.row(i), matN3.row(i).cross(vec3));
  95. }
  96. void test_geo_orthomethods()
  97. {
  98. for(int i = 0; i < g_repeat; i++) {
  99. CALL_SUBTEST_1( orthomethods_3<float>() );
  100. CALL_SUBTEST_2( orthomethods_3<double>() );
  101. CALL_SUBTEST_4( orthomethods_3<std::complex<double> >() );
  102. CALL_SUBTEST_1( (orthomethods<float,2>()) );
  103. CALL_SUBTEST_2( (orthomethods<double,2>()) );
  104. CALL_SUBTEST_1( (orthomethods<float,3>()) );
  105. CALL_SUBTEST_2( (orthomethods<double,3>()) );
  106. CALL_SUBTEST_3( (orthomethods<float,7>()) );
  107. CALL_SUBTEST_4( (orthomethods<std::complex<double>,8>()) );
  108. CALL_SUBTEST_5( (orthomethods<float,Dynamic>(36)) );
  109. CALL_SUBTEST_6( (orthomethods<double,Dynamic>(35)) );
  110. }
  111. }