source: pacpussensors/trunk/Vislab/lib3dv/eigen/test/sparse_solver.h@ 136

Last change on this file since 136 was 136, checked in by ldecherf, 7 years ago

Doc

File size: 11.5 KB
Line 
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2011 Gael Guennebaud <g.gael@free.fr>
5//
6// This Source Code Form is subject to the terms of the Mozilla
7// Public License v. 2.0. If a copy of the MPL was not distributed
8// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10#include "sparse.h"
11#include <Eigen/SparseCore>
12
13template<typename Solver, typename Rhs, typename DenseMat, typename DenseRhs>
14void check_sparse_solving(Solver& solver, const typename Solver::MatrixType& A, const Rhs& b, const DenseMat& dA, const DenseRhs& db)
15{
16 typedef typename Solver::MatrixType Mat;
17 typedef typename Mat::Scalar Scalar;
18
19 DenseRhs refX = dA.lu().solve(db);
20 {
21 Rhs x(b.rows(), b.cols());
22 Rhs oldb = b;
23
24 solver.compute(A);
25 if (solver.info() != Success)
26 {
27 std::cerr << "sparse solver testing: factorization failed (check_sparse_solving)\n";
28 exit(0);
29 return;
30 }
31 x = solver.solve(b);
32 if (solver.info() != Success)
33 {
34 std::cerr << "sparse solver testing: solving failed\n";
35 return;
36 }
37 VERIFY(oldb.isApprox(b) && "sparse solver testing: the rhs should not be modified!");
38
39 VERIFY(x.isApprox(refX,test_precision<Scalar>()));
40 x.setZero();
41 // test the analyze/factorize API
42 solver.analyzePattern(A);
43 solver.factorize(A);
44 if (solver.info() != Success)
45 {
46 std::cerr << "sparse solver testing: factorization failed (check_sparse_solving)\n";
47 exit(0);
48 return;
49 }
50 x = solver.solve(b);
51 if (solver.info() != Success)
52 {
53 std::cerr << "sparse solver testing: solving failed\n";
54 return;
55 }
56 VERIFY(oldb.isApprox(b) && "sparse solver testing: the rhs should not be modified!");
57
58 VERIFY(x.isApprox(refX,test_precision<Scalar>()));
59 }
60
61 // test dense Block as the result and rhs:
62 {
63 DenseRhs x(db.rows(), db.cols());
64 DenseRhs oldb(db);
65 x.setZero();
66 x.block(0,0,x.rows(),x.cols()) = solver.solve(db.block(0,0,db.rows(),db.cols()));
67 VERIFY(oldb.isApprox(db) && "sparse solver testing: the rhs should not be modified!");
68 VERIFY(x.isApprox(refX,test_precision<Scalar>()));
69 }
70
71 // if not too large, do some extra check:
72 if(A.rows()<2000)
73 {
74
75 // test expression as input
76 {
77 solver.compute(0.5*(A+A));
78 Rhs x = solver.solve(b);
79 VERIFY(x.isApprox(refX,test_precision<Scalar>()));
80
81 Solver solver2(0.5*(A+A));
82 Rhs x2 = solver2.solve(b);
83 VERIFY(x2.isApprox(refX,test_precision<Scalar>()));
84 }
85 }
86}
87
88template<typename Solver, typename Rhs>
89void check_sparse_solving_real_cases(Solver& solver, const typename Solver::MatrixType& A, const Rhs& b, const Rhs& refX)
90{
91 typedef typename Solver::MatrixType Mat;
92 typedef typename Mat::Scalar Scalar;
93 typedef typename Mat::RealScalar RealScalar;
94
95 Rhs x(b.rows(), b.cols());
96
97 solver.compute(A);
98 if (solver.info() != Success)
99 {
100 std::cerr << "sparse solver testing: factorization failed (check_sparse_solving_real_cases)\n";
101 exit(0);
102 return;
103 }
104 x = solver.solve(b);
105 if (solver.info() != Success)
106 {
107 std::cerr << "sparse solver testing: solving failed\n";
108 return;
109 }
110
111 RealScalar res_error;
112 // Compute the norm of the relative error
113 if(refX.size() != 0)
114 res_error = (refX - x).norm()/refX.norm();
115 else
116 {
117 // Compute the relative residual norm
118 res_error = (b - A * x).norm()/b.norm();
119 }
120 if (res_error > test_precision<Scalar>() ){
121 std::cerr << "Test " << g_test_stack.back() << " failed in "EI_PP_MAKE_STRING(__FILE__)
122 << " (" << EI_PP_MAKE_STRING(__LINE__) << ")" << std::endl << std::endl;
123 abort();
124 }
125
126}
127template<typename Solver, typename DenseMat>
128void check_sparse_determinant(Solver& solver, const typename Solver::MatrixType& A, const DenseMat& dA)
129{
130 typedef typename Solver::MatrixType Mat;
131 typedef typename Mat::Scalar Scalar;
132
133 solver.compute(A);
134 if (solver.info() != Success)
135 {
136 std::cerr << "sparse solver testing: factorization failed (check_sparse_determinant)\n";
137 return;
138 }
139
140 Scalar refDet = dA.determinant();
141 VERIFY_IS_APPROX(refDet,solver.determinant());
142}
143template<typename Solver, typename DenseMat>
144void check_sparse_abs_determinant(Solver& solver, const typename Solver::MatrixType& A, const DenseMat& dA)
145{
146 using std::abs;
147 typedef typename Solver::MatrixType Mat;
148 typedef typename Mat::Scalar Scalar;
149
150 solver.compute(A);
151 if (solver.info() != Success)
152 {
153 std::cerr << "sparse solver testing: factorization failed (check_sparse_abs_determinant)\n";
154 return;
155 }
156
157 Scalar refDet = abs(dA.determinant());
158 VERIFY_IS_APPROX(refDet,solver.absDeterminant());
159}
160
161template<typename Solver, typename DenseMat>
162int generate_sparse_spd_problem(Solver& , typename Solver::MatrixType& A, typename Solver::MatrixType& halfA, DenseMat& dA, int maxSize = 300)
163{
164 typedef typename Solver::MatrixType Mat;
165 typedef typename Mat::Scalar Scalar;
166 typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
167
168 int size = internal::random<int>(1,maxSize);
169 double density = (std::max)(8./(size*size), 0.01);
170
171 Mat M(size, size);
172 DenseMatrix dM(size, size);
173
174 initSparse<Scalar>(density, dM, M, ForceNonZeroDiag);
175
176 A = M * M.adjoint();
177 dA = dM * dM.adjoint();
178
179 halfA.resize(size,size);
180 if(Solver::UpLo==(Lower|Upper))
181 halfA = A;
182 else
183 halfA.template selfadjointView<Solver::UpLo>().rankUpdate(M);
184
185 return size;
186}
187
188
189#ifdef TEST_REAL_CASES
190template<typename Scalar>
191inline std::string get_matrixfolder()
192{
193 std::string mat_folder = TEST_REAL_CASES;
194 if( internal::is_same<Scalar, std::complex<float> >::value || internal::is_same<Scalar, std::complex<double> >::value )
195 mat_folder = mat_folder + static_cast<std::string>("/complex/");
196 else
197 mat_folder = mat_folder + static_cast<std::string>("/real/");
198 return mat_folder;
199}
200#endif
201
202template<typename Solver> void check_sparse_spd_solving(Solver& solver)
203{
204 typedef typename Solver::MatrixType Mat;
205 typedef typename Mat::Scalar Scalar;
206 typedef SparseMatrix<Scalar,ColMajor> SpMat;
207 typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
208 typedef Matrix<Scalar,Dynamic,1> DenseVector;
209
210 // generate the problem
211 Mat A, halfA;
212 DenseMatrix dA;
213 for (int i = 0; i < g_repeat; i++) {
214 int size = generate_sparse_spd_problem(solver, A, halfA, dA);
215
216 // generate the right hand sides
217 int rhsCols = internal::random<int>(1,16);
218 double density = (std::max)(8./(size*rhsCols), 0.1);
219 SpMat B(size,rhsCols);
220 DenseVector b = DenseVector::Random(size);
221 DenseMatrix dB(size,rhsCols);
222 initSparse<Scalar>(density, dB, B, ForceNonZeroDiag);
223
224 check_sparse_solving(solver, A, b, dA, b);
225 check_sparse_solving(solver, halfA, b, dA, b);
226 check_sparse_solving(solver, A, dB, dA, dB);
227 check_sparse_solving(solver, halfA, dB, dA, dB);
228 check_sparse_solving(solver, A, B, dA, dB);
229 check_sparse_solving(solver, halfA, B, dA, dB);
230
231 // check only once
232 if(i==0)
233 {
234 b = DenseVector::Zero(size);
235 check_sparse_solving(solver, A, b, dA, b);
236 }
237 }
238
239 // First, get the folder
240#ifdef TEST_REAL_CASES
241 if (internal::is_same<Scalar, float>::value
242 || internal::is_same<Scalar, std::complex<float> >::value)
243 return ;
244
245 std::string mat_folder = get_matrixfolder<Scalar>();
246 MatrixMarketIterator<Scalar> it(mat_folder);
247 for (; it; ++it)
248 {
249 if (it.sym() == SPD){
250 Mat halfA;
251 PermutationMatrix<Dynamic, Dynamic, Index> pnull;
252 halfA.template selfadjointView<Solver::UpLo>() = it.matrix().template triangularView<Eigen::Lower>().twistedBy(pnull);
253
254 std::cout<< " ==== SOLVING WITH MATRIX " << it.matname() << " ==== \n";
255 check_sparse_solving_real_cases(solver, it.matrix(), it.rhs(), it.refX());
256 check_sparse_solving_real_cases(solver, halfA, it.rhs(), it.refX());
257 }
258 }
259#endif
260}
261
262template<typename Solver> void check_sparse_spd_determinant(Solver& solver)
263{
264 typedef typename Solver::MatrixType Mat;
265 typedef typename Mat::Scalar Scalar;
266 typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
267
268 // generate the problem
269 Mat A, halfA;
270 DenseMatrix dA;
271 generate_sparse_spd_problem(solver, A, halfA, dA, 30);
272
273 for (int i = 0; i < g_repeat; i++) {
274 check_sparse_determinant(solver, A, dA);
275 check_sparse_determinant(solver, halfA, dA );
276 }
277}
278
279template<typename Solver, typename DenseMat>
280int generate_sparse_square_problem(Solver&, typename Solver::MatrixType& A, DenseMat& dA, int maxSize = 300)
281{
282 typedef typename Solver::MatrixType Mat;
283 typedef typename Mat::Scalar Scalar;
284
285 int size = internal::random<int>(1,maxSize);
286 double density = (std::max)(8./(size*size), 0.01);
287
288 A.resize(size,size);
289 dA.resize(size,size);
290
291 initSparse<Scalar>(density, dA, A, ForceNonZeroDiag);
292
293 return size;
294}
295
296
297struct prune_column {
298 int m_col;
299 prune_column(int col) : m_col(col) {}
300 template<class Scalar>
301 bool operator()(int, int col, const Scalar&) const {
302 return col != m_col;
303 }
304};
305
306template<typename Solver> void check_sparse_square_solving(Solver& solver, bool checkDeficient = false)
307{
308 typedef typename Solver::MatrixType Mat;
309 typedef typename Mat::Scalar Scalar;
310 typedef SparseMatrix<Scalar,ColMajor> SpMat;
311 typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
312 typedef Matrix<Scalar,Dynamic,1> DenseVector;
313
314 int rhsCols = internal::random<int>(1,16);
315
316 Mat A;
317 DenseMatrix dA;
318 for (int i = 0; i < g_repeat; i++) {
319 int size = generate_sparse_square_problem(solver, A, dA);
320
321 A.makeCompressed();
322 DenseVector b = DenseVector::Random(size);
323 DenseMatrix dB(size,rhsCols);
324 SpMat B(size,rhsCols);
325 double density = (std::max)(8./(size*rhsCols), 0.1);
326 initSparse<Scalar>(density, dB, B, ForceNonZeroDiag);
327 B.makeCompressed();
328 check_sparse_solving(solver, A, b, dA, b);
329 check_sparse_solving(solver, A, dB, dA, dB);
330 check_sparse_solving(solver, A, B, dA, dB);
331
332 // check only once
333 if(i==0)
334 {
335 b = DenseVector::Zero(size);
336 check_sparse_solving(solver, A, b, dA, b);
337 }
338 // regression test for Bug 792 (structurally rank deficient matrices):
339 if(checkDeficient && size>1) {
340 int col = internal::random<int>(0,size-1);
341 A.prune(prune_column(col));
342 solver.compute(A);
343 VERIFY_IS_EQUAL(solver.info(), NumericalIssue);
344 }
345 }
346
347 // First, get the folder
348#ifdef TEST_REAL_CASES
349 if (internal::is_same<Scalar, float>::value
350 || internal::is_same<Scalar, std::complex<float> >::value)
351 return ;
352
353 std::string mat_folder = get_matrixfolder<Scalar>();
354 MatrixMarketIterator<Scalar> it(mat_folder);
355 for (; it; ++it)
356 {
357 std::cout<< " ==== SOLVING WITH MATRIX " << it.matname() << " ==== \n";
358 check_sparse_solving_real_cases(solver, it.matrix(), it.rhs(), it.refX());
359 }
360#endif
361
362}
363
364template<typename Solver> void check_sparse_square_determinant(Solver& solver)
365{
366 typedef typename Solver::MatrixType Mat;
367 typedef typename Mat::Scalar Scalar;
368 typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
369
370 // generate the problem
371 Mat A;
372 DenseMatrix dA;
373 generate_sparse_square_problem(solver, A, dA, 30);
374 A.makeCompressed();
375 for (int i = 0; i < g_repeat; i++) {
376 check_sparse_determinant(solver, A, dA);
377 }
378}
379
380template<typename Solver> void check_sparse_square_abs_determinant(Solver& solver)
381{
382 typedef typename Solver::MatrixType Mat;
383 typedef typename Mat::Scalar Scalar;
384 typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
385
386 // generate the problem
387 Mat A;
388 DenseMatrix dA;
389 generate_sparse_square_problem(solver, A, dA, 30);
390 A.makeCompressed();
391 for (int i = 0; i < g_repeat; i++) {
392 check_sparse_abs_determinant(solver, A, dA);
393 }
394}
395
Note: See TracBrowser for help on using the repository browser.