28#ifndef EWOMS_FLASH_BOUNDARY_RATE_VECTOR_HH
29#define EWOMS_FLASH_BOUNDARY_RATE_VECTOR_HH
31#include <opm/material/common/MathToolbox.hpp>
32#include <opm/material/common/Valgrind.hpp>
47template <
class TypeTag>
57 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
58 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
59 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
60 enum { conti0EqIdx = Indices::conti0EqIdx };
61 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
64 using Toolbox = MathToolbox<Evaluation>;
87 template <
class Context,
class Flu
idState>
91 const FluidState& fluidState)
93 ExtensiveQuantities extQuants;
94 extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState);
95 const auto& insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
96 const unsigned focusDofIdx = context.focusDofIndex();
97 const unsigned interiorDofIdx = context.interiorScvIndex(bfIdx, timeIdx);
102 (*this) = Evaluation(0.0);
103 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
105 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
106 if (focusDofIdx == interiorDofIdx) {
107 density = fluidState.density(phaseIdx);
110 density = getValue(fluidState.density(phaseIdx));
113 else if (focusDofIdx == interiorDofIdx) {
114 density = insideIntQuants.fluidState().density(phaseIdx);
117 density = getValue(insideIntQuants.fluidState().density(phaseIdx));
120 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
122 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
123 if (focusDofIdx == interiorDofIdx) {
124 molarity = fluidState.molarity(phaseIdx, compIdx);
127 molarity = getValue(fluidState.molarity(phaseIdx, compIdx));
130 else if (focusDofIdx == interiorDofIdx) {
131 molarity = insideIntQuants.fluidState().molarity(phaseIdx, compIdx);
134 molarity = getValue(insideIntQuants.fluidState().molarity(phaseIdx, compIdx));
139 (*this)[conti0EqIdx + compIdx] += extQuants.volumeFlux(phaseIdx)*molarity;
143 Evaluation specificEnthalpy;
144 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
145 if (focusDofIdx == interiorDofIdx) {
146 specificEnthalpy = fluidState.enthalpy(phaseIdx);
149 specificEnthalpy = getValue(fluidState.enthalpy(phaseIdx));
152 else if (focusDofIdx == interiorDofIdx) {
153 specificEnthalpy = insideIntQuants.fluidState().enthalpy(phaseIdx);
156 specificEnthalpy = getValue(insideIntQuants.fluidState().enthalpy(phaseIdx));
159 Evaluation enthalpyRate = density*extQuants.volumeFlux(phaseIdx)*specificEnthalpy;
160 EnergyModule::addToEnthalpyRate(*
this, enthalpyRate);
165 EnergyModule::addToEnthalpyRate(*
this, EnergyModule::thermalConductionRate(extQuants));
168 for (
unsigned i = 0; i < numEq; ++i) {
169 Valgrind::CheckDefined((*
this)[i]);
177 template <
class Context,
class Flu
idState>
181 const FluidState& fluidState)
183 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
186 std::for_each(this->begin(), this->end(),
187 [](
auto& val) { val = Toolbox::min(0.0, val); });
193 template <
class Context,
class Flu
idState>
197 const FluidState& fluidState)
199 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
202 std::for_each(this->begin(), this->end(),
203 [](
auto& val) { val = Toolbox::max(0.0, val); });
210 { (*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 model which is based on...
Definition: flashboundaryratevector.hh:49
FlashBoundaryRateVector(const FlashBoundaryRateVector &value)=default
FlashBoundaryRateVector(const Evaluation &value)
Definition: flashboundaryratevector.hh:73
void setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition: flashboundaryratevector.hh:209
FlashBoundaryRateVector & operator=(const FlashBoundaryRateVector &value)=default
void setFreeFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify a free-flow boundary.
Definition: flashboundaryratevector.hh:88
FlashBoundaryRateVector()=default
void setInFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition: flashboundaryratevector.hh:178
void setOutFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an outflow boundary.
Definition: flashboundaryratevector.hh:194
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
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