28#ifndef EWOMS_RICHARDS_BOUNDARY_RATE_VECTOR_HH
29#define EWOMS_RICHARDS_BOUNDARY_RATE_VECTOR_HH
31#include <opm/material/common/Valgrind.hpp>
32#include <opm/material/constraintsolvers/NcpFlash.hpp>
43template <
class TypeTag>
53 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
54 enum { contiEqIdx = Indices::contiEqIdx };
55 enum { liquidPhaseIdx = getPropValue<TypeTag, Properties::LiquidPhaseIndex>() };
57 using Toolbox = Opm::MathToolbox<Evaluation>;
81 template <
class Context,
class Flu
idState>
82 void setFreeFlow(
const Context& context,
unsigned bfIdx,
unsigned timeIdx,
const FluidState& fluidState)
84 ExtensiveQuantities extQuants;
85 extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState);
86 const auto& insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
87 unsigned focusDofIdx = context.focusDofIndex();
88 unsigned interiorDofIdx = context.interiorScvIndex(bfIdx, timeIdx);
93 (*this) = Evaluation(0.0);
95 unsigned phaseIdx = liquidPhaseIdx;
97 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
98 if (focusDofIdx == interiorDofIdx)
99 density = fluidState.density(phaseIdx);
101 density = Opm::getValue(fluidState.density(phaseIdx));
103 else if (focusDofIdx == interiorDofIdx)
104 density = insideIntQuants.fluidState().density(phaseIdx);
106 density = Opm::getValue(insideIntQuants.fluidState().density(phaseIdx));
110 (*this)[contiEqIdx] += extQuants.volumeFlux(phaseIdx) * density;
113 for (
unsigned i = 0; i < numEq; ++i) {
114 Opm::Valgrind::CheckDefined((*
this)[i]);
116 Opm::Valgrind::CheckDefined(*
this);
123 template <
class Context,
class Flu
idState>
127 const FluidState& fluidState)
129 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
132 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
133 Evaluation& val = this->operator[](eqIdx);
134 val = Toolbox::min(0.0, val);
141 template <
class Context,
class Flu
idState>
145 const FluidState& fluidState)
147 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
150 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
151 Evaluation& val = this->operator[](eqIdx);
152 val = Toolbox::max(0.0, val);
160 { (*this) = Evaluation(0.0); }
Implements a boundary vector for the fully implicit Richards model.
Definition: richardsboundaryratevector.hh:45
RichardsBoundaryRateVector(const RichardsBoundaryRateVector &value)=default
void setFreeFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify a free-flow boundary.
Definition: richardsboundaryratevector.hh:82
void setOutFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an outflow boundary.
Definition: richardsboundaryratevector.hh:142
void setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition: richardsboundaryratevector.hh:159
void setInFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition: richardsboundaryratevector.hh:124
RichardsBoundaryRateVector()
Definition: richardsboundaryratevector.hh:60
RichardsBoundaryRateVector & operator=(const RichardsBoundaryRateVector &value)=default
RichardsBoundaryRateVector(const Evaluation &value)
Definition: richardsboundaryratevector.hh:67
Definition: blackoilboundaryratevector.hh:37
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:235