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>
58 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
59 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
60 enum { enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>() };
61 enum { enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>() };
62 enum { enableEnergy = (getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::FullyImplicitThermal) };
63 enum { contiEnergyEqIdx = Indices::contiEnergyEqIdx };
64 enum { enableFoam = getPropValue<TypeTag, Properties::EnableFoam>() };
65 enum { enableMICP = Indices::enableMICP };
67 static constexpr bool blackoilConserveSurfaceVolume =
68 getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>();
70 static constexpr EnergyModules energyModuleType = getPropValue<TypeTag, Properties::EnergyModuleType>();
94 template <
class Context,
class Flu
idState>
98 const FluidState& fluidState)
100 ExtensiveQuantities extQuants;
101 extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState);
102 const auto& insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
103 const unsigned focusDofIdx = context.focusDofIndex();
104 const unsigned interiorDofIdx = context.interiorScvIndex(bfIdx, timeIdx);
110 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
111 if (!FluidSystem::phaseIsActive(phaseIdx)) {
114 const auto& pBoundary = fluidState.pressure(phaseIdx);
115 const Evaluation& pInside = insideIntQuants.fluidState().pressure(phaseIdx);
120 if (pBoundary < pInside) {
122 LocalResidual::template evalPhaseFluxes_<Evaluation>(tmp,
124 insideIntQuants.pvtRegionIndex(),
126 insideIntQuants.fluidState());
128 else if (pBoundary > pInside) {
129 using RhsEval = std::conditional_t<std::is_same_v<typename FluidState::Scalar, Evaluation>,
132 LocalResidual::template evalPhaseFluxes_<RhsEval>(tmp,
134 insideIntQuants.pvtRegionIndex(),
139 for (
unsigned i = 0; i < tmp.size(); ++i) {
140 (*this)[i] += tmp[i];
144 if constexpr (enableEnergy) {
146 Evaluation specificEnthalpy;
147 if (pBoundary > pInside) {
148 if (focusDofIdx == interiorDofIdx) {
149 density = fluidState.density(phaseIdx);
150 specificEnthalpy = fluidState.enthalpy(phaseIdx);
153 density = getValue(fluidState.density(phaseIdx));
154 specificEnthalpy = getValue(fluidState.enthalpy(phaseIdx));
157 else if (focusDofIdx == interiorDofIdx) {
158 density = insideIntQuants.fluidState().density(phaseIdx);
159 specificEnthalpy = insideIntQuants.fluidState().enthalpy(phaseIdx);
162 density = getValue(insideIntQuants.fluidState().density(phaseIdx));
163 specificEnthalpy = getValue(insideIntQuants.fluidState().enthalpy(phaseIdx));
166 const Evaluation enthalpyRate = density * extQuants.volumeFlux(phaseIdx) * specificEnthalpy;
168 getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>());
172 if constexpr (enableSolvent) {
173 (*this)[Indices::contiSolventEqIdx] = extQuants.solventVolumeFlux();
174 if (blackoilConserveSurfaceVolume) {
175 (*this)[Indices::contiSolventEqIdx] *= insideIntQuants.solventInverseFormationVolumeFactor();
178 (*this)[Indices::contiSolventEqIdx] *= insideIntQuants.solventDensity();
182 if constexpr (enablePolymer) {
183 (*this)[Indices::contiPolymerEqIdx] = extQuants.volumeFlux(FluidSystem::waterPhaseIdx) *
184 insideIntQuants.polymerConcentration();
187 if constexpr (enableMICP) {
188 (*this)[Indices::contiMicrobialEqIdx] = extQuants.volumeFlux(FluidSystem::waterPhaseIdx) *
189 insideIntQuants.microbialConcentration();
190 (*this)[Indices::contiOxygenEqIdx] = extQuants.volumeFlux(FluidSystem::waterPhaseIdx) *
191 insideIntQuants.oxygenConcentration();
192 (*this)[Indices::contiUreaEqIdx] = extQuants.volumeFlux(FluidSystem::waterPhaseIdx) *
193 insideIntQuants.ureaConcentration();
195 (*this)[Indices::contiUreaEqIdx] *= getPropValue<TypeTag, Properties::BlackOilUreaScalingFactor>();
199 LocalResidual::adaptMassConservationQuantities_(*
this, insideIntQuants.pvtRegionIndex());
202 if constexpr (enableEnergy) {
204 getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>());
208 for (
unsigned i = 0; i < numEq; ++i) {
209 Valgrind::CheckDefined((*
this)[i]);
211 Valgrind::CheckDefined(*
this);
218 template <
class Context,
class Flu
idState>
222 const FluidState& fluidState)
224 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
228 std::for_each(this->begin(), this->end(),
229 [](
auto& val) { val = std::min(Scalar(0), val); });
235 template <
class Context,
class Flu
idState>
239 const FluidState& fluidState)
241 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
245 std::for_each(this->begin(), this->end(),
246 [](
auto& val) { val = std::max(Scalar(0), val); });
253 { (*this) = Scalar(0); }
262 template <
class Context,
class Flu
idState>
264 [[maybe_unused]]
unsigned bfIdx,
265 [[maybe_unused]]
unsigned timeIdx,
266 [[maybe_unused]]
const FluidState& boundaryFluidState)
272 if constexpr (enableEnergy) {
273 ExtensiveQuantities extQuants;
274 extQuants.updateBoundary(context, bfIdx, timeIdx, boundaryFluidState);
276 (*this)[contiEnergyEqIdx] += extQuants.energyFlux();
279 for (
unsigned i = 0; i < numEq; ++i) {
280 Valgrind::CheckDefined((*
this)[i]);
282 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:82
void setOutFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an outflow boundary.
Definition: blackoilboundaryratevector.hh:236
BlackOilBoundaryRateVector & operator=(const BlackOilBoundaryRateVector &value)=default
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:95
BlackOilBoundaryRateVector()=default
Default constructor.
void setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition: blackoilboundaryratevector.hh:252
void setInFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition: blackoilboundaryratevector.hh:219
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:263
Contains the high level supplements required to extend the black oil model by energy.
Definition: blackoilenergymodules.hh:60
static void addToEnthalpyRate(RateVector &flux, const Evaluation &hRate)
Definition: blackoilenergymodules.hh:253
Definition: blackoilbioeffectsmodules.hh:43
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