Skip to content
Snippets Groups Projects
Commit 52260791 authored by Oliver Sander's avatar Oliver Sander
Browse files

Merge branch 'fix-method-printsparsematrix' into 'master'

Various fixes and tests for the printSparseMatrix method

See merge request !558
parents 338462bb 965bbada
No related branches found
No related tags found
1 merge request!558Various fixes and tests for the printSparseMatrix method
Pipeline #68273 passed
......@@ -246,6 +246,36 @@ namespace Dune {
s.precision(oldprec);
}
namespace Impl
{
template<class B>
void printInnerMatrixElement(std::ostream& s,
const B& innerMatrixElement,
int innerrow, int innercol)
{
s<<innerMatrixElement<<" ";
}
template<class B, int n>
void printInnerMatrixElement(std::ostream& s,
const ScaledIdentityMatrix<B,n> innerMatrixElement,
int innerrow, int innercol)
{
if (innerrow == innercol)
s<<innerMatrixElement.scalar()<<" ";
else
s<<"-";
}
template<class B, int n, int m>
void printInnerMatrixElement(std::ostream& s,
const FieldMatrix<B,n,m> innerMatrixElement,
int innerrow, int innercol)
{
s<<innerMatrixElement[innerrow][innercol]<<" ";
}
}
/**
* \brief Prints a BCRSMatrix with fixed sized blocks.
*
......@@ -290,8 +320,8 @@ namespace Dune {
typedef typename Matrix::ConstRowIterator Row;
int n = InnerMatrixType::rows;
int m = InnerMatrixType::cols;
constexpr int n = std::decay_t<decltype(Impl::asMatrix(std::declval<InnerMatrixType>()))>::rows;
constexpr int m = std::decay_t<decltype(Impl::asMatrix(std::declval<InnerMatrixType>()))>::cols;
for(Row row=mat.begin(); row != mat.end(); ++row) {
int skipcols=0;
bool reachedEnd=false;
......@@ -325,7 +355,7 @@ namespace Dune {
}
for(int innercol=0; innercol < m; ++innercol) {
s.width(9);
printInnerMatrixElement(s,*col,innerrow,innercol);
Impl::printInnerMatrixElement(s,*col,innerrow,innercol);
}
s<<"|";
......@@ -346,25 +376,6 @@ namespace Dune {
s.precision(oldprec);
}
template<class B, int n>
void printInnerMatrixElement(std::ostream& s,
const ScaledIdentityMatrix<B,n> innerMatrixElement,
int innerrow, int innercol)
{
if (innerrow == innercol)
s<<innerMatrixElement.scalar()<<" ";
else
s<<"-";
}
template<class B, int n, int m, class A>
void printInnerMatrixElement(std::ostream& s,
const FieldMatrix<B,n,m> innerMatrixElement,
int innerrow, int innercol)
{
s<<innerMatrixElement[innerrow][innercol]<<" ";
}
namespace
{
template<typename T>
......
......@@ -89,22 +89,30 @@ int main(int argc, char** argv)
v3.push_back(2);
Dune::writeVectorToMatlabHelper(v3, std::cout);
// Test the printmatrix method
// Test the printmatrix and printSparseMatrix methods
// BCRSMatrix
{
Dune::BCRSMatrix<double> matrix;
setupLaplacian(matrix, 3);
Dune::printmatrix(std::cout, matrix, "BCRSMatrix<double>", "--");
Dune::printSparseMatrix(std::cout, matrix, "BCRSMatrix<double>", "--");
}
{
Dune::BCRSMatrix<Dune::FieldMatrix<double,1,1> > matrix;
setupLaplacian(matrix, 3);
Dune::printmatrix(std::cout, matrix, "BCRSMatrix<FieldMatrix<double,1,1> >", "--");
Dune::printSparseMatrix(std::cout, matrix, "BCRSMatrix<FieldMatrix<double,1,1> >", "--");
}
{
Dune::BCRSMatrix<Dune::ScaledIdentityMatrix<double,1> > matrix;
setupLaplacian(matrix, 3);
Dune::printSparseMatrix(std::cout, matrix, "BCRSMatrix<ScaledIdentityMatrix<double,1> >", "--");
}
{
Dune::BCRSMatrix<Dune::FieldMatrix<double,2,3> > matrix;
setupLaplacian(matrix, 3);
Dune::printmatrix(std::cout, matrix, "BCRSMatrix<FieldMatrix<double,2,3> >", "--");
Dune::printSparseMatrix(std::cout, matrix, "BCRSMatrix<FieldMatrix<double,2,3> >", "--");
}
// Matrix
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment