ncpmodel.hh File Reference
#include <opm/material/densead/Math.hpp>
#include "ncpproperties.hh"
#include "ncplocalresidual.hh"
#include "ncpextensivequantities.hh"
#include "ncpprimaryvariables.hh"
#include "ncpboundaryratevector.hh"
#include "ncpratevector.hh"
#include "ncpintensivequantities.hh"
#include "ncpnewtonmethod.hh"
#include "ncpindices.hh"
#include <opm/common/Exceptions.hpp>
#include <opm/models/common/multiphasebasemodel.hh>
#include <opm/models/common/energymodule.hh>
#include <opm/models/common/diffusionmodule.hh>
#include <opm/models/io/vtkcompositionmodule.hh>
#include <opm/models/io/vtkenergymodule.hh>
#include <opm/models/io/vtkdiffusionmodule.hh>
#include <opm/material/common/Valgrind.hpp>
#include <dune/common/fvector.hh>
#include <sstream>
#include <string>
#include <vector>
#include <array>
Include dependency graph for ncpmodel.hh:

Go to the source code of this file.

Classes

struct  Opm::Properties::TTag::NcpModel
 Define the type tag for the compositional NCP model. More...
 
struct  Opm::Properties::LocalResidual< TypeTag, TTag::NcpModel >
 Use the Ncp local jacobian operator for the compositional NCP model. More...
 
struct  Opm::Properties::NewtonMethod< TypeTag, TTag::NcpModel >
 Use the Ncp specific newton method for the compositional NCP model. More...
 
struct  Opm::Properties::Model< TypeTag, TTag::NcpModel >
 the Model property More...
 
struct  Opm::Properties::BaseProblem< TypeTag, TTag::NcpModel >
 The type of the base base class for actual problems. More...
 
struct  Opm::Properties::EnableEnergy< TypeTag, TTag::NcpModel >
 Disable the energy equation by default. More...
 
struct  Opm::Properties::EnableDiffusion< TypeTag, TTag::NcpModel >
 disable diffusion by default More...
 
struct  Opm::Properties::RateVector< TypeTag, TTag::NcpModel >
 the RateVector property More...
 
struct  Opm::Properties::BoundaryRateVector< TypeTag, TTag::NcpModel >
 the BoundaryRateVector property More...
 
struct  Opm::Properties::PrimaryVariables< TypeTag, TTag::NcpModel >
 the PrimaryVariables property More...
 
struct  Opm::Properties::IntensiveQuantities< TypeTag, TTag::NcpModel >
 the IntensiveQuantities property More...
 
struct  Opm::Properties::ExtensiveQuantities< TypeTag, TTag::NcpModel >
 the ExtensiveQuantities property More...
 
struct  Opm::Properties::Indices< TypeTag, TTag::NcpModel >
 The indices required by the compositional NCP model. More...
 
struct  Opm::Properties::NcpPressureBaseWeight< TypeTag, TTag::NcpModel >
 The unmodified weight for the pressure primary variable. More...
 
struct  Opm::Properties::NcpSaturationsBaseWeight< TypeTag, TTag::NcpModel >
 The weight for the saturation primary variables. More...
 
struct  Opm::Properties::NcpFugacitiesBaseWeight< TypeTag, TTag::NcpModel >
 The unmodified weight for the fugacity primary variables. More...
 
class  Opm::NcpModel< TypeTag >
 A compositional multi-phase model based on non-linear complementarity functions. More...
 

Namespaces

namespace  Opm
 
namespace  Opm::Properties
 
namespace  Opm::Properties::TTag
 The generic type tag for problems using the immiscible multi-phase model.
 

Detailed Description

A compositional multi-phase model based on non-linear complementarity functions.

This model implements a $M$-phase flow of a fluid mixture composed of $N$ chemical species. The phases are denoted by lower index $\alpha \in \{ 1, \dots, M \}$. All fluid phases are mixtures of $N \geq M - 1$ chemical species which are denoted by the upper index $\kappa \in \{ 1, \dots, N \} $.

By default, the standard multi-phase Darcy approach is used to determine the velocity, i.e.

\[
     \mathbf{v}_\alpha = - \frac{k_{r\alpha}}{\mu_\alpha} \mathbf{K}
     \left(\mathbf{grad}\, p_\alpha - \varrho_{\alpha} \mathbf{g} \right) \;,
   \]

although the actual approach which is used can be specified via the FluxModule property. For example, the velocity model can by changed to the Forchheimer approach by

template<class TypeTag>
struct FluxModule<TypeTag, TTag::MyProblemTypeTag> { using type = ForchheimerFluxModule<TypeTag>; };

The core of the model is the conservation mass of each component by means of the equation

\[
   \sum_\alpha \frac{\partial\;\phi c_\alpha^\kappa S_\alpha }{\partial t}
   - \sum_\alpha \mathrm{div} \left\{ c_\alpha^\kappa \mathbf{v}_\alpha \right\}
   - q^\kappa = 0 \;.
   \]

For the missing $M$ model assumptions, the model uses non-linear complementarity functions. These are based on the observation that if a fluid phase is not present, the sum of the mole fractions of this fluid phase is smaller than $1$, i.e.

\[ \forall \alpha: S_\alpha = 0 \implies \sum_\kappa
   x_\alpha^\kappa \leq 1 \]

Also, if a fluid phase may be present at a given spatial location its saturation must be non-negative:

\[ \forall \alpha: \sum_\kappa x_\alpha^\kappa = 1 \implies S_\alpha \geq 0
 *\]

Since at any given spatial location, a phase is always either present or not present, one of the strict equalities on the right hand side is always true, i.e.

\[
   \forall \alpha: S_\alpha \left( \sum_\kappa x_\alpha^\kappa - 1 \right) = 0
   \]

always holds.

These three equations constitute a non-linear complementarity problem, which can be solved using so-called non-linear complementarity functions $\Phi(a, b)$. Such functions have the property

\[\Phi(a,b) = 0 \iff a \geq0 \land b \geq0  \land a \cdot b = 0 \]

Several non-linear complementarity functions have been suggested, e.g. the Fischer-Burmeister function

\[ \Phi(a,b) = a + b - \sqrt{a^2 + b^2} \;. \]

This model uses

\[ \Phi(a,b) = \min \{a,  b \}\;, \]

because of its piecewise linearity.

The model assumes local thermodynamic equilibrium and uses the following primary variables:

  • The pressure of the first phase $p_1$
  • The component fugacities $f^1, \dots, f^{N}$
  • The saturations of the first $M-1$ phases $S_1, \dots, S_{M-1}$
  • Temperature $T$ if the energy equation is enabled