Use std::vector instead of DiagonalMatrix.

This is because DiagonalMatrix lacks a swap() method.
This commit is contained in:
Atgeirr Flø Rasmussen
2015-08-27 09:33:30 +02:00
committed by babrodtk
parent 87f677af02
commit 3c905845f9
2 changed files with 448 additions and 433 deletions

View File

@@ -28,6 +28,8 @@
#include <opm/core/utility/platform_dependent/reenable_warnings.h>
#include <opm/autodiff/fastSparseProduct.hpp>
#include <vector>
namespace Opm
{
@@ -69,7 +71,7 @@ namespace Opm
: type_(D),
rows_(d.rows()),
cols_(d.cols()),
d_(d)
d_(d.diagonal().array().data(), d.diagonal().array().data() + d.rows())
{
}
@@ -85,6 +87,17 @@ namespace Opm
void swap(AutoDiffMatrix& other)
{
std::swap(type_, other.type_);
std::swap(rows_, other.rows_);
std::swap(cols_, other.cols_);
d_.swap(other.d_);
s_.swap(other.s_);
}
AutoDiffMatrix operator+(const AutoDiffMatrix& rhs) const
{
assert(rows_ == rhs.rows_);
@@ -207,13 +220,15 @@ namespace Opm
{
AutoDiffMatrix retval(*this);
retval.type_ = D;
retval.d_ = rhs * Eigen::VectorXd::Ones(rows_).asDiagonal();
retval.d_.assign(rows_, rhs);
return retval;
}
case D:
{
AutoDiffMatrix retval(*this);
retval.d_ = rhs * retval.d_;
for (double& elem : retval.d_) {
elem *= rhs;
}
return retval;
}
case S:
@@ -239,7 +254,7 @@ namespace Opm
case I:
return rhs;
case D:
return d_ * rhs;
return Eigen::Map<const Eigen::VectorXd>(d_.data(), rows_) * rhs;
case S:
return s_ * rhs;
}
@@ -258,17 +273,18 @@ namespace Opm
retval.type_ = D;
retval.rows_ = lhs.rows_;
retval.cols_ = rhs.cols_;
retval.d_ = Eigen::VectorXd::Constant(lhs.rows_, 2.0).asDiagonal();
retval.d_.assign(lhs.rows_, 2.0);
return retval;
}
static AutoDiffMatrix sumDI(const AutoDiffMatrix& lhs, const AutoDiffMatrix& rhs)
{
static_cast<void>(rhs); // Silence release-mode warning.
assert(lhs.type_ == D);
assert(rhs.type_ == I);
AutoDiffMatrix retval = lhs;
for (int r = 0; r < lhs.rows_; ++r) {
retval.d_.diagonal()(r) += 1.0;
retval.d_[r] += 1.0;
}
return retval;
}
@@ -279,7 +295,7 @@ namespace Opm
assert(rhs.type_ == D);
AutoDiffMatrix retval = lhs;
for (int r = 0; r < lhs.rows_; ++r) {
retval.d_.diagonal()(r) += rhs.d_.diagonal()(r);
retval.d_[r] += rhs.d_[r];
}
return retval;
}
@@ -302,7 +318,7 @@ 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_);
retval.type_ = S;
retval.rows_ = lhs.rows_;
retval.cols_ = rhs.cols_;
@@ -334,7 +350,10 @@ namespace Opm
retval.type_ = D;
retval.rows_ = lhs.rows_;
retval.cols_ = rhs.cols_;
retval.d_ = (lhs.d_.diagonal().array() * rhs.d_.diagonal().array()).matrix().asDiagonal();
retval.d_.resize(lhs.rows_);
for (int r = 0; r < lhs.rows_; ++r) {
retval.d_[r] = lhs.d_[r] * rhs.d_[r];
}
return retval;
}
@@ -343,13 +362,10 @@ namespace Opm
assert(lhs.type_ == D);
assert(rhs.type_ == S);
AutoDiffMatrix retval;
// 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_);
fastDiagSparseProduct(lhs.d_, rhs.s_, retval.s_);
// retval.s_ = lhs.d_ * rhs.s_;
return retval;
}
@@ -358,13 +374,10 @@ namespace Opm
assert(lhs.type_ == S);
assert(rhs.type_ == D);
AutoDiffMatrix retval;
// 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_);
fastSparseDiagProduct(lhs.s_, rhs.d_, retval.s_);
// retval.s_ = lhs.s_ * rhs.d_;
return retval;
}
@@ -392,7 +405,7 @@ namespace Opm
s = spdiag(Eigen::VectorXd::Ones(rows_));
return;
case D:
s = spdiag(d_.diagonal());
s = spdiag(d_);
return;
case S:
s = s_;
@@ -434,7 +447,7 @@ namespace Opm
case I:
return (row == col) ? 1.0 : 0.0;
case D:
return (row == col) ? d_.diagonal()(row) : 0.0;
return (row == col) ? d_[row] : 0.0;
case S:
return s_.coeff(row, col);
}
@@ -445,7 +458,7 @@ namespace Opm
MatrixType type_;
int rows_;
int cols_;
Eigen::DiagonalMatrix<double, Eigen::Dynamic> d_;
std::vector<double> d_;
Eigen::SparseMatrix<double> s_;
template <class V>

View File

@@ -142,7 +142,8 @@ void fastSparseProduct(const Lhs& lhs, const Rhs& rhs, ResultType& res)
inline void fastDiagSparseProduct(const Eigen::DiagonalMatrix<double, Eigen::Dynamic>& lhs,
inline void fastDiagSparseProduct(// const Eigen::DiagonalMatrix<double, Eigen::Dynamic>& lhs,
const std::vector<double>& lhs,
const Eigen::SparseMatrix<double>& rhs,
Eigen::SparseMatrix<double>& res)
{
@@ -153,7 +154,7 @@ inline void fastDiagSparseProduct(const Eigen::DiagonalMatrix<double, Eigen::Dyn
for (int col = 0; col < n; ++col) {
typedef Eigen::SparseMatrix<double>::InnerIterator It;
for (It it(res, col); it; ++it) {
it.valueRef() *= lhs.diagonal()(it.row());
it.valueRef() *= lhs[it.row()]; // lhs.diagonal()(it.row());
}
}
}
@@ -162,7 +163,8 @@ inline void fastDiagSparseProduct(const Eigen::DiagonalMatrix<double, Eigen::Dyn
inline void fastSparseDiagProduct(const Eigen::SparseMatrix<double>& lhs,
const Eigen::DiagonalMatrix<double, Eigen::Dynamic>& rhs,
// const Eigen::DiagonalMatrix<double, Eigen::Dynamic>& rhs,
const std::vector<double>& rhs,
Eigen::SparseMatrix<double>& res)
{
res = lhs;
@@ -172,7 +174,7 @@ inline void fastSparseDiagProduct(const Eigen::SparseMatrix<double>& lhs,
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);
it.valueRef() *= rhs[col]; // rhs.diagonal()(col);
}
}
}