My Project
Loading...
Searching...
No Matches
Opm::NcpModel< TypeTag > Class Template Reference

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

#include <ncpmodel.hh>

Inheritance diagram for Opm::NcpModel< TypeTag >:
Opm::MultiPhaseBaseModel< TypeTag >

Public Member Functions

 NcpModel (Simulator &simulator)
 
void finishInit ()
 Apply the initial conditions to the model.
 
void adaptGrid ()
 
std::string primaryVarName (unsigned pvIdx) const
 Given an primary variable index, return a human readable name.
 
std::string eqName (unsigned eqIdx) const
 Given an equation index, return a human readable name.
 
void updateBegin ()
 Called by the update() method before it tries to apply the newton method.
 
void updatePVWeights (const ElementContext &elemCtx) const
 Update the weights of all primary variables within an element given the complete set of intensive quantities.
 
Scalar primaryVarWeight (unsigned globalDofIdx, unsigned pvIdx) const
 Returns the relative weight of a primary variable for calculating relative errors.
 
Scalar eqWeight (unsigned globalDofIdx, unsigned eqIdx) const
 Returns the relative weight of an equation.
 
Scalar minActivityCoeff (unsigned globalDofIdx, unsigned compIdx) const
 Returns the smallest activity coefficient of a component for the most current solution at a vertex.
 
void registerOutputModules_ ()
 
- Public Member Functions inherited from Opm::MultiPhaseBaseModel< TypeTag >
 MultiPhaseBaseModel (Simulator &simulator)
 
bool phaseIsConsidered (unsigned) const
 Returns true iff a fluid phase is used by the model.
 
void globalPhaseStorage (EqVector &storage, unsigned phaseIdx)
 Compute the total storage inside one phase of all conservation quantities.
 
void registerOutputModules_ ()
 

Static Public Member Functions

static void registerParameters ()
 Register all run-time parameters for the immiscible model.
 
static std::string name ()
 
- Static Public Member Functions inherited from Opm::MultiPhaseBaseModel< TypeTag >
static void registerParameters ()
 Register all run-time parameters for the immiscible model.
 

Public Attributes

Scalar referencePressure_
 
std::vector< ComponentVector > minActivityCoeff_
 

Detailed Description

template<class TypeTag>
class Opm::NcpModel< TypeTag >

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>; };
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:242
Specifies a flux module which uses the Forchheimer relation.
Definition forchheimerfluxmodule.hh:62

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

Member Function Documentation

◆ eqName()

template<class TypeTag >
std::string Opm::NcpModel< TypeTag >::eqName ( unsigned  eqIdx) const
inline

Given an equation index, return a human readable name.

Parameters
eqIdxThe index of the conservation equation of interest.

◆ eqWeight()

template<class TypeTag >
Scalar Opm::NcpModel< TypeTag >::eqWeight ( unsigned  globalDofIdx,
unsigned  eqIdx 
) const
inline

Returns the relative weight of an equation.

Parameters
globalVertexIdxThe global index of the vertex
eqIdxThe index of the equation

◆ finishInit()

template<class TypeTag >
void Opm::NcpModel< TypeTag >::finishInit ( )
inline

Apply the initial conditions to the model.

◆ minActivityCoeff()

template<class TypeTag >
Scalar Opm::NcpModel< TypeTag >::minActivityCoeff ( unsigned  globalDofIdx,
unsigned  compIdx 
) const
inline

Returns the smallest activity coefficient of a component for the most current solution at a vertex.

Parameters
globalDofIdxThe global index of the vertex (i.e. finite volume) of interest.
compIdxThe index of the component of interest.

◆ name()

template<class TypeTag >
static std::string Opm::NcpModel< TypeTag >::name ( )
inlinestatic

◆ primaryVarName()

template<class TypeTag >
std::string Opm::NcpModel< TypeTag >::primaryVarName ( unsigned  pvIdx) const
inline

Given an primary variable index, return a human readable name.

Parameters
pvIdxThe index of the primary variable of interest.

◆ primaryVarWeight()

template<class TypeTag >
Scalar Opm::NcpModel< TypeTag >::primaryVarWeight ( unsigned  globalDofIdx,
unsigned  pvIdx 
) const
inline

Returns the relative weight of a primary variable for calculating relative errors.

Parameters
globalDofIdxThe global index of the degree of freedom
pvIdxThe index of the primary variable

◆ updateBegin()

template<class TypeTag >
void Opm::NcpModel< TypeTag >::updateBegin ( )
inline

Called by the update() method before it tries to apply the newton method.

This is primary a hook which the actual model can overload.

◆ updatePVWeights()

template<class TypeTag >
void Opm::NcpModel< TypeTag >::updatePVWeights ( const ElementContext &  elemCtx) const
inline

Update the weights of all primary variables within an element given the complete set of intensive quantities.


The documentation for this class was generated from the following file: