3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_IMMISCIBLE_LOCAL_RESIDUAL_HH
26#define DUMUX_IMMISCIBLE_LOCAL_RESIDUAL_HH
27
30
31namespace Dumux {
32
38template<class TypeTag>
39class ImmiscibleLocalResidual : public GetPropType<TypeTag, Properties::BaseLocalResidual>
40{
46 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
48 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
49 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
50 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
51 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
53 using Element = typename GridView::template Codim<0>::Entity;
55
57 static constexpr int numPhases = ModelTraits::numFluidPhases();
58 static constexpr int conti0EqIdx = ModelTraits::Indices::conti0EqIdx;
59
60public:
61 using ParentType::ParentType;
62
75 NumEqVector computeStorage(const Problem& problem,
76 const SubControlVolume& scv,
77 const VolumeVariables& volVars) const
78 {
79 // partial time derivative of the phase mass
80 NumEqVector storage;
81 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
82 {
83 auto eqIdx = conti0EqIdx + phaseIdx;
84 storage[eqIdx] = volVars.porosity()
85 * volVars.density(phaseIdx)
86 * volVars.saturation(phaseIdx);
87
89 EnergyLocalResidual::fluidPhaseStorage(storage, scv, volVars, phaseIdx);
90 }
91
93 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
94
95 return storage;
96 }
97
98
109 NumEqVector computeFlux(const Problem& problem,
110 const Element& element,
111 const FVElementGeometry& fvGeometry,
112 const ElementVolumeVariables& elemVolVars,
113 const SubControlVolumeFace& scvf,
114 const ElementFluxVariablesCache& elemFluxVarsCache) const
115 {
116 FluxVariables fluxVars;
117 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
118
119 NumEqVector flux;
120 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
121 {
122 // the physical quantities for which we perform upwinding
123 auto upwindTerm = [phaseIdx](const auto& volVars)
124 { return volVars.density(phaseIdx)*volVars.mobility(phaseIdx); };
125
126 auto eqIdx = conti0EqIdx + phaseIdx;
127 flux[eqIdx] = fluxVars.advectiveFlux(phaseIdx, upwindTerm);
128
130 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
131 }
132
134 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
135
136 return flux;
137 }
138};
139
140} // end namespace Dumux
141
142#endif
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
Definition: adapt.hh:29
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.