3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
cvfelocalresidual.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_CVFE_LOCAL_RESIDUAL_HH
26#define DUMUX_CVFE_LOCAL_RESIDUAL_HH
27
28#include <dune/common/std/type_traits.hh>
29#include <dune/geometry/type.hh>
30#include <dune/istl/matrix.hh>
31
36
37namespace Dumux::Detail {
38
39template<class Imp, class P, class G, class S, class V>
41 std::declval<Imp>().computeStorage(
42 std::declval<P>(), std::declval<G>(), std::declval<S>(), std::declval<V>(), true
43 )
44);
45
46template<class Imp, class P, class G, class S, class V>
47constexpr inline bool hasTimeInfoInterfaceCVFE()
48{ return Dune::Std::is_detected<TimeInfoInterfaceCVFEDetector, Imp, P, G, S, V>::value; }
49
50template<class Imp>
52 std::declval<Imp>().isOverlapping()
53);
54
55template<class Imp>
56constexpr inline bool hasScvfIsOverlapping()
57{ return Dune::Std::is_detected<SCVFIsOverlappingDetector, Imp>::value; }
58
59} // end namespace Dumux::Detail
60
61
62namespace Dumux {
63
70template<class TypeTag>
71class CVFELocalResidual : public FVLocalResidual<TypeTag>
72{
78 using GridView = typename GridGeometry::GridView;
79 using Element = typename GridView::template Codim<0>::Entity;
81 using FVElementGeometry = typename GridGeometry::LocalView;
83 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
84 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
85 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
86 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
87 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
89 using NumEqVector = Dumux::NumEqVector<PrimaryVariables>;
90 using Extrusion = Extrusion_t<GridGeometry>;
91
92public:
94 using ParentType::ParentType;
95
98 const Problem& problem,
99 const Element& element,
100 const FVElementGeometry& fvGeometry,
101 const ElementVolumeVariables& elemVolVars,
102 const ElementBoundaryTypes& elemBcTypes,
103 const ElementFluxVariablesCache& elemFluxVarsCache,
104 const SubControlVolumeFace& scvf) const
105 {
106 const auto flux = evalFlux(problem, element, fvGeometry, elemVolVars, elemBcTypes, elemFluxVarsCache, scvf);
107 if (!scvf.boundary())
108 {
109 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
110 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
111 residual[insideScv.localDofIndex()] += flux;
112
113 // for control-volume finite element schemes with overlapping control volumes
114 if constexpr (Detail::hasScvfIsOverlapping<SubControlVolumeFace>())
115 {
116 if (!scvf.isOverlapping())
117 residual[outsideScv.localDofIndex()] -= flux;
118 }
119 else
120 residual[outsideScv.localDofIndex()] -= flux;
121 }
122 else
123 {
124 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
125 residual[insideScv.localDofIndex()] += flux;
126 }
127 }
128
130 NumEqVector evalFlux(const Problem& problem,
131 const Element& element,
132 const FVElementGeometry& fvGeometry,
133 const ElementVolumeVariables& elemVolVars,
134 const ElementBoundaryTypes& elemBcTypes,
135 const ElementFluxVariablesCache& elemFluxVarsCache,
136 const SubControlVolumeFace& scvf) const
137 {
138 NumEqVector flux(0.0);
139
140 // inner faces
141 if (!scvf.boundary())
142 flux += this->asImp().computeFlux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
143
144 // boundary faces
145 else
146 {
147 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
148 const auto& bcTypes = elemBcTypes.get(fvGeometry, scv);
149
150 // Treat Neumann and Robin ("solution dependent Neumann") boundary conditions.
151 // For Dirichlet there is no addition to the residual here but they
152 // are enforced strongly by replacing the residual entry afterwards.
153 if (bcTypes.hasNeumann())
154 {
155 auto neumannFluxes = problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
156
157 // multiply neumann fluxes with the area and the extrusion factor
158 neumannFluxes *= Extrusion::area(fvGeometry, scvf)*elemVolVars[scv].extrusionFactor();
159
160 // only add fluxes to equations for which Neumann is set
161 for (int eqIdx = 0; eqIdx < NumEqVector::dimension; ++eqIdx)
162 if (bcTypes.isNeumann(eqIdx))
163 flux[eqIdx] += neumannFluxes[eqIdx];
164 }
165 }
166
167 return flux;
168 }
169
187 const Problem& problem,
188 const Element& element,
189 const FVElementGeometry& fvGeometry,
190 const ElementVolumeVariables& prevElemVolVars,
191 const ElementVolumeVariables& curElemVolVars,
192 const SubControlVolume& scv) const
193 {
194 const auto& curVolVars = curElemVolVars[scv];
195 const auto& prevVolVars = prevElemVolVars[scv];
196
197 // Compute storage with the model specific storage residual
198 // This addresses issues #792/#940 in ad-hoc way by additionally providing crude time level information (previous or current)
199 // to the low-level interfaces if this is supported by the LocalResidual implementation
200 NumEqVector prevStorage = computeStorageImpl_(problem, fvGeometry, scv, prevVolVars, /*previous time level?*/true);
201 NumEqVector storage = computeStorageImpl_(problem, fvGeometry, scv, curVolVars, /*previous time level?*/false);
202
203 prevStorage *= prevVolVars.extrusionFactor();
204 storage *= curVolVars.extrusionFactor();
205
206 storage -= prevStorage;
207 storage *= Extrusion::volume(fvGeometry, scv);
208 storage /= this->timeLoop().timeStepSize();
209
210 residual[scv.localDofIndex()] += storage;
211 }
212
213private:
214 NumEqVector computeStorageImpl_(const Problem& problem,
215 const FVElementGeometry& fvGeometry,
216 const SubControlVolume& scv,
217 const VolumeVariables& volVars,
218 [[maybe_unused]] bool isPreviousTimeStep) const
219 {
220 if constexpr (Detail::hasTimeInfoInterfaceCVFE<Implementation, Problem, FVElementGeometry, SubControlVolume, VolumeVariables>())
221 return this->asImp().computeStorage(problem, fvGeometry, scv, volVars, isPreviousTimeStep);
222 else
223 return this->asImp().computeStorage(problem, scv, volVars);
224 }
225};
226
227} // end namespace Dumux
228
229#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
Distance implementation details.
Definition: cvfelocalresidual.hh:37
constexpr bool hasScvfIsOverlapping()
Definition: cvfelocalresidual.hh:56
constexpr bool hasTimeInfoInterfaceCVFE()
Definition: cvfelocalresidual.hh:47
decltype(std::declval< Imp >().computeStorage(std::declval< P >(), std::declval< G >(), std::declval< S >(), std::declval< V >(), true)) TimeInfoInterfaceCVFEDetector
Definition: cvfelocalresidual.hh:44
decltype(std::declval< Imp >().isOverlapping()) SCVFIsOverlappingDetector
Definition: cvfelocalresidual.hh:53
The element-wise residual for control-volume finite element schemes.
Definition: cvfelocalresidual.hh:72
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: cvfelocalresidual.hh:130
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: cvfelocalresidual.hh:186
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: cvfelocalresidual.hh:97
typename ParentType::ElementResidualVector ElementResidualVector
Definition: cvfelocalresidual.hh:93
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.