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>
44template <
class TypeTag>
56 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
57 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
58 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
59 enum { enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>() };
60 enum { enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>() };
61 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
62 enum { conti0EqIdx = Indices::conti0EqIdx };
63 enum { contiEnergyEqIdx = Indices::contiEnergyEqIdx };
64 enum { enableFoam = getPropValue<TypeTag, Properties::EnableFoam>() };
65 enum { enableMICP = getPropValue<TypeTag, Properties::EnableMICP>() };
67 static constexpr bool blackoilConserveSurfaceVolume = getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>();
93 template <
class Context,
class Flu
idState>
97 const FluidState& fluidState)
99 ExtensiveQuantities extQuants;
100 extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState);
101 const auto& insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
102 unsigned focusDofIdx = context.focusDofIndex();
103 unsigned interiorDofIdx = context.interiorScvIndex(bfIdx, timeIdx);
109 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
110 if (!FluidSystem::phaseIsActive(phaseIdx)) {
113 const auto& pBoundary = fluidState.pressure(phaseIdx);
114 const Evaluation& pInside = insideIntQuants.fluidState().pressure(phaseIdx);
119 if (pBoundary < pInside)
121 LocalResidual::template evalPhaseFluxes_<Evaluation>(tmp,
123 insideIntQuants.pvtRegionIndex(),
125 insideIntQuants.fluidState());
126 else if (pBoundary > pInside) {
127 using RhsEval =
typename std::conditional<std::is_same<typename FluidState::Scalar, Evaluation>::value,
128 Evaluation, Scalar>::type;
130 LocalResidual::template evalPhaseFluxes_<RhsEval>(tmp,
132 insideIntQuants.pvtRegionIndex(),
137 for (
unsigned i = 0; i < tmp.size(); ++i)
138 (*
this)[i] += tmp[i];
141 if constexpr (enableEnergy) {
143 Evaluation specificEnthalpy;
144 if (pBoundary > pInside) {
145 if (focusDofIdx == interiorDofIdx) {
146 density = fluidState.density(phaseIdx);
147 specificEnthalpy = fluidState.enthalpy(phaseIdx);
150 density = getValue(fluidState.density(phaseIdx));
151 specificEnthalpy = getValue(fluidState.enthalpy(phaseIdx));
154 else if (focusDofIdx == interiorDofIdx) {
155 density = insideIntQuants.fluidState().density(phaseIdx);
156 specificEnthalpy = insideIntQuants.fluidState().enthalpy(phaseIdx);
159 density = getValue(insideIntQuants.fluidState().density(phaseIdx));
160 specificEnthalpy = getValue(insideIntQuants.fluidState().enthalpy(phaseIdx));
163 Evaluation enthalpyRate = density*extQuants.volumeFlux(phaseIdx)*specificEnthalpy;
168 if constexpr (enableSolvent) {
169 (*this)[Indices::contiSolventEqIdx] = extQuants.solventVolumeFlux();
170 if (blackoilConserveSurfaceVolume)
171 (*this)[Indices::contiSolventEqIdx] *= insideIntQuants.solventInverseFormationVolumeFactor();
173 (*
this)[Indices::contiSolventEqIdx] *= insideIntQuants.solventDensity();
177 if constexpr (enablePolymer) {
178 (*this)[Indices::contiPolymerEqIdx] = extQuants.volumeFlux(FluidSystem::waterPhaseIdx) * insideIntQuants.polymerConcentration();
181 if constexpr (enableMICP) {
182 (*this)[Indices::contiMicrobialEqIdx] = extQuants.volumeFlux(FluidSystem::waterPhaseIdx) * insideIntQuants.microbialConcentration();
183 (*this)[Indices::contiOxygenEqIdx] = extQuants.volumeFlux(FluidSystem::waterPhaseIdx) * insideIntQuants.oxygenConcentration();
184 (*this)[Indices::contiUreaEqIdx] = extQuants.volumeFlux(FluidSystem::waterPhaseIdx) * insideIntQuants.ureaConcentration();
188 LocalResidual::adaptMassConservationQuantities_(*
this, insideIntQuants.pvtRegionIndex());
191 if constexpr (enableEnergy)
195 for (
unsigned i = 0; i < numEq; ++i) {
196 Valgrind::CheckDefined((*
this)[i]);
198 Valgrind::CheckDefined(*
this);
205 template <
class Context,
class Flu
idState>
209 const FluidState& fluidState)
211 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
215 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
216 Scalar& val = this->operator[](eqIdx);
217 val = std::min<Scalar>(0.0, val);
224 template <
class Context,
class Flu
idState>
228 const FluidState& fluidState)
230 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
234 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
235 Scalar& val = this->operator[](eqIdx);
236 val = std::max( Scalar(0), val);
244 { (*this) = Scalar(0); }
253 template <
class Context,
class Flu
idState>
255 [[maybe_unused]]
unsigned bfIdx,
256 [[maybe_unused]]
unsigned timeIdx,
257 [[maybe_unused]]
const FluidState& boundaryFluidState)
263 if constexpr (enableEnergy) {
264 ExtensiveQuantities extQuants;
265 extQuants.updateBoundary(context, bfIdx, timeIdx, boundaryFluidState);
267 (*this)[contiEnergyEqIdx] += extQuants.energyFlux();
270 for (
unsigned i = 0; i < numEq; ++i)
271 Valgrind::CheckDefined((*
this)[i]);
272 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:46
BlackOilBoundaryRateVector()
Default constructor.
Definition: blackoilboundaryratevector.hh:75
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:225
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:94
void setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition: blackoilboundaryratevector.hh:243
void setInFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition: blackoilboundaryratevector.hh:206
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:254
Contains the high level supplements required to extend the black oil model by energy.
Definition: blackoilenergymodules.hh:52
static void addToEnthalpyRate(RateVector &flux, const Evaluation &hRate)
Definition: blackoilenergymodules.hh:237
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