5 #ifndef __pinocchio_autodiff_casadi_hpp__
6 #define __pinocchio_autodiff_casadi_hpp__
8 #define PINOCCHIO_WITH_CASADI_SUPPORT
12 #include <casadi/casadi.hpp>
14 #include <Eigen/Sparse>
25 struct constant_pi<::casadi::SX> : constant_pi<double>
29 template<
typename Scalar>
30 struct constant_pi<::casadi::Matrix<Scalar>> : constant_pi<Scalar>
41 inline bool operator||(
const bool x, const ::casadi::Matrix<SXElem> & )
49 template<
typename Scalar>
57 return ::casadi::Matrix<Scalar>(Base::template precision<degree>());
68 template<
typename Scalar>
69 struct cast_impl<::casadi::Matrix<Scalar>,
Scalar>
71 #if EIGEN_VERSION_AT_LEAST(3, 2, 90)
74 static inline Scalar run(const ::casadi::Matrix<Scalar> & x)
80 #if EIGEN_VERSION_AT_LEAST(3, 2, 90) && !EIGEN_VERSION_AT_LEAST(3, 2, 93)
81 template<
typename Scalar,
bool IsInteger>
82 struct significant_decimals_default_impl<::casadi::Matrix<Scalar>, IsInteger>
84 static inline int run()
86 return std::numeric_limits<Scalar>::digits10;
97 template<
typename Scalar>
98 struct NumTraits<::casadi::Matrix<Scalar>>
100 using Real = ::casadi::Matrix<Scalar>;
114 RequireInitialization = 1,
123 return ::casadi::Matrix<Scalar>(std::numeric_limits<double>::epsilon());
128 return ::casadi::Matrix<Scalar>(NumTraits<double>::dummy_precision());
143 return std::numeric_limits<double>::digits10;
153 template<
typename MT,
typename Scalar>
154 inline void copy(::casadi::Matrix<Scalar>
const & src, Eigen::MatrixBase<MT> & dst)
156 Eigen::DenseIndex
const m = src.size1();
157 Eigen::DenseIndex
const n = src.size2();
161 for (Eigen::DenseIndex
i = 0;
i <
m; ++
i)
162 for (Eigen::DenseIndex j = 0; j <
n; ++j)
163 dst(
i, j) = src(
i, j);
167 template<
typename TensorDerived,
typename Scalar>
168 inline void copy(::casadi::Matrix<Scalar>
const & src, Eigen::TensorBase<TensorDerived> & dst_)
170 TensorDerived & dst =
static_cast<TensorDerived &
>(dst_);
179 dst(0,
i, j) = src(
i, j);
183 template<
typename TensorDerived,
typename Scalar>
184 inline void copy(Eigen::TensorBase<TensorDerived>
const & _src, ::casadi::Matrix<Scalar> & dst)
186 const TensorDerived & src =
static_cast<const TensorDerived &
>(_src);
195 dst(
i, j) = src(0,
i, j);
199 template<
typename MT,
typename Scalar>
200 inline void copy(Eigen::MatrixBase<MT>
const & src, ::casadi::Matrix<Scalar> & dst)
202 Eigen::DenseIndex
const m = src.rows();
203 Eigen::DenseIndex
const n = src.cols();
207 for (Eigen::DenseIndex
i = 0;
i <
m; ++
i)
208 for (Eigen::DenseIndex j = 0; j <
n; ++j)
209 dst(
i, j) = src(
i, j);
213 template<
typename MT,
typename Scalar>
214 inline void copy(Eigen::SparseMatrixBase<MT>
const & src, ::casadi::Matrix<Scalar> & dst)
216 Eigen::DenseIndex
const m = src.rows();
217 Eigen::DenseIndex
const n = src.cols();
225 dst(
i, j) = src.derived().coeff(
i, j);
229 template<
typename MatrixDerived>
230 inline void sym(
const Eigen::MatrixBase<MatrixDerived> & eig_mat, std::string
const &
name)
235 for (Eigen::DenseIndex
i = 0;
i < eig_mat.rows(); ++
i)
236 for (Eigen::DenseIndex j = 0; j < eig_mat.cols(); ++j)
237 eig_mat_(
i, j) =
SX::sym(
name +
"_" + std::to_string(
i) +
"_" + std::to_string(j));
250 template<
typename Scalar>
251 struct return_type_min<::casadi::Matrix<Scalar>, ::casadi::Matrix<Scalar>>
253 typedef ::casadi::Matrix<Scalar>
type;
256 template<
typename Scalar,
typename T>
257 struct return_type_min<::casadi::Matrix<Scalar>,
T>
259 typedef ::casadi::Matrix<Scalar>
type;
262 template<
typename Scalar,
typename T>
263 struct return_type_min<
T, ::casadi::Matrix<Scalar>>
265 typedef ::casadi::Matrix<Scalar>
type;
268 template<
typename Scalar>
269 struct call_min<::casadi::Matrix<Scalar>, ::casadi::Matrix<Scalar>>
271 static inline ::casadi::Matrix<Scalar>
272 run(const ::casadi::Matrix<Scalar> &
a, const ::casadi::Matrix<Scalar> &
b)
278 template<
typename S1,
typename S2>
279 struct call_min<::casadi::Matrix<S1>, S2>
282 static inline ::casadi::Matrix<S1>
run(const ::casadi::Matrix<S1> &
a,
const S2 &
b)
288 template<
typename S1,
typename S2>
289 struct call_min<S1, ::casadi::Matrix<S2>>
292 static inline ::casadi::Matrix<S2>
run(
const S1 &
a, const ::casadi::Matrix<S2> &
b)
298 template<
typename Scalar>
299 struct return_type_max<::casadi::Matrix<Scalar>, ::casadi::Matrix<Scalar>>
301 typedef ::casadi::Matrix<Scalar>
type;
304 template<
typename Scalar,
typename T>
305 struct return_type_max<::casadi::Matrix<Scalar>,
T>
307 typedef ::casadi::Matrix<Scalar>
type;
310 template<
typename Scalar,
typename T>
311 struct return_type_max<
T, ::casadi::Matrix<Scalar>>
313 typedef ::casadi::Matrix<Scalar>
type;
316 template<
typename Scalar>
317 struct call_max<::casadi::Matrix<Scalar>, ::casadi::Matrix<Scalar>>
319 static inline ::casadi::Matrix<Scalar>
320 run(const ::casadi::Matrix<Scalar> &
a, const ::casadi::Matrix<Scalar> &
b)
326 template<
typename S1,
typename S2>
327 struct call_max<::casadi::Matrix<S1>, S2>
330 static inline ::casadi::Matrix<S1>
run(const ::casadi::Matrix<S1> &
a,
const S2 &
b)
336 template<
typename S1,
typename S2>
337 struct call_max<S1, ::casadi::Matrix<S2>>
340 static inline ::casadi::Matrix<S2>
run(
const S1 &
a, const ::casadi::Matrix<S2> &
b)
355 template<
typename Scalar>
356 struct ScalarBinaryOpTraits<
357 ::casadi::Matrix<Scalar>,
358 ::casadi::Matrix<Scalar>,
359 internal::scalar_max_op<::casadi::Matrix<Scalar>, ::casadi::Matrix<Scalar>>>
364 template<
typename S1,
typename S2>
365 struct ScalarBinaryOpTraits<
366 ::casadi::Matrix<S1>,
368 internal::scalar_max_op<::casadi::Matrix<S1>, S2>>
373 template<
typename S1,
typename S2>
374 struct ScalarBinaryOpTraits<
376 ::casadi::Matrix<S2>,
377 internal::scalar_max_op<S1, ::casadi::Matrix<S2>>>
385 template<
typename Scalar>
386 struct scalar_max_op<::casadi::Matrix<Scalar>, ::casadi::Matrix<Scalar>>
387 : binary_op_base<::casadi::Matrix<Scalar>, ::casadi::Matrix<Scalar>>
392 typename ScalarBinaryOpTraits<LhsScalar, RhsScalar, scalar_max_op>::ReturnType
result_type;
394 EIGEN_EMPTY_STRUCT_CTOR(scalar_max_op)
407 template<
typename S1,
typename S2>
408 struct scalar_max_op<S1, ::casadi::Matrix<S2>> : binary_op_base<S1, ::casadi::Matrix<S2>>
413 typename ScalarBinaryOpTraits<LhsScalar, RhsScalar, scalar_max_op>::ReturnType
result_type;
415 EIGEN_EMPTY_STRUCT_CTOR(scalar_max_op)
423 template<
typename S1,
typename S2>
424 struct scalar_max_op<::casadi::Matrix<S1>, S2> : binary_op_base<::casadi::Matrix<S1>, S2>
429 typename ScalarBinaryOpTraits<LhsScalar, RhsScalar, scalar_max_op>::ReturnType
result_type;
431 EIGEN_EMPTY_STRUCT_CTOR(scalar_max_op)
447 #endif // #ifndef __pinocchio_autodiff_casadi_hpp__