28#ifndef EWOMS_NCP_BOUNDARY_RATE_VECTOR_HH
29#define EWOMS_NCP_BOUNDARY_RATE_VECTOR_HH
34#include <opm/material/common/Valgrind.hpp>
44template <
class TypeTag>
54 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
55 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
56 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
57 enum { conti0EqIdx = Indices::conti0EqIdx };
58 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
62 using Toolbox = Opm::MathToolbox<Evaluation>;
85 template <
class Context,
class Flu
idState>
89 const FluidState& fluidState)
91 ExtensiveQuantities extQuants;
92 extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState);
93 const auto& insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
94 unsigned focusDofIdx = context.focusDofIndex();
95 unsigned interiorDofIdx = context.interiorScvIndex(bfIdx, timeIdx);
100 (*this) = Evaluation(0.0);
101 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
103 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
104 if (focusDofIdx == interiorDofIdx)
105 density = fluidState.density(phaseIdx);
107 density = Opm::getValue(fluidState.density(phaseIdx));
109 else if (focusDofIdx == interiorDofIdx)
110 density = insideIntQuants.fluidState().density(phaseIdx);
112 density = Opm::getValue(insideIntQuants.fluidState().density(phaseIdx));
114 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
116 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
117 if (focusDofIdx == interiorDofIdx)
118 molarity = fluidState.molarity(phaseIdx, compIdx);
120 molarity = Opm::getValue(fluidState.molarity(phaseIdx, compIdx));
122 else if (focusDofIdx == interiorDofIdx)
123 molarity = insideIntQuants.fluidState().molarity(phaseIdx, compIdx);
125 molarity = Opm::getValue(insideIntQuants.fluidState().molarity(phaseIdx, compIdx));
129 (*this)[conti0EqIdx + compIdx] += extQuants.volumeFlux(phaseIdx)*molarity;
133 Evaluation specificEnthalpy;
134 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
135 if (focusDofIdx == interiorDofIdx)
136 specificEnthalpy = fluidState.enthalpy(phaseIdx);
138 specificEnthalpy = Opm::getValue(fluidState.enthalpy(phaseIdx));
140 else if (focusDofIdx == interiorDofIdx)
141 specificEnthalpy = insideIntQuants.fluidState().enthalpy(phaseIdx);
143 specificEnthalpy = Opm::getValue(insideIntQuants.fluidState().enthalpy(phaseIdx));
145 Evaluation enthalpyRate = density*extQuants.volumeFlux(phaseIdx)*specificEnthalpy;
146 EnergyModule::addToEnthalpyRate(*
this, enthalpyRate);
151 EnergyModule::addToEnthalpyRate(*
this, EnergyModule::thermalConductionRate(extQuants));
154 for (
unsigned i = 0; i < numEq; ++i) {
155 Opm::Valgrind::CheckDefined((*
this)[i]);
163 template <
class Context,
class Flu
idState>
167 const FluidState& fluidState)
169 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
172 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
173 this->operator[](eqIdx) = Toolbox::min(0.0, this->
operator[](eqIdx));
180 template <
class Context,
class Flu
idState>
184 const FluidState& fluidState)
186 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
189 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
190 this->operator[](eqIdx) = Toolbox::max(0.0, this->
operator[](eqIdx));
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 NCP model.
Definition: ncpboundaryratevector.hh:46
NcpBoundaryRateVector(const Evaluation &value)
Definition: ncpboundaryratevector.hh:72
void setInFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition: ncpboundaryratevector.hh:164
void setFreeFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify a free-flow boundary.
Definition: ncpboundaryratevector.hh:86
void setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition: ncpboundaryratevector.hh:197
NcpBoundaryRateVector & operator=(const NcpBoundaryRateVector &value)=default
void setOutFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an outflow boundary.
Definition: ncpboundaryratevector.hh:181
NcpBoundaryRateVector(const NcpBoundaryRateVector &value)=default
NcpBoundaryRateVector()
Definition: ncpboundaryratevector.hh:65
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
Declares the properties required for the NCP compositional multi-phase model.