AutoDiffMatrix: use fastSparse{Add,Substract} when sparsisty patterns match.

This commit is contained in:
Robert Kloefkorn 2016-02-12 16:00:24 +01:00
parent 704603e2b2
commit d5b6566e06
2 changed files with 104 additions and 2 deletions

View File

@ -264,7 +264,12 @@ namespace Opm
AutoDiffMatrix& operator+=(const AutoDiffMatrix& rhs)
{
*this = *this + rhs;
if( type_ == Sparse && rhs.type_ == Sparse )
{
fastSparseAdd( sparse_, rhs.sparse_ );
}
else
*this = *this + rhs;
return *this;
}
@ -275,7 +280,12 @@ namespace Opm
AutoDiffMatrix& operator-=(const AutoDiffMatrix& rhs)
{
*this = *this + (rhs * -1.0);
if( type_ == Sparse && rhs.type_ == Sparse )
{
fastSparseSubstract( sparse_, rhs.sparse_ );
}
else
*this = *this + (rhs * -1.0);
return *this;
}

View File

@ -178,6 +178,98 @@ inline void fastSparseDiagProduct(const Eigen::SparseMatrix<double>& lhs,
}
namespace detail {
struct AddTo
{
template <class T>
inline void operator()( T& a, const T& b ) const
{
a += b;
}
};
struct SubstractFrom
{
template <class T>
inline void operator()( T& a, const T& b ) const
{
a -= b;
}
};
} // end namespace detail
// this function substracts two sparse matrices
// if the sparsity pattern is the same a faster add/substract is performed
template<typename Lhs, typename Rhs>
void fastSparseBinaryOp(Lhs& lhs, const Rhs& rhs, const bool add )
{
// if non zeros don't match, matrices don't have the same sparsity pattern
bool equal = lhs.nonZeros() == rhs.nonZeros();
typedef typename Eigen::internal::remove_all<Lhs>::type::Scalar Scalar;
typedef typename Eigen::internal::remove_all<Lhs>::type::Index Index;
const Index nnz = lhs.nonZeros();
const Index* rhsJ = rhs.innerIndexPtr();
const Index* lhsJ = lhs.innerIndexPtr();
// check complete sparsity pattern
if( equal )
{
for(Index i=0; i<nnz; ++i )
{
equal &= ( lhsJ[ i ] == rhsJ[ i ] );
}
}
// if sparsity pattern does not match,
// fallback to default implementation
if( ! equal )
{
if( add ) {
lhs += rhs;
}
else {
lhs -= rhs;
}
return;
}
// fast add/substract using only the data pointers
const Scalar* rhsV = rhs.valuePtr();
Scalar* lhsV = lhs.valuePtr();
if( add )
{
for(Index i=0; i<nnz; ++i )
{
lhsV[ i ] += rhsV[ i ];
}
}
else
{
for(Index i=0; i<nnz; ++i )
{
lhsV[ i ] -= rhsV[ i ];
}
}
}
template<typename Lhs, typename Rhs>
void fastSparseAdd(Lhs& lhs, const Rhs& rhs )
{
fastSparseBinaryOp( lhs, rhs, true ); //detail::AddTo() );
}
template<typename Lhs, typename Rhs>
void fastSparseSubstract(Lhs& lhs, const Rhs& rhs )
{
fastSparseBinaryOp( lhs, rhs, false ); //detail::SubstractFrom() );
}
} // end namespace Opm