28#ifndef EWOMS_FLASH_LOCAL_RESIDUAL_HH
29#define EWOMS_FLASH_LOCAL_RESIDUAL_HH
31#include <opm/material/common/MathToolbox.hpp>
32#include <opm/material/common/Valgrind.hpp>
45template <
class TypeTag>
55 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
56 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
57 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
58 enum { conti0EqIdx = Indices::conti0EqIdx };
60 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
63 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
66 using Toolbox = Opm::MathToolbox<Evaluation>;
72 template <
class LhsEval>
74 const ElementContext& elemCtx,
77 unsigned phaseIdx)
const
79 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
80 const auto& fs = intQuants.fluidState();
83 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
84 unsigned eqIdx = conti0EqIdx + compIdx;
86 Toolbox::template decay<LhsEval>(fs.molarity(phaseIdx, compIdx)) *
87 Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx)) *
88 Toolbox::template decay<LhsEval>(intQuants.porosity());
91 EnergyModule::addPhaseStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx), phaseIdx);
97 template <
class LhsEval>
99 const ElementContext& elemCtx,
101 unsigned timeIdx)
const
104 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
108 EnergyModule::addSolidEnergyStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx));
115 const ElementContext& elemCtx,
117 unsigned timeIdx)
const
121 Valgrind::CheckDefined(flux);
124 Valgrind::CheckDefined(flux);
131 const ElementContext& elemCtx,
133 unsigned timeIdx)
const
135 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
137 const unsigned focusDofIdx = elemCtx.focusDofIndex();
138 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
140 const unsigned upIdx =
static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
141 const IntensiveQuantities& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
146 if (upIdx == focusDofIdx) {
148 up.fluidState().molarDensity(phaseIdx) *
149 extQuants.volumeFlux(phaseIdx);
151 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
152 flux[conti0EqIdx + compIdx] +=
153 tmp * up.fluidState().moleFraction(phaseIdx, compIdx);
158 Toolbox::value(up.fluidState().molarDensity(phaseIdx)) *
159 extQuants.volumeFlux(phaseIdx);
161 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
162 flux[conti0EqIdx + compIdx] +=
163 tmp * Toolbox::value(up.fluidState().moleFraction(phaseIdx, compIdx));
168 EnergyModule::addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
175 const ElementContext& elemCtx,
177 unsigned timeIdx)
const
179 DiffusionModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
180 EnergyModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
187 const ElementContext& elemCtx,
189 unsigned timeIdx)
const
191 Valgrind::SetUndefined(source);
192 elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
193 Valgrind::CheckDefined(source);
Provides the auxiliary methods required for consideration of the diffusion equation.
Definition: diffusionmodule.hh:51
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:54
Calculates the local residual of the compositional multi-phase model based on flash calculations.
Definition: flash/flashlocalresidual.hh:47
void computeStorage(Dune::FieldVector< LhsEval, numEq > &storage, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Evaluate the amount all conservation quantities (e.g. phase mass) within a finite sub-control volume.
Definition: flash/flashlocalresidual.hh:98
void addDiffusiveFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Adds the diffusive flux at a given flux integration point.
Definition: flash/flashlocalresidual.hh:174
void addPhaseStorage(Dune::FieldVector< LhsEval, numEq > &storage, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx, unsigned phaseIdx) const
Adds the amount all conservation quantities (e.g. phase mass) within a single fluid phase.
Definition: flash/flashlocalresidual.hh:73
void computeSource(RateVector &source, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Calculate the source term of the equation.
Definition: flash/flashlocalresidual.hh:186
void addAdvectiveFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Add the advective mass flux at a given flux integration point.
Definition: flash/flashlocalresidual.hh:130
void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Evaluates the total mass flux of all conservation quantities over a face of a sub-control volume.
Definition: flash/flashlocalresidual.hh:114
Classes required for molecular diffusion.
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
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