From 9cf77ce1d80cf17aa79c5da95b578ee2a4490152 Mon Sep 17 00:00:00 2001 From: Desire NUENTSA Date: Mon, 12 Nov 2012 15:20:37 +0100 Subject: Add support for Sparse QR factorization --- test/spqr_support.cpp | 62 +++++++++++++++++++++++++++++++++++++++++++++++++++ 1 file changed, 62 insertions(+) create mode 100644 test/spqr_support.cpp (limited to 'test/spqr_support.cpp') diff --git a/test/spqr_support.cpp b/test/spqr_support.cpp new file mode 100644 index 000000000..9ff4f9cba --- /dev/null +++ b/test/spqr_support.cpp @@ -0,0 +1,62 @@ +// This file is part of Eigen, a lightweight C++ template library +// for linear algebra. +// +// Copyright (C) 2012 Desire Nuentsa Wakam +// +// This Source Code Form is subject to the terms of the Mozilla +// Public License v. 2.0. If a copy of the MPL was not distributed +#include "sparse.h" +#include + + +template +int generate_sparse_rectangular_problem(MatrixType& A, DenseMat& dA, int maxRows = 300, int maxCols = 300) +{ + eigen_assert(maxRows >= maxCols); + typedef typename MatrixType::Scalar Scalar; + int rows = internal::random(1,maxRows); + int cols = internal::random(1,maxCols); + double density = (std::max)(8./(rows*cols), 0.01); + + A.resize(rows,rows); + dA.resize(rows,rows); + initSparse(density, dA, A,ForceNonZeroDiag); + A.makeCompressed(); + return rows; +} + +template void test_spqr_scalar() +{ + typedef SparseMatrix MatrixType; + MatrixType A; + Matrix dA; + typedef Matrix DenseVector; + DenseVector refX,x,b; + SPQR solver; + generate_sparse_rectangular_problem(A,dA); + + int n = A.cols(); + b = DenseVector::Random(n); + solver.compute(A); + if (solver.info() != Success) + { + std::cerr << "sparse QR factorization failed\n"; + exit(0); + return; + } + solver._solve(b,x); + if (solver.info() != Success) + { + std::cerr << "sparse QR factorization failed\n"; + exit(0); + return; + } + //Compare with a dense solver + refX = dA.colPivHouseholderQr().solve(b); + VERIFY(x.isApprox(refX,test_precision())); +} +void test_spqr_support() +{ + CALL_SUBTEST_1(test_spqr_scalar()); + CALL_SUBTEST_2(test_spqr_scalar >()); +} \ No newline at end of file -- cgit v1.2.3