3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
fclocalresidual.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_FACECENTERED_LOCAL_RESIDUAL_HH
26#define DUMUX_FACECENTERED_LOCAL_RESIDUAL_HH
27
28#include <dune/geometry/type.hh>
29#include <dune/istl/matrix.hh>
30
35
36namespace Dumux {
37
44template<class TypeTag>
46{
51 using GridView = typename GridGeometry::GridView;
52 using Element = typename GridView::template Codim<0>::Entity;
54 using FVElementGeometry = typename GridGeometry::LocalView;
55 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
56 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
57 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
58 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
60 using Extrusion = Extrusion_t<GridGeometry>;
61
62public:
64 using ParentType::ParentType;
65
68 const Problem& problem,
69 const Element& element,
70 const FVElementGeometry& fvGeometry,
71 const ElementVolumeVariables& elemVolVars,
72 const ElementBoundaryTypes& elemBcTypes,
73 const ElementFluxVariablesCache& elemFluxVarsCache,
74 const SubControlVolumeFace& scvf) const
75 {
76 const auto flux = evalFlux(problem, element, fvGeometry, elemVolVars, elemBcTypes, elemFluxVarsCache, scvf);
77 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
78 residual[insideScv.localDofIndex()] += flux;
79 }
80
82 NumEqVector evalFlux(const Problem& problem,
83 const Element& element,
84 const FVElementGeometry& fvGeometry,
85 const ElementVolumeVariables& elemVolVars,
86 const ElementBoundaryTypes& elemBcTypes,
87 const ElementFluxVariablesCache& elemFluxVarsCache,
88 const SubControlVolumeFace& scvf) const
89 {
90 if (elemBcTypes.hasDirichlet())
91 return this->asImp().maybeHandleDirichletBoundary(problem, element, fvGeometry, elemVolVars, elemBcTypes, elemFluxVarsCache, scvf);
92
93 if (elemBcTypes.hasNeumann())
94 return this->asImp().maybeHandleNeumannBoundary(problem, element, fvGeometry, elemVolVars, elemBcTypes, elemFluxVarsCache, scvf);
95
96 return this->asImp().computeFlux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache, elemBcTypes);
97 }
98
112 NumEqVector computeSource(const Problem& problem,
113 const Element& element,
114 const FVElementGeometry& fvGeometry,
115 const ElementVolumeVariables& elemVolVars,
116 const SubControlVolume& scv) const
117 {
118 NumEqVector source(0.0);
119
120 // add contributions from volume flux sources
121 source += problem.source(element, fvGeometry, elemVolVars, scv)[scv.dofAxis()];
122
123 // add contribution from possible point sources
124 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv)[scv.dofAxis()];
125
126 return source;
127 }
128
130
147 const Problem& problem,
148 const Element& element,
149 const FVElementGeometry& fvGeometry,
150 const ElementVolumeVariables& prevElemVolVars,
151 const ElementVolumeVariables& curElemVolVars,
152 const SubControlVolume& scv) const
153 {
154 const auto& curVolVars = curElemVolVars[scv];
155 const auto& prevVolVars = prevElemVolVars[scv];
156
157 // mass balance within the element. this is the
158 // \f$\frac{m}{\partial t}\f$ term if using implicit or explicit
159 // euler as time discretization.
160 //
161
163 NumEqVector prevStorage = this->asImp().computeStorage(problem, scv, prevVolVars, true/*isPreviousStorage*/);
164 NumEqVector storage = this->asImp().computeStorage(problem, scv, curVolVars, false/*isPreviousStorage*/);
165
166 prevStorage *= prevVolVars.extrusionFactor();
167 storage *= curVolVars.extrusionFactor();
168
169 storage -= prevStorage;
170 storage *= Extrusion::volume(fvGeometry, scv);
171 storage /= this->timeLoop().timeStepSize();
172
173 residual[scv.localDofIndex()] += storage;
174 }
175};
176
177} // end namespace Dumux
178
179#endif
The element-wise residual for finite volume schemes.
A helper to deduce a vector with the same size as numbers of equations.
Helper classes to compute the integration elements.
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:171
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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
The element-wise residual for the box scheme.
Definition: fclocalresidual.hh:46
void evalStorage(ElementResidualVector &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &prevElemVolVars, const ElementVolumeVariables &curElemVolVars, const SubControlVolume &scv) const
Compute the storage local residual, i.e. the deviation of the storage term from zero for instationary...
Definition: fclocalresidual.hh:146
typename ParentType::ElementResidualVector ElementResidualVector
Definition: fclocalresidual.hh:63
void evalFlux(ElementResidualVector &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementBoundaryTypes &elemBcTypes, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
evaluate flux residuals for one sub control volume face and add to residual
Definition: fclocalresidual.hh:67
NumEqVector evalFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementBoundaryTypes &elemBcTypes, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
evaluate flux residuals for one sub control volume face
Definition: fclocalresidual.hh:82
NumEqVector computeSource(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const
Calculate the source term of the equation.
Definition: fclocalresidual.hh:112
The element-wise residual for finite volume schemes.
Definition: fvlocalresidual.hh:47
Implementation & asImp()
Definition: fvlocalresidual.hh:499
const Problem & problem() const
the problem
Definition: fvlocalresidual.hh:484
ReservedBlockVector< NumEqVector, FVElementGeometry::maxNumElementScvs > ElementResidualVector
the container storing all element residuals
Definition: fvlocalresidual.hh:69
const TimeLoop & timeLoop() const
Definition: fvlocalresidual.hh:489
ElementResidualVector evalStorage(const Problem &problem, const Element &element, const GridGeometry &gridGeometry, const GridVariables &gridVariables, const SolutionVector &sol) const
Compute the storage term for the current solution.
Definition: fvlocalresidual.hh:98
virtual Scalar timeStepSize() const =0
Returns the suggested time step length .
Declares all properties used in Dumux.