5 #ifndef __pinocchio_math_tridiagonal_matrix_hpp__
6 #define __pinocchio_math_tridiagonal_matrix_hpp__
13 template<
typename Scalar,
int Options = 0>
16 template<
typename _Scalar,
int _Options>
24 typedef Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic, Options>
PlainMatrixType;
27 template<
typename Tr
idiagonalSymmetricMatrix,
typename MatrixDerived>
30 template<
typename Tr
idiagonalSymmetricMatrix,
typename MatrixDerived>
33 :
public traits<TridiagonalSymmetricMatrix>
37 template<
typename MatrixDerived,
typename Tr
idiagonalSymmetricMatrix>
40 template<
typename MatrixDerived,
typename Tr
idiagonalSymmetricMatrix>
43 :
public traits<TridiagonalSymmetricMatrix>
47 template<
typename Tr
idiagonalSymmetricMatrix>
50 template<
typename Tr
idiagonalSymmetricMatrix>
52 :
public traits<TridiagonalSymmetricMatrix>
56 template<
typename Tr
idiagonalSymmetricMatrixInverse,
typename MatrixDerived>
59 template<
typename Tr
idiagonalSymmetricMatrixInverse,
typename MatrixDerived>
62 MatrixDerived>> :
public traits<TridiagonalSymmetricMatrixInverse>
72 template<
typename Scalar,
int Options>
74 :
public traits<typename pinocchio::traits<
75 pinocchio::TridiagonalSymmetricMatrixTpl<Scalar, Options>>::PlainMatrixType>
85 template<
typename Tr
idiagonalSymmetricMatrix,
typename MatrixDerived>
87 TridiagonalSymmetricMatrix,
90 typename pinocchio::traits<pinocchio::TridiagonalSymmetricMatrixApplyOnTheRightReturnType<
91 TridiagonalSymmetricMatrix,
92 MatrixDerived>>::PlainMatrixType>
95 TridiagonalSymmetricMatrix,
105 template<
typename MatrixDerived,
typename Tr
idiagonalSymmetricMatrix>
108 TridiagonalSymmetricMatrix>>
110 typename pinocchio::traits<pinocchio::TridiagonalSymmetricMatrixApplyOnTheLeftReturnType<
112 TridiagonalSymmetricMatrix>>::PlainMatrixType>
116 TridiagonalSymmetricMatrix>>
125 template<
typename Tr
idiagonalSymmetricMatrix>
126 struct traits<
pinocchio::TridiagonalSymmetricMatrixInverse<TridiagonalSymmetricMatrix>>
127 :
public traits<TridiagonalSymmetricMatrix>
131 template<
typename Tr
idiagonalSymmetricMatrixInverse,
typename MatrixDerived>
133 TridiagonalSymmetricMatrixInverse,
135 :
public traits<typename pinocchio::traits<
136 pinocchio::TridiagonalSymmetricMatrixInverseApplyOnTheRightReturnType<
137 TridiagonalSymmetricMatrixInverse,
138 MatrixDerived>>::PlainMatrixType>
142 TridiagonalSymmetricMatrixInverse,
160 template<
typename _Scalar,
int _Options>
161 struct TridiagonalSymmetricMatrixTpl
162 :
public Eigen::ReturnByValue<TridiagonalSymmetricMatrixTpl<_Scalar, _Options>>
164 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
181 assert(
size > 0 &&
"size should be greater than 0.");
193 return !(*
this == other);
235 bool isZero(
const Scalar prec = Eigen::NumTraits<Scalar>::dummy_precision())
const
251 template<
typename VectorCoeffType>
252 void setDiagonal(
const Eigen::MatrixBase<VectorCoeffType> & diagonal_coefficients)
256 VectorCoeffType::IsVectorAtCompileTime,
257 "VectorCoeffType should be a vector type at compile time");
277 template<
typename ResultType>
278 inline void evalTo(ResultType & result)
const
281 result.template diagonal<1>() =
subDiagonal().conjugate();
286 template<
typename MatrixDerived>
291 return ReturnType(*
this,
mat.derived());
294 template<
typename MatrixDerived>
299 return ReturnType(
mat.derived(), *
this);
302 template<
typename MatrixDerived>
304 operator*(
const Eigen::MatrixBase<MatrixDerived> & mat)
const
315 template<
typename LhsMatrixType,
typename S,
int O>
322 return rhs.applyOnTheLeft(lhs);
325 template<
typename Tr
idiagonalSymmetricMatrix,
typename RhsMatrixType>
326 struct TridiagonalSymmetricMatrixApplyOnTheRightReturnType
327 :
public Eigen::ReturnByValue<TridiagonalSymmetricMatrixApplyOnTheRightReturnType<
328 TridiagonalSymmetricMatrix,
335 const TridiagonalSymmetricMatrix & lhs,
const RhsMatrixType &
rhs)
341 template<
typename ResultType>
342 inline void evalTo(ResultType & result)
const
350 const Eigen::DenseIndex reduced_size =
cols() - 1;
352 result.noalias() =
m_lhs.diagonal().asDiagonal() *
m_rhs;
354 result.topRows(reduced_size).noalias() +=
355 m_lhs.subDiagonal().conjugate().asDiagonal() *
m_rhs.bottomRows(reduced_size);
357 result.bottomRows(reduced_size).noalias() +=
358 m_lhs.subDiagonal().asDiagonal() *
m_rhs.topRows(reduced_size);
371 const TridiagonalSymmetricMatrix &
m_lhs;
375 template<
typename LhsMatrixType,
typename Tr
idiagonalSymmetricMatrix>
377 :
public Eigen::ReturnByValue<
378 TridiagonalSymmetricMatrixApplyOnTheLeftReturnType<LhsMatrixType, TridiagonalSymmetricMatrix>>
384 const LhsMatrixType & lhs,
const TridiagonalSymmetricMatrix &
rhs)
390 template<
typename ResultType>
391 inline void evalTo(ResultType & result)
const
399 const Eigen::DenseIndex reduced_size =
cols() - 1;
401 result.noalias() =
m_lhs *
m_rhs.diagonal().asDiagonal();
403 result.rightCols(reduced_size).noalias() +=
404 m_lhs.leftCols(reduced_size) *
m_rhs.subDiagonal().conjugate().asDiagonal();
406 result.leftCols(reduced_size).noalias() +=
407 m_lhs.rightCols(reduced_size) *
m_rhs.subDiagonal().asDiagonal();
421 const TridiagonalSymmetricMatrix &
m_rhs;
424 template<
typename _Tr
idiagonalSymmetricMatrix>
426 :
public Eigen::ReturnByValue<TridiagonalSymmetricMatrixInverse<_TridiagonalSymmetricMatrix>>
454 template<
typename MatrixDerived>
460 return ReturnType(*
this,
mat.derived());
463 template<
typename MatrixDerived>
465 operator*(
const Eigen::MatrixBase<MatrixDerived> & mat)
const
470 template<
typename ResultType>
471 inline void evalTo(ResultType & result)
const
484 result.setIdentity();
485 for (Eigen::DenseIndex
i = 1;
i <
m_size; ++
i)
487 result.row(
i).head(
i) -=
w[
i - 1] * result.row(
i - 1).head(
i);
492 for (Eigen::DenseIndex
i =
m_size - 2;
i >= 0; --
i)
494 result.row(
i) -=
c[
i] * result.row(
i + 1);
495 result.row(
i) /=
b[
i];
509 template<
typename T,
typename MatrixDerived>
520 for (Eigen::DenseIndex
i = 1;
i <
m_size; ++
i)
522 w.coeffRef(
i - 1) /=
b[
i - 1];
533 template<
typename Tr
idiagonalSymmetricMatrixInverse,
typename RhsMatrixType>
535 :
public Eigen::ReturnByValue<TridiagonalSymmetricMatrixInverseApplyOnTheRightReturnType<
536 TridiagonalSymmetricMatrixInverse,
549 template<
typename ResultType>
550 inline void evalTo(ResultType & result)
const
565 for (Eigen::DenseIndex
i = 1;
i <
size; ++
i)
567 result.row(
i) -=
w[
i - 1] * result.row(
i - 1);
572 for (Eigen::DenseIndex
i =
size - 2;
i >= 0; --
i)
574 result.row(
i) -=
c[
i] * result.row(
i + 1);
575 result.row(
i) /=
b[
i];
594 #endif // #ifndef __pinocchio_math_tridiagonal_matrix_hpp__