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

[!267] Implement printmatrix for scalar entries

Merge branch 'implement-printmatrix-for-scalar-entries' into 'master'

See merge request [!267]

  [!267]: Nonecore/dune-istl/merge_requests/267
parents fdcfc577 142b88f2
No related branches found
No related tags found
1 merge request!267Implement printmatrix for scalar entries
Pipeline #15580 failed
......@@ -151,6 +151,31 @@ namespace Dune {
}
}
/**
* \brief Print one row of a matrix, specialization for number types
*
* \code
* #include <dune/istl/io.hh>
* \endcode
*/
template<class K>
void print_row (std::ostream& s, const K& value,
typename FieldMatrix<K,1,1>::size_type I,
typename FieldMatrix<K,1,1>::size_type J,
typename FieldMatrix<K,1,1>::size_type therow,
int width, int precision,
typename std::enable_if_t<Dune::IsNumber<K>::value>* sfinae = nullptr)
{
DUNE_UNUSED_PARAMETER(I);
DUNE_UNUSED_PARAMETER(J);
DUNE_UNUSED_PARAMETER(therow);
DUNE_UNUSED_PARAMETER(precision);
s << " "; // space in front of each entry
s.width(width); // set width for each entry anew
s << value;
}
/**
* \brief Print one row of a matrix
*
......@@ -161,7 +186,8 @@ namespace Dune {
template<class M>
void print_row (std::ostream& s, const M& A, typename M::size_type I,
typename M::size_type J, typename M::size_type therow,
int width, int precision)
int width, int precision,
typename std::enable_if_t<!Dune::IsNumber<M>::value>* sfinae = nullptr)
{
typename M::size_type i0=I;
for (typename M::size_type i=0; i<A.N(); i++)
......@@ -190,60 +216,6 @@ namespace Dune {
}
}
/**
* \brief Print one row of a matrix, specialization for FieldMatrix
*
* \code
* #include <dune/istl/io.hh>
* \endcode
*/
template<class K, int n, int m>
void print_row (std::ostream& s, const FieldMatrix<K,n,m>& A,
typename FieldMatrix<K,n,m>::size_type I,
typename FieldMatrix<K,n,m>::size_type J,
typename FieldMatrix<K,n,m>::size_type therow, int width,
int precision)
{
DUNE_UNUSED_PARAMETER(J);
DUNE_UNUSED_PARAMETER(precision);
typedef typename FieldMatrix<K,n,m>::size_type size_type;
for (size_type i=0; i<n; i++)
if (I+i==therow)
for (int j=0; j<m; j++)
{
s << " "; // space in front of each entry
s.width(width); // set width for each entry anew
s << A[i][j]; // yeah, the number !
}
}
/**
* \brief Print one row of a matrix, specialization for FieldMatrix<K,1,1>
*
* \code
* #include <dune/istl/io.hh>
* \endcode
*/
template<class K>
void print_row (std::ostream& s, const FieldMatrix<K,1,1>& A,
typename FieldMatrix<K,1,1>::size_type I,
typename FieldMatrix<K,1,1>::size_type J,
typename FieldMatrix<K,1,1>::size_type therow,
int width, int precision)
{
DUNE_UNUSED_PARAMETER(J);
DUNE_UNUSED_PARAMETER(precision);
if (I==therow)
{
s << " "; // space in front of each entry
s.width(width); // set width for each entry anew
s << static_cast<K>(A); // yeah, the number !
}
}
/**
* \brief Print a generic block matrix
*
......@@ -502,12 +474,7 @@ namespace Dune {
outStream.precision(oldPrecision);
}
// Write vector entries to a stream
// TODO: The writeVectorToMatlab method does not actually need this helper anymore:
// As there is no template specialization involved, the code may as well
// be simply folded into writeVectorToMatlab. On the other hand,
// writeVectorToMatlabHelper writes the vector content into a stream,
// which may be useful to some people.
// Recursively write vector entries to a stream
template<class V>
void writeVectorToMatlabHelper (const V& v, std::ostream& stream)
{
......
......@@ -83,4 +83,21 @@ int main(int argc, char** argv)
v3.push_back(1);
v3.push_back(2);
Dune::writeVectorToMatlabHelper(v3, std::cout);
// Test the printmatrix method
{
Dune::BCRSMatrix<double> matrix;
setupLaplacian(matrix, 3);
Dune::printmatrix(std::cout, matrix, "matrix", "--");
}
{
Dune::BCRSMatrix<Dune::FieldMatrix<double,1,1> > matrix;
setupLaplacian(matrix, 3);
Dune::printmatrix(std::cout, matrix, "matrix", "--");
}
{
Dune::BCRSMatrix<Dune::FieldMatrix<double,2,3> > matrix;
setupLaplacian(matrix, 3);
Dune::printmatrix(std::cout, matrix, "matrix", "--");
}
}
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment