autodiff.cpp 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2009 Gael Guennebaud <g.gael@free.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 <unsupported/Eigen/AutoDiff>
  11. template<typename Scalar>
  12. EIGEN_DONT_INLINE Scalar foo(const Scalar& x, const Scalar& y)
  13. {
  14. using namespace std;
  15. // return x+std::sin(y);
  16. EIGEN_ASM_COMMENT("mybegin");
  17. // pow(float, int) promotes to pow(double, double)
  18. return x*2 - 1 + static_cast<Scalar>(pow(1+x,2)) + 2*sqrt(y*y+0) - 4 * sin(0+x) + 2 * cos(y+0) - exp(Scalar(-0.5)*x*x+0);
  19. //return x+2*y*x;//x*2 -std::pow(x,2);//(2*y/x);// - y*2;
  20. EIGEN_ASM_COMMENT("myend");
  21. }
  22. template<typename Vector>
  23. EIGEN_DONT_INLINE typename Vector::Scalar foo(const Vector& p)
  24. {
  25. typedef typename Vector::Scalar Scalar;
  26. return (p-Vector(Scalar(-1),Scalar(1.))).norm() + (p.array() * p.array()).sum() + p.dot(p);
  27. }
  28. template<typename _Scalar, int NX=Dynamic, int NY=Dynamic>
  29. struct TestFunc1
  30. {
  31. typedef _Scalar Scalar;
  32. enum {
  33. InputsAtCompileTime = NX,
  34. ValuesAtCompileTime = NY
  35. };
  36. typedef Matrix<Scalar,InputsAtCompileTime,1> InputType;
  37. typedef Matrix<Scalar,ValuesAtCompileTime,1> ValueType;
  38. typedef Matrix<Scalar,ValuesAtCompileTime,InputsAtCompileTime> JacobianType;
  39. int m_inputs, m_values;
  40. TestFunc1() : m_inputs(InputsAtCompileTime), m_values(ValuesAtCompileTime) {}
  41. TestFunc1(int inputs_, int values_) : m_inputs(inputs_), m_values(values_) {}
  42. int inputs() const { return m_inputs; }
  43. int values() const { return m_values; }
  44. template<typename T>
  45. void operator() (const Matrix<T,InputsAtCompileTime,1>& x, Matrix<T,ValuesAtCompileTime,1>* _v) const
  46. {
  47. Matrix<T,ValuesAtCompileTime,1>& v = *_v;
  48. v[0] = 2 * x[0] * x[0] + x[0] * x[1];
  49. v[1] = 3 * x[1] * x[0] + 0.5 * x[1] * x[1];
  50. if(inputs()>2)
  51. {
  52. v[0] += 0.5 * x[2];
  53. v[1] += x[2];
  54. }
  55. if(values()>2)
  56. {
  57. v[2] = 3 * x[1] * x[0] * x[0];
  58. }
  59. if (inputs()>2 && values()>2)
  60. v[2] *= x[2];
  61. }
  62. void operator() (const InputType& x, ValueType* v, JacobianType* _j) const
  63. {
  64. (*this)(x, v);
  65. if(_j)
  66. {
  67. JacobianType& j = *_j;
  68. j(0,0) = 4 * x[0] + x[1];
  69. j(1,0) = 3 * x[1];
  70. j(0,1) = x[0];
  71. j(1,1) = 3 * x[0] + 2 * 0.5 * x[1];
  72. if (inputs()>2)
  73. {
  74. j(0,2) = 0.5;
  75. j(1,2) = 1;
  76. }
  77. if(values()>2)
  78. {
  79. j(2,0) = 3 * x[1] * 2 * x[0];
  80. j(2,1) = 3 * x[0] * x[0];
  81. }
  82. if (inputs()>2 && values()>2)
  83. {
  84. j(2,0) *= x[2];
  85. j(2,1) *= x[2];
  86. j(2,2) = 3 * x[1] * x[0] * x[0];
  87. j(2,2) = 3 * x[1] * x[0] * x[0];
  88. }
  89. }
  90. }
  91. };
  92. #if EIGEN_HAS_VARIADIC_TEMPLATES
  93. /* Test functor for the C++11 features. */
  94. template <typename Scalar>
  95. struct integratorFunctor
  96. {
  97. typedef Matrix<Scalar, 2, 1> InputType;
  98. typedef Matrix<Scalar, 2, 1> ValueType;
  99. /*
  100. * Implementation starts here.
  101. */
  102. integratorFunctor(const Scalar gain) : _gain(gain) {}
  103. integratorFunctor(const integratorFunctor& f) : _gain(f._gain) {}
  104. const Scalar _gain;
  105. template <typename T1, typename T2>
  106. void operator() (const T1 &input, T2 *output, const Scalar dt) const
  107. {
  108. T2 &o = *output;
  109. /* Integrator to test the AD. */
  110. o[0] = input[0] + input[1] * dt * _gain;
  111. o[1] = input[1] * _gain;
  112. }
  113. /* Only needed for the test */
  114. template <typename T1, typename T2, typename T3>
  115. void operator() (const T1 &input, T2 *output, T3 *jacobian, const Scalar dt) const
  116. {
  117. T2 &o = *output;
  118. /* Integrator to test the AD. */
  119. o[0] = input[0] + input[1] * dt * _gain;
  120. o[1] = input[1] * _gain;
  121. if (jacobian)
  122. {
  123. T3 &j = *jacobian;
  124. j(0, 0) = 1;
  125. j(0, 1) = dt * _gain;
  126. j(1, 0) = 0;
  127. j(1, 1) = _gain;
  128. }
  129. }
  130. };
  131. template<typename Func> void forward_jacobian_cpp11(const Func& f)
  132. {
  133. typedef typename Func::ValueType::Scalar Scalar;
  134. typedef typename Func::ValueType ValueType;
  135. typedef typename Func::InputType InputType;
  136. typedef typename AutoDiffJacobian<Func>::JacobianType JacobianType;
  137. InputType x = InputType::Random(InputType::RowsAtCompileTime);
  138. ValueType y, yref;
  139. JacobianType j, jref;
  140. const Scalar dt = internal::random<double>();
  141. jref.setZero();
  142. yref.setZero();
  143. f(x, &yref, &jref, dt);
  144. //std::cerr << "y, yref, jref: " << "\n";
  145. //std::cerr << y.transpose() << "\n\n";
  146. //std::cerr << yref << "\n\n";
  147. //std::cerr << jref << "\n\n";
  148. AutoDiffJacobian<Func> autoj(f);
  149. autoj(x, &y, &j, dt);
  150. //std::cerr << "y j (via autodiff): " << "\n";
  151. //std::cerr << y.transpose() << "\n\n";
  152. //std::cerr << j << "\n\n";
  153. VERIFY_IS_APPROX(y, yref);
  154. VERIFY_IS_APPROX(j, jref);
  155. }
  156. #endif
  157. template<typename Func> void forward_jacobian(const Func& f)
  158. {
  159. typename Func::InputType x = Func::InputType::Random(f.inputs());
  160. typename Func::ValueType y(f.values()), yref(f.values());
  161. typename Func::JacobianType j(f.values(),f.inputs()), jref(f.values(),f.inputs());
  162. jref.setZero();
  163. yref.setZero();
  164. f(x,&yref,&jref);
  165. // std::cerr << y.transpose() << "\n\n";;
  166. // std::cerr << j << "\n\n";;
  167. j.setZero();
  168. y.setZero();
  169. AutoDiffJacobian<Func> autoj(f);
  170. autoj(x, &y, &j);
  171. // std::cerr << y.transpose() << "\n\n";;
  172. // std::cerr << j << "\n\n";;
  173. VERIFY_IS_APPROX(y, yref);
  174. VERIFY_IS_APPROX(j, jref);
  175. }
  176. // TODO also check actual derivatives!
  177. template <int>
  178. void test_autodiff_scalar()
  179. {
  180. Vector2f p = Vector2f::Random();
  181. typedef AutoDiffScalar<Vector2f> AD;
  182. AD ax(p.x(),Vector2f::UnitX());
  183. AD ay(p.y(),Vector2f::UnitY());
  184. AD res = foo<AD>(ax,ay);
  185. VERIFY_IS_APPROX(res.value(), foo(p.x(),p.y()));
  186. }
  187. // TODO also check actual derivatives!
  188. template <int>
  189. void test_autodiff_vector()
  190. {
  191. Vector2f p = Vector2f::Random();
  192. typedef AutoDiffScalar<Vector2f> AD;
  193. typedef Matrix<AD,2,1> VectorAD;
  194. VectorAD ap = p.cast<AD>();
  195. ap.x().derivatives() = Vector2f::UnitX();
  196. ap.y().derivatives() = Vector2f::UnitY();
  197. AD res = foo<VectorAD>(ap);
  198. VERIFY_IS_APPROX(res.value(), foo(p));
  199. }
  200. template <int>
  201. void test_autodiff_jacobian()
  202. {
  203. CALL_SUBTEST(( forward_jacobian(TestFunc1<double,2,2>()) ));
  204. CALL_SUBTEST(( forward_jacobian(TestFunc1<double,2,3>()) ));
  205. CALL_SUBTEST(( forward_jacobian(TestFunc1<double,3,2>()) ));
  206. CALL_SUBTEST(( forward_jacobian(TestFunc1<double,3,3>()) ));
  207. CALL_SUBTEST(( forward_jacobian(TestFunc1<double>(3,3)) ));
  208. #if EIGEN_HAS_VARIADIC_TEMPLATES
  209. CALL_SUBTEST(( forward_jacobian_cpp11(integratorFunctor<double>(10)) ));
  210. #endif
  211. }
  212. template <int>
  213. void test_autodiff_hessian()
  214. {
  215. typedef AutoDiffScalar<VectorXd> AD;
  216. typedef Matrix<AD,Eigen::Dynamic,1> VectorAD;
  217. typedef AutoDiffScalar<VectorAD> ADD;
  218. typedef Matrix<ADD,Eigen::Dynamic,1> VectorADD;
  219. VectorADD x(2);
  220. double s1 = internal::random<double>(), s2 = internal::random<double>(), s3 = internal::random<double>(), s4 = internal::random<double>();
  221. x(0).value()=s1;
  222. x(1).value()=s2;
  223. //set unit vectors for the derivative directions (partial derivatives of the input vector)
  224. x(0).derivatives().resize(2);
  225. x(0).derivatives().setZero();
  226. x(0).derivatives()(0)= 1;
  227. x(1).derivatives().resize(2);
  228. x(1).derivatives().setZero();
  229. x(1).derivatives()(1)=1;
  230. //repeat partial derivatives for the inner AutoDiffScalar
  231. x(0).value().derivatives() = VectorXd::Unit(2,0);
  232. x(1).value().derivatives() = VectorXd::Unit(2,1);
  233. //set the hessian matrix to zero
  234. for(int idx=0; idx<2; idx++) {
  235. x(0).derivatives()(idx).derivatives() = VectorXd::Zero(2);
  236. x(1).derivatives()(idx).derivatives() = VectorXd::Zero(2);
  237. }
  238. ADD y = sin(AD(s3)*x(0) + AD(s4)*x(1));
  239. VERIFY_IS_APPROX(y.value().derivatives()(0), y.derivatives()(0).value());
  240. VERIFY_IS_APPROX(y.value().derivatives()(1), y.derivatives()(1).value());
  241. VERIFY_IS_APPROX(y.value().derivatives()(0), s3*std::cos(s1*s3+s2*s4));
  242. VERIFY_IS_APPROX(y.value().derivatives()(1), s4*std::cos(s1*s3+s2*s4));
  243. VERIFY_IS_APPROX(y.derivatives()(0).derivatives(), -std::sin(s1*s3+s2*s4)*Vector2d(s3*s3,s4*s3));
  244. VERIFY_IS_APPROX(y.derivatives()(1).derivatives(), -std::sin(s1*s3+s2*s4)*Vector2d(s3*s4,s4*s4));
  245. ADD z = x(0)*x(1);
  246. VERIFY_IS_APPROX(z.derivatives()(0).derivatives(), Vector2d(0,1));
  247. VERIFY_IS_APPROX(z.derivatives()(1).derivatives(), Vector2d(1,0));
  248. }
  249. double bug_1222() {
  250. typedef Eigen::AutoDiffScalar<Eigen::Vector3d> AD;
  251. const double _cv1_3 = 1.0;
  252. const AD chi_3 = 1.0;
  253. // this line did not work, because operator+ returns ADS<DerType&>, which then cannot be converted to ADS<DerType>
  254. const AD denom = chi_3 + _cv1_3;
  255. return denom.value();
  256. }
  257. #ifdef EIGEN_TEST_PART_5
  258. double bug_1223() {
  259. using std::min;
  260. typedef Eigen::AutoDiffScalar<Eigen::Vector3d> AD;
  261. const double _cv1_3 = 1.0;
  262. const AD chi_3 = 1.0;
  263. const AD denom = 1.0;
  264. // failed because implementation of min attempts to construct ADS<DerType&> via constructor AutoDiffScalar(const Real& value)
  265. // without initializing m_derivatives (which is a reference in this case)
  266. #define EIGEN_TEST_SPACE
  267. const AD t = min EIGEN_TEST_SPACE (denom / chi_3, 1.0);
  268. const AD t2 = min EIGEN_TEST_SPACE (denom / (chi_3 * _cv1_3), 1.0);
  269. return t.value() + t2.value();
  270. }
  271. // regression test for some compilation issues with specializations of ScalarBinaryOpTraits
  272. void bug_1260() {
  273. Matrix4d A = Matrix4d::Ones();
  274. Vector4d v = Vector4d::Ones();
  275. A*v;
  276. }
  277. // check a compilation issue with numext::max
  278. double bug_1261() {
  279. typedef AutoDiffScalar<Matrix2d> AD;
  280. typedef Matrix<AD,2,1> VectorAD;
  281. VectorAD v(0.,0.);
  282. const AD maxVal = v.maxCoeff();
  283. const AD minVal = v.minCoeff();
  284. return maxVal.value() + minVal.value();
  285. }
  286. double bug_1264() {
  287. typedef AutoDiffScalar<Vector2d> AD;
  288. const AD s = 0.;
  289. const Matrix<AD, 3, 1> v1(0.,0.,0.);
  290. const Matrix<AD, 3, 1> v2 = (s + 3.0) * v1;
  291. return v2(0).value();
  292. }
  293. // check with expressions on constants
  294. double bug_1281() {
  295. int n = 2;
  296. typedef AutoDiffScalar<VectorXd> AD;
  297. const AD c = 1.;
  298. AD x0(2,n,0);
  299. AD y1 = (AD(c)+AD(c))*x0;
  300. y1 = x0 * (AD(c)+AD(c));
  301. AD y2 = (-AD(c))+x0;
  302. y2 = x0+(-AD(c));
  303. AD y3 = (AD(c)*(-AD(c))+AD(c))*x0;
  304. y3 = x0 * (AD(c)*(-AD(c))+AD(c));
  305. return (y1+y2+y3).value();
  306. }
  307. #endif
  308. EIGEN_DECLARE_TEST(autodiff)
  309. {
  310. for(int i = 0; i < g_repeat; i++) {
  311. CALL_SUBTEST_1( test_autodiff_scalar<1>() );
  312. CALL_SUBTEST_2( test_autodiff_vector<1>() );
  313. CALL_SUBTEST_3( test_autodiff_jacobian<1>() );
  314. CALL_SUBTEST_4( test_autodiff_hessian<1>() );
  315. }
  316. CALL_SUBTEST_5( bug_1222() );
  317. CALL_SUBTEST_5( bug_1223() );
  318. CALL_SUBTEST_5( bug_1260() );
  319. CALL_SUBTEST_5( bug_1261() );
  320. CALL_SUBTEST_5( bug_1281() );
  321. }