Skip to content
Snippets Groups Projects
Commit 2693915e authored by Elias Pipping's avatar Elias Pipping
Browse files

Incorporate addToDiagonal from StaticMatrix

parent 273fa630
No related branches found
No related tags found
No related merge requests found
#install headers
install(FILES
addtodiagonal.hh
axpy.hh
componentwisematrixmap.hh
genericvectortools.hh
......
#ifndef DUNE_MATRIX_VECTOR_ADDTODIAGONAL_HH
#define DUNE_MATRIX_VECTOR_ADDTODIAGONAL_HH
#include <dune/istl/scaledidmatrix.hh>
namespace Dune {
namespace MatrixVector {
template <class Matrix>
static void addToDiagonal(Matrix& x, const typename Matrix::field_type a) {
for (int i = 0; i < Matrix::rows; ++i)
x[i][i] += a;
}
/*
the line
template <typename FieldType, int n> static void
addToDiagonal(Dune::ScaledIdentityMatrix<FieldType,n>& x, const
FieldType a)
should NOT be used as it may lead to ambiguities in which case the
general method above will be used. an example is the line taken
from massassembler.hh (line 83 at the time):
StaticMatrix::addToDiagonal(localMatrix[i][i], values[i] * zi);
where values[i] is of type FieldVector<FieldType,1> and zi is a
double. This call then does not exactly fit this specialization
(without the implicit cast of FV<FT,1>) and hence some wild
template voodoo decides which of the templates is to be taken - in
this case with a gcc 4.4.5 that was the one above leading to a
wrong factor n in the diagonal value
*/
template <typename FieldType, int n>
static void addToDiagonal(
Dune::ScaledIdentityMatrix<FieldType, n>& x,
const typename Dune::ScaledIdentityMatrix<FieldType, n>::field_type a) {
x.scalar() += a;
}
}
}
#endif
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