28 #ifndef EWOMS_FLASH_BOUNDARY_RATE_VECTOR_HH 29 #define EWOMS_FLASH_BOUNDARY_RATE_VECTOR_HH 34 #include <opm/common/Valgrind.hpp> 45 template <
class TypeTag>
48 typedef typename GET_PROP_TYPE(TypeTag, RateVector) ParentType;
49 typedef typename GET_PROP_TYPE(TypeTag, ExtensiveQuantities) ExtensiveQuantities;
50 typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
51 typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
52 typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
53 typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
58 enum { conti0EqIdx = Indices::conti0EqIdx };
62 typedef Opm::MathToolbox<Evaluation> Toolbox;
85 template <
class Context,
class Flu
idState>
89 const FluidState& fluidState)
91 typename FluidSystem::template ParameterCache<typename FluidState::Scalar> paramCache;
92 paramCache.updateAll(fluidState);
94 ExtensiveQuantities extQuants;
95 extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState, paramCache);
96 const auto& insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
101 (*this) = Evaluation(0.0);
102 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
103 Evaluation meanMBoundary = 0;
104 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
106 fluidState.moleFraction(phaseIdx, compIdx)*FluidSystem::molarMass(compIdx);
109 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx))
110 density = FluidSystem::density(fluidState, paramCache, phaseIdx);
112 density = insideIntQuants.fluidState().density(phaseIdx);
114 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
116 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx))
117 molarity = fluidState.moleFraction(phaseIdx, compIdx)*density/meanMBoundary;
119 molarity = insideIntQuants.fluidState().molarity(phaseIdx, compIdx);
123 (*this)[conti0EqIdx + compIdx] += extQuants.volumeFlux(phaseIdx)*molarity;
127 Evaluation specificEnthalpy;
128 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx))
129 specificEnthalpy = FluidSystem::enthalpy(fluidState, paramCache, phaseIdx);
131 specificEnthalpy = insideIntQuants.fluidState().enthalpy(phaseIdx);
134 Evaluation enthalpyRate = density*extQuants.volumeFlux(phaseIdx)*specificEnthalpy;
135 EnergyModule::addToEnthalpyRate(*
this, enthalpyRate);
139 EnergyModule::addToEnthalpyRate(*
this, EnergyModule::heatConductionRate(extQuants));
142 for (
unsigned i = 0; i < numEq; ++i) {
143 Opm::Valgrind::CheckDefined((*
this)[i]);
151 template <
class Context,
class Flu
idState>
155 const FluidState& fluidState)
157 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
160 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
161 Evaluation& val = this->operator[](eqIdx);
162 val = Toolbox::min(0.0, val);
169 template <
class Context,
class Flu
idState>
173 const FluidState& fluidState)
175 this->
setFreeFlow(context, bfIdx, timeIdx, fluidState);
178 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
179 Evaluation& val = this->operator[](eqIdx);
180 val = Toolbox::max(0.0, val);
188 { (*this) = Evaluation(0.0); }
void setFreeFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify a free-flow boundary.
Definition: flashboundaryratevector.hh:86
FlashBoundaryRateVector(const Evaluation &value)
Definition: flashboundaryratevector.hh:72
Definition: baseauxiliarymodule.hh:37
Provides the auxiliary methods required for consideration of the energy equation. ...
Definition: energymodule.hh:59
#define GET_PROP_VALUE(TypeTag, PropTagName)
Access the value attribute of a property for a type tag.
Definition: propertysystem.hh:469
Contains the classes required to consider energy as a conservation quantity in a multi-phase module...
#define GET_PROP_TYPE(TypeTag, PropTagName)
Access the type attribute of a property for a type tag.
Definition: propertysystem.hh:486
void setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition: flashboundaryratevector.hh:187
void setInFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition: flashboundaryratevector.hh:152
Declares the properties required by the compositional multi-phase model based on flash calculations...
void setOutFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an outflow boundary.
Definition: flashboundaryratevector.hh:170
Implements a boundary vector for the fully implicit compositional multi-phase model which is based on...
Definition: flashboundaryratevector.hh:46