New implementation of AutoDiffMatrix, some tests.

Compiles and tests successfully, but test coverage very
limited. New approach based on relatively primitive
run-time switching instead of trying to use inheritance.
This commit is contained in:
Atgeirr Flø Rasmussen 2015-08-24 13:55:16 +02:00 committed by babrodtk
parent 6a5a48e728
commit 47e7dbe943
2 changed files with 381 additions and 240 deletions

View File

@ -1,5 +1,5 @@
/* /*
Copyright 2014 SINTEF ICT, Applied Mathematics. Copyright 2014, 2015 SINTEF ICT, Applied Mathematics.
This file is part of the Open Porous Media project (OPM). This file is part of the Open Porous Media project (OPM).
@ -31,217 +31,308 @@
namespace Opm namespace Opm
{ {
/// Implementation details for class AutoDiffMatrix. class AutoDiffMatrix
namespace AutoDiffMatrixDetail
{ {
public:
AutoDiffMatrix()
: type_(Z),
rows_(0),
cols_(0)
{
}
class Zero;
class Identity; enum CreationType { ZeroMatrix, IdentityMatrix };
class Diagonal;
class Sparse;
typedef std::shared_ptr<Zero> ZeroMat; AutoDiffMatrix(const CreationType t, const int rows)
typedef std::shared_ptr<Identity> IdentityMat; : type_(t == ZeroMatrix ? Z : I),
typedef std::shared_ptr<Diagonal> DiagonalMat; rows_(rows),
typedef std::shared_ptr<Sparse> SparseMat; cols_(rows)
{
}
explicit AutoDiffMatrix(const Eigen::DiagonalMatrix<double, Eigen::Dynamic>& d)
: type_(D),
rows_(d.rows()),
cols_(d.cols()),
d_(d)
{
}
explicit AutoDiffMatrix(const Eigen::SparseMatrix<double>& s)
: type_(S),
rows_(s.rows()),
cols_(s.cols()),
s_(s)
{
}
AutoDiffMatrix operator+(const AutoDiffMatrix& rhs) const
{
switch (type_) {
case Z:
return rhs;
case I:
switch (rhs.type_) {
case Z:
return *this;
case I:
return sumII(*this, rhs);
case D:
return rhs + (*this);
case S:
return rhs + (*this);
}
case D:
switch (rhs.type_) {
case Z:
return *this;
case I:
return sumDI(*this, rhs);
case D:
return sumDD(*this, rhs);
case S:
return rhs + (*this);
}
case S:
switch (rhs.type_) {
case Z:
return *this;
case I:
return sumSI(*this, rhs);
case D:
return sumSD(*this, rhs);
case S:
return sumSS(*this, rhs);
}
}
}
AutoDiffMatrix operator*(const AutoDiffMatrix& rhs) const
{
switch (type_) {
case Z:
return *this;
case I:
switch (rhs.type_) {
case Z:
return rhs;
case I:
return rhs;
case D:
return rhs;
case S:
return rhs;
}
case D:
switch (rhs.type_) {
case Z:
return rhs;
case I:
return *this;
case D:
return prodDD(*this, rhs);
case S:
return prodDS(*this, rhs);
}
case S:
switch (rhs.type_) {
case Z:
return rhs;
case I:
return *this;
case D:
return prodSD(*this, rhs);
case S:
return prodSS(*this, rhs);
}
}
}
static AutoDiffMatrix sumII(const AutoDiffMatrix& lhs, const AutoDiffMatrix& rhs)
{
assert(lhs.type_ == I);
assert(rhs.type_ == I);
AutoDiffMatrix retval;
retval.type_ = D;
retval.rows_ = lhs.rows_;
retval.cols_ = rhs.cols_;
retval.d_ = Eigen::VectorXd::Constant(lhs.rows_, 2.0).asDiagonal();
return retval;
}
class Interface static AutoDiffMatrix sumDI(const AutoDiffMatrix& lhs, const AutoDiffMatrix& rhs)
{ {
public: assert(lhs.type_ == D);
typedef std::shared_ptr<Interface> Mat; assert(rhs.type_ == I);
AutoDiffMatrix retval = lhs;
for (int r = 0; r < lhs.rows_; ++r) {
retval.d_.diagonal()(r) += 1.0;
}
return retval;
}
virtual ~Interface() static AutoDiffMatrix sumDD(const AutoDiffMatrix& lhs, const AutoDiffMatrix& rhs)
{ {
} assert(lhs.type_ == D);
assert(rhs.type_ == D);
AutoDiffMatrix retval = lhs;
for (int r = 0; r < lhs.rows_; ++r) {
retval.d_.diagonal()(r) += rhs.d_.diagonal()(r);
}
return retval;
}
virtual Mat operator+(const Mat& rhs) = 0; static AutoDiffMatrix sumSI(const AutoDiffMatrix& lhs, const AutoDiffMatrix& rhs)
virtual Mat addIdentity(const IdentityMat& rhs) = 0; {
virtual Mat addDiagonal(const DiagonalMat& rhs) = 0; assert(lhs.type_ == S);
virtual Mat addSparse(const SparsMate& rhs) = 0; assert(rhs.type_ == I);
AutoDiffMatrix retval;
Eigen::SparseMatrix<double> ident = spdiag(Eigen::VectorXd::Ones(lhs.rows_));
retval.type_ = S;
retval.rows_ = lhs.rows_;
retval.cols_ = rhs.cols_;
retval.s_ = lhs.s_ + ident;
return retval;
}
virtual Mat operator*(const Mat& rhs) = 0; static AutoDiffMatrix sumSD(const AutoDiffMatrix& lhs, const AutoDiffMatrix& rhs)
virtual Mat leftMulDiagonal(const DiagonalMat& rhs) = 0; {
virtual Mat leftMulSparse(const SparseMat& rhs) = 0; 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_;
retval.s_ = lhs.s_ + diag;
return retval;
}
static AutoDiffMatrix sumSS(const AutoDiffMatrix& lhs, const AutoDiffMatrix& rhs)
{
assert(lhs.type_ == S);
assert(rhs.type_ == S);
AutoDiffMatrix retval;
retval.type_ = S;
retval.rows_ = lhs.rows_;
retval.cols_ = rhs.cols_;
retval.s_ = lhs.s_ + rhs.s_;
return retval;
}
static AutoDiffMatrix prodDD(const AutoDiffMatrix& lhs, const AutoDiffMatrix& rhs)
{
assert(lhs.type_ == D);
assert(rhs.type_ == D);
AutoDiffMatrix retval = lhs;
for (int r = 0; r < lhs.rows_; ++r) {
retval.d_.diagonal().array() *= rhs.d_.diagonal().array();
}
return retval;
}
class Zero : public Interface static AutoDiffMatrix prodDS(const AutoDiffMatrix& lhs, const AutoDiffMatrix& rhs)
{ {
public: assert(lhs.type_ == S);
virtual Mat operator+(const Mat& rhs) assert(rhs.type_ == D);
{ AutoDiffMatrix retval;
return rhs; Eigen::SparseMatrix<double> diag = spdiag(rhs.d_.diagonal());
} retval.type_ = S;
retval.rows_ = lhs.rows_;
retval.cols_ = rhs.cols_;
retval.s_ = lhs.s_ * diag;
return retval;
}
virtual Mat addIdentity(const IdentityMat& rhs) static AutoDiffMatrix prodSD(const AutoDiffMatrix& lhs, const AutoDiffMatrix& rhs)
{ {
return rhs; 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_;
retval.s_ = diag * lhs.s_;
return retval;
}
virtual Mat addDiagonal(const DiagonalMat& rhs) static AutoDiffMatrix prodSS(const AutoDiffMatrix& lhs, const AutoDiffMatrix& rhs)
{ {
return rhs; assert(lhs.type_ == S);
} assert(rhs.type_ == S);
AutoDiffMatrix retval;
virtual Mat addSparse(const SparseMat& rhs) retval.type_ = S;
{ retval.rows_ = lhs.rows_;
return rhs; retval.cols_ = rhs.cols_;
} retval.s_ = lhs.s_ * rhs.s_;
return retval;
}
virtual Mat operator*(const Mat& rhs)
{
return std::make_shared<Zero>();
}
virtual Mat leftMulDiagonal(const DiagonalMat& rhs)
{
return std::make_shared<Zero>();
}
virtual Mat leftMulSparse(const SparseMat& rhs)
{
return std::make_shared<Zero>();
}
};
void toSparse(Eigen::SparseMatrix<double>& s) const
{
switch (type_) {
case Z:
s = Eigen::SparseMatrix<double>(rows_, cols_);
return;
case I:
s = spdiag(Eigen::VectorXd::Ones(rows_));
return;
case D:
s = spdiag(d_.diagonal());
return;
case S:
s = s_;
return;
}
}
private:
enum MatrixType { Z, I, D, S };
MatrixType type_;
int rows_;
int cols_;
Eigen::DiagonalMatrix<double, Eigen::Dynamic> d_;
Eigen::SparseMatrix<double> s_;
class Identity : public Interface template <class V>
{ static inline
public: Eigen::SparseMatrix<double>
virtual Mat operator+(const Mat& rhs) spdiag(const V& d)
{ {
return rhs->addIdentity(*this); typedef Eigen::SparseMatrix<double> M;
} const int n = d.size();
M mat(n, n);
mat.reserve(Eigen::ArrayXi::Ones(n, 1));
for (M::Index i = 0; i < n; ++i) {
mat.insert(i, i) = d[i];
}
virtual Mat addIdentity(const IdentityMat& rhs) return mat;
{ }
return rhs;
}
virtual Mat addDiagonal(const DiagonalMat& rhs) };
{
return rhs;
}
virtual Mat addSparse(const SparseMat& rhs)
{
return rhs;
}
virtual Mat operator*(const Mat& rhs)
{
return std::make_shared<Zero>();
}
virtual Mat leftMulDiagonal(const DiagonalMat& rhs)
{
return std::make_shared<Zero>();
}
virtual Mat leftMulSparse(const SparseMat& rhs)
{
return std::make_shared<Zero>();
}
};
class Diagonal : public Interface
{
public:
virtual Mat operator+(const Mat& rhs)
{
return (*rhs) + (*this);
}
operator+(const IdentityMat& rhs)
{
// TODO return Diagnonal(...);
}
operator+(const DiagonalMat& rhs)
{
// TODO return Diagonal(...);
}
operator+(const SparseMat& rhs)
{
// TODO return Sparse(...);
}
virtual Mat operator*(const Mat& rhs)
{
return (*rhs) * (;
}
Mat operator*(const IdentityMat& rhs)
{
return *this;
}
Mat operator*(const DiagonalMat& rhs)
{
// TODO return Diagonal(...);
}
Mat operator*(const SparseMat& rhs)
{
// TODO return Sparse(...);
}
};
class Sparse : public Interface
{
virtual Mat operator+(const Mat& rhs)
{
return (*rhs) + (*this);
}
operator+(const IdentityMat& rhs)
{
// TODO return Sparse(...);
}
operator+(const DiagonalMat& rhs)
{
// TODO return Sparse(...);
}
operator+(const SparseMat& rhs)
{
// TODO return Sparse(...);
}
virtual Mat operator*(const Mat& rhs)
{
return rhs;
}
Mat operator*(const IdentityMat& rhs)
{
return *this;
}
Mat operator*(const DiagonalMat& rhs)
{
// TODO return Sparse(...);
}
Mat operator*(const SparseMat& rhs)
{
// TODO return Sparse(...);
}
};
} // namespace AutoDiffMatrixDetail
} // namespace Opm } // namespace Opm

View File

@ -26,100 +26,150 @@
#define BOOST_TEST_MODULE AutoDiffMatrixTest #define BOOST_TEST_MODULE AutoDiffMatrixTest
#include <opm/autodiff/AutoDiffMatrix.hpp> #include <opm/autodiff/AutoDiffMatrix.hpp>
#include <opm/autodiff/AutoDiffHelpers.hpp>
#include <boost/test/unit_test.hpp> #include <boost/test/unit_test.hpp>
using namespace Opm::AutoDiffMatrix;
using std::make_shared;
typedef Eigen::SparseMatrix<double> Sp; typedef Eigen::SparseMatrix<double> Sp;
typedef Opm::AutoDiffMatrix Mat;
using namespace Opm;
namespace { bool
template <typename Scalar> operator ==(const Eigen::SparseMatrix<double>& A,
bool const Eigen::SparseMatrix<double>& B)
operator ==(const Eigen::SparseMatrix<Scalar>& A, {
const Eigen::SparseMatrix<Scalar>& B) // Two SparseMatrices are equal if
{ // 0) They have the same ordering (enforced by equal types)
// Two SparseMatrices are equal if // 1) They have the same outer and inner dimensions
// 0) They have the same ordering (enforced by equal types) // 2) They have the same number of non-zero elements
// 1) They have the same outer and inner dimensions // 3) They have the same sparsity structure
// 2) They have the same number of non-zero elements // 4) The non-zero elements are equal
// 3) They have the same sparsity structure
// 4) The non-zero elements are equal
// 1) Outer and inner dimensions // 1) Outer and inner dimensions
bool eq = (A.outerSize() == B.outerSize()); bool eq = (A.outerSize() == B.outerSize());
eq = eq && (A.innerSize() == B.innerSize()); eq = eq && (A.innerSize() == B.innerSize());
// 2) Equal number of non-zero elements // 2) Equal number of non-zero elements
eq = eq && (A.nonZeros() == B.nonZeros()); eq = eq && (A.nonZeros() == B.nonZeros());
for (typename Eigen::SparseMatrix<Scalar>::Index for (typename Eigen::SparseMatrix<double>::Index
k0 = 0, kend = A.outerSize(); eq && (k0 < kend); ++k0) { k0 = 0, kend = A.outerSize(); eq && (k0 < kend); ++k0) {
for (typename Eigen::SparseMatrix<Scalar>::InnerIterator for (typename Eigen::SparseMatrix<double>::InnerIterator
iA(A, k0), iB(B, k0); eq && (iA && iB); ++iA, ++iB) { iA(A, k0), iB(B, k0); eq && (iA && iB); ++iA, ++iB) {
// 3) Sparsity structure // 3) Sparsity structure
eq = (iA.row() == iB.row()) && (iA.col() == iB.col()); eq = (iA.row() == iB.row()) && (iA.col() == iB.col());
// 4) Equal non-zero elements // 4) Equal non-zero elements
eq = eq && (iA.value() == iB.value()); eq = eq && (iA.value() == iB.value());
} }
}
return eq;
// Note: Investigate implementing this operator as
// return A.cwiseNotEqual(B).count() == 0;
} }
}
return eq;
// Note: Investigate implementing this operator as
// return A.cwiseNotEqual(B).count() == 0;
}
BOOST_AUTO_TEST_CASE(Initialization) BOOST_AUTO_TEST_CASE(Initialization)
{ {
// Setup. // Setup.
Mat z = make_shared<Zero>(3,3); Mat z = Mat(AutoDiffMatrix::ZeroMatrix, 3);
Mat i = make_shared<Identity>(3); Mat i = Mat(AutoDiffMatrix::IdentityMatrix, 3);
Eigen::Array<double, Eigen::Dynamic> d1(3); Eigen::Array<double, Eigen::Dynamic, 1> d1(3);
d1 << 0.2, 1.2, 13.4; d1 << 0.2, 1.2, 13.4;
Mat d = make_shared<Diagonal>(d1); Mat d = Mat(d1.matrix().asDiagonal());
Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic> s1(3,2); Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic> s1(3,2);
s1 << s1 <<
1.0, 0.0, 2.0, 1.0, 0.0, 2.0,
0.0, 1.0, 0.0; 0.0, 1.0, 0.0;
Sp s2(s1); Sp s2(s1.sparseView());
Mat s = make_shared<Sparse>(s2); Mat s = Mat(s2);
} }
BOOST_AUTO_TEST_CASE(EigenConversion) BOOST_AUTO_TEST_CASE(EigenConversion)
{ {
// Setup // Setup.
Mat z = make_shared<Zero>(3,3); Mat z = Mat(AutoDiffMatrix::ZeroMatrix, 3);
Mat i = make_shared<Identity>(3); Mat i = Mat(AutoDiffMatrix::IdentityMatrix, 3);
Eigen::Array<double, Eigen::Dynamic> d1(3); Eigen::Array<double, Eigen::Dynamic, 1> d1(3);
d1 << 0.2, 1.2, 13.4; d1 << 0.2, 1.2, 13.4;
Mat d = make_shared<Diagonal>(d1); Mat d = Mat(d1.matrix().asDiagonal());
Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic> s1(3,2); Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic> s1(3,2);
s1 << s1 <<
1.0, 0.0, 2.0, 1.0, 0.0, 2.0,
0.0, 1.0, 0.0; 0.0, 1.0, 0.0;
Mat s = make_shared<Sparse>(Sp(s1)); Sp s2(s1.sparseView());
Mat s = Mat(s2);
// Convert to Eigen::SparseMatrix // Convert to Eigen::SparseMatrix
Sp x; Sp x;
z->toSparse(x); z.toSparse(x);
BOOST_CHECK_EQUAL(x, Sp(3,3)); Sp z1(3,3);
i->toSparse(x); BOOST_CHECK(x == z1);
Sp i1(Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic>::Identity(3,3)); i.toSparse(x);
BOOST_CHECK_EQUAL(x, i1); Sp i1(Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic>::Identity(3,3).sparseView());
d->toSparse(x); BOOST_CHECK(x == i1);
BOOST_CHECK_EQUAL(x, Sp(d1.matrix().asDiagonal())); d.toSparse(x);
s->toSparse(x); Sp d2 = spdiag(d1);
BOOST_CHECK_EQUAL(x, Sp(s1)); BOOST_CHECK(x == d2);
s.toSparse(x);
BOOST_CHECK(x == s2);
}
BOOST_AUTO_TEST_CASE(AdditionOps)
{
// Setup.
Mat z = Mat(AutoDiffMatrix::ZeroMatrix, 3);
Sp zs(3,3);
Mat i = Mat(AutoDiffMatrix::IdentityMatrix, 3);
Sp is(Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic>::Identity(3,3).sparseView());
Eigen::Array<double, Eigen::Dynamic, 1> d1(3);
d1 << 0.2, 1.2, 13.4;
Mat d = Mat(d1.matrix().asDiagonal());
Sp ds = spdiag(d1);
Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic> s1(3,3);
s1 <<
1.0, 0.0, 2.0,
0.0, 1.0, 0.0,
0.0, 0.0, 2.0;
Sp ss(s1.sparseView());
Mat s = Mat(ss);
// Convert to Eigen::SparseMatrix
Sp x;
z.toSparse(x);
BOOST_CHECK(x == zs);
i.toSparse(x);
BOOST_CHECK(x == is);
d.toSparse(x);
BOOST_CHECK(x == ds);
s.toSparse(x);
BOOST_CHECK(x == ss);
// Adding zero.
auto zpz = z + z;
zpz.toSparse(x);
BOOST_CHECK(x == zs);
auto ipz = i + z;
ipz.toSparse(x);
BOOST_CHECK(x == is);
auto dpz = d + z;
dpz.toSparse(x);
BOOST_CHECK(x == ds);
auto spz = s + z;
spz.toSparse(x);
BOOST_CHECK(x == ss);
} }