2011-10-11 17:32:26 +08:00
|
|
|
// This file is part of Eigen, a lightweight C++ template library
|
|
|
|
// for linear algebra.
|
|
|
|
//
|
|
|
|
// Copyright (C) 2011 Gael Guennebaud <g.gael@free.fr>
|
|
|
|
//
|
2012-07-14 02:42:47 +08:00
|
|
|
// 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
|
|
|
|
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
|
2011-10-11 17:32:26 +08:00
|
|
|
|
|
|
|
#include "sparse.h"
|
2011-11-12 21:11:27 +08:00
|
|
|
#include <Eigen/SparseCore>
|
2011-10-11 17:32:26 +08:00
|
|
|
|
|
|
|
template<typename Solver, typename Rhs, typename DenseMat, typename DenseRhs>
|
|
|
|
void check_sparse_solving(Solver& solver, const typename Solver::MatrixType& A, const Rhs& b, const DenseMat& dA, const DenseRhs& db)
|
|
|
|
{
|
|
|
|
typedef typename Solver::MatrixType Mat;
|
|
|
|
typedef typename Mat::Scalar Scalar;
|
2014-12-05 05:48:53 +08:00
|
|
|
typedef typename Mat::StorageIndex StorageIndex;
|
2011-10-11 17:32:26 +08:00
|
|
|
|
|
|
|
DenseRhs refX = dA.lu().solve(db);
|
2013-06-24 01:11:32 +08:00
|
|
|
{
|
|
|
|
Rhs x(b.rows(), b.cols());
|
|
|
|
Rhs oldb = b;
|
2011-10-11 17:32:26 +08:00
|
|
|
|
2013-06-24 01:11:32 +08:00
|
|
|
solver.compute(A);
|
|
|
|
if (solver.info() != Success)
|
|
|
|
{
|
|
|
|
std::cerr << "sparse solver testing: factorization failed (check_sparse_solving)\n";
|
|
|
|
exit(0);
|
|
|
|
return;
|
|
|
|
}
|
|
|
|
x = solver.solve(b);
|
|
|
|
if (solver.info() != Success)
|
|
|
|
{
|
2015-02-09 18:41:25 +08:00
|
|
|
std::cerr << "sparse solver testing: solving failed (" << typeid(Solver).name() << ")\n";
|
2013-06-24 01:11:32 +08:00
|
|
|
return;
|
|
|
|
}
|
|
|
|
VERIFY(oldb.isApprox(b) && "sparse solver testing: the rhs should not be modified!");
|
|
|
|
VERIFY(x.isApprox(refX,test_precision<Scalar>()));
|
2014-10-20 22:46:47 +08:00
|
|
|
|
2013-06-24 01:11:32 +08:00
|
|
|
x.setZero();
|
|
|
|
// test the analyze/factorize API
|
|
|
|
solver.analyzePattern(A);
|
|
|
|
solver.factorize(A);
|
|
|
|
if (solver.info() != Success)
|
|
|
|
{
|
|
|
|
std::cerr << "sparse solver testing: factorization failed (check_sparse_solving)\n";
|
|
|
|
exit(0);
|
|
|
|
return;
|
|
|
|
}
|
|
|
|
x = solver.solve(b);
|
|
|
|
if (solver.info() != Success)
|
|
|
|
{
|
|
|
|
std::cerr << "sparse solver testing: solving failed\n";
|
|
|
|
return;
|
|
|
|
}
|
|
|
|
VERIFY(oldb.isApprox(b) && "sparse solver testing: the rhs should not be modified!");
|
|
|
|
VERIFY(x.isApprox(refX,test_precision<Scalar>()));
|
2014-10-20 22:46:47 +08:00
|
|
|
|
|
|
|
|
|
|
|
x.setZero();
|
|
|
|
// test with Map
|
2014-12-05 05:48:53 +08:00
|
|
|
MappedSparseMatrix<Scalar,Mat::Options,StorageIndex> Am(A.rows(), A.cols(), A.nonZeros(), const_cast<StorageIndex*>(A.outerIndexPtr()), const_cast<StorageIndex*>(A.innerIndexPtr()), const_cast<Scalar*>(A.valuePtr()));
|
2014-10-20 22:46:47 +08:00
|
|
|
solver.compute(Am);
|
|
|
|
if (solver.info() != Success)
|
|
|
|
{
|
|
|
|
std::cerr << "sparse solver testing: factorization failed (check_sparse_solving)\n";
|
|
|
|
exit(0);
|
|
|
|
return;
|
|
|
|
}
|
|
|
|
DenseRhs dx(refX);
|
|
|
|
dx.setZero();
|
|
|
|
Map<DenseRhs> xm(dx.data(), dx.rows(), dx.cols());
|
|
|
|
Map<const DenseRhs> bm(db.data(), db.rows(), db.cols());
|
|
|
|
xm = solver.solve(bm);
|
|
|
|
if (solver.info() != Success)
|
|
|
|
{
|
2015-02-09 18:41:25 +08:00
|
|
|
std::cerr << "sparse solver testing: solving with a Map failed\n";
|
|
|
|
exit(0);
|
2014-10-20 22:46:47 +08:00
|
|
|
return;
|
|
|
|
}
|
|
|
|
VERIFY(oldb.isApprox(bm) && "sparse solver testing: the rhs should not be modified!");
|
|
|
|
VERIFY(xm.isApprox(refX,test_precision<Scalar>()));
|
2012-02-27 20:22:38 +08:00
|
|
|
}
|
2012-06-06 15:40:01 +08:00
|
|
|
|
2015-02-17 00:05:10 +08:00
|
|
|
// test initialization ctor
|
|
|
|
{
|
|
|
|
Rhs x(b.rows(), b.cols());
|
|
|
|
Solver solver2(A);
|
|
|
|
VERIFY(solver2.info() == Success);
|
|
|
|
x = solver2.solve(b);
|
|
|
|
VERIFY(x.isApprox(refX,test_precision<Scalar>()));
|
|
|
|
}
|
|
|
|
|
2013-06-24 01:11:32 +08:00
|
|
|
// test dense Block as the result and rhs:
|
2012-06-06 15:40:01 +08:00
|
|
|
{
|
|
|
|
DenseRhs x(db.rows(), db.cols());
|
2013-06-24 01:11:32 +08:00
|
|
|
DenseRhs oldb(db);
|
2012-06-06 15:40:01 +08:00
|
|
|
x.setZero();
|
2013-06-24 01:11:32 +08:00
|
|
|
x.block(0,0,x.rows(),x.cols()) = solver.solve(db.block(0,0,db.rows(),db.cols()));
|
|
|
|
VERIFY(oldb.isApprox(db) && "sparse solver testing: the rhs should not be modified!");
|
2012-06-06 15:40:01 +08:00
|
|
|
VERIFY(x.isApprox(refX,test_precision<Scalar>()));
|
|
|
|
}
|
2014-10-06 17:41:50 +08:00
|
|
|
|
|
|
|
// test uncompressed inputs
|
|
|
|
{
|
|
|
|
Mat A2 = A;
|
2014-12-05 05:48:53 +08:00
|
|
|
A2.reserve((ArrayXf::Random(A.outerSize())+2).template cast<typename Mat::StorageIndex>().eval());
|
2014-10-06 17:41:50 +08:00
|
|
|
solver.compute(A2);
|
|
|
|
Rhs x = solver.solve(b);
|
|
|
|
VERIFY(x.isApprox(refX,test_precision<Scalar>()));
|
|
|
|
}
|
2011-10-11 17:32:26 +08:00
|
|
|
}
|
|
|
|
|
2012-03-29 20:32:54 +08:00
|
|
|
template<typename Solver, typename Rhs>
|
|
|
|
void check_sparse_solving_real_cases(Solver& solver, const typename Solver::MatrixType& A, const Rhs& b, const Rhs& refX)
|
|
|
|
{
|
|
|
|
typedef typename Solver::MatrixType Mat;
|
|
|
|
typedef typename Mat::Scalar Scalar;
|
|
|
|
typedef typename Mat::RealScalar RealScalar;
|
|
|
|
|
|
|
|
Rhs x(b.rows(), b.cols());
|
|
|
|
|
|
|
|
solver.compute(A);
|
|
|
|
if (solver.info() != Success)
|
|
|
|
{
|
|
|
|
std::cerr << "sparse solver testing: factorization failed (check_sparse_solving_real_cases)\n";
|
|
|
|
exit(0);
|
|
|
|
return;
|
|
|
|
}
|
|
|
|
x = solver.solve(b);
|
|
|
|
if (solver.info() != Success)
|
|
|
|
{
|
|
|
|
std::cerr << "sparse solver testing: solving failed\n";
|
|
|
|
return;
|
|
|
|
}
|
|
|
|
|
|
|
|
RealScalar res_error;
|
|
|
|
// Compute the norm of the relative error
|
|
|
|
if(refX.size() != 0)
|
|
|
|
res_error = (refX - x).norm()/refX.norm();
|
|
|
|
else
|
|
|
|
{
|
|
|
|
// Compute the relative residual norm
|
|
|
|
res_error = (b - A * x).norm()/b.norm();
|
|
|
|
}
|
|
|
|
if (res_error > test_precision<Scalar>() ){
|
|
|
|
std::cerr << "Test " << g_test_stack.back() << " failed in "EI_PP_MAKE_STRING(__FILE__)
|
|
|
|
<< " (" << EI_PP_MAKE_STRING(__LINE__) << ")" << std::endl << std::endl;
|
|
|
|
abort();
|
|
|
|
}
|
|
|
|
|
|
|
|
}
|
2011-10-11 17:32:26 +08:00
|
|
|
template<typename Solver, typename DenseMat>
|
|
|
|
void check_sparse_determinant(Solver& solver, const typename Solver::MatrixType& A, const DenseMat& dA)
|
|
|
|
{
|
|
|
|
typedef typename Solver::MatrixType Mat;
|
|
|
|
typedef typename Mat::Scalar Scalar;
|
|
|
|
|
|
|
|
solver.compute(A);
|
|
|
|
if (solver.info() != Success)
|
|
|
|
{
|
2011-11-12 21:11:27 +08:00
|
|
|
std::cerr << "sparse solver testing: factorization failed (check_sparse_determinant)\n";
|
2011-10-11 17:32:26 +08:00
|
|
|
return;
|
|
|
|
}
|
|
|
|
|
|
|
|
Scalar refDet = dA.determinant();
|
|
|
|
VERIFY_IS_APPROX(refDet,solver.determinant());
|
|
|
|
}
|
2014-10-17 22:52:56 +08:00
|
|
|
template<typename Solver, typename DenseMat>
|
|
|
|
void check_sparse_abs_determinant(Solver& solver, const typename Solver::MatrixType& A, const DenseMat& dA)
|
|
|
|
{
|
|
|
|
using std::abs;
|
|
|
|
typedef typename Solver::MatrixType Mat;
|
|
|
|
typedef typename Mat::Scalar Scalar;
|
|
|
|
|
|
|
|
solver.compute(A);
|
|
|
|
if (solver.info() != Success)
|
|
|
|
{
|
|
|
|
std::cerr << "sparse solver testing: factorization failed (check_sparse_abs_determinant)\n";
|
|
|
|
return;
|
|
|
|
}
|
2011-10-11 17:32:26 +08:00
|
|
|
|
2014-10-17 22:52:56 +08:00
|
|
|
Scalar refDet = abs(dA.determinant());
|
|
|
|
VERIFY_IS_APPROX(refDet,solver.absDeterminant());
|
|
|
|
}
|
2011-10-11 17:32:26 +08:00
|
|
|
|
|
|
|
template<typename Solver, typename DenseMat>
|
|
|
|
int generate_sparse_spd_problem(Solver& , typename Solver::MatrixType& A, typename Solver::MatrixType& halfA, DenseMat& dA, int maxSize = 300)
|
|
|
|
{
|
|
|
|
typedef typename Solver::MatrixType Mat;
|
|
|
|
typedef typename Mat::Scalar Scalar;
|
|
|
|
typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
|
|
|
|
|
|
|
|
int size = internal::random<int>(1,maxSize);
|
|
|
|
double density = (std::max)(8./(size*size), 0.01);
|
|
|
|
|
|
|
|
Mat M(size, size);
|
|
|
|
DenseMatrix dM(size, size);
|
|
|
|
|
|
|
|
initSparse<Scalar>(density, dM, M, ForceNonZeroDiag);
|
|
|
|
|
|
|
|
A = M * M.adjoint();
|
|
|
|
dA = dM * dM.adjoint();
|
|
|
|
|
|
|
|
halfA.resize(size,size);
|
2015-02-11 01:57:41 +08:00
|
|
|
if(Solver::UpLo==(Lower|Upper))
|
|
|
|
halfA = A;
|
|
|
|
else
|
|
|
|
halfA.template selfadjointView<Solver::UpLo>().rankUpdate(M);
|
2011-10-11 17:32:26 +08:00
|
|
|
|
|
|
|
return size;
|
|
|
|
}
|
|
|
|
|
2012-03-29 20:53:42 +08:00
|
|
|
|
|
|
|
#ifdef TEST_REAL_CASES
|
|
|
|
template<typename Scalar>
|
|
|
|
inline std::string get_matrixfolder()
|
|
|
|
{
|
|
|
|
std::string mat_folder = TEST_REAL_CASES;
|
|
|
|
if( internal::is_same<Scalar, std::complex<float> >::value || internal::is_same<Scalar, std::complex<double> >::value )
|
2012-08-03 19:05:45 +08:00
|
|
|
mat_folder = mat_folder + static_cast<std::string>("/complex/");
|
2012-03-29 20:53:42 +08:00
|
|
|
else
|
2012-08-03 19:05:45 +08:00
|
|
|
mat_folder = mat_folder + static_cast<std::string>("/real/");
|
2012-03-29 20:53:42 +08:00
|
|
|
return mat_folder;
|
|
|
|
}
|
|
|
|
#endif
|
|
|
|
|
2011-10-11 17:32:26 +08:00
|
|
|
template<typename Solver> void check_sparse_spd_solving(Solver& solver)
|
|
|
|
{
|
|
|
|
typedef typename Solver::MatrixType Mat;
|
|
|
|
typedef typename Mat::Scalar Scalar;
|
2012-02-04 21:20:56 +08:00
|
|
|
typedef SparseMatrix<Scalar,ColMajor> SpMat;
|
2011-10-11 17:32:26 +08:00
|
|
|
typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
|
|
|
|
typedef Matrix<Scalar,Dynamic,1> DenseVector;
|
|
|
|
|
|
|
|
// generate the problem
|
|
|
|
Mat A, halfA;
|
|
|
|
DenseMatrix dA;
|
2012-03-29 20:32:54 +08:00
|
|
|
for (int i = 0; i < g_repeat; i++) {
|
2013-03-20 23:15:18 +08:00
|
|
|
int size = generate_sparse_spd_problem(solver, A, halfA, dA);
|
|
|
|
|
|
|
|
// generate the right hand sides
|
|
|
|
int rhsCols = internal::random<int>(1,16);
|
|
|
|
double density = (std::max)(8./(size*rhsCols), 0.1);
|
|
|
|
SpMat B(size,rhsCols);
|
|
|
|
DenseVector b = DenseVector::Random(size);
|
|
|
|
DenseMatrix dB(size,rhsCols);
|
|
|
|
initSparse<Scalar>(density, dB, B, ForceNonZeroDiag);
|
|
|
|
|
2012-03-29 20:32:54 +08:00
|
|
|
check_sparse_solving(solver, A, b, dA, b);
|
|
|
|
check_sparse_solving(solver, halfA, b, dA, b);
|
|
|
|
check_sparse_solving(solver, A, dB, dA, dB);
|
|
|
|
check_sparse_solving(solver, halfA, dB, dA, dB);
|
|
|
|
check_sparse_solving(solver, A, B, dA, dB);
|
|
|
|
check_sparse_solving(solver, halfA, B, dA, dB);
|
2013-03-20 23:15:18 +08:00
|
|
|
|
|
|
|
// check only once
|
|
|
|
if(i==0)
|
|
|
|
{
|
|
|
|
b = DenseVector::Zero(size);
|
|
|
|
check_sparse_solving(solver, A, b, dA, b);
|
|
|
|
}
|
2012-03-29 20:32:54 +08:00
|
|
|
}
|
2013-03-20 23:15:18 +08:00
|
|
|
|
2012-03-29 20:32:54 +08:00
|
|
|
// First, get the folder
|
2012-03-29 20:53:42 +08:00
|
|
|
#ifdef TEST_REAL_CASES
|
2012-03-29 20:32:54 +08:00
|
|
|
if (internal::is_same<Scalar, float>::value
|
|
|
|
|| internal::is_same<Scalar, std::complex<float> >::value)
|
|
|
|
return ;
|
|
|
|
|
|
|
|
std::string mat_folder = get_matrixfolder<Scalar>();
|
|
|
|
MatrixMarketIterator<Scalar> it(mat_folder);
|
|
|
|
for (; it; ++it)
|
|
|
|
{
|
|
|
|
if (it.sym() == SPD){
|
|
|
|
Mat halfA;
|
|
|
|
PermutationMatrix<Dynamic, Dynamic, Index> pnull;
|
|
|
|
halfA.template selfadjointView<Solver::UpLo>() = it.matrix().template triangularView<Eigen::Lower>().twistedBy(pnull);
|
|
|
|
|
|
|
|
std::cout<< " ==== SOLVING WITH MATRIX " << it.matname() << " ==== \n";
|
|
|
|
check_sparse_solving_real_cases(solver, it.matrix(), it.rhs(), it.refX());
|
|
|
|
check_sparse_solving_real_cases(solver, halfA, it.rhs(), it.refX());
|
|
|
|
}
|
|
|
|
}
|
|
|
|
#endif
|
2011-10-11 17:32:26 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
template<typename Solver> void check_sparse_spd_determinant(Solver& solver)
|
|
|
|
{
|
|
|
|
typedef typename Solver::MatrixType Mat;
|
|
|
|
typedef typename Mat::Scalar Scalar;
|
|
|
|
typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
|
|
|
|
|
|
|
|
// generate the problem
|
|
|
|
Mat A, halfA;
|
|
|
|
DenseMatrix dA;
|
|
|
|
generate_sparse_spd_problem(solver, A, halfA, dA, 30);
|
2012-03-29 20:32:54 +08:00
|
|
|
|
|
|
|
for (int i = 0; i < g_repeat; i++) {
|
|
|
|
check_sparse_determinant(solver, A, dA);
|
|
|
|
check_sparse_determinant(solver, halfA, dA );
|
|
|
|
}
|
2011-10-11 17:32:26 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
template<typename Solver, typename DenseMat>
|
2015-02-17 02:09:48 +08:00
|
|
|
int generate_sparse_square_problem(Solver&, typename Solver::MatrixType& A, DenseMat& dA, int maxSize = 300, int options = ForceNonZeroDiag)
|
2011-10-11 17:32:26 +08:00
|
|
|
{
|
|
|
|
typedef typename Solver::MatrixType Mat;
|
|
|
|
typedef typename Mat::Scalar Scalar;
|
|
|
|
|
|
|
|
int size = internal::random<int>(1,maxSize);
|
|
|
|
double density = (std::max)(8./(size*size), 0.01);
|
|
|
|
|
|
|
|
A.resize(size,size);
|
|
|
|
dA.resize(size,size);
|
|
|
|
|
2015-02-17 02:09:48 +08:00
|
|
|
initSparse<Scalar>(density, dA, A, options);
|
2011-10-11 17:32:26 +08:00
|
|
|
|
|
|
|
return size;
|
|
|
|
}
|
|
|
|
|
|
|
|
template<typename Solver> void check_sparse_square_solving(Solver& solver)
|
|
|
|
{
|
|
|
|
typedef typename Solver::MatrixType Mat;
|
|
|
|
typedef typename Mat::Scalar Scalar;
|
2013-01-26 01:17:17 +08:00
|
|
|
typedef SparseMatrix<Scalar,ColMajor> SpMat;
|
2011-10-11 17:32:26 +08:00
|
|
|
typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
|
|
|
|
typedef Matrix<Scalar,Dynamic,1> DenseVector;
|
|
|
|
|
|
|
|
int rhsCols = internal::random<int>(1,16);
|
|
|
|
|
|
|
|
Mat A;
|
|
|
|
DenseMatrix dA;
|
2012-03-29 20:32:54 +08:00
|
|
|
for (int i = 0; i < g_repeat; i++) {
|
2013-03-20 23:15:18 +08:00
|
|
|
int size = generate_sparse_square_problem(solver, A, dA);
|
|
|
|
|
|
|
|
A.makeCompressed();
|
|
|
|
DenseVector b = DenseVector::Random(size);
|
|
|
|
DenseMatrix dB(size,rhsCols);
|
|
|
|
SpMat B(size,rhsCols);
|
|
|
|
double density = (std::max)(8./(size*rhsCols), 0.1);
|
|
|
|
initSparse<Scalar>(density, dB, B, ForceNonZeroDiag);
|
|
|
|
B.makeCompressed();
|
2012-03-29 20:32:54 +08:00
|
|
|
check_sparse_solving(solver, A, b, dA, b);
|
|
|
|
check_sparse_solving(solver, A, dB, dA, dB);
|
2013-01-26 01:17:17 +08:00
|
|
|
check_sparse_solving(solver, A, B, dA, dB);
|
2013-03-20 23:15:18 +08:00
|
|
|
|
|
|
|
// check only once
|
|
|
|
if(i==0)
|
|
|
|
{
|
|
|
|
b = DenseVector::Zero(size);
|
|
|
|
check_sparse_solving(solver, A, b, dA, b);
|
|
|
|
}
|
2012-03-29 20:32:54 +08:00
|
|
|
}
|
2013-03-20 23:15:18 +08:00
|
|
|
|
2012-03-29 20:32:54 +08:00
|
|
|
// First, get the folder
|
2012-03-29 20:53:42 +08:00
|
|
|
#ifdef TEST_REAL_CASES
|
2012-03-29 20:32:54 +08:00
|
|
|
if (internal::is_same<Scalar, float>::value
|
|
|
|
|| internal::is_same<Scalar, std::complex<float> >::value)
|
|
|
|
return ;
|
|
|
|
|
|
|
|
std::string mat_folder = get_matrixfolder<Scalar>();
|
|
|
|
MatrixMarketIterator<Scalar> it(mat_folder);
|
|
|
|
for (; it; ++it)
|
|
|
|
{
|
|
|
|
std::cout<< " ==== SOLVING WITH MATRIX " << it.matname() << " ==== \n";
|
|
|
|
check_sparse_solving_real_cases(solver, it.matrix(), it.rhs(), it.refX());
|
|
|
|
}
|
|
|
|
#endif
|
2011-10-11 17:32:26 +08:00
|
|
|
|
|
|
|
}
|
|
|
|
|
|
|
|
template<typename Solver> void check_sparse_square_determinant(Solver& solver)
|
|
|
|
{
|
|
|
|
typedef typename Solver::MatrixType Mat;
|
|
|
|
typedef typename Mat::Scalar Scalar;
|
|
|
|
typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
|
2015-02-17 02:09:48 +08:00
|
|
|
|
2012-03-29 20:32:54 +08:00
|
|
|
for (int i = 0; i < g_repeat; i++) {
|
2015-02-17 02:09:48 +08:00
|
|
|
// generate the problem
|
|
|
|
Mat A;
|
|
|
|
DenseMatrix dA;
|
|
|
|
|
|
|
|
int size = internal::random<int>(1,30);
|
|
|
|
dA.setRandom(size,size);
|
|
|
|
|
|
|
|
dA = (dA.array().abs()<0.3).select(0,dA);
|
|
|
|
dA.diagonal() = (dA.diagonal().array()==0).select(1,dA.diagonal());
|
|
|
|
A = dA.sparseView();
|
|
|
|
A.makeCompressed();
|
|
|
|
|
2012-03-29 20:32:54 +08:00
|
|
|
check_sparse_determinant(solver, A, dA);
|
|
|
|
}
|
2011-10-11 17:32:26 +08:00
|
|
|
}
|
2014-10-17 22:52:56 +08:00
|
|
|
|
|
|
|
template<typename Solver> void check_sparse_square_abs_determinant(Solver& solver)
|
|
|
|
{
|
|
|
|
typedef typename Solver::MatrixType Mat;
|
|
|
|
typedef typename Mat::Scalar Scalar;
|
|
|
|
typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
|
|
|
|
|
|
|
|
for (int i = 0; i < g_repeat; i++) {
|
2015-02-17 02:09:48 +08:00
|
|
|
// generate the problem
|
|
|
|
Mat A;
|
|
|
|
DenseMatrix dA;
|
|
|
|
generate_sparse_square_problem(solver, A, dA, 30);
|
|
|
|
A.makeCompressed();
|
2014-10-17 22:52:56 +08:00
|
|
|
check_sparse_abs_determinant(solver, A, dA);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|