26 #ifndef EWOMS_RICHARDS_BOUNDARY_RATE_VECTOR_HH
27 #define EWOMS_RICHARDS_BOUNDARY_RATE_VECTOR_HH
29 #include <opm/material/common/Valgrind.hpp>
30 #include <opm/material/constraintsolvers/NcpFlash.hpp>
41 template <
class TypeTag>
44 typedef typename GET_PROP_TYPE(TypeTag, RateVector) ParentType;
45 typedef typename GET_PROP_TYPE(TypeTag, ExtensiveQuantities) ExtensiveQuantities;
46 typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
47 typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
48 typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
49 typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
52 enum { contiEqIdx = Indices::contiEqIdx };
53 enum { liquidPhaseIdx =
GET_PROP_VALUE(TypeTag, LiquidPhaseIndex) };
55 typedef Opm::MathToolbox<Evaluation> Toolbox;
80 template <
class Context,
class Flu
idState>
81 void setFreeFlow(
const Context &context,
int bfIdx,
int timeIdx,
const FluidState &fluidState)
83 typename FluidSystem::ParameterCache paramCache;
84 paramCache.updateAll(fluidState);
86 ExtensiveQuantities extQuants;
87 extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState, paramCache);
88 const auto &insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
93 (*this) = Toolbox::createConstant(0.0);
95 int phaseIdx = liquidPhaseIdx;
97 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx))
98 density = FluidSystem::density(fluidState, paramCache, phaseIdx);
100 density = insideIntQuants.fluidState().density(phaseIdx);
104 (*this)[contiEqIdx] += extQuants.volumeFlux(phaseIdx) * density;
107 for (
int i = 0; i < numEq; ++i) {
108 Valgrind::CheckDefined((*
this)[i]);
110 Valgrind::CheckDefined(*
this);
117 template <
class Context,
class Flu
idState>
118 void setInFlow(
const Context &context,
int bfIdx,
int timeIdx,
119 const FluidState &fluidState)
121 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
124 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
125 Evaluation& val = this->operator[](eqIdx);
126 val = Toolbox::min(0.0, val);
133 template <
class Context,
class Flu
idState>
134 void setOutFlow(
const Context &context,
int bfIdx,
int timeIdx,
135 const FluidState &fluidState)
137 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
140 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
141 Evaluation& val = this->operator[](eqIdx);
142 val = Toolbox::max(0.0, val);
150 { (*this) = Toolbox::createConstant(0.0); }
RichardsBoundaryRateVector(const RichardsBoundaryRateVector &value)
Definition: richardsboundaryratevector.hh:73
void setInFlow(const Context &context, int bfIdx, int timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition: richardsboundaryratevector.hh:118
void setFreeFlow(const Context &context, int bfIdx, int timeIdx, const FluidState &fluidState)
Specify a free-flow boundary.
Definition: richardsboundaryratevector.hh:81
Implements a boundary vector for the fully implicit Richards model.
Definition: richardsboundaryratevector.hh:42
#define GET_PROP_VALUE(TypeTag, PropTagName)
Access the value attribute of a property for a type tag.
Definition: propertysystem.hh:468
#define GET_PROP_TYPE(TypeTag, PropTagName)
Access the type attribute of a property for a type tag.
Definition: propertysystem.hh:485
Intensive quantities required by the Richards model.
Definition: baseauxiliarymodule.hh:35
RichardsBoundaryRateVector()
Definition: richardsboundaryratevector.hh:58
RichardsBoundaryRateVector(const Evaluation &value)
Definition: richardsboundaryratevector.hh:65
void setOutFlow(const Context &context, int bfIdx, int timeIdx, const FluidState &fluidState)
Specify an outflow boundary.
Definition: richardsboundaryratevector.hh:134
void setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition: richardsboundaryratevector.hh:149