diff --git a/Makefile.am b/Makefile.am
index 305ed0bd..6cf45cd0 100644
--- a/Makefile.am
+++ b/Makefile.am
@@ -29,6 +29,7 @@ opm/core/fluid/blackoil/SinglePvtLiveOil.cpp \
opm/core/fluid/blackoil/SinglePvtInterface.cpp \
opm/core/fluid/BlackoilPropertiesBasic.cpp \
opm/core/fluid/BlackoilPropertiesFromDeck.cpp \
+opm/core/fluid/IncompPropertiesBasic.cpp \
opm/core/fluid/PvtPropertiesBasic.cpp \
opm/core/fluid/RockBasic.cpp \
opm/core/fluid/RockFromDeck.cpp \
@@ -97,6 +98,8 @@ opm/core/fluid/BlackoilFluid.hpp \
opm/core/fluid/BlackoilPropertiesInterface.hpp \
opm/core/fluid/BlackoilPropertiesBasic.hpp \
opm/core/fluid/BlackoilPropertiesFromDeck.hpp \
+opm/core/fluid/IncompPropertiesInterface.hpp \
+opm/core/fluid/IncompPropertiesBasic.hpp \
opm/core/fluid/PvtPropertiesBasic.hpp \
opm/core/fluid/RockBasic.hpp \
opm/core/fluid/RockFromDeck.hpp \
diff --git a/opm/core/fluid/IncompPropertiesBasic.cpp b/opm/core/fluid/IncompPropertiesBasic.cpp
new file mode 100644
index 00000000..2d0bbbf6
--- /dev/null
+++ b/opm/core/fluid/IncompPropertiesBasic.cpp
@@ -0,0 +1,141 @@
+/*
+ Copyright 2012 SINTEF ICT, Applied Mathematics.
+
+ This file is part of the Open Porous Media project (OPM).
+
+ OPM is free software: you can redistribute it and/or modify
+ it under the terms of the GNU General Public License as published by
+ the Free Software Foundation, either version 3 of the License, or
+ (at your option) any later version.
+
+ OPM is distributed in the hope that it will be useful,
+ but WITHOUT ANY WARRANTY; without even the implied warranty of
+ MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
+ GNU General Public License for more details.
+
+ You should have received a copy of the GNU General Public License
+ along with OPM. If not, see .
+*/
+
+
+
+#include
+#include
+#include
+#include
+
+namespace Opm
+{
+
+ IncompPropertiesBasic::IncompPropertiesBasic(const Dune::parameter::ParameterGroup& param,
+ const int dim,
+ const int num_cells)
+ {
+ double poro = param.getDefault("porosity", 1.0);
+ using namespace Dune::unit;
+ using namespace Dune::prefix;
+ double perm = param.getDefault("permeability", 100)*milli*darcy;
+ rock_.init(dim, num_cells, poro, perm);
+ pvt_.init(param);
+ satprops_.init(param);
+ if (pvt_.numPhases() != satprops_.numPhases()) {
+ THROW("IncompPropertiesBasic::IncompPropertiesBasic() - Inconsistent number of phases in pvt data ("
+ << pvt_.numPhases() << ") and saturation-dependent function data (" << satprops_.numPhases() << ").");
+ }
+ viscosity_.resize(pvt_.numPhases());
+ pvt_.mu(1, 0, 0, &viscosity_[0]);
+ }
+
+ IncompPropertiesBasic::~IncompPropertiesBasic()
+ {
+ }
+
+
+ /// \return D, the number of spatial dimensions.
+ int IncompPropertiesBasic::numDimensions() const
+ {
+ return rock_.numDimensions();
+ }
+
+ /// \return N, the number of cells.
+ int IncompPropertiesBasic::numCells() const
+ {
+ return rock_.numCells();
+ }
+
+ /// \return Array of N porosity values.
+ const double* IncompPropertiesBasic::porosity() const
+ {
+ return rock_.porosity();
+ }
+
+ /// \return Array of ND^2 permeability values.
+ /// The D^2 permeability values for a cell are organized as a matrix,
+ /// which is symmetric (so ordering does not matter).
+ const double* IncompPropertiesBasic::permeability() const
+ {
+ return rock_.permeability();
+ }
+
+
+ // ---- Fluid interface ----
+
+ /// \return P, the number of phases (also the number of components).
+ int IncompPropertiesBasic::numPhases() const
+ {
+ return pvt_.numPhases();
+ }
+
+ /// \return Array of P viscosity values.
+ const double* IncompPropertiesBasic::viscosity() const
+ {
+ return &viscosity_[0];
+ }
+
+ /// \return Array of P density values.
+ const double* IncompPropertiesBasic::density() const
+ {
+ return pvt_.surfaceDensities();
+ }
+
+ /// \param[in] n Number of data points.
+ /// \param[in] s Array of nP saturation values.
+ /// \param[in] cells Array of n cell indices to be associated with the s values.
+ /// \param[out] kr Array of nP relperm values, array must be valid before calling.
+ /// \param[out] dkrds If non-null: array of nP^2 relperm derivative values,
+ /// array must be valid before calling.
+ /// The P^2 derivative matrix is
+ /// m_{ij} = \frac{dkr_i}{ds^j},
+ /// and is output in Fortran order (m_00 m_10 m_20 m01 ...)
+ void IncompPropertiesBasic::relperm(const int n,
+ const double* s,
+ const int* /*cells*/,
+ double* kr,
+ double* dkrds) const
+ {
+ satprops_.relperm(n, s, kr, dkrds);
+ }
+
+
+ /// \param[in] n Number of data points.
+ /// \param[in] s Array of nP saturation values.
+ /// \param[in] cells Array of n cell indices to be associated with the s values.
+ /// \param[out] pc Array of nP capillary pressure values, array must be valid before calling.
+ /// \param[out] dpcds If non-null: array of nP^2 derivative values,
+ /// array must be valid before calling.
+ /// The P^2 derivative matrix is
+ /// m_{ij} = \frac{dpc_i}{ds^j},
+ /// and is output in Fortran order (m_00 m_10 m_20 m01 ...)
+ void IncompPropertiesBasic::capPress(const int n,
+ const double* s,
+ const int* /*cells*/,
+ double* pc,
+ double* dpcds) const
+ {
+ satprops_.relperm(n, s, pc, dpcds);
+ }
+
+
+
+} // namespace Opm
+
diff --git a/opm/core/fluid/IncompPropertiesBasic.hpp b/opm/core/fluid/IncompPropertiesBasic.hpp
new file mode 100644
index 00000000..951dcf18
--- /dev/null
+++ b/opm/core/fluid/IncompPropertiesBasic.hpp
@@ -0,0 +1,129 @@
+/*
+ Copyright 2012 SINTEF ICT, Applied Mathematics.
+
+ This file is part of the Open Porous Media project (OPM).
+
+ OPM is free software: you can redistribute it and/or modify
+ it under the terms of the GNU General Public License as published by
+ the Free Software Foundation, either version 3 of the License, or
+ (at your option) any later version.
+
+ OPM is distributed in the hope that it will be useful,
+ but WITHOUT ANY WARRANTY; without even the implied warranty of
+ MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
+ GNU General Public License for more details.
+
+ You should have received a copy of the GNU General Public License
+ along with OPM. If not, see .
+*/
+
+#ifndef OPM_INCOMPPROPERTIESBASIC_HEADER_INCLUDED
+#define OPM_INCOMPPROPERTIESBASIC_HEADER_INCLUDED
+
+#include
+#include
+#include
+#include
+
+namespace Opm
+{
+
+ /// Concrete class implementing the incompressible property
+ /// interface, reading all necessary input from parameters.
+ ///
+ /// Supports variable number of spatial dimensions, called D.
+ /// Supports variable number of phases, called P.
+ /// In general, when arguments call for n values of some vector or
+ /// matrix property, such as saturation, they shall always be
+ /// ordered cellwise:
+ /// [s^1_0 s^2_0 s^3_0 s^1_1 s^2_2 ... ]
+ /// in which s^i_j denotes saturation of phase i in cell j.
+ class IncompPropertiesBasic : public IncompPropertiesInterface
+ {
+ public:
+ /// Construct from parameters.
+ /// The following parameters are accepted (defaults):
+ /// num_phases (2) Must be 1 or 2.
+ /// relperm_func ("Linear") Must be "Constant", "Linear" or "Quadratic".
+ /// rho1 [rho2, rho3] (1.0e3) Density in kg/m^3
+ /// mu1 [mu2, mu3] (1.0) Viscosity in cP
+ /// porosity (1.0) Porosity
+ /// permeability (100.0) Permeability in mD
+ IncompPropertiesBasic(const Dune::parameter::ParameterGroup& param,
+ const int dim,
+ const int num_cells);
+
+ /// Destructor.
+ virtual ~IncompPropertiesBasic();
+
+ // ---- Rock interface ----
+
+ /// \return D, the number of spatial dimensions.
+ virtual int numDimensions() const;
+
+ /// \return N, the number of cells.
+ virtual int numCells() const;
+
+ /// \return Array of N porosity values.
+ virtual const double* porosity() const;
+
+ /// \return Array of ND^2 permeability values.
+ /// The D^2 permeability values for a cell are organized as a matrix,
+ /// which is symmetric (so ordering does not matter).
+ virtual const double* permeability() const;
+
+
+ // ---- Fluid interface ----
+
+ /// \return P, the number of phases (also the number of components).
+ virtual int numPhases() const;
+
+ /// \return Array of P viscosity values.
+ virtual const double* viscosity() const;
+
+ /// \return Array of P density values.
+ virtual const double* density() const;
+
+ /// \param[in] n Number of data points.
+ /// \param[in] s Array of nP saturation values.
+ /// \param[in] cells Array of n cell indices to be associated with the s values.
+ /// \param[out] kr Array of nP relperm values, array must be valid before calling.
+ /// \param[out] dkrds If non-null: array of nP^2 relperm derivative values,
+ /// array must be valid before calling.
+ /// The P^2 derivative matrix is
+ /// m_{ij} = \frac{dkr_i}{ds^j},
+ /// and is output in Fortran order (m_00 m_10 m_20 m01 ...)
+ virtual void relperm(const int n,
+ const double* s,
+ const int* cells,
+ double* kr,
+ double* dkrds) const;
+
+
+ /// \param[in] n Number of data points.
+ /// \param[in] s Array of nP saturation values.
+ /// \param[in] cells Array of n cell indices to be associated with the s values.
+ /// \param[out] pc Array of nP capillary pressure values, array must be valid before calling.
+ /// \param[out] dpcds If non-null: array of nP^2 derivative values,
+ /// array must be valid before calling.
+ /// The P^2 derivative matrix is
+ /// m_{ij} = \frac{dpc_i}{ds^j},
+ /// and is output in Fortran order (m_00 m_10 m_20 m01 ...)
+ virtual void capPress(const int n,
+ const double* s,
+ const int* cells,
+ double* pc,
+ double* dpcds) const;
+ private:
+ RockBasic rock_;
+ PvtPropertiesBasic pvt_;
+ SaturationPropsBasic satprops_;
+ std::vector viscosity_;
+ };
+
+
+
+} // namespace Opm
+
+
+#endif // OPM_INCOMPPROPERTIESBASIC_HEADER_INCLUDED
diff --git a/opm/core/fluid/IncompPropertiesInterface.hpp b/opm/core/fluid/IncompPropertiesInterface.hpp
new file mode 100644
index 00000000..5ac84fbb
--- /dev/null
+++ b/opm/core/fluid/IncompPropertiesInterface.hpp
@@ -0,0 +1,105 @@
+/*
+ Copyright 2012 SINTEF ICT, Applied Mathematics.
+
+ This file is part of the Open Porous Media project (OPM).
+
+ OPM is free software: you can redistribute it and/or modify
+ it under the terms of the GNU General Public License as published by
+ the Free Software Foundation, either version 3 of the License, or
+ (at your option) any later version.
+
+ OPM is distributed in the hope that it will be useful,
+ but WITHOUT ANY WARRANTY; without even the implied warranty of
+ MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
+ GNU General Public License for more details.
+
+ You should have received a copy of the GNU General Public License
+ along with OPM. If not, see .
+*/
+
+#ifndef OPM_INCOMPPROPERTIESINTERFACE_HEADER_INCLUDED
+#define OPM_INCOMPPROPERTIESINTERFACE_HEADER_INCLUDED
+
+namespace Opm
+{
+
+ /// Abstract base class for incompressible fluid and reservoir properties.
+ ///
+ /// Supports variable number of spatial dimensions, called D.
+ /// Supports variable number of phases, called P.
+ /// In general, when arguments call for n values of some vector or
+ /// matrix property, such as saturation, they shall always be
+ /// ordered cellwise:
+ /// [s^1_0 s^2_0 s^3_0 s^1_1 s^2_2 ... ]
+ /// in which s^i_j denotes saturation of phase i in cell j.
+ class IncompPropertiesInterface
+ {
+ public:
+ virtual ~IncompPropertiesInterface() {}
+
+ // ---- Rock interface ----
+
+ /// \return D, the number of spatial dimensions.
+ virtual int numDimensions() const = 0;
+
+ /// \return N, the number of cells.
+ virtual int numCells() const = 0;
+
+ /// \return Array of N porosity values.
+ virtual const double* porosity() const = 0;
+
+ /// \return Array of ND^2 permeability values.
+ /// The D^2 permeability values for a cell are organized as a matrix,
+ /// which is symmetric (so ordering does not matter).
+ virtual const double* permeability() const = 0;
+
+
+ // ---- Fluid interface ----
+
+ /// \return P, the number of phases (also the number of components).
+ virtual int numPhases() const = 0;
+
+ /// \return Array of P viscosity values.
+ virtual const double* viscosity() const = 0;
+
+ /// \return Array of P density values.
+ virtual const double* density() const = 0;
+
+ /// \param[in] n Number of data points.
+ /// \param[in] s Array of nP saturation values.
+ /// \param[in] cells Array of n cell indices to be associated with the s values.
+ /// \param[out] kr Array of nP relperm values, array must be valid before calling.
+ /// \param[out] dkrds If non-null: array of nP^2 relperm derivative values,
+ /// array must be valid before calling.
+ /// The P^2 derivative matrix is
+ /// m_{ij} = \frac{dkr_i}{ds^j},
+ /// and is output in Fortran order (m_00 m_10 m_20 m01 ...)
+ virtual void relperm(const int n,
+ const double* s,
+ const int* cells,
+ double* kr,
+ double* dkrds) const = 0;
+
+
+ /// \param[in] n Number of data points.
+ /// \param[in] s Array of nP saturation values.
+ /// \param[in] cells Array of n cell indices to be associated with the s values.
+ /// \param[out] pc Array of nP capillary pressure values, array must be valid before calling.
+ /// \param[out] dpcds If non-null: array of nP^2 derivative values,
+ /// array must be valid before calling.
+ /// The P^2 derivative matrix is
+ /// m_{ij} = \frac{dpc_i}{ds^j},
+ /// and is output in Fortran order (m_00 m_10 m_20 m01 ...)
+ virtual void capPress(const int n,
+ const double* s,
+ const int* cells,
+ double* pc,
+ double* dpcds) const = 0;
+ };
+
+
+
+} // namespace Opm
+
+
+#endif // OPM_INCOMPPROPERTIESINTERFACE_HEADER_INCLUDED