14#ifndef DUMUX_3P3C_LOCAL_RESIDUAL_HH
15#define DUMUX_3P3C_LOCAL_RESIDUAL_HH
28template<
class TypeTag>
35 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
36 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
42 using Element =
typename GridView::template Codim<0>::Entity;
52 contiWEqIdx = Indices::conti0EqIdx + FluidSystem::wPhaseIdx,
53 contiNEqIdx = Indices::conti0EqIdx + FluidSystem::nPhaseIdx,
54 contiGEqIdx = Indices::conti0EqIdx + FluidSystem::gPhaseIdx,
56 wPhaseIdx = FluidSystem::wPhaseIdx,
57 nPhaseIdx = FluidSystem::nPhaseIdx,
58 gPhaseIdx = FluidSystem::gPhaseIdx,
60 wCompIdx = FluidSystem::wCompIdx,
61 nCompIdx = FluidSystem::nCompIdx,
62 gCompIdx = FluidSystem::gCompIdx
67 using ParentType::ParentType;
81 const SubControlVolume& scv,
82 const VolumeVariables& volVars)
const
84 NumEqVector storage(0.0);
87 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
89 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
91 auto eqIdx = Indices::conti0EqIdx + compIdx;
92 storage[eqIdx] += volVars.porosity()
93 * volVars.saturation(phaseIdx)
94 * volVars.molarDensity(phaseIdx)
95 * volVars.moleFraction(phaseIdx, compIdx);
99 EnergyLocalResidual::fluidPhaseStorage(storage, problem, scv, volVars, phaseIdx);
103 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
120 const Element& element,
121 const FVElementGeometry& fvGeometry,
122 const ElementVolumeVariables& elemVolVars,
123 const SubControlVolumeFace& scvf,
124 const ElementFluxVariablesCache& elemFluxVarsCache)
const
126 FluxVariables fluxVars;
127 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
128 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
131 NumEqVector flux(0.0);
134 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
136 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
138 auto upwindTerm = [phaseIdx, compIdx](
const VolumeVariables& volVars)
139 {
return volVars.molarDensity(phaseIdx)*volVars.moleFraction(phaseIdx, compIdx)*volVars.mobility(phaseIdx); };
142 auto eqIdx = Indices::conti0EqIdx + compIdx;
143 flux[eqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
147 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
151 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
154 const auto diffusionFluxesWPhase = fluxVars.molecularDiffusionFlux(wPhaseIdx);
155 Scalar jGW = diffusionFluxesWPhase[gCompIdx];
156 Scalar jNW = diffusionFluxesWPhase[nCompIdx];
157 Scalar jWW = -(jGW+jNW);
162 jGW /= FluidSystem::molarMass(gCompIdx);
163 jNW /= FluidSystem::molarMass(nCompIdx);
164 jWW /= FluidSystem::molarMass(wCompIdx);
167 const auto diffusionFluxesGPhase = fluxVars.molecularDiffusionFlux(gPhaseIdx);
168 Scalar jWG = diffusionFluxesGPhase[wCompIdx];
169 Scalar jNG = diffusionFluxesGPhase[nCompIdx];
170 Scalar jGG = -(jWG+jNG);
175 jWG /= FluidSystem::molarMass(wCompIdx);
176 jNG /= FluidSystem::molarMass(nCompIdx);
177 jGG /= FluidSystem::molarMass(gCompIdx);
181 const Scalar jWN = 0.0;
182 const Scalar jGN = 0.0;
183 const Scalar jNN = 0.0;
185 flux[contiWEqIdx] += jWW+jWG+jWN;
186 flux[contiNEqIdx] += jNW+jNG+jNN;
187 flux[contiGEqIdx] += jGW+jGG+jGN;
Element-wise calculation of the Jacobian matrix for problems using the three-phase three-component fu...
Definition: porousmediumflow/3p3c/localresidual.hh:30
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:80
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:119
Defines all properties used in Dumux.
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:34
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
A helper to deduce a vector with the same size as numbers of equations.