forward_adolc.cpp
Go to the documentation of this file.
00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2008 Gael Guennebaud <g.gael@free.fr>
00005 //
00006 // This Source Code Form is subject to the terms of the Mozilla
00007 // Public License v. 2.0. If a copy of the MPL was not distributed
00008 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
00009 
00010 #include "main.h"
00011 #include <Eigen/Dense>
00012 
00013 #define NUMBER_DIRECTIONS 16
00014 #include <unsupported/Eigen/AdolcForward>
00015 
00016 int adtl::ADOLC_numDir;
00017 
00018 template<typename Vector>
00019 EIGEN_DONT_INLINE typename Vector::Scalar foo(const Vector& p)
00020 {
00021   typedef typename Vector::Scalar Scalar;
00022   return (p-Vector(Scalar(-1),Scalar(1.))).norm() + (p.array().sqrt().abs() * p.array().sin()).sum() + p.dot(p);
00023 }
00024 
00025 template<typename _Scalar, int NX=Dynamic, int NY=Dynamic>
00026 struct TestFunc1
00027 {
00028   typedef _Scalar Scalar;
00029   enum {
00030     InputsAtCompileTime = NX,
00031     ValuesAtCompileTime = NY
00032   };
00033   typedef Matrix<Scalar,InputsAtCompileTime,1> InputType;
00034   typedef Matrix<Scalar,ValuesAtCompileTime,1> ValueType;
00035   typedef Matrix<Scalar,ValuesAtCompileTime,InputsAtCompileTime> JacobianType;
00036 
00037   int m_inputs, m_values;
00038 
00039   TestFunc1() : m_inputs(InputsAtCompileTime), m_values(ValuesAtCompileTime) {}
00040   TestFunc1(int inputs, int values) : m_inputs(inputs), m_values(values) {}
00041 
00042   int inputs() const { return m_inputs; }
00043   int values() const { return m_values; }
00044 
00045   template<typename T>
00046   void operator() (const Matrix<T,InputsAtCompileTime,1>& x, Matrix<T,ValuesAtCompileTime,1>* _v) const
00047   {
00048     Matrix<T,ValuesAtCompileTime,1>& v = *_v;
00049 
00050     v[0] = 2 * x[0] * x[0] + x[0] * x[1];
00051     v[1] = 3 * x[1] * x[0] + 0.5 * x[1] * x[1];
00052     if(inputs()>2)
00053     {
00054       v[0] += 0.5 * x[2];
00055       v[1] += x[2];
00056     }
00057     if(values()>2)
00058     {
00059       v[2] = 3 * x[1] * x[0] * x[0];
00060     }
00061     if (inputs()>2 && values()>2)
00062       v[2] *= x[2];
00063   }
00064 
00065   void operator() (const InputType& x, ValueType* v, JacobianType* _j) const
00066   {
00067     (*this)(x, v);
00068 
00069     if(_j)
00070     {
00071       JacobianType& j = *_j;
00072 
00073       j(0,0) = 4 * x[0] + x[1];
00074       j(1,0) = 3 * x[1];
00075 
00076       j(0,1) = x[0];
00077       j(1,1) = 3 * x[0] + 2 * 0.5 * x[1];
00078 
00079       if (inputs()>2)
00080       {
00081         j(0,2) = 0.5;
00082         j(1,2) = 1;
00083       }
00084       if(values()>2)
00085       {
00086         j(2,0) = 3 * x[1] * 2 * x[0];
00087         j(2,1) = 3 * x[0] * x[0];
00088       }
00089       if (inputs()>2 && values()>2)
00090       {
00091         j(2,0) *= x[2];
00092         j(2,1) *= x[2];
00093 
00094         j(2,2) = 3 * x[1] * x[0] * x[0];
00095         j(2,2) = 3 * x[1] * x[0] * x[0];
00096       }
00097     }
00098   }
00099 };
00100 
00101 template<typename Func> void adolc_forward_jacobian(const Func& f)
00102 {
00103     typename Func::InputType x = Func::InputType::Random(f.inputs());
00104     typename Func::ValueType y(f.values()), yref(f.values());
00105     typename Func::JacobianType j(f.values(),f.inputs()), jref(f.values(),f.inputs());
00106 
00107     jref.setZero();
00108     yref.setZero();
00109     f(x,&yref,&jref);
00110 //     std::cerr << y.transpose() << "\n\n";;
00111 //     std::cerr << j << "\n\n";;
00112 
00113     j.setZero();
00114     y.setZero();
00115     AdolcForwardJacobian<Func> autoj(f);
00116     autoj(x, &y, &j);
00117 //     std::cerr << y.transpose() << "\n\n";;
00118 //     std::cerr << j << "\n\n";;
00119 
00120     VERIFY_IS_APPROX(y, yref);
00121     VERIFY_IS_APPROX(j, jref);
00122 }
00123 
00124 void test_forward_adolc()
00125 {
00126   adtl::ADOLC_numDir = NUMBER_DIRECTIONS;
00127 
00128   for(int i = 0; i < g_repeat; i++) {
00129     CALL_SUBTEST(( adolc_forward_jacobian(TestFunc1<double,2,2>()) ));
00130     CALL_SUBTEST(( adolc_forward_jacobian(TestFunc1<double,2,3>()) ));
00131     CALL_SUBTEST(( adolc_forward_jacobian(TestFunc1<double,3,2>()) ));
00132     CALL_SUBTEST(( adolc_forward_jacobian(TestFunc1<double,3,3>()) ));
00133     CALL_SUBTEST(( adolc_forward_jacobian(TestFunc1<double>(3,3)) ));
00134   }
00135 
00136   {
00137     // simple instanciation tests
00138     Matrix<adtl::adouble,2,1> x;
00139     foo(x);
00140     Matrix<adtl::adouble,Dynamic,Dynamic> A(4,4);;
00141     A.selfadjointView<Lower>().eigenvalues();
00142   }
00143 }


shape_reconstruction
Author(s): Roberto Martín-Martín
autogenerated on Sat Jun 8 2019 18:31:02