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 { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
61 enum { enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>() };
62 enum { enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>() };
63 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
64 enum { conti0EqIdx = Indices::conti0EqIdx };
65 enum { contiEnergyEqIdx = Indices::contiEnergyEqIdx };
66 enum { enableFoam = getPropValue<TypeTag, Properties::EnableFoam>() };
67 enum { enableMICP = getPropValue<TypeTag, Properties::EnableMICP>() };
69 static constexpr bool blackoilConserveSurfaceVolume =
70 getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>();
95 template <
class Context,
class Flu
idState>
99 const FluidState& fluidState)
101 ExtensiveQuantities extQuants;
102 extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState);
103 const auto& insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
104 const unsigned focusDofIdx = context.focusDofIndex();
105 const unsigned interiorDofIdx = context.interiorScvIndex(bfIdx, timeIdx);
111 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
112 if (!FluidSystem::phaseIsActive(phaseIdx)) {
115 const auto& pBoundary = fluidState.pressure(phaseIdx);
116 const Evaluation& pInside = insideIntQuants.fluidState().pressure(phaseIdx);
121 if (pBoundary < pInside) {
123 LocalResidual::template evalPhaseFluxes_<Evaluation>(tmp,
125 insideIntQuants.pvtRegionIndex(),
127 insideIntQuants.fluidState());
129 else if (pBoundary > pInside) {
130 using RhsEval = std::conditional_t<std::is_same_v<typename FluidState::Scalar, Evaluation>,
133 LocalResidual::template evalPhaseFluxes_<RhsEval>(tmp,
135 insideIntQuants.pvtRegionIndex(),
140 for (
unsigned i = 0; i < tmp.size(); ++i) {
141 (*this)[i] += tmp[i];
145 if constexpr (enableEnergy) {
147 Evaluation specificEnthalpy;
148 if (pBoundary > pInside) {
149 if (focusDofIdx == interiorDofIdx) {
150 density = fluidState.density(phaseIdx);
151 specificEnthalpy = fluidState.enthalpy(phaseIdx);
154 density = getValue(fluidState.density(phaseIdx));
155 specificEnthalpy = getValue(fluidState.enthalpy(phaseIdx));
158 else if (focusDofIdx == interiorDofIdx) {
159 density = insideIntQuants.fluidState().density(phaseIdx);
160 specificEnthalpy = insideIntQuants.fluidState().enthalpy(phaseIdx);
163 density = getValue(insideIntQuants.fluidState().density(phaseIdx));
164 specificEnthalpy = getValue(insideIntQuants.fluidState().enthalpy(phaseIdx));
167 const Evaluation enthalpyRate = density * extQuants.volumeFlux(phaseIdx) * specificEnthalpy;
169 getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>());
173 if constexpr (enableSolvent) {
174 (*this)[Indices::contiSolventEqIdx] = extQuants.solventVolumeFlux();
175 if (blackoilConserveSurfaceVolume) {
176 (*this)[Indices::contiSolventEqIdx] *= insideIntQuants.solventInverseFormationVolumeFactor();
179 (*this)[Indices::contiSolventEqIdx] *= insideIntQuants.solventDensity();
183 if constexpr (enablePolymer) {
184 (*this)[Indices::contiPolymerEqIdx] = extQuants.volumeFlux(FluidSystem::waterPhaseIdx) *
185 insideIntQuants.polymerConcentration();
188 if constexpr (enableMICP) {
189 (*this)[Indices::contiMicrobialEqIdx] = extQuants.volumeFlux(FluidSystem::waterPhaseIdx) *
190 insideIntQuants.microbialConcentration();
191 (*this)[Indices::contiOxygenEqIdx] = extQuants.volumeFlux(FluidSystem::waterPhaseIdx) *
192 insideIntQuants.oxygenConcentration();
193 (*this)[Indices::contiUreaEqIdx] = extQuants.volumeFlux(FluidSystem::waterPhaseIdx) *
194 insideIntQuants.ureaConcentration();
196 (*this)[Indices::contiUreaEqIdx] *= getPropValue<TypeTag, Properties::BlackOilUreaScalingFactor>();
200 LocalResidual::adaptMassConservationQuantities_(*
this, insideIntQuants.pvtRegionIndex());
203 if constexpr (enableEnergy) {
205 getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>());
209 for (
unsigned i = 0; i < numEq; ++i) {
210 Valgrind::CheckDefined((*
this)[i]);
212 Valgrind::CheckDefined(*
this);
219 template <
class Context,
class Flu
idState>
223 const FluidState& fluidState)
225 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
229 std::for_each(this->begin(), this->end(),
230 [](
auto& val) { val = std::min(Scalar(0), val); });
236 template <
class Context,
class Flu
idState>
240 const FluidState& fluidState)
242 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
246 std::for_each(this->begin(), this->end(),
247 [](
auto& val) { val = std::max(Scalar(0), val); });
254 { (*this) = Scalar(0); }
263 template <
class Context,
class Flu
idState>
265 [[maybe_unused]]
unsigned bfIdx,
266 [[maybe_unused]]
unsigned timeIdx,
267 [[maybe_unused]]
const FluidState& boundaryFluidState)
273 if constexpr (enableEnergy) {
274 ExtensiveQuantities extQuants;
275 extQuants.updateBoundary(context, bfIdx, timeIdx, boundaryFluidState);
277 (*this)[contiEnergyEqIdx] += extQuants.energyFlux();
280 for (
unsigned i = 0; i < numEq; ++i) {
281 Valgrind::CheckDefined((*
this)[i]);
283 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:83
void setOutFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an outflow boundary.
Definition: blackoilboundaryratevector.hh:237
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:96
BlackOilBoundaryRateVector()=default
Default constructor.
void setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition: blackoilboundaryratevector.hh:253
void setInFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition: blackoilboundaryratevector.hh:220
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:264
Contains the high level supplements required to extend the black oil model by energy.
Definition: blackoilenergymodules.hh:58
static void addToEnthalpyRate(RateVector &flux, const Evaluation &hRate)
Definition: blackoilenergymodules.hh:251
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