28#ifndef EWOMS_PVS_BOUNDARY_RATE_VECTOR_HH
29#define EWOMS_PVS_BOUNDARY_RATE_VECTOR_HH
31#include <opm/material/common/MathToolbox.hpp>
32#include <opm/material/common/Valgrind.hpp>
49template <
class TypeTag>
59 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
60 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
61 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
62 enum { conti0EqIdx = Indices::conti0EqIdx };
63 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
66 using Toolbox = MathToolbox<Evaluation>;
89 template <
class Context,
class Flu
idState>
90 void setFreeFlow(
const Context& context,
unsigned bfIdx,
unsigned timeIdx,
const FluidState& fluidState)
92 ExtensiveQuantities extQuants;
93 extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState);
94 const auto& insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
95 const unsigned focusDofIdx = context.focusDofIndex();
96 const unsigned interiorDofIdx = context.interiorScvIndex(bfIdx, timeIdx);
101 (*this) = Evaluation(0.0);
102 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
104 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
105 if (focusDofIdx == interiorDofIdx) {
106 density = fluidState.density(phaseIdx);
109 density = getValue(fluidState.density(phaseIdx));
112 else if (focusDofIdx == interiorDofIdx) {
113 density = insideIntQuants.fluidState().density(phaseIdx);
116 density = getValue(insideIntQuants.fluidState().density(phaseIdx));
119 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
121 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
122 if (focusDofIdx == interiorDofIdx) {
123 molarity = fluidState.molarity(phaseIdx, compIdx);
126 molarity = getValue(fluidState.molarity(phaseIdx, compIdx));
129 else if (focusDofIdx == interiorDofIdx) {
130 molarity = insideIntQuants.fluidState().molarity(phaseIdx, compIdx);
133 molarity = getValue(insideIntQuants.fluidState().molarity(phaseIdx, compIdx));
138 (*this)[conti0EqIdx + compIdx] += extQuants.volumeFlux(phaseIdx) * molarity;
141 if constexpr (enableEnergy) {
142 Evaluation specificEnthalpy;
143 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
144 if (focusDofIdx == interiorDofIdx) {
145 specificEnthalpy = fluidState.enthalpy(phaseIdx);
148 specificEnthalpy = getValue(fluidState.enthalpy(phaseIdx));
151 else if (focusDofIdx == interiorDofIdx) {
152 specificEnthalpy = insideIntQuants.fluidState().enthalpy(phaseIdx);
155 specificEnthalpy = getValue(insideIntQuants.fluidState().enthalpy(phaseIdx));
158 const Evaluation enthalpyRate = density * extQuants.volumeFlux(phaseIdx) * specificEnthalpy;
159 EnergyModule::addToEnthalpyRate(*
this, enthalpyRate);
163 if constexpr (enableEnergy) {
165 EnergyModule::addToEnthalpyRate(*
this, EnergyModule::thermalConductionRate(extQuants));
169 for (
unsigned i = 0; i < numEq; ++i) {
170 Opm::Valgrind::CheckDefined((*
this)[i]);
178 template <
class Context,
class Flu
idState>
182 const FluidState& fluidState)
184 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
187 std::for_each(this->begin(), this->end(),
188 [](
auto& val) { 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 std::for_each(this->begin(), this->end(),
204 [](
auto& val) { val = Toolbox::max(0.0, val); });
211 { (*this) = Evaluation(0.0); }
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:54
Implements a rate vector on the boundary for the fully implicit compositional multi-phase primary var...
Definition: pvsboundaryratevector.hh:51
PvsBoundaryRateVector & operator=(const PvsBoundaryRateVector &value)=default
PvsBoundaryRateVector(const PvsBoundaryRateVector &value)=default
void setOutFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an outflow boundary.
Definition: pvsboundaryratevector.hh:195
void setInFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition: pvsboundaryratevector.hh:179
PvsBoundaryRateVector()=default
PvsBoundaryRateVector(const Evaluation &value)
Definition: pvsboundaryratevector.hh:75
void setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition: pvsboundaryratevector.hh:210
void setFreeFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify a free-flow boundary.
Definition: pvsboundaryratevector.hh:90
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