pinocchio  3.7.0
A fast and flexible implementation of Rigid Body Dynamics algorithms and their analytical derivatives
 
Loading...
Searching...
No Matches
eigenvalues.hpp
1//
2// Copyright (c) 2022-2024 INRIA
3//
4
5#ifndef __pinocchio_math_eigenvalues_hpp__
6#define __pinocchio_math_eigenvalues_hpp__
7
8#include "pinocchio/math/fwd.hpp"
9#include <Eigen/Core>
10
11namespace pinocchio
12{
13
16 template<typename _Vector>
18 {
19 typedef typename PINOCCHIO_EIGEN_PLAIN_TYPE(_Vector) Vector;
20 typedef typename Vector::Scalar Scalar;
21
22 explicit PowerIterationAlgoTpl(
23 const Eigen::DenseIndex size, const int max_it = 10, const Scalar rel_tol = Scalar(1e-8))
24 : principal_eigen_vector(size)
25 , lowest_eigen_vector(size)
26 , max_it(max_it)
27 , it(0)
28 , rel_tol(rel_tol)
29 , eigen_vector_prev(size)
30 {
31 reset();
32 }
33
34 template<typename MatrixLike>
35 void run(const MatrixLike & mat)
36 {
37 PINOCCHIO_CHECK_INPUT_ARGUMENT(max_it >= 1);
38 PINOCCHIO_CHECK_INPUT_ARGUMENT((check_expression_if_real<Scalar, true>(rel_tol > 0)));
39 Scalar eigenvalue_est = principal_eigen_vector.norm();
40
41 for (it = 0; it < max_it; ++it)
42 {
44 principal_eigen_vector /= eigenvalue_est;
45 eigen_vector_prev = principal_eigen_vector;
46 principal_eigen_vector.noalias() = mat * eigen_vector_prev;
47
48 eigenvalue_est = principal_eigen_vector.norm();
49
50 convergence_criteria = math::fabs(eigenvalue_est_prev - eigenvalue_est);
52 convergence_criteria
53 <= rel_tol * math::max(math::fabs(eigenvalue_est_prev), math::fabs(eigenvalue_est))))
54 break;
55 }
56
57 largest_eigen_value = eigenvalue_est;
58 }
59
60 template<typename MatrixLike, typename VectorLike>
61 void run(const MatrixLike & mat, const Eigen::PlainObjectBase<VectorLike> & eigenvector_est)
62 {
63 principal_eigen_vector = eigenvector_est;
64 run(mat);
65 }
66
67 template<typename MatrixLike>
68 void lowest(const MatrixLike & mat, const bool compute_largest = true)
69 {
70 PINOCCHIO_CHECK_INPUT_ARGUMENT(max_it >= 1);
71 PINOCCHIO_CHECK_INPUT_ARGUMENT((check_expression_if_real<Scalar, true>(rel_tol > 0)));
72
74 run(mat);
75
76 Scalar eigenvalue_est = lowest_eigen_vector.norm();
77
78 for (it = 0; it < max_it; ++it)
79 {
81 lowest_eigen_vector /= eigenvalue_est;
82 eigen_vector_prev = lowest_eigen_vector;
83 lowest_eigen_vector.noalias() = mat * eigen_vector_prev;
84 lowest_eigen_vector -= largest_eigen_value * eigen_vector_prev;
85
86 eigenvalue_est = lowest_eigen_vector.norm();
87
88 convergence_criteria = math::fabs(eigenvalue_est_prev - eigenvalue_est);
90 convergence_criteria
91 <= rel_tol * math::max(math::fabs(eigenvalue_est_prev), math::fabs(eigenvalue_est))))
92 break;
93 }
94
95 lowest_eigen_value = largest_eigen_value - eigenvalue_est;
96 }
97
98 template<typename MatrixLike, typename VectorLike>
99 void lowest(
100 const MatrixLike & mat,
101 const Eigen::PlainObjectBase<VectorLike> & largest_eigenvector_est,
102 const Eigen::PlainObjectBase<VectorLike> & lowest_eigenvector_est,
103 const bool compute_largest = true)
104 {
105 principal_eigen_vector = largest_eigenvector_est;
106 lowest_eigen_vector = lowest_eigenvector_est;
107 lowest(mat, compute_largest);
108 }
109
110 void reset()
111 {
112 const Scalar normalized_value = Scalar(1) / math::sqrt(Scalar(principal_eigen_vector.size()));
113 principal_eigen_vector.fill(normalized_value);
114 lowest_eigen_vector.fill(normalized_value);
115
116 largest_eigen_value = std::numeric_limits<Scalar>::min();
117 lowest_eigen_value = std::numeric_limits<Scalar>::max();
118 }
119
120 Vector principal_eigen_vector;
121 Vector lowest_eigen_vector;
122 Scalar largest_eigen_value;
123 Scalar lowest_eigen_value;
124 int max_it;
125 int it;
126 Scalar rel_tol;
127 Scalar convergence_criteria;
128
129 protected:
130 Vector eigen_vector_prev;
131 }; // struct PowerIterationAlgoTpl
132
137 template<typename MatrixLike, typename VectorLike>
139 const MatrixLike & mat,
140 const Eigen::PlainObjectBase<VectorLike> & _eigenvector_est,
141 const int max_it = 10,
142 const typename MatrixLike::Scalar rel_tol = 1e-8)
143 {
144 PowerIterationAlgoTpl<VectorLike> algo(mat.rows(), max_it, rel_tol);
145 algo.run(mat, _eigenvector_est.derived());
146 _eigenvector_est.const_cast_derived() = algo.principal_eigen_vector;
147 }
148
152 template<typename MatrixLike>
153 Eigen::Matrix<typename MatrixLike::Scalar, MatrixLike::RowsAtCompileTime, 1>
155 const MatrixLike & mat, const int max_it = 10, const typename MatrixLike::Scalar rel_tol = 1e-8)
156 {
157 typedef Eigen::Matrix<typename MatrixLike::Scalar, MatrixLike::RowsAtCompileTime, 1> Vector;
158 PowerIterationAlgoTpl<Vector> algo(mat.rows(), max_it, rel_tol);
159 algo.run(mat);
160 return algo.principal_eigen_vector;
161 }
162
167 template<typename MatrixLike, typename VectorLike1, typename VectorLike2>
169 const MatrixLike & mat,
170 const Eigen::PlainObjectBase<VectorLike1> & largest_eigenvector_est,
171 const Eigen::PlainObjectBase<VectorLike2> & lowest_eigenvector_est,
172 const bool compute_largest = true,
173 const int max_it = 10,
174 const typename MatrixLike::Scalar rel_tol = 1e-8)
175 {
176 PowerIterationAlgoTpl<VectorLike1> algo(mat.rows(), max_it, rel_tol);
177 algo.lowest(
179 largest_eigenvector_est.const_cast_derived() = algo.principal_eigen_vector;
180 lowest_eigenvector_est.const_cast_derived() = algo.lowest_eigen_vector;
181 }
182
186 template<typename MatrixLike>
187 Eigen::Matrix<typename MatrixLike::Scalar, MatrixLike::RowsAtCompileTime, 1>
189 const MatrixLike & mat,
190 const bool compute_largest = true,
191 const int max_it = 10,
192 const typename MatrixLike::Scalar rel_tol = 1e-8)
193 {
194 typedef Eigen::Matrix<typename MatrixLike::Scalar, MatrixLike::RowsAtCompileTime, 1> Vector;
195 PowerIterationAlgoTpl<Vector> algo(mat.rows(), max_it, rel_tol);
196 algo.lowest(mat, compute_largest);
197 return algo.lowest_eigen_vector;
198 }
199
204 template<typename VectorLike>
205 typename VectorLike::Scalar
206 retrieveLargestEigenvalue(const Eigen::MatrixBase<VectorLike> & eigenvector)
207 {
208 return eigenvector.norm();
209 }
210} // namespace pinocchio
211
212#endif // #ifndef __pinocchio_math_eigenvalues_hpp__
Main pinocchio namespace.
Definition treeview.dox:11
VectorLike::Scalar retrieveLargestEigenvalue(const Eigen::MatrixBase< VectorLike > &eigenvector)
Compute the lagest eigenvalue of a given matrix. This is taking the eigenvector computed by the funct...
void computeLargestEigenvector(const MatrixLike &mat, const Eigen::PlainObjectBase< VectorLike > &_eigenvector_est, const int max_it=10, const typename MatrixLike::Scalar rel_tol=1e-8)
Compute the lagest eigenvector of a given matrix according to a given eigenvector estimate.
void computeLowestEigenvector(const MatrixLike &mat, const Eigen::PlainObjectBase< VectorLike1 > &largest_eigenvector_est, const Eigen::PlainObjectBase< VectorLike2 > &lowest_eigenvector_est, const bool compute_largest=true, const int max_it=10, const typename MatrixLike::Scalar rel_tol=1e-8)
Compute the lagest eigenvector of a given matrix according to a given eigenvector estimate.
Compute the largest eigenvalues and the associated principle eigenvector via power iteration.