28#ifndef EWOMS_NCP_LOCAL_RESIDUAL_HH
29#define EWOMS_NCP_LOCAL_RESIDUAL_HH
31#include <dune/istl/bvector.hh>
33#include <opm/material/common/MathToolbox.hpp>
34#include <opm/material/common/Valgrind.hpp>
51template <
class TypeTag>
63 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
64 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
65 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
66 enum { ncp0EqIdx = Indices::ncp0EqIdx };
67 enum { conti0EqIdx = Indices::conti0EqIdx };
69 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
72 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
75 using EvalEqVector = Dune::FieldVector<Evaluation, numEq>;
76 using ElemEvalEqVector = Dune::BlockVector<EvalEqVector>;
77 using Toolbox = MathToolbox<Evaluation>;
83 template <
class LhsEval>
85 const ElementContext& elemCtx,
88 unsigned phaseIdx)
const
90 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
91 const auto& fluidState = intQuants.fluidState();
94 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
95 const unsigned eqIdx = conti0EqIdx + compIdx;
97 Toolbox::template decay<LhsEval>(fluidState.molarity(phaseIdx, compIdx)) *
98 Toolbox::template decay<LhsEval>(fluidState.saturation(phaseIdx)) *
99 Toolbox::template decay<LhsEval>(intQuants.porosity());
102 EnergyModule::addPhaseStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx), phaseIdx);
108 template <
class LhsEval>
110 const ElementContext& elemCtx,
112 unsigned timeIdx)
const
115 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
119 EnergyModule::addSolidEnergyStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx));
126 const ElementContext& elemCtx,
128 unsigned timeIdx)
const
132 Valgrind::CheckDefined(flux);
135 Valgrind::CheckDefined(flux);
142 const ElementContext& elemCtx,
144 unsigned timeIdx)
const
146 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
148 unsigned focusDofIdx = elemCtx.focusDofIndex();
149 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
152 const unsigned upIdx =
static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
153 const IntensiveQuantities& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
158 if (upIdx == focusDofIdx) {
159 const Evaluation tmp =
160 up.fluidState().molarDensity(phaseIdx) *
161 extQuants.volumeFlux(phaseIdx);
163 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
164 flux[conti0EqIdx + compIdx] +=
165 tmp * up.fluidState().moleFraction(phaseIdx, compIdx);
169 const Evaluation tmp =
170 Toolbox::value(up.fluidState().molarDensity(phaseIdx)) *
171 extQuants.volumeFlux(phaseIdx);
173 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
174 flux[conti0EqIdx + compIdx] +=
175 tmp * Toolbox::value(up.fluidState().moleFraction(phaseIdx, compIdx));
180 EnergyModule::addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
187 const ElementContext& elemCtx,
189 unsigned timeIdx)
const
191 DiffusionModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
192 EnergyModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
202 const ElementContext& elemCtx,
204 unsigned timeIdx)
const
206 Valgrind::SetUndefined(source);
207 elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
208 Valgrind::CheckDefined(source);
211 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
212 source[ncp0EqIdx + phaseIdx] =
phaseNcp(elemCtx, dofIdx, timeIdx, phaseIdx);
219 template <
class LhsEval = Evaluation>
223 unsigned phaseIdx)
const
225 const auto& fluidState = elemCtx.intensiveQuantities(dofIdx, timeIdx).fluidState();
226 using FluidState = std::remove_const_t<std::remove_reference_t<
decltype(fluidState)>>;
228 using LhsToolbox = MathToolbox<LhsEval>;
230 const LhsEval& a = phaseNotPresentIneq_<FluidState, LhsEval>(fluidState, phaseIdx);
231 const LhsEval& b = phasePresentIneq_<FluidState, LhsEval>(fluidState, phaseIdx);
232 return LhsToolbox::min(a, b);
240 template <
class Flu
idState,
class LhsEval>
241 LhsEval phasePresentIneq_(
const FluidState& fluidState,
unsigned phaseIdx)
const
243 using FsToolbox = MathToolbox<typename FluidState::Scalar>;
244 return FsToolbox::template decay<LhsEval>(fluidState.saturation(phaseIdx));
251 template <
class Flu
idState,
class LhsEval>
252 LhsEval phaseNotPresentIneq_(
const FluidState& fluidState,
unsigned phaseIdx)
const
254 using FsToolbox = MathToolbox<typename FluidState::Scalar>;
258 for (
unsigned i = 0; i < numComponents; ++i) {
259 a -= FsToolbox::template decay<LhsEval>(fluidState.moleFraction(phaseIdx, i));
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
Details needed to calculate the local residual in the compositional multi-phase NCP-model .
Definition: ncplocalresidual.hh:53
LhsEval phaseNcp(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx, unsigned phaseIdx) const
Returns the value of the NCP-function for a phase.
Definition: ncplocalresidual.hh:220
void addAdvectiveFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Add the advective mass flux at a given flux integration point.
Definition: ncplocalresidual.hh:141
void addDiffusiveFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Adds the diffusive flux at a given flux integration point.
Definition: ncplocalresidual.hh:186
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: ncplocalresidual.hh:125
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: ncplocalresidual.hh:84
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: ncplocalresidual.hh:109
void computeSource(RateVector &source, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Calculate the source term of the equation.
Definition: ncplocalresidual.hh:201
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