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>
48template <
class TypeTag>
58 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
59 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
60 enum { conti0EqIdx = Indices::conti0EqIdx };
61 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
63 using Toolbox = MathToolbox<Evaluation>;
99 template <
class Context,
class Flu
idState>
100 void setFreeFlow(
const Context& context,
unsigned bfIdx,
unsigned timeIdx,
const FluidState& fluidState)
102 ExtensiveQuantities extQuants;
103 extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState);
104 const auto& insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
105 const unsigned focusDofIdx = context.focusDofIndex();
106 const unsigned interiorDofIdx = context.interiorScvIndex(bfIdx, timeIdx);
111 (*this) = Evaluation(0.0);
112 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
113 const auto& pBoundary = fluidState.pressure(phaseIdx);
114 const Evaluation& pInside = insideIntQuants.fluidState().pressure(phaseIdx);
118 if (pBoundary > pInside) {
119 if (focusDofIdx == interiorDofIdx) {
120 density = fluidState.density(phaseIdx);
123 density = getValue(fluidState.density(phaseIdx));
126 else if (focusDofIdx == interiorDofIdx) {
127 density = insideIntQuants.fluidState().density(phaseIdx);
130 density = getValue(insideIntQuants.fluidState().density(phaseIdx));
133 Valgrind::CheckDefined(density);
134 Valgrind::CheckDefined(extQuants.volumeFlux(phaseIdx));
136 (*this)[conti0EqIdx + phaseIdx] += extQuants.volumeFlux(phaseIdx) * density;
139 if constexpr (enableEnergy) {
140 Evaluation specificEnthalpy;
141 if (pBoundary > pInside) {
142 if (focusDofIdx == interiorDofIdx) {
143 specificEnthalpy = fluidState.enthalpy(phaseIdx);
146 specificEnthalpy = getValue(fluidState.enthalpy(phaseIdx));
149 else if (focusDofIdx == interiorDofIdx) {
150 specificEnthalpy = insideIntQuants.fluidState().enthalpy(phaseIdx);
153 specificEnthalpy = getValue(insideIntQuants.fluidState().enthalpy(phaseIdx));
156 const Evaluation enthalpyRate = density * extQuants.volumeFlux(phaseIdx) * specificEnthalpy;
157 EnergyModule::addToEnthalpyRate(*
this, enthalpyRate);
162 EnergyModule::addToEnthalpyRate(*
this, EnergyModule::thermalConductionRate(extQuants));
165 for (
unsigned i = 0; i < numEq; ++i) {
166 Valgrind::CheckDefined((*
this)[i]);
168 Valgrind::CheckDefined(*
this);
181 template <
class Context,
class Flu
idState>
185 const FluidState& fluidState)
187 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
190 std::for_each(this->begin(), this->end(),
191 [](
auto& val) { val = Toolbox::min(0.0, val); });
203 template <
class Context,
class Flu
idState>
207 const FluidState& fluidState)
209 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
212 std::for_each(this->begin(), this->end(),
213 [](
auto& val) { val = Toolbox::max(0.0, val); });
220 { (*this) = Evaluation(0.0); }
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:54
Implements a boundary vector for the fully implicit multi-phase model which assumes immiscibility.
Definition: immiscibleboundaryratevector.hh:50
void setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition: immiscibleboundaryratevector.hh:219
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:100
ImmiscibleBoundaryRateVector(const Evaluation &value)
Constructor that assigns all entries to a scalar value.
Definition: immiscibleboundaryratevector.hh:75
void setOutFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an outflow boundary.
Definition: immiscibleboundaryratevector.hh:204
void setInFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition: immiscibleboundaryratevector.hh:182
ImmiscibleBoundaryRateVector()=default
ImmiscibleBoundaryRateVector(const ImmiscibleBoundaryRateVector &value)=default
Copy constructor.
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
Declare the properties used by the infrastructure code of the finite volume discretizations.
Defines the common properties required by the porous medium multi-phase models.
Definition: blackoilboundaryratevector.hh:39
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:233