25#ifndef DUMUX_IMMISCIBLE_LOCAL_RESIDUAL_HH
26#define DUMUX_IMMISCIBLE_LOCAL_RESIDUAL_HH
38template<
class TypeTag>
50 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
51 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
53 using Element =
typename GridView::template Codim<0>::Entity;
57 static constexpr int numPhases = ModelTraits::numFluidPhases();
58 static constexpr int conti0EqIdx = ModelTraits::Indices::conti0EqIdx;
61 using ParentType::ParentType;
76 const SubControlVolume& scv,
77 const VolumeVariables& volVars)
const
81 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
83 auto eqIdx = conti0EqIdx + phaseIdx;
84 storage[eqIdx] = volVars.porosity()
85 * volVars.density(phaseIdx)
86 * volVars.saturation(phaseIdx);
89 EnergyLocalResidual::fluidPhaseStorage(storage, scv, volVars, phaseIdx);
93 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
110 const Element& element,
111 const FVElementGeometry& fvGeometry,
112 const ElementVolumeVariables& elemVolVars,
113 const SubControlVolumeFace& scvf,
114 const ElementFluxVariablesCache& elemFluxVarsCache)
const
116 FluxVariables fluxVars;
117 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
120 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
123 auto upwindTerm = [phaseIdx](
const auto& volVars)
124 {
return volVars.density(phaseIdx)*volVars.mobility(phaseIdx); };
126 auto eqIdx = conti0EqIdx + phaseIdx;
127 flux[eqIdx] = fluxVars.advectiveFlux(phaseIdx, upwindTerm);
130 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
134 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
A helper to deduce a vector with the same size as numbers of equations.
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition: numeqvector.hh:46
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
Element-wise calculation of the residual for problems using the n-phase immiscible fully implicit mod...
Definition: porousmediumflow/immiscible/localresidual.hh:40
NumEqVector computeFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache) const
Evaluatex the mass flux over a face of a sub control volume.
Definition: porousmediumflow/immiscible/localresidual.hh:109
NumEqVector computeStorage(const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars) const
Evaluatex the rate of change of all conservation quantites (e.g. phase mass) within a sub-control vol...
Definition: porousmediumflow/immiscible/localresidual.hh:75
Declares all properties used in Dumux.