26#ifndef DUMUX_3P3C_LOCAL_RESIDUAL_HH
27#define DUMUX_3P3C_LOCAL_RESIDUAL_HH
40template<
class TypeTag>
47 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
48 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
54 using Element =
typename GridView::template Codim<0>::Entity;
64 contiWEqIdx = Indices::conti0EqIdx + FluidSystem::wPhaseIdx,
65 contiNEqIdx = Indices::conti0EqIdx + FluidSystem::nPhaseIdx,
66 contiGEqIdx = Indices::conti0EqIdx + FluidSystem::gPhaseIdx,
68 wPhaseIdx = FluidSystem::wPhaseIdx,
69 nPhaseIdx = FluidSystem::nPhaseIdx,
70 gPhaseIdx = FluidSystem::gPhaseIdx,
72 wCompIdx = FluidSystem::wCompIdx,
73 nCompIdx = FluidSystem::nCompIdx,
74 gCompIdx = FluidSystem::gCompIdx
79 using ParentType::ParentType;
93 const SubControlVolume& scv,
94 const VolumeVariables& volVars)
const
96 NumEqVector storage(0.0);
99 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
101 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
103 auto eqIdx = Indices::conti0EqIdx + compIdx;
104 storage[eqIdx] += volVars.porosity()
105 * volVars.saturation(phaseIdx)
106 * volVars.molarDensity(phaseIdx)
107 * volVars.moleFraction(phaseIdx, compIdx);
111 EnergyLocalResidual::fluidPhaseStorage(storage, scv, volVars, phaseIdx);
115 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
132 const Element& element,
133 const FVElementGeometry& fvGeometry,
134 const ElementVolumeVariables& elemVolVars,
135 const SubControlVolumeFace& scvf,
136 const ElementFluxVariablesCache& elemFluxVarsCache)
const
138 FluxVariables fluxVars;
139 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
140 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
143 NumEqVector flux(0.0);
146 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
148 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
150 auto upwindTerm = [phaseIdx, compIdx](
const VolumeVariables& volVars)
151 {
return volVars.molarDensity(phaseIdx)*volVars.moleFraction(phaseIdx, compIdx)*volVars.mobility(phaseIdx); };
154 auto eqIdx = Indices::conti0EqIdx + compIdx;
155 flux[eqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
159 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
163 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
166 const auto diffusionFluxesWPhase = fluxVars.molecularDiffusionFlux(wPhaseIdx);
167 Scalar jGW = diffusionFluxesWPhase[gCompIdx];
168 Scalar jNW = diffusionFluxesWPhase[nCompIdx];
169 Scalar jWW = -(jGW+jNW);
174 jGW /= FluidSystem::molarMass(gCompIdx);
175 jNW /= FluidSystem::molarMass(nCompIdx);
176 jWW /= FluidSystem::molarMass(wCompIdx);
179 const auto diffusionFluxesGPhase = fluxVars.molecularDiffusionFlux(gPhaseIdx);
180 Scalar jWG = diffusionFluxesGPhase[wCompIdx];
181 Scalar jNG = diffusionFluxesGPhase[nCompIdx];
182 Scalar jGG = -(jWG+jNG);
187 jWG /= FluidSystem::molarMass(wCompIdx);
188 jNG /= FluidSystem::molarMass(nCompIdx);
189 jGG /= FluidSystem::molarMass(gCompIdx);
193 const Scalar jWN = 0.0;
194 const Scalar jGN = 0.0;
195 const Scalar jNN = 0.0;
197 flux[contiWEqIdx] += jWW+jWG+jWN;
198 flux[contiNEqIdx] += jNW+jNG+jNN;
199 flux[contiGEqIdx] += jGW+jGG+jGN;
A helper to deduce a vector with the same size as numbers of equations.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
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 Jacobian matrix for problems using the three-phase three-component fu...
Definition: porousmediumflow/3p3c/localresidual.hh:42
NumEqVector computeStorage(const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars) const
Evaluates the amount of all conservation quantities (e.g. phase mass) within a sub-control volume.
Definition: porousmediumflow/3p3c/localresidual.hh:92
NumEqVector computeFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache) const
Evaluates the total flux of all conservation quantities over a face of a sub-control volume.
Definition: porousmediumflow/3p3c/localresidual.hh:131
Declares all properties used in Dumux.