28#ifndef EWOMS_NCP_BOUNDARY_RATE_VECTOR_HH
29#define EWOMS_NCP_BOUNDARY_RATE_VECTOR_HH
31#include <opm/material/common/MathToolbox.hpp>
32#include <opm/material/common/Valgrind.hpp>
48template <
class TypeTag>
58 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
59 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
60 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
61 enum { conti0EqIdx = Indices::conti0EqIdx };
62 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
66 using Toolbox = MathToolbox<Evaluation>;
89 template <
class Context,
class Flu
idState>
93 const FluidState& fluidState)
95 ExtensiveQuantities extQuants;
96 extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState);
97 const auto& insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
98 const unsigned focusDofIdx = context.focusDofIndex();
99 const unsigned interiorDofIdx = context.interiorScvIndex(bfIdx, timeIdx);
104 (*this) = Evaluation(0.0);
105 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
107 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
108 if (focusDofIdx == interiorDofIdx) {
109 density = fluidState.density(phaseIdx);
112 density = getValue(fluidState.density(phaseIdx));
115 else if (focusDofIdx == interiorDofIdx) {
116 density = insideIntQuants.fluidState().density(phaseIdx);
119 density = getValue(insideIntQuants.fluidState().density(phaseIdx));
122 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
124 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
125 if (focusDofIdx == interiorDofIdx) {
126 molarity = fluidState.molarity(phaseIdx, compIdx);
129 molarity = getValue(fluidState.molarity(phaseIdx, compIdx));
132 else if (focusDofIdx == interiorDofIdx) {
133 molarity = insideIntQuants.fluidState().molarity(phaseIdx, compIdx);
136 molarity = getValue(insideIntQuants.fluidState().molarity(phaseIdx, compIdx));
141 (*this)[conti0EqIdx + compIdx] += extQuants.volumeFlux(phaseIdx) * molarity;
145 Evaluation specificEnthalpy;
146 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
147 if (focusDofIdx == interiorDofIdx) {
148 specificEnthalpy = fluidState.enthalpy(phaseIdx);
151 specificEnthalpy = getValue(fluidState.enthalpy(phaseIdx));
154 else if (focusDofIdx == interiorDofIdx) {
155 specificEnthalpy = insideIntQuants.fluidState().enthalpy(phaseIdx);
158 specificEnthalpy = getValue(insideIntQuants.fluidState().enthalpy(phaseIdx));
161 const Evaluation enthalpyRate = density * extQuants.volumeFlux(phaseIdx) * specificEnthalpy;
162 EnergyModule::addToEnthalpyRate(*
this, enthalpyRate);
167 EnergyModule::addToEnthalpyRate(*
this, EnergyModule::thermalConductionRate(extQuants));
170 for (
unsigned i = 0; i < numEq; ++i) {
171 Valgrind::CheckDefined((*
this)[i]);
179 template <
class Context,
class Flu
idState>
183 const FluidState& fluidState)
185 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
188 std::for_each(this->begin(), this->end(),
189 [](
auto& val) { val = Toolbox::min(Scalar{0.0}, val); });
195 template <
class Context,
class Flu
idState>
199 const FluidState& fluidState)
201 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
204 std::for_each(this->begin(), this->end(),
205 [](
auto& val) { val = Toolbox::max(Scalar{0.0}, val); });
212 { (*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 compositional multi-phase NCP model.
Definition: ncpboundaryratevector.hh:50
NcpBoundaryRateVector(const Evaluation &value)
Definition: ncpboundaryratevector.hh:75
void setInFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition: ncpboundaryratevector.hh:180
void setFreeFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify a free-flow boundary.
Definition: ncpboundaryratevector.hh:90
void setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition: ncpboundaryratevector.hh:211
NcpBoundaryRateVector()=default
NcpBoundaryRateVector & operator=(const NcpBoundaryRateVector &value)=default
void setOutFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an outflow boundary.
Definition: ncpboundaryratevector.hh:196
NcpBoundaryRateVector(const NcpBoundaryRateVector &value)=default
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