Skip to content
Snippets Groups Projects
Commit ec3db135 authored by Carsten Gräser's avatar Carsten Gräser Committed by Martin Nolte
Browse files

[test] Add check for BCRSMatrix::nonzeroe()

parent 572d972c
No related branches found
No related tags found
1 merge request!93[test] Add check for BCRSMatrix::nonzeroes()
Pipeline #
......@@ -4,8 +4,20 @@
#include <dune/istl/bcrsmatrix.hh>
#include <dune/common/fmatrix.hh>
#include <dune/common/exceptions.hh>
#include <dune/common/test/testsuite.hh>
#include <dune/istl/io.hh>
template<class M>
std::size_t computeNNZ(M&& matrix)
{
std::size_t nnz = 0;
for (auto&& row : matrix)
for (auto&& entry : row)
++nnz;
return nnz;
}
template<class M>
struct Builder
{
......@@ -15,11 +27,14 @@ struct Builder
}
};
template<class B, class A>
struct Builder<Dune::BCRSMatrix<B,A> >
{
void randomBuild(int rows, int cols)
Dune::TestSuite randomBuild(int rows, int cols)
{
Dune::TestSuite testSuite("BCRSMatrix::random");
int maxNZCols = 15; // maximal number of nonzeros per row
{
Dune::BCRSMatrix<B,A> matrix( rows, cols, Dune::BCRSMatrix<B,A>::random );
......@@ -43,6 +58,10 @@ struct Builder<Dune::BCRSMatrix<B,A> >
Dune::printmatrix(std::cout, matrix, "random", "row");
Dune::BCRSMatrix<B,A> matrix1(matrix);
Dune::printmatrix(std::cout, matrix, "random", "row");
auto nnz = computeNNZ(matrix);
testSuite.check(nnz == matrix.nonzeroes(), "nnz after building")
<< "BCRSMatrix::nonzeroes() returns " << matrix.nonzeroes() << " while there are " << nnz << " stored entries.";
}
/*{
......@@ -63,11 +82,11 @@ struct Builder<Dune::BCRSMatrix<B,A> >
Dune::printmatrix(std::cout, matrix, "random", "row");
}*/
return testSuite;
}
void rowWiseBuild(Dune::BCRSMatrix<B,A>& matrix, int /* rows */, int cols)
void rowWiseBuild(Dune::TestSuite& testSuite, Dune::BCRSMatrix<B,A>& matrix, int /* rows */, int cols)
{
for(typename Dune::BCRSMatrix<B,A>::CreateIterator ci=matrix.createbegin(), cend=matrix.createend();
ci!=cend; ++ci)
{
......@@ -80,26 +99,49 @@ struct Builder<Dune::BCRSMatrix<B,A> >
ci.insert(i+1);
}
{
auto nnz = computeNNZ(matrix);
testSuite.check(nnz == matrix.nonzeroes(), "nnz after building")
<< "BCRSMatrix::nonzeroes() returns " << matrix.nonzeroes() << " while there are " << nnz << " stored entries.";
}
Dune::printmatrix(std::cout, matrix, "row_wise", "row");
// test copy ctor
Dune::BCRSMatrix<B,A> matrix1(matrix);
Dune::printmatrix(std::cout, matrix1, "row_wise", "row");
{
auto nnz = computeNNZ(matrix1);
testSuite.check(nnz == matrix1.nonzeroes(), "nnz after constructed from")
<< "BCRSMatrix::nonzeroes() returns " << matrix1.nonzeroes() << " while there are " << nnz << " stored entries.";
}
// test copy assignment
Dune::BCRSMatrix<B,A> matrix2;
matrix2 = matrix;
Dune::printmatrix(std::cout, matrix2, "row_wise", "row");
{
auto nnz = computeNNZ(matrix2);
testSuite.check(nnz == matrix2.nonzeroes(), "nnz after assigned from")
<< "BCRSMatrix::nonzeroes() returns " << matrix2.nonzeroes() << " while there are " << nnz << " stored entries.";
}
}
void rowWiseBuild(int rows, int cols)
Dune::TestSuite rowWiseBuild(int rows, int cols)
{
Dune::TestSuite testSuite("BCRSMatrix::row_wise without nnz");
Dune::BCRSMatrix<B,A> matrix( rows, cols, Dune::BCRSMatrix<B,A>::row_wise );
rowWiseBuild(matrix, rows, cols);
rowWiseBuild(testSuite, matrix, rows, cols);
return testSuite;
}
void rowWiseBuild(int rows, int cols, int nnz)
Dune::TestSuite rowWiseBuild(int rows, int cols, int nnz)
{
Dune::TestSuite testSuite("BCRSMatrix::row_wise with nnz");
Dune::BCRSMatrix<B,A> matrix( rows, cols, nnz, Dune::BCRSMatrix<B,A>::row_wise );
rowWiseBuild(matrix, rows, cols);
rowWiseBuild(testSuite, matrix, rows, cols);
return testSuite;
}
};
......@@ -116,11 +158,17 @@ void testDoubleSetSize()
int main()
{
try{
Dune::TestSuite testSuite;
Builder<Dune::BCRSMatrix<Dune::FieldMatrix<double,1,1> > > builder;
builder.randomBuild(5,4);
builder.rowWiseBuild(5,4,13);
builder.rowWiseBuild(5,4);
testSuite.subTest(builder.randomBuild(5,4));
testSuite.subTest(builder.rowWiseBuild(5,4,13));
testSuite.subTest(builder.rowWiseBuild(5,4));
testDoubleSetSize();
return testSuite.exit();
}catch(Dune::Exception e) {
std::cerr << e<<std::endl;
return 1;
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment