bugfix, Mismatched free() / delete / delete [] in DuneMatrix due to use of new instead

of the matrix internal allocators.

This fix also avoid the copying of the BCRSMatrix by providing a contructor that creates
the DuneMatrix for a given Eigen SparseMatrix.
This commit is contained in:
Robert Kloefkorn 2014-10-09 14:06:02 +02:00
parent ed75a02ac0
commit 9f58ad5476
2 changed files with 31 additions and 24 deletions

View File

@ -1,5 +1,6 @@
/* /*
Copyright 2014 SINTEF ICT, Applied Mathematics. Copyright 2014 SINTEF ICT, Applied Mathematics.
Copyright 2014 IRIS AS
This file is part of the Open Porous Media project (OPM). This file is part of the Open Porous Media project (OPM).
@ -26,6 +27,9 @@
#include <opm/core/utility/platform_dependent/disable_warnings.h> #include <opm/core/utility/platform_dependent/disable_warnings.h>
#include <Eigen/Eigen>
#include <Eigen/Sparse>
#include <dune/common/fmatrix.hh> #include <dune/common/fmatrix.hh>
#include <dune/common/version.hh> #include <dune/common/version.hh>
@ -43,9 +47,26 @@ namespace Opm
{ {
public: public:
DuneMatrix(const int rows, const int cols, const int* ia, const int* ja, const double* sa) DuneMatrix(const int rows, const int cols, const int* ia, const int* ja, const double* sa)
// : build_mode(unknown), ready(built), n(rows), m(cols), nnz(ia[rows]), {
// allocationSize(nnz), r(0), a(0), // create BCRSMatrix from given CSR storage
// avg(0), overflowsize(-1.0) init( rows, cols, ia, ja, sa );
}
/// \brief create an ISTL BCRSMatrix from a Eigen::SparseMatrix
DuneMatrix( const Eigen::SparseMatrix<double, Eigen::RowMajor>& matrix )
{
// Create ISTL matrix.
const int rows = matrix.rows();
const int cols = matrix.cols();
const int* ia = matrix.outerIndexPtr();
const int* ja = matrix.innerIndexPtr();
const double* sa = matrix.valuePtr();
// create BCRSMatrix from Eigen matrix
init( rows, cols, ia, ja, sa );
}
protected:
void init(const int rows, const int cols, const int* ia, const int* ja, const double* sa)
{ {
typedef Dune::BCRSMatrix< Dune::FieldMatrix<double, 1, 1> > Super; typedef Dune::BCRSMatrix< Dune::FieldMatrix<double, 1, 1> > Super;
typedef Super::block_type block_type; typedef Super::block_type block_type;
@ -61,12 +82,14 @@ namespace Opm
this->overflowsize = -1.0; this->overflowsize = -1.0;
#endif #endif
this->a = new block_type[this->nnz]; // make sure to use the allocators of this matrix
// because the same allocators are used to deallocate the data
this->a = this->allocator_.allocate(this->nnz);
static_assert(sizeof(block_type) == sizeof(double), "This constructor requires a block type that is the same as a double."); static_assert(sizeof(block_type) == sizeof(double), "This constructor requires a block type that is the same as a double.");
std::copy(sa, sa + this->nnz, reinterpret_cast<double*>(this->a)); std::copy(sa, sa + this->nnz, reinterpret_cast<double*>(this->a));
this->j.reset(new Super::size_type[this->nnz]); this->j.reset(this->sizeAllocator_.allocate(this->nnz));
std::copy(ja, ja +this-> nnz, this->j.get()); std::copy(ja, ja +this-> nnz, this->j.get());
this->r = new Super::row_type[rows]; this->r = rowAllocator_.allocate(rows);
for (int row = 0; row < rows; ++row) { for (int row = 0; row < rows; ++row) {
this->r[row].set(ia[row+1] - ia[row], this->a + ia[row], this->j.get() + ia[row]); this->r[row].set(ia[row+1] - ia[row], this->a + ia[row], this->j.get() + ia[row]);
} }

View File

@ -168,10 +168,10 @@ namespace Opm
SolutionVector dx(SolutionVector::Zero(b.size())); SolutionVector dx(SolutionVector::Zero(b.size()));
// Create ISTL matrix. // Create ISTL matrix.
Mat istlA = makeIstlMatrix(A); DuneMatrix istlA( A );
// Create ISTL matrix for elliptic part. // Create ISTL matrix for elliptic part.
Mat istlAe = makeIstlMatrix(A.topLeftCorner(nc, nc)); DuneMatrix istlAe( A.topLeftCorner(nc, nc) );
// Construct operator, scalar product and vectors needed. // Construct operator, scalar product and vectors needed.
typedef Dune::MatrixAdapter<Mat,Vector,Vector> Operator; typedef Dune::MatrixAdapter<Mat,Vector,Vector> Operator;
@ -457,22 +457,6 @@ namespace Opm
b = L * total_residual.value().matrix(); b = L * total_residual.value().matrix();
} }
Mat makeIstlMatrix(const Eigen::SparseMatrix<double, Eigen::RowMajor>& matrix)
{
// Create ISTL matrix.
const int size = matrix.rows();
const int* ia = matrix.outerIndexPtr();
const int* ja = matrix.innerIndexPtr();
const double* sa = matrix.valuePtr();
return Opm::DuneMatrix(size, size, ia, ja, sa);
}
} // anonymous namespace } // anonymous namespace