mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-20 11:48:25 -06:00
Experiment with fast diagonal-sparse products.
This commit is contained in:
parent
1f32594f79
commit
11a33b3017
@ -343,11 +343,13 @@ namespace Opm
|
||||
assert(lhs.type_ == D);
|
||||
assert(rhs.type_ == S);
|
||||
AutoDiffMatrix retval;
|
||||
Eigen::SparseMatrix<double> diag = spdiag(lhs.d_.diagonal());
|
||||
// Eigen::SparseMatrix<double> diag = spdiag(lhs.d_.diagonal());
|
||||
retval.type_ = S;
|
||||
retval.rows_ = lhs.rows_;
|
||||
retval.cols_ = rhs.cols_;
|
||||
fastSparseProduct(diag, rhs.s_, retval.s_);
|
||||
// fastSparseProduct(diag, rhs.s_, retval.s_);
|
||||
fastDiagSparseProduct(lhs.d_, rhs.s_, retval.s_);
|
||||
// retval.s_ = lhs.d_ * rhs.s_;
|
||||
return retval;
|
||||
}
|
||||
|
||||
@ -356,11 +358,13 @@ namespace Opm
|
||||
assert(lhs.type_ == S);
|
||||
assert(rhs.type_ == D);
|
||||
AutoDiffMatrix retval;
|
||||
Eigen::SparseMatrix<double> diag = spdiag(rhs.d_.diagonal());
|
||||
// Eigen::SparseMatrix<double> diag = spdiag(rhs.d_.diagonal());
|
||||
retval.type_ = S;
|
||||
retval.rows_ = lhs.rows_;
|
||||
retval.cols_ = rhs.cols_;
|
||||
fastSparseProduct(lhs.s_, diag, retval.s_);
|
||||
// fastSparseProduct(lhs.s_, diag, retval.s_);
|
||||
fastSparseDiagProduct(lhs.s_, rhs.d_, retval.s_);
|
||||
// retval.s_ = lhs.s_ * rhs.d_;
|
||||
return retval;
|
||||
}
|
||||
|
||||
|
@ -141,6 +141,44 @@ void fastSparseProduct(const Lhs& lhs, const Rhs& rhs, ResultType& res)
|
||||
|
||||
|
||||
|
||||
|
||||
inline void fastDiagSparseProduct(const Eigen::DiagonalMatrix<double, Eigen::Dynamic>& lhs,
|
||||
const Eigen::SparseMatrix<double>& rhs,
|
||||
Eigen::SparseMatrix<double>& res)
|
||||
{
|
||||
res = rhs;
|
||||
|
||||
// Multiply rows by diagonal lhs.
|
||||
int n = res.cols();
|
||||
for (int col = 0; col < n; ++col) {
|
||||
typedef Eigen::SparseMatrix<double>::InnerIterator It;
|
||||
for (It it(res, col); it; ++it) {
|
||||
it.valueRef() *= rhs.diagonal()(it.row());
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
|
||||
|
||||
inline void fastSparseDiagProduct(const Eigen::SparseMatrix<double>& lhs,
|
||||
const Eigen::DiagonalMatrix<double, Eigen::Dynamic>& rhs,
|
||||
Eigen::SparseMatrix<double>& res)
|
||||
{
|
||||
res = lhs;
|
||||
|
||||
// Multiply columns by diagonal rhs.
|
||||
int n = res.cols();
|
||||
for (int col = 0; col < n; ++col) {
|
||||
typedef Eigen::SparseMatrix<double>::InnerIterator It;
|
||||
for (It it(res, col); it; ++it) {
|
||||
it.valueRef() *= rhs.diagonal()(col);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
|
||||
} // end namespace Opm
|
||||
|
||||
#endif // OPM_FASTSPARSEPRODUCT_HEADER_INCLUDED
|
||||
|
Loading…
Reference in New Issue
Block a user