28#ifndef EWOMS_FLASH_BOUNDARY_RATE_VECTOR_HH
29#define EWOMS_FLASH_BOUNDARY_RATE_VECTOR_HH
32#include <opm/material/common/Valgrind.hpp>
43template <
class TypeTag>
53 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
54 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
55 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
56 enum { conti0EqIdx = Indices::conti0EqIdx };
57 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
60 using Toolbox = Opm::MathToolbox<Evaluation>;
83 template <
class Context,
class Flu
idState>
87 const FluidState& fluidState)
89 ExtensiveQuantities extQuants;
90 extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState);
91 const auto& insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
92 unsigned focusDofIdx = context.focusDofIndex();
93 unsigned interiorDofIdx = context.interiorScvIndex(bfIdx, timeIdx);
98 (*this) = Evaluation(0.0);
99 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
101 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
102 if (focusDofIdx == interiorDofIdx)
103 density = fluidState.density(phaseIdx);
105 density = Opm::getValue(fluidState.density(phaseIdx));
107 else if (focusDofIdx == interiorDofIdx)
108 density = insideIntQuants.fluidState().density(phaseIdx);
110 density = Opm::getValue(insideIntQuants.fluidState().density(phaseIdx));
112 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
114 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
115 if (focusDofIdx == interiorDofIdx)
116 molarity = fluidState.molarity(phaseIdx, compIdx);
118 molarity = Opm::getValue(fluidState.molarity(phaseIdx, compIdx));
120 else if (focusDofIdx == interiorDofIdx)
121 molarity = insideIntQuants.fluidState().molarity(phaseIdx, compIdx);
123 molarity = Opm::getValue(insideIntQuants.fluidState().molarity(phaseIdx, compIdx));
127 (*this)[conti0EqIdx + compIdx] += extQuants.volumeFlux(phaseIdx)*molarity;
131 Evaluation specificEnthalpy;
132 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
133 if (focusDofIdx == interiorDofIdx)
134 specificEnthalpy = fluidState.enthalpy(phaseIdx);
136 specificEnthalpy = Opm::getValue(fluidState.enthalpy(phaseIdx));
138 else if (focusDofIdx == interiorDofIdx)
139 specificEnthalpy = insideIntQuants.fluidState().enthalpy(phaseIdx);
141 specificEnthalpy = Opm::getValue(insideIntQuants.fluidState().enthalpy(phaseIdx));
143 Evaluation enthalpyRate = density*extQuants.volumeFlux(phaseIdx)*specificEnthalpy;
144 EnergyModule::addToEnthalpyRate(*
this, enthalpyRate);
149 EnergyModule::addToEnthalpyRate(*
this, EnergyModule::thermalConductionRate(extQuants));
152 for (
unsigned i = 0; i < numEq; ++i) {
153 Opm::Valgrind::CheckDefined((*
this)[i]);
161 template <
class Context,
class Flu
idState>
165 const FluidState& fluidState)
167 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
170 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
171 Evaluation& val = this->operator[](eqIdx);
172 val = Toolbox::min(0.0, val);
179 template <
class Context,
class Flu
idState>
183 const FluidState& fluidState)
185 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
188 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
189 Evaluation& val = this->operator[](eqIdx);
190 val = Toolbox::max(0.0, val);
198 { (*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 compositional multi-phase model which is based on...
Definition: flashboundaryratevector.hh:45
FlashBoundaryRateVector(const FlashBoundaryRateVector &value)=default
FlashBoundaryRateVector(const Evaluation &value)
Definition: flashboundaryratevector.hh:70
void setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition: flashboundaryratevector.hh:197
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:84
void setInFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition: flashboundaryratevector.hh:162
FlashBoundaryRateVector()
Definition: flashboundaryratevector.hh:63
void setOutFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an outflow boundary.
Definition: flashboundaryratevector.hh:180
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
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