v2.0.0
Loading...
Searching...
No Matches
linalg.cpp
Go to the documentation of this file.
1//=============================================================================================================
34
35//=============================================================================================================
36// INCLUDES
37//=============================================================================================================
38
39#include "linalg.h"
40
41//=============================================================================================================
42// EIGEN INCLUDES
43//=============================================================================================================
44
45#include <Eigen/Eigen>
46
47//=============================================================================================================
48// QT INCLUDES
49//=============================================================================================================
50
51#include <QDebug>
52
53//=============================================================================================================
54// USED NAMESPACES
55//=============================================================================================================
56
57using namespace UTILSLIB;
58using namespace Eigen;
59
60//=============================================================================================================
61// DEFINE MEMBER METHODS
62//=============================================================================================================
63
64VectorXd Linalg::combine_xyz(const VectorXd& vec)
65{
66 if (vec.size() % 3 != 0)
67 {
68 qWarning("Linalg::combine_xyz: Input must be a row or column vector with 3N components.");
69 return VectorXd();
70 }
71
72 MatrixXd tmp = MatrixXd(vec.transpose());
73 SparseMatrix<double> s = make_block_diag(tmp, 3);
74
75 SparseMatrix<double> sC = s * s.transpose();
76 VectorXd comb(sC.rows());
77
78 for (qint32 i = 0; i < sC.rows(); ++i)
79 comb[i] = sC.coeff(i, i);
80
81 return comb;
82}
83
84//=============================================================================================================
85
86double Linalg::getConditionNumber(const MatrixXd& A,
87 VectorXd &s)
88{
89 JacobiSVD<MatrixXd> svd(A);
90 s = svd.singularValues();
91
92 double c = s.maxCoeff()/s.minCoeff();
93
94 return c;
95}
96
97//=============================================================================================================
98
99double Linalg::getConditionSlope(const MatrixXd& A,
100 VectorXd &s)
101{
102 JacobiSVD<MatrixXd> svd(A);
103 s = svd.singularValues();
104
105 double c = s.maxCoeff()/s.mean();
106
107 return c;
108}
109
110//=============================================================================================================
111
112void Linalg::get_whitener(MatrixXd &A,
113 bool pca,
114 QString ch_type,
115 VectorXd &eig,
116 MatrixXd &eigvec)
117{
118 SelfAdjointEigenSolver<MatrixXd> t_eigenSolver(A);
119
120 eig = t_eigenSolver.eigenvalues();
121 eigvec = t_eigenSolver.eigenvectors().transpose();
122
123 Linalg::sort<double>(eig, eigvec, false);
124 qint32 rnk = Linalg::rank(A);
125
126 for(qint32 i = 0; i < eig.size()-rnk; ++i)
127 eig(i) = 0;
128
129 printf("Setting small %s eigenvalues to zero.\n", ch_type.toUtf8().constData());
130 if (!pca)
131 printf("Not doing PCA for %s\n", ch_type.toUtf8().constData());
132 else
133 {
134 printf("Doing PCA for %s.",ch_type.toUtf8().constData());
135 eigvec = eigvec.block(eigvec.rows()-rnk, 0, rnk, eigvec.cols());
136 }
137}
138
139//=============================================================================================================
140
141void Linalg::get_whitener(MatrixXd &A,
142 bool pca,
143 const std::string& ch_type,
144 VectorXd &eig,
145 MatrixXd &eigvec)
146{
147 SelfAdjointEigenSolver<MatrixXd> t_eigenSolver(A);
148
149 eig = t_eigenSolver.eigenvalues();
150 eigvec = t_eigenSolver.eigenvectors().transpose();
151
152 Linalg::sort<double>(eig, eigvec, false);
153 qint32 rnk = Linalg::rank(A);
154
155 for(qint32 i = 0; i < eig.size()-rnk; ++i)
156 eig(i) = 0;
157
158 printf("Setting small %s eigenvalues to zero.\n", ch_type.c_str());
159 if (!pca)
160 printf("Not doing PCA for %s\n", ch_type.c_str());
161 else
162 {
163 printf("Doing PCA for %s.",ch_type.c_str());
164 eigvec = eigvec.block(eigvec.rows()-rnk, 0, rnk, eigvec.cols());
165 }
166}
167
168//=============================================================================================================
169
170VectorXi Linalg::intersect(const VectorXi &v1,
171 const VectorXi &v2,
172 VectorXi &idx_sel)
173{
174 std::vector<int> tmp;
175
176 std::vector< std::pair<int,int> > t_vecIntIdxValue;
177
178 for(qint32 i = 0; i < v1.size(); ++i)
179 tmp.push_back(v1[i]);
180
181 std::vector<int>::iterator it;
182 for(qint32 i = 0; i < v2.size(); ++i)
183 {
184 it = std::search(tmp.begin(), tmp.end(), &v2[i], &v2[i]+1);
185 if(it != tmp.end())
186 t_vecIntIdxValue.push_back(std::pair<int,int>(v2[i], it-tmp.begin()));
187 }
188
189 std::sort(t_vecIntIdxValue.begin(), t_vecIntIdxValue.end(), Linalg::compareIdxValuePairSmallerThan<int>);
190
191 VectorXi p_res(t_vecIntIdxValue.size());
192 idx_sel = VectorXi(t_vecIntIdxValue.size());
193
194 for(quint32 i = 0; i < t_vecIntIdxValue.size(); ++i)
195 {
196 p_res[i] = t_vecIntIdxValue[i].first;
197 idx_sel[i] = t_vecIntIdxValue[i].second;
198 }
199
200 return p_res;
201}
202
203//=============================================================================================================
204
205SparseMatrix<double> Linalg::make_block_diag(const MatrixXd &A,
206 qint32 n)
207{
208 qint32 ma = A.rows();
209 qint32 na = A.cols();
210 float bdn = static_cast<float>(na) / n;
211
212 if (bdn - std::floor(bdn))
213 {
214 qWarning("Linalg::make_block_diag: Width of matrix must be an even multiple of n.");
215 return SparseMatrix<double>();
216 }
217
218 typedef Eigen::Triplet<double> T;
219 std::vector<T> tripletList;
220 tripletList.reserve(static_cast<size_t>(bdn * ma * n));
221
222 for (qint32 i = 0; i < static_cast<qint32>(bdn); ++i)
223 {
224 qint32 current_col = i * n;
225 qint32 current_row = i * ma;
226
227 for (qint32 r = 0; r < ma; ++r)
228 for (qint32 c = 0; c < n; ++c)
229 tripletList.push_back(T(r + current_row, c + current_col, A(r, c + current_col)));
230 }
231
232 SparseMatrix<double> bd(static_cast<int>(std::floor(ma * bdn + 0.5f)), na);
233 bd.setFromTriplets(tripletList.begin(), tripletList.end());
234
235 return bd;
236}
237
238//=============================================================================================================
239
240qint32 Linalg::rank(const MatrixXd& A,
241 double tol)
242{
243 JacobiSVD<MatrixXd> t_svdA(A);
244 VectorXd s = t_svdA.singularValues();
245 double t_dMax = s.maxCoeff();
246 t_dMax *= tol;
247 qint32 sum = 0;
248 for(qint32 i = 0; i < s.size(); ++i)
249 sum += s[i] > t_dMax ? 1 : 0;
250 return sum;
251}
Linalg class declaration.
Shared utilities (I/O helpers, spectral analysis, layout management, warp algorithms).
static bool compareIdxValuePairSmallerThan(const std::pair< int, T > &lhs, const std::pair< int, T > &rhs)
Definition linalg.h:361
static Eigen::VectorXd combine_xyz(const Eigen::VectorXd &vec)
Definition linalg.cpp:64
static Eigen::VectorXi sort(Eigen::Matrix< T, Eigen::Dynamic, 1 > &v, bool desc=true)
Definition linalg.h:284
static qint32 rank(const Eigen::MatrixXd &A, double tol=1e-8)
Definition linalg.cpp:240
static void get_whitener(Eigen::MatrixXd &A, bool pca, QString ch_type, Eigen::VectorXd &eig, Eigen::MatrixXd &eigvec)
static Eigen::VectorXi intersect(const Eigen::VectorXi &v1, const Eigen::VectorXi &v2, Eigen::VectorXi &idx_sel)
Definition linalg.cpp:170
static double getConditionNumber(const Eigen::MatrixXd &A, Eigen::VectorXd &s)
Definition linalg.cpp:86
static double getConditionSlope(const Eigen::MatrixXd &A, Eigen::VectorXd &s)
Definition linalg.cpp:99
static Eigen::SparseMatrix< double > make_block_diag(const Eigen::MatrixXd &A, qint32 n)
Definition linalg.cpp:205