Distinguish fluid from geology properties.

Most of class ImpesTPFAAD deals (or will deal) with fluid properties
with only a token nod to the geology--mostly in the form of derived
quantities such as transmissibility and pore volumes.  There is no
need to conflate the roles of fluid interface and geology interface
for this purpose.
This commit is contained in:
Bård Skaflestad
2013-05-08 09:39:48 +02:00
parent 1b3e43068f
commit 3ff5347ac9
2 changed files with 61 additions and 15 deletions

View File

@@ -49,17 +49,17 @@ namespace {
namespace Opm {
class LinearSolverInterface;
template <typename Scalar, class BOProps>
template <typename Scalar, class BOFluid>
class PressureDependentFluidData {
public:
typedef AutoDiff::ForwardBlock<Scalar> ADB;
PressureDependentFluidData(const int nc,
const BOProps& props)
const BOFluid& fluid)
: nc_ (nc)
, np_ (props.numPhases())
, np_ (fluid.numPhases())
, cells_(buildAllCells(nc))
, props_(props)
, fluid_(fluid)
, A_ (nc_, np_ * np_)
, dA_ (nc_, np_ * np_)
, mu_ (nc_, np_ )
@@ -76,7 +76,7 @@ namespace Opm {
const std::vector<double>& s = state.saturation();
double* dkrds = 0; // Ignore rel-perm derivatives
props_.relperm(nc_, & s[0], & cells_[0],
fluid_.relperm(nc_, & s[0], & cells_[0],
kr_.data(), dkrds);
}
@@ -86,10 +86,10 @@ namespace Opm {
const std::vector<double>& p = state.pressure();
const std::vector<double>& z = state.surfacevol();
props_.matrix (nc_, & p[0], & z[0], & cells_[0],
fluid_.matrix (nc_, & p[0], & z[0], & cells_[0],
A_ .data(), dA_ .data());
props_.viscosity(nc_, & p[0], & z[0], & cells_[0],
fluid_.viscosity(nc_, & p[0], & z[0], & cells_[0],
mu_.data(), dmu_.data());
}
@@ -132,7 +132,7 @@ namespace Opm {
const int np_;
const std::vector<int> cells_;
const BOProps& props_;
const BOFluid& fluid_;
typedef Eigen::Array<Scalar,
Eigen::Dynamic,
@@ -152,14 +152,15 @@ namespace Opm {
const typename ADB::V one_ ;
};
template <class BOProps>
template <class BOFluid, class GeoProps>
class ImpesTPFAAD {
public:
ImpesTPFAAD(const UnstructuredGrid& grid ,
const BOProps& props)
const BOFluid& fluid,
const GeoProps& geo )
: grid_ (grid)
, pdepfdata_(grid.number_of_cells, props)
, geo_ (geo)
, pdepfdata_(grid.number_of_cells, fluid)
{
}
@@ -168,9 +169,11 @@ namespace Opm {
ImpesTPFAAD(const ImpesTPFAAD& rhs);
ImpesTPFAAD& operator=(const ImpesTPFAAD& rhs);
typedef PressureDependentFluidData<double,BOProps> PDepFData;
typedef PressureDependentFluidData<double, BOFluid> PDepFData;
typedef typename PDepFData::ADB ADB;
const UnstructuredGrid& grid_;
const GeoProps& geo_ ;
PDepFData pdepfdata_;
};
} // namespace Opm

View File

@@ -25,8 +25,45 @@
#include <opm/core/props/BlackoilPropertiesBasic.hpp>
#include <opm/core/pressure/tpfa/trans_tpfa.h>
#include <opm/core/utility/parameters/ParameterGroup.hpp>
#include <algorithm>
namespace {
template <class Geology, class Vector>
class DerivedGeology {
public:
typedef Vector V;
DerivedGeology(const UnstructuredGrid& grid,
const Geology& geo)
: pvol_ (grid.number_of_cells)
, trans_(grid.number_of_faces)
{
// Pore volume
const typename Vector::Index nc = grid.number_of_cells;
std::transform(grid.cell_volumes, grid.cell_volumes + nc,
geo.porosity(), pvol_.data(),
std::multiplies<double>());
// Transmissibility
Vector htrans(grid.cell_facepos[nc]);
UnstructuredGrid* ug = const_cast<UnstructuredGrid*>(& grid);
tpfa_htrans_compute(ug, geo.permeability(), htrans.data());
tpfa_trans_compute (ug, htrans.data() , trans_.data());
}
const Vector& poreVolume() const { return pvol_ ; }
const Vector& transmissibility() const { return trans_; }
private:
Vector pvol_ ;
Vector trans_;
};
}
int
main(int argc, char* argv[])
{
@@ -37,8 +74,14 @@ main(int argc, char* argv[])
const int nc = g->number_of_cells;
const Opm::BlackoilPropertiesBasic props(param, 2, nc);
typedef Opm::ImpesTPFAAD<Opm::BlackoilPropertiesInterface> PSolver;
PSolver ps(*g, props);
typedef AutoDiff::ForwardBlock<double> ADB;
typedef Opm::BlackoilPropertiesInterface Geology;
typedef DerivedGeology<Geology, ADB::V> GeoProps;
typedef Opm::BlackoilPropertiesInterface BOFluid;
typedef Opm::ImpesTPFAAD<BOFluid, GeoProps> PSolver;
GeoProps geo(*g, props);
PSolver ps (*g, props, geo);
return 0;
}