21 #include <boost/test/unit_test.hpp>
22 #include <boost/utility/binary.hpp>
25 #define DLL_EXT ".dll"
30 BOOST_AUTO_TEST_SUITE(BOOST_TEST_MODULE)
35 using CppAD::NearEqual;
38 typedef AD<Scalar> ADScalar;
48 model.lowerPositionLimit.head<3>().
fill(-1.);
49 model.upperPositionLimit.head<3>().
fill(1.);
52 ADModel ad_model =
model.cast<ADScalar>();
53 ADData ad_data(ad_model);
56 typedef Model::ConfigVectorType ConfigVectorType;
57 typedef Model::TangentVectorType TangentVectorType;
58 ConfigVectorType
q(
model.nq);
61 TangentVectorType
v(TangentVectorType::Random(
model.nv));
62 TangentVectorType
a(TangentVectorType::Random(
model.nv));
64 typedef ADModel::ConfigVectorType ADConfigVectorType;
65 typedef ADModel::TangentVectorType ADTangentVectorType;
67 ADConfigVectorType ad_q =
q.cast<ADScalar>();
68 ADTangentVectorType ad_v =
v.cast<ADScalar>();
69 ADTangentVectorType ad_a =
a.cast<ADScalar>();
71 typedef Eigen::Matrix<ADScalar, Eigen::Dynamic, 1> VectorXAD;
73 data.M.triangularView<Eigen::StrictlyLower>() =
74 data.M.transpose().triangularView<Eigen::StrictlyLower>();
79 CppAD::Independent(ad_a);
83 Eigen::Map<ADData::TangentVectorType>(
Y.data(),
model.nv, 1) = ad_data.tau;
85 CppAD::ADFun<Scalar> ad_fun(ad_a,
Y);
88 Eigen::Map<Data::TangentVectorType>(
x.data(),
model.nv, 1) =
a;
90 CPPAD_TESTVECTOR(
Scalar)
tau = ad_fun.Forward(0,
x);
91 BOOST_CHECK(Eigen::Map<Data::TangentVectorType>(
tau.data(),
model.nv, 1).isApprox(
data.tau));
94 Data::MatrixXs M = Eigen::Map<PINOCCHIO_EIGEN_PLAIN_ROW_MAJOR_TYPE(Data::MatrixXs)>(
96 BOOST_CHECK(
M.isApprox(
data.M));
99 ADTangentVectorType ad_tau =
tau.cast<ADScalar>();
102 data.Minv.triangularView<Eigen::StrictlyLower>() =
103 data.Minv.transpose().triangularView<Eigen::StrictlyLower>();
107 CppAD::Independent(ad_tau);
111 Eigen::Map<ADData::TangentVectorType>(
Y.data(),
model.nv, 1) = ad_data.ddq;
113 CppAD::ADFun<Scalar> ad_fun(ad_tau,
Y);
116 Eigen::Map<Data::TangentVectorType>(
x.data(),
model.nv, 1) =
tau;
118 CPPAD_TESTVECTOR(
Scalar) ddq = ad_fun.Forward(0,
x);
119 BOOST_CHECK(Eigen::Map<Data::TangentVectorType>(ddq.data(),
model.nv, 1).isApprox(
a));
121 CPPAD_TESTVECTOR(
Scalar) dddq_da = ad_fun.Jacobian(
x);
122 Data::MatrixXs Minv = Eigen::Map<PINOCCHIO_EIGEN_PLAIN_ROW_MAJOR_TYPE(Data::MatrixXs)>(
124 BOOST_CHECK(Minv.isApprox(
data.Minv));
131 using CppAD::NearEqual;
134 typedef AD<Scalar> ADScalar;
145 model.lowerPositionLimit.head<3>().
fill(-1.);
146 model.upperPositionLimit.head<3>().
fill(1.);
149 ADModel ad_model =
model.cast<ADScalar>();
150 ADData ad_data(ad_model);
153 typedef Model::ConfigVectorType ConfigVectorType;
154 typedef Model::TangentVectorType TangentVectorType;
155 ConfigVectorType
q(
model.nq);
158 TangentVectorType
v(TangentVectorType::Random(
model.nv));
159 TangentVectorType
a(TangentVectorType::Random(
model.nv));
161 typedef ADModel::ConfigVectorType ADConfigVectorType;
162 typedef ADModel::TangentVectorType ADTangentVectorType;
164 ADConfigVectorType ad_q =
q.cast<ADScalar>();
165 ADTangentVectorType ad_v =
v.cast<ADScalar>();
166 ADTangentVectorType ad_a =
a.cast<ADScalar>();
189 typedef Eigen::Matrix<ADScalar, Eigen::Dynamic, 1> VectorXAD;
192 CppAD::Independent(ad_v);
198 Eigen::DenseIndex current_id = 0;
199 for (Eigen::DenseIndex k = 0; k < 3; ++k)
201 Y[current_id + k + Motion::LINEAR] = v_local.linear()[k];
202 Y[current_id + k + Motion::ANGULAR] = v_local.angular()[k];
206 for (Eigen::DenseIndex k = 0; k < 3; ++k)
208 Y[current_id + k + Motion::LINEAR] = v_global.linear()[k];
209 Y[current_id + k + Motion::ANGULAR] = v_global.angular()[k];
213 for (Eigen::DenseIndex k = 0; k < 3; ++k)
215 Y[current_id + k + Motion::LINEAR] = a_local.linear()[k];
216 Y[current_id + k + Motion::ANGULAR] = a_local.angular()[k];
220 CppAD::ADFun<Scalar> vjoint(ad_v,
Y);
223 for (Eigen::DenseIndex k = 0; k <
model.nv; ++k)
228 CPPAD_TESTVECTOR(
Scalar)
y = vjoint.Forward(0,
x);
234 .isApprox(
Motion(Eigen::Map<Motion::Vector6>(y_ptr))));
239 CPPAD_TESTVECTOR(
Scalar) dY_dv = vjoint.Jacobian(
x);
241 Scalar * dY_dv_ptr = dY_dv.data();
243 Eigen::Map<PINOCCHIO_EIGEN_PLAIN_ROW_MAJOR_TYPE(Data::Matrix6x)>(dY_dv_ptr, 6,
model.nv);
244 dY_dv_ptr += ad_J_local.size();
246 Eigen::Map<PINOCCHIO_EIGEN_PLAIN_ROW_MAJOR_TYPE(Data::Matrix6x)>(dY_dv_ptr, 6,
model.nv);
247 dY_dv_ptr += ad_J_global.size();
249 BOOST_CHECK(ad_J_local.isApprox(J_local));
250 BOOST_CHECK(ad_J_global.isApprox(J_global));
254 CppAD::Independent(ad_a);
258 Eigen::DenseIndex current_id = 0;
259 for (Eigen::DenseIndex k = 0; k < 3; ++k)
261 Y[current_id + k + Motion::LINEAR] = v_local.linear()[k];
262 Y[current_id + k + Motion::ANGULAR] = v_local.angular()[k];
266 for (Eigen::DenseIndex k = 0; k < 3; ++k)
268 Y[current_id + k + Motion::LINEAR] = a_local.linear()[k];
269 Y[current_id + k + Motion::ANGULAR] = a_local.angular()[k];
273 CppAD::ADFun<Scalar> ajoint(ad_a,
Y);
276 for (Eigen::DenseIndex k = 0; k <
model.nv; ++k)
281 CPPAD_TESTVECTOR(
Scalar)
y = ajoint.Forward(0,
x);
286 CPPAD_TESTVECTOR(
Scalar) dY_da = ajoint.Jacobian(
x);
288 Scalar * dY_da_ptr = dY_da.data();
290 Eigen::Map<PINOCCHIO_EIGEN_PLAIN_ROW_MAJOR_TYPE(Data::Matrix6x)>(dY_da_ptr, 6,
model.nv);
291 dY_da_ptr += ad_dv_da.size();
293 Eigen::Map<PINOCCHIO_EIGEN_PLAIN_ROW_MAJOR_TYPE(Data::Matrix6x)>(dY_da_ptr, 6,
model.nv);
294 dY_da_ptr += ad_J_local.size();
296 BOOST_CHECK(ad_dv_da.isZero());
297 BOOST_CHECK(ad_J_local.isApprox(J_local));
304 using CppAD::NearEqual;
307 typedef AD<Scalar> ADScalar;
317 model.lowerPositionLimit.head<3>().
fill(-1.);
318 model.upperPositionLimit.head<3>().
fill(1.);
321 ADModel ad_model =
model.cast<ADScalar>();
322 ADData ad_data(ad_model);
325 typedef Model::ConfigVectorType ConfigVectorType;
326 typedef Model::TangentVectorType TangentVectorType;
327 ConfigVectorType
q(
model.nq);
330 TangentVectorType
v(TangentVectorType::Random(
model.nv));
331 TangentVectorType
a(TangentVectorType::Random(
model.nv));
333 typedef ADModel::ConfigVectorType ADConfigVectorType;
334 typedef ADModel::TangentVectorType ADTangentVectorType;
336 ADConfigVectorType ad_q =
q.cast<ADScalar>();
337 ADTangentVectorType ad_v =
v.cast<ADScalar>();
338 ADTangentVectorType ad_a =
a.cast<ADScalar>();
340 typedef Eigen::Matrix<ADScalar, Eigen::Dynamic, 1> VectorXAD;
341 typedef Eigen::Matrix<ADScalar, Eigen::Dynamic, Eigen::Dynamic> MatrixXAD;
343 data.M.triangularView<Eigen::StrictlyLower>() =
344 data.M.transpose().triangularView<Eigen::StrictlyLower>();
348 CppAD::Independent(ad_a);
352 Eigen::Map<ADData::TangentVectorType>(
Y.data(),
model.nv, 1) = ad_data.tau;
354 CppAD::ADFun<Scalar>
f(ad_a,
Y);
355 CppAD::ADFun<ADScalar, Scalar> af =
f.base2ad();
357 CppAD::Independent(ad_a);
358 MatrixXAD
dtau_da = af.Jacobian(ad_a);
359 VectorXAD dtau_da_vector(
model.nv *
model.nv);
361 CppAD::ADFun<double> ad_fun(ad_a, dtau_da_vector);
363 ad_fun.function_name_set(
"ad_fun");
366 std::string c_type =
"double";
367 std::string csrc_file =
"jit_JSIM.c";
369 ofs.open(csrc_file, std::ofstream::out);
370 ad_fun.to_csrc(ofs, c_type);
374 std::string dll_file =
"./jit_JSIM" DLL_EXT;
375 CPPAD_TESTVECTOR(std::string) csrc_files(1);
376 csrc_files[0] = csrc_file;
377 std::map<std::string, std::string> options;
378 std::string err_msg = CppAD::create_dll_lib(dll_file, csrc_files, options);
381 std::cerr <<
"jit_JSIM: err_msg = " << err_msg <<
"\n";
385 CppAD::link_dll_lib dll_linker(dll_file, err_msg);
388 std::cerr <<
"jit_JSIM: err_msg = " << err_msg <<
"\n";
391 std::string function_name =
"cppad_jit_ad_fun";
392 void * void_ptr = dll_linker(function_name, err_msg);
395 std::cerr <<
"jit_JSIM: err_msg = " << err_msg <<
"\n";
398 using CppAD::jit_double;
399 jit_double ad_fun_ptr =
reinterpret_cast<jit_double
>(void_ptr);
402 Eigen::Map<Data::TangentVectorType>(
x.data(),
model.nv, 1) =
a;
404 size_t compare_change = 0,
nx = (size_t)
model.nv,
405 ndtau_da_ = (
size_t)
model.nv * (size_t)
model.nv;
406 std::vector<double> dtau_da_jit(ndtau_da_);
407 ad_fun_ptr(
nx,
x.data(), ndtau_da_, dtau_da_jit.data(), &compare_change);
409 Data::MatrixXs M = Eigen::Map<PINOCCHIO_EIGEN_PLAIN_ROW_MAJOR_TYPE(Data::MatrixXs)>(
411 BOOST_CHECK(
M.isApprox(
data.M));
415 BOOST_AUTO_TEST_SUITE_END()