28#ifndef EWOMS_PVS_LOCAL_RESIDUAL_HH
29#define EWOMS_PVS_LOCAL_RESIDUAL_HH
31#include <opm/material/common/MathToolbox.hpp>
32#include <opm/material/common/Valgrind.hpp>
47template <
class TypeTag>
57 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
58 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
59 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
60 enum { conti0EqIdx = Indices::conti0EqIdx };
62 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
65 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
68 using Toolbox = MathToolbox<Evaluation>;
74 template <
class LhsEval>
76 const ElementContext& elemCtx,
79 unsigned phaseIdx)
const
81 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
82 const auto& fs = intQuants.fluidState();
85 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
86 const unsigned eqIdx = conti0EqIdx + compIdx;
88 Toolbox::template decay<LhsEval>(fs.molarity(phaseIdx, compIdx)) *
89 Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx)) *
90 Toolbox::template decay<LhsEval>(intQuants.porosity());
93 EnergyModule::addPhaseStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx), phaseIdx);
99 template <
class LhsEval>
101 const ElementContext& elemCtx,
103 unsigned timeIdx)
const
106 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
110 EnergyModule::addSolidEnergyStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx));
117 const ElementContext& elemCtx,
119 unsigned timeIdx)
const
123 Valgrind::CheckDefined(flux);
126 Valgrind::CheckDefined(flux);
133 const ElementContext& elemCtx,
135 unsigned timeIdx)
const
137 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
139 const unsigned focusDofIdx = elemCtx.focusDofIndex();
140 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
143 const unsigned upIdx =
static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
144 const IntensiveQuantities& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
149 if (upIdx == focusDofIdx) {
150 const Evaluation tmp =
151 up.fluidState().molarDensity(phaseIdx) *
152 extQuants.volumeFlux(phaseIdx);
154 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
155 flux[conti0EqIdx + compIdx] +=
156 tmp * up.fluidState().moleFraction(phaseIdx, compIdx);
160 const Evaluation tmp =
161 Toolbox::value(up.fluidState().molarDensity(phaseIdx)) *
162 extQuants.volumeFlux(phaseIdx);
164 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
165 flux[conti0EqIdx + compIdx] +=
166 tmp * Toolbox::value(up.fluidState().moleFraction(phaseIdx, compIdx));
171 EnergyModule::addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
178 const ElementContext& elemCtx,
180 unsigned timeIdx)
const
182 DiffusionModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
183 EnergyModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
190 const ElementContext& elemCtx,
192 unsigned timeIdx)
const
194 Valgrind::SetUndefined(source);
195 elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
196 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
Element-wise calculation of the local residual for the compositional multi-phase primary variable swi...
Definition: pvslocalresidual.hh:49
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:75
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:100
void computeSource(RateVector &source, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Calculate the source term of the equation.
Definition: pvslocalresidual.hh:189
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:177
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:132
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:116
Classes required for molecular diffusion.
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
Declare the properties used by the infrastructure code of the finite volume discretizations.
Defines the common properties required by the porous medium multi-phase models.
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