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

[!268] Implement MatrixDimension for dense Matrices

Merge branch 'implement-printmatrix-matrix-double' into 'master'

See merge request [!268]

  [!268]: Nonecore/dune-istl/merge_requests/268
parents 9e485731 bb9f42a1
No related branches found
No related tags found
1 merge request!268Implement MatrixDimension for dense Matrices
Pipeline #15597 passed
......@@ -219,6 +219,40 @@ namespace Dune
}
};
// Default implementation for scalar types
template<typename B, typename TA>
struct MatrixDimension<Matrix<B,TA> >
{
using block_type = typename Matrix<B,TA>::block_type;
using size_type = typename Matrix<B,TA>::size_type;
static size_type rowdim (const Matrix<B,TA>& A, size_type i)
{
return MatrixDimension<block_type>::rowdim(A[i][0]);
}
static size_type coldim (const Matrix<B,TA>& A, size_type c)
{
return MatrixDimension<block_type>::coldim(A[0][c]);
}
static size_type rowdim (const Matrix<B,TA>& A)
{
size_type nn=0;
for (size_type i=0; i<A.N(); i++)
nn += rowdim(A,i);
return nn;
}
static size_type coldim (const Matrix<B,TA>& A)
{
size_type nn=0;
for (size_type i=0; i<A.M(); i++)
nn += coldim(A,i);
return nn;
}
};
template<typename B, typename TA>
struct MatrixDimension<BCRSMatrix<B,TA> >
......
......@@ -5,6 +5,7 @@
#include <dune/common/diagonalmatrix.hh>
#include <dune/istl/scaledidmatrix.hh>
#include <dune/istl/bcrsmatrix.hh>
#include <dune/istl/matrix.hh>
#include <dune/istl/io.hh>
#include "laplacian.hh"
......@@ -85,19 +86,37 @@ int main(int argc, char** argv)
Dune::writeVectorToMatlabHelper(v3, std::cout);
// Test the printmatrix method
// BCRSMatrix
{
Dune::BCRSMatrix<double> matrix;
setupLaplacian(matrix, 3);
Dune::printmatrix(std::cout, matrix, "matrix", "--");
Dune::printmatrix(std::cout, matrix, "BCRSMatrix<double>", "--");
}
{
Dune::BCRSMatrix<Dune::FieldMatrix<double,1,1> > matrix;
setupLaplacian(matrix, 3);
Dune::printmatrix(std::cout, matrix, "matrix", "--");
Dune::printmatrix(std::cout, matrix, "BCRSMatrix<FieldMatrix<double,1,1> >", "--");
}
{
Dune::BCRSMatrix<Dune::FieldMatrix<double,2,3> > matrix;
setupLaplacian(matrix, 3);
Dune::printmatrix(std::cout, matrix, "matrix", "--");
Dune::printmatrix(std::cout, matrix, "BCRSMatrix<FieldMatrix<double,2,3> >", "--");
}
// Matrix
{
Dune::Matrix<double> matrix(3,3);
matrix = 0;
Dune::printmatrix(std::cout, matrix, "Matrix<double>", "--");
}
{
Dune::Matrix<Dune::FieldMatrix<double,1,1> > matrix(3,3);
matrix = 0;
Dune::printmatrix(std::cout, matrix, "Matrix<FieldMatrix<double,1,1> >", "--");
}
{
Dune::Matrix<Dune::FieldMatrix<double,2,3> > matrix(3,3);
matrix = 0;
Dune::printmatrix(std::cout, matrix, "Matrix<FieldMatrix<double,2,3> >", "--");
}
}
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment