28#ifndef EWOMS_BLACK_OIL_BOUNDARY_RATE_VECTOR_HH
29#define EWOMS_BLACK_OIL_BOUNDARY_RATE_VECTOR_HH
31#include <opm/material/common/Valgrind.hpp>
32#include <opm/material/constraintsolvers/NcpFlash.hpp>
46template <
class TypeTag>
63 enum { contiEnergyEqIdx = Indices::contiEnergyEqIdx };
65 enum { enableMICP = Indices::enableMICP };
67 static constexpr bool blackoilConserveSurfaceVolume =
93 template <
class Context,
class Flu
idState>
97 const FluidState& fluidState)
102 const unsigned focusDofIdx = context.focusDofIndex();
110 if (!FluidSystem::phaseIsActive(
phaseIdx)) {
128 using RhsEval = std::conditional_t<std::is_same_v<typename FluidState::Scalar, Evaluation>,
138 for (
unsigned i = 0; i < tmp.size(); ++i) {
139 (*this)[i] += tmp[i];
143 if constexpr (enableEnergy) {
148 density = fluidState.density(
phaseIdx);
171 if constexpr (enableSolvent) {
172 (*this)[Indices::contiSolventEqIdx] =
extQuants.solventVolumeFlux();
173 if (blackoilConserveSurfaceVolume) {
174 (*this)[Indices::contiSolventEqIdx] *=
insideIntQuants.solventInverseFormationVolumeFactor();
177 (*this)[Indices::contiSolventEqIdx] *=
insideIntQuants.solventDensity();
181 if constexpr (enablePolymer) {
182 (*this)[Indices::contiPolymerEqIdx] =
extQuants.volumeFlux(FluidSystem::waterPhaseIdx) *
186 if constexpr (enableMICP) {
187 (*this)[Indices::contiMicrobialEqIdx] =
extQuants.volumeFlux(FluidSystem::waterPhaseIdx) *
189 (*this)[Indices::contiOxygenEqIdx] =
extQuants.volumeFlux(FluidSystem::waterPhaseIdx) *
191 (*this)[Indices::contiUreaEqIdx] =
extQuants.volumeFlux(FluidSystem::waterPhaseIdx) *
198 LocalResidual::adaptMassConservationQuantities_(*
this,
insideIntQuants.pvtRegionIndex());
201 if constexpr (enableEnergy) {
202 EnergyModule::addToEnthalpyRate(*
this,
extQuants.energyFlux() *
207 for (
unsigned i = 0; i < numEq; ++i) {
208 Valgrind::CheckDefined((*
this)[i]);
210 Valgrind::CheckDefined(*
this);
217 template <
class Context,
class Flu
idState>
221 const FluidState& fluidState)
227 std::for_each(this->begin(), this->end(),
228 [](
auto&
val) {
val = std::min(Scalar(0),
val); });
234 template <
class Context,
class Flu
idState>
238 const FluidState& fluidState)
244 std::for_each(this->begin(), this->end(),
245 [](
auto&
val) {
val = std::max(Scalar(0),
val); });
252 { (*this) = Scalar(0); }
261 template <
class Context,
class Flu
idState>
271 if constexpr (enableEnergy) {
275 (*this)[contiEnergyEqIdx] +=
extQuants.energyFlux();
278 for (
unsigned i = 0; i < numEq; ++i) {
279 Valgrind::CheckDefined((*
this)[i]);
281 Valgrind::CheckDefined(*
this);
Contains the classes required to extend the black-oil model by energy.
Implements a boundary vector for the fully implicit black-oil model.
Definition blackoilboundaryratevector.hh:48
BlackOilBoundaryRateVector(Scalar value)
Definition blackoilboundaryratevector.hh:81
void setOutFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an outflow boundary.
Definition blackoilboundaryratevector.hh:235
BlackOilBoundaryRateVector(const BlackOilBoundaryRateVector &value)=default
Copy constructor.
void setFreeFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify a free-flow boundary.
Definition blackoilboundaryratevector.hh:94
BlackOilBoundaryRateVector()=default
Default constructor.
void setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition blackoilboundaryratevector.hh:251
void setInFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition blackoilboundaryratevector.hh:218
void setThermalFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &boundaryFluidState)
an energy flux that corresponds to the thermal conduction from
Definition blackoilboundaryratevector.hh:262
Contains the high level supplements required to extend the black oil model by energy.
Definition blackoilenergymodules.hh:58
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
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