28#ifndef EWOMS_PVS_LOCAL_RESIDUAL_HH
29#define EWOMS_PVS_LOCAL_RESIDUAL_HH
36#include <opm/material/common/Valgrind.hpp>
46template <
class TypeTag>
56 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
57 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
58 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
59 enum { conti0EqIdx = Indices::conti0EqIdx };
61 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
64 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
67 using Toolbox = Opm::MathToolbox<Evaluation>;
73 template <
class LhsEval>
75 const ElementContext& elemCtx,
78 unsigned phaseIdx)
const
80 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
81 const auto& fs = intQuants.fluidState();
84 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
85 unsigned eqIdx = conti0EqIdx + compIdx;
87 Toolbox::template decay<LhsEval>(fs.molarity(phaseIdx, compIdx))
88 * Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx))
89 * Toolbox::template decay<LhsEval>(intQuants.porosity());
92 EnergyModule::addPhaseStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx), phaseIdx);
98 template <
class LhsEval>
100 const ElementContext& elemCtx,
102 unsigned timeIdx)
const
105 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
108 EnergyModule::addSolidEnergyStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx));
115 const ElementContext& elemCtx,
117 unsigned timeIdx)
const
121 Opm::Valgrind::CheckDefined(flux);
124 Opm::Valgrind::CheckDefined(flux);
131 const ElementContext& elemCtx,
133 unsigned timeIdx)
const
135 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
137 unsigned focusDofIdx = elemCtx.focusDofIndex();
138 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
141 unsigned upIdx =
static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
142 const IntensiveQuantities& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
147 if (upIdx == focusDofIdx) {
149 up.fluidState().molarDensity(phaseIdx)
150 * extQuants.volumeFlux(phaseIdx);
152 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
153 flux[conti0EqIdx + compIdx] +=
154 tmp*up.fluidState().moleFraction(phaseIdx, compIdx);
159 Toolbox::value(up.fluidState().molarDensity(phaseIdx))
160 * extQuants.volumeFlux(phaseIdx);
162 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
163 flux[conti0EqIdx + compIdx] +=
164 tmp*Toolbox::value(up.fluidState().moleFraction(phaseIdx, compIdx));
169 EnergyModule::addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
176 const ElementContext& elemCtx,
178 unsigned timeIdx)
const
180 DiffusionModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
181 EnergyModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
188 const ElementContext& elemCtx,
190 unsigned timeIdx)
const
192 Opm::Valgrind::SetUndefined(source);
193 elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
194 Opm::Valgrind::CheckDefined(source);
Provides the auxiliary methods required for consideration of the diffusion equation.
Definition: diffusionmodule.hh:48
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:50
Element-wise calculation of the local residual for the compositional multi-phase primary variable swi...
Definition: pvslocalresidual.hh:48
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: pvslocalresidual.hh:74
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: pvslocalresidual.hh:99
void computeSource(RateVector &source, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Calculate the source term of the equation.
Definition: pvslocalresidual.hh:187
void addDiffusiveFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Adds the diffusive flux at a given flux integration point.
Definition: pvslocalresidual.hh:175
void addAdvectiveFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Add the advective mass flux at a given flux integration point.
Definition: pvslocalresidual.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: pvslocalresidual.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: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 compositional multi-phase primary variable switching model.