28#ifndef OPM_PTFLASH_LOCAL_RESIDUAL_HH
29#define OPM_PTFLASH_LOCAL_RESIDUAL_HH
31#include <opm/material/common/MathToolbox.hpp>
32#include <opm/material/common/Valgrind.hpp>
45template <
class TypeTag>
46class FlashLocalResidual:
public GetPropType<TypeTag, Properties::DiscLocalResidual>
48 using Evaluation = GetPropType<TypeTag, Properties::Evaluation>;
49 using EqVector = GetPropType<TypeTag, Properties::EqVector>;
50 using RateVector = GetPropType<TypeTag, Properties::RateVector>;
51 using Indices = GetPropType<TypeTag, Properties::Indices>;
52 using IntensiveQuantities = GetPropType<TypeTag, Properties::IntensiveQuantities>;
53 using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
54 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
56 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
57 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
58 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
59 enum { water0Idx = Indices::water0Idx };
60 enum { conti0EqIdx = Indices::conti0EqIdx };
62 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
64 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
67 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
70 static constexpr bool waterEnabled = Indices::waterEnabled;
72 using Toolbox = MathToolbox<Evaluation>;
78 template <
class LhsEval>
80 const ElementContext& elemCtx,
83 unsigned phaseIdx)
const
85 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
86 const auto& fs = intQuants.fluidState();
89 if (waterEnabled && phaseIdx ==
static_cast<unsigned int>(waterPhaseIdx)) {
90 const unsigned eqIdx = conti0EqIdx + numComponents;
92 Toolbox::template decay<LhsEval>(fs.density(phaseIdx)) *
93 Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx)) *
94 Toolbox::template decay<LhsEval>(intQuants.porosity());
98 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
99 const unsigned eqIdx = conti0EqIdx + compIdx;
101 Toolbox::template decay<LhsEval>(fs.massFraction(phaseIdx, compIdx)) *
102 Toolbox::template decay<LhsEval>(fs.density(phaseIdx)) *
103 Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx)) *
104 Toolbox::template decay<LhsEval>(intQuants.porosity());
108 EnergyModule::addPhaseStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx), phaseIdx);
114 template <
class LhsEval>
116 const ElementContext& elemCtx,
118 unsigned timeIdx)
const
121 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
125 EnergyModule::addSolidEnergyStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx));
132 const ElementContext& elemCtx,
134 unsigned timeIdx)
const
138 Valgrind::CheckDefined(flux);
141 Valgrind::CheckDefined(flux);
148 const ElementContext& elemCtx,
150 unsigned timeIdx)
const
152 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
154 const unsigned focusDofIdx = elemCtx.focusDofIndex();
155 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
157 const auto upIdx =
static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
158 const IntensiveQuantities& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
163 if (upIdx == focusDofIdx) {
164 const Evaluation tmp =
165 up.fluidState().density(phaseIdx) *
166 extQuants.volumeFlux(phaseIdx);
168 if (waterEnabled && phaseIdx ==
static_cast<unsigned int>(waterPhaseIdx)) {
169 const unsigned eqIdx = conti0EqIdx + numComponents;
173 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
174 flux[conti0EqIdx + compIdx] +=
175 tmp * up.fluidState().massFraction(phaseIdx, compIdx);
180 const Evaluation tmp =
181 Toolbox::value(up.fluidState().density(phaseIdx)) *
182 extQuants.volumeFlux(phaseIdx);
184 if (waterEnabled && phaseIdx ==
static_cast<unsigned int>(waterPhaseIdx)) {
185 const unsigned eqIdx = conti0EqIdx + numComponents;
189 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
190 flux[conti0EqIdx + compIdx] +=
191 tmp * Toolbox::value(up.fluidState().massFraction(phaseIdx, compIdx));
197 EnergyModule::addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
204 const ElementContext& elemCtx,
206 unsigned timeIdx)
const
208 DiffusionModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
209 EnergyModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
216 const ElementContext& elemCtx,
218 unsigned timeIdx)
const
220 Valgrind::SetUndefined(source);
221 elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
222 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
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: ptflash/flashlocalresidual.hh:115
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: ptflash/flashlocalresidual.hh:79
void computeSource(RateVector &source, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Calculate the source term of the equation.
Definition: ptflash/flashlocalresidual.hh:215
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: ptflash/flashlocalresidual.hh:131
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