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

Implement printSparseMatrix for scalar inner types

And test it!
parent 6a999a4e
No related branches found
No related tags found
1 merge request!558Various fixes and tests for the printSparseMatrix method
......@@ -248,6 +248,14 @@ namespace Dune {
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,
......@@ -312,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;
......
......@@ -95,6 +95,7 @@ int main(int argc, char** argv)
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;
......
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