version 3.10-dev
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// 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_FACECENTERED_LOCAL_RESIDUAL_HH
14#define DUMUX_FACECENTERED_LOCAL_RESIDUAL_HH
15
16#include <dune/geometry/type.hh>
17#include <dune/istl/matrix.hh>
18
23
24namespace Dumux {
25
32template<class TypeTag>
34{
39 using GridView = typename GridGeometry::GridView;
40 using Element = typename GridView::template Codim<0>::Entity;
42 using FVElementGeometry = typename GridGeometry::LocalView;
43 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
44 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
45 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
46 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
48 using Extrusion = Extrusion_t<GridGeometry>;
49
50public:
52 using ParentType::ParentType;
53
59 const Problem& problem,
60 const Element& element,
61 const FVElementGeometry& fvGeometry,
62 const ElementVolumeVariables& elemVolVars,
63 const ElementBoundaryTypes& elemBcTypes,
64 const ElementFluxVariablesCache& elemFluxVarsCache,
65 const SubControlVolumeFace& scvf) const
66 {
67 const auto flux = evalFlux(problem, element, fvGeometry, elemVolVars, elemBcTypes, elemFluxVarsCache, scvf);
68 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
69 residual[insideScv.localDofIndex()] += flux;
70 }
71
73 NumEqVector evalFlux(const Problem& problem,
74 const Element& element,
75 const FVElementGeometry& fvGeometry,
76 const ElementVolumeVariables& elemVolVars,
77 const ElementBoundaryTypes& elemBcTypes,
78 const ElementFluxVariablesCache& elemFluxVarsCache,
79 const SubControlVolumeFace& scvf) const
80 {
81 if (elemBcTypes.hasDirichlet())
82 return this->asImp().maybeHandleDirichletBoundary(problem, element, fvGeometry, elemVolVars, elemBcTypes, elemFluxVarsCache, scvf);
83
84 if (elemBcTypes.hasNeumann())
85 return this->asImp().maybeHandleNeumannBoundary(problem, element, fvGeometry, elemVolVars, elemBcTypes, elemFluxVarsCache, scvf);
86
87 return this->asImp().computeFlux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache, elemBcTypes);
88 }
89
103 NumEqVector computeSource(const Problem& problem,
104 const Element& element,
105 const FVElementGeometry& fvGeometry,
106 const ElementVolumeVariables& elemVolVars,
107 const SubControlVolume& scv) const
108 {
109 NumEqVector source(0.0);
110
111 // add contributions from volume flux sources
112 source += problem.source(element, fvGeometry, elemVolVars, scv)[scv.dofAxis()];
113
114 // add contribution from possible point sources
115 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv)[scv.dofAxis()];
116
117 return source;
118 }
119
121
138 const Problem& problem,
139 const Element& element,
140 const FVElementGeometry& fvGeometry,
141 const ElementVolumeVariables& prevElemVolVars,
142 const ElementVolumeVariables& curElemVolVars,
143 const SubControlVolume& scv) const
144 {
145 const auto& curVolVars = curElemVolVars[scv];
146 const auto& prevVolVars = prevElemVolVars[scv];
147
148 // mass balance within the element. this is the
149 // \f$\frac{m}{\partial t}\f$ term if using implicit or explicit
150 // euler as time discretization.
151 //
152
154 NumEqVector prevStorage = this->asImp().computeStorage(problem, scv, prevVolVars, true/*isPreviousStorage*/);
155 NumEqVector storage = this->asImp().computeStorage(problem, scv, curVolVars, false/*isPreviousStorage*/);
156
157 prevStorage *= prevVolVars.extrusionFactor();
158 storage *= curVolVars.extrusionFactor();
159
160 storage -= prevStorage;
161 storage *= Extrusion::volume(fvGeometry, scv);
162 storage /= this->timeLoop().timeStepSize();
163
164 residual[scv.localDofIndex()] += storage;
165 }
166};
167
168} // end namespace Dumux
169
170#endif
The element-wise residual for finite volume schemes.
Definition: fvlocalresidual.hh:35
Implementation & asImp()
Definition: fvlocalresidual.hh:488
const Problem & problem() const
the problem
Definition: fvlocalresidual.hh:473
ReservedBlockVector< NumEqVector, FVElementGeometry::maxNumElementScvs > ElementResidualVector
the container storing all element residuals
Definition: fvlocalresidual.hh:57
const TimeLoop & timeLoop() const
Definition: fvlocalresidual.hh:478
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:86
The element-wise residual for the box scheme.
Definition: fclocalresidual.hh:34
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:137
typename ParentType::ElementResidualVector ElementResidualVector
Definition: fclocalresidual.hh:51
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:58
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:73
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:103
virtual Scalar timeStepSize() const =0
Returns the suggested time step length .
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
The element-wise residual for finite volume schemes.
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
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:159
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
A helper to deduce a vector with the same size as numbers of equations.