28#ifndef EWOMS_IMMISCIBLE_BOUNDARY_RATE_VECTOR_HH
29#define EWOMS_IMMISCIBLE_BOUNDARY_RATE_VECTOR_HH
31#include <opm/material/common/Valgrind.hpp>
32#include <opm/material/constraintsolvers/NcpFlash.hpp>
44template <
class TypeTag>
54 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
55 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
56 enum { conti0EqIdx = Indices::conti0EqIdx };
57 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
59 using Toolbox = Opm::MathToolbox<Evaluation>;
97 template <
class Context,
class Flu
idState>
98 void setFreeFlow(
const Context& context,
unsigned bfIdx,
unsigned timeIdx,
const FluidState& fluidState)
100 ExtensiveQuantities extQuants;
101 extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState);
102 const auto& insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
103 unsigned focusDofIdx = context.focusDofIndex();
104 unsigned interiorDofIdx = context.interiorScvIndex(bfIdx, timeIdx);
109 (*this) = Evaluation(0.0);
110 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
111 const auto& pBoundary = fluidState.pressure(phaseIdx);
112 const Evaluation& pInside = insideIntQuants.fluidState().pressure(phaseIdx);
116 if (pBoundary > pInside) {
117 if (focusDofIdx == interiorDofIdx)
118 density = fluidState.density(phaseIdx);
120 density = Opm::getValue(fluidState.density(phaseIdx));
122 else if (focusDofIdx == interiorDofIdx)
123 density = insideIntQuants.fluidState().density(phaseIdx);
125 density = Opm::getValue(insideIntQuants.fluidState().density(phaseIdx));
127 Opm::Valgrind::CheckDefined(density);
128 Opm::Valgrind::CheckDefined(extQuants.volumeFlux(phaseIdx));
130 (*this)[conti0EqIdx + phaseIdx] += extQuants.volumeFlux(phaseIdx)*density;
134 Evaluation specificEnthalpy;
135 if (pBoundary > pInside) {
136 if (focusDofIdx == interiorDofIdx)
137 specificEnthalpy = fluidState.enthalpy(phaseIdx);
139 specificEnthalpy = Opm::getValue(fluidState.enthalpy(phaseIdx));
141 else if (focusDofIdx == interiorDofIdx)
142 specificEnthalpy = insideIntQuants.fluidState().enthalpy(phaseIdx);
144 specificEnthalpy = Opm::getValue(insideIntQuants.fluidState().enthalpy(phaseIdx));
146 Evaluation enthalpyRate = density*extQuants.volumeFlux(phaseIdx)*specificEnthalpy;
147 EnergyModule::addToEnthalpyRate(*
this, enthalpyRate);
152 EnergyModule::addToEnthalpyRate(*
this, EnergyModule::thermalConductionRate(extQuants));
155 for (
unsigned i = 0; i < numEq; ++i)
156 Opm::Valgrind::CheckDefined((*
this)[i]);
157 Opm::Valgrind::CheckDefined(*
this);
170 template <
class Context,
class Flu
idState>
174 const FluidState& fluidState)
176 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
179 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
180 Evaluation& val = this->operator[](eqIdx);
181 val = Toolbox::min(0.0, val);
194 template <
class Context,
class Flu
idState>
198 const FluidState& fluidState)
200 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
203 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
204 Evaluation& val = this->operator[](eqIdx);
205 val = Toolbox::max(0.0, val);
213 { (*this) = Evaluation(0.0); }
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:50
Implements a boundary vector for the fully implicit multi-phase model which assumes immiscibility.
Definition: immiscibleboundaryratevector.hh:46
void setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition: immiscibleboundaryratevector.hh:212
ImmiscibleBoundaryRateVector & operator=(const ImmiscibleBoundaryRateVector &value)=default
void setFreeFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify a free-flow boundary.
Definition: immiscibleboundaryratevector.hh:98
ImmiscibleBoundaryRateVector(const Evaluation &value)
Constructor that assigns all entries to a scalar value.
Definition: immiscibleboundaryratevector.hh:73
void setOutFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an outflow boundary.
Definition: immiscibleboundaryratevector.hh:195
ImmiscibleBoundaryRateVector()
Definition: immiscibleboundaryratevector.hh:63
void setInFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition: immiscibleboundaryratevector.hh:171
ImmiscibleBoundaryRateVector(const ImmiscibleBoundaryRateVector &value)=default
Copy constructor.
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