version 3.8
porousmediumflow/immiscible/localresidual.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_IMMISCIBLE_LOCAL_RESIDUAL_HH
14#define DUMUX_IMMISCIBLE_LOCAL_RESIDUAL_HH
15
18
19namespace Dumux {
20
26template<class TypeTag>
27class ImmiscibleLocalResidual : public GetPropType<TypeTag, Properties::BaseLocalResidual>
28{
34 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
36 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
37 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
38 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
39 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
41 using Element = typename GridView::template Codim<0>::Entity;
43
45 static constexpr int numPhases = ModelTraits::numFluidPhases();
46 static constexpr int conti0EqIdx = ModelTraits::Indices::conti0EqIdx;
47
48public:
49 using ParentType::ParentType;
50
63 NumEqVector computeStorage(const Problem& problem,
64 const SubControlVolume& scv,
65 const VolumeVariables& volVars) const
66 {
67 // partial time derivative of the phase mass
68 NumEqVector storage;
69 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
70 {
71 auto eqIdx = conti0EqIdx + phaseIdx;
72 storage[eqIdx] = volVars.porosity()
73 * volVars.density(phaseIdx)
74 * volVars.saturation(phaseIdx);
75
77 EnergyLocalResidual::fluidPhaseStorage(storage, scv, volVars, phaseIdx);
78 }
79
81 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
82
83 return storage;
84 }
85
86
97 NumEqVector computeFlux(const Problem& problem,
98 const Element& element,
99 const FVElementGeometry& fvGeometry,
100 const ElementVolumeVariables& elemVolVars,
101 const SubControlVolumeFace& scvf,
102 const ElementFluxVariablesCache& elemFluxVarsCache) const
103 {
104 FluxVariables fluxVars;
105 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
106
107 NumEqVector flux;
108 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
109 {
110 // the physical quantities for which we perform upwinding
111 auto upwindTerm = [phaseIdx](const auto& volVars)
112 { return volVars.density(phaseIdx)*volVars.mobility(phaseIdx); };
113
114 auto eqIdx = conti0EqIdx + phaseIdx;
115 flux[eqIdx] = fluxVars.advectiveFlux(phaseIdx, upwindTerm);
116
118 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
119 }
120
122 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
123
124 return flux;
125 }
126};
127
128} // end namespace Dumux
129
130#endif
Element-wise calculation of the residual for problems using the n-phase immiscible fully implicit mod...
Definition: porousmediumflow/immiscible/localresidual.hh:28
NumEqVector computeFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache) const
Evaluate the mass flux over a face of a sub control volume.
Definition: porousmediumflow/immiscible/localresidual.hh:97
NumEqVector computeStorage(const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars) const
Evaluate the rate of change of all conservation quantites (e.g. phase mass) within a sub-control volu...
Definition: porousmediumflow/immiscible/localresidual.hh:63
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
Definition: adapt.hh:17
A helper to deduce a vector with the same size as numbers of equations.