version 3.11-dev
Loading...
Searching...
No Matches
assembly/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-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_BASE_LOCAL_RESIDUAL_HH
13#define DUMUX_BASE_LOCAL_RESIDUAL_HH
14
15#include <cassert>
16
17#include <dune/common/exceptions.hh>
18
19#include <dumux/common/typetraits/localdofs_.hh>
28
29namespace Dumux::Experimental {
30
37template<class TypeTag>
39{
44 using Element = typename GridView::template Codim<0>::Entity;
45 using ElementDiscretization = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
48 using Extrusion = Extrusion_t<GridGeometry>;
50 using GridVariablesCache = typename GridVariables::GridVariablesCache;
51 using ElementVariables = typename GridVariablesCache::LocalView;
52 using TimeLoop = TimeLoopBase<Scalar>;
53
54public:
57
59 LocalResidual(const Problem* problem,
60 const TimeLoop* timeLoop = nullptr)
61 : problem_(problem)
62 , timeLoop_(timeLoop)
63 {}
64
69 // \{
70
81 ElementResidualVector evalStorage(const Element& element,
82 const ElementDiscretization& elemDisc,
83 const ElementVariables& prevElemVars,
84 const ElementVariables& curElemVars) const
85 {
86 assert(!this->isStationary() && "no time loop set for storage term evaluation");
87
88 // initialize the residual vector for all local dofs in this element
89 ElementResidualVector residual(Dumux::Detail::LocalDofs::numLocalDofs(elemDisc));
90
91 // evaluate the volume terms (storage + source terms)
92 // forward to the local residual specialized for the discretization methods
94 {
95 for (const auto& scv : scvs(elemDisc))
96 this->asImp().evalStorage(residual, this->problem(), element, elemDisc, prevElemVars, curElemVars, scv);
97 }
98
99 // allow for additional contributions (e.g. hybrid CVFE / FE schemes)
100 this->asImp().addToElementStorageResidual(residual, this->problem(), element, elemDisc, prevElemVars, curElemVars);
101
102 return residual;
103 }
104
114 const ElementDiscretization& elemDisc,
115 const ElementVariables& elemVars) const
116 {
117 // initialize the residual vector for all local dofs in this element
118 ElementResidualVector residual(Dumux::Detail::LocalDofs::numLocalDofs(elemDisc));
119
121 {
122 // evaluate the volume terms (storage + source terms)
123 // forward to the local residual specialized for the discretization methods
124 for (const auto& scv : scvs(elemDisc))
125 this->asImp().evalSource(residual, this->problem(), element, elemDisc, elemVars, scv);
126
127 // forward to the local residual specialized for the discretization methods
128 // TODO: Provide scvfs range only over interior scvfs
129 for (auto&& scvf : scvfs(elemDisc))
130 if(!scvf.boundary())
131 this->asImp().evalFlux(residual, this->problem(), element, elemDisc, elemVars, scvf);
132 }
133
134 // allow for additional contributions (e.g. hybrid CVFE / FE schemes)
135 this->asImp().addToElementFluxAndSourceResidual(residual, this->problem(), element, elemDisc, elemVars);
136
137 // add boundary flux contributions
138 this->asImp().addBoundaryFluxIntegral(residual, this->problem(), elemDisc, elemVars);
139
140 return residual;
141 }
142
145 const Problem& problem,
146 const Element& element,
147 const ElementDiscretization& elemDisc,
148 const ElementVariables& prevElemVars,
149 const ElementVariables& curElemVars) const
150 {}
151
154 const Problem& problem,
155 const Element& element,
156 const ElementDiscretization& elemDisc,
157 const ElementVariables& curElemVars) const
158 {}
159
162 const Problem& problem,
163 const ElementDiscretization& elemDisc,
164 const ElementVariables& elemVars) const
165 {
166 if(!elemDisc.hasBoundaryFaces())
167 return;
168
169 for (const auto& boundaryFace : boundaryFaces(elemDisc))
170 {
171 const auto& bcTypes = problem.boundaryTypes(elemDisc, boundaryFace);
172 if(!bcTypes.hasFluxBoundary())
173 continue;
174
175 problem.addFEBoundaryFluxIntegral(residual, elemDisc, elemVars, boundaryFace, bcTypes);
176
177 for(const auto& scvf : scvfs(elemDisc, boundaryFace))
178 problem.addFVBoundaryFluxIntegral(residual, elemDisc, elemVars, scvf, bcTypes);
179 }
180 }
181
182 // \}
183
188 // \{
189
200 template<class SubControlVolume>
201 NumEqVector storageIntegral(const ElementDiscretization& elemDisc,
202 const ElementVariables& elemVars,
203 const SubControlVolume& scv,
204 bool isPreviousTimeLevel) const
205 {
206 DUNE_THROW(Dune::NotImplemented, "This model does not implement a storageIntegral method!");
207 }
208
219 template<class SubControlVolume>
220 NumEqVector sourceIntegral(const ElementDiscretization& elemDisc,
221 const ElementVariables& elemVars,
222 const SubControlVolume& scv) const
223 {
224 NumEqVector source(0.0);
225
226 const auto& problem = this->asImp().problem();
227 for (const auto& qpData : CVFE::quadratureRule(elemDisc, scv))
228 source += qpData.weight() * problem.source(elemDisc, elemVars, qpData.ipData());
229
230 source *= elemVars[scv].extrusionFactor();
231
232 // add contribution from possible point sources
233 const auto& pointSources = problem.pointSources();
234 if (!pointSources.empty())
235 for (const auto& context : pointSources.contexts(elemDisc, scv))
236 {
237 auto psValues = pointSources.eval(elemDisc, elemVars, context);
238 source += psValues;
239 }
240
241 return source;
242 }
243
254 template<class SubControlVolumeFace>
255 NumEqVector fluxIntegral(const ElementDiscretization& elemDisc,
256 const ElementVariables& elemVars,
257 const SubControlVolumeFace& scvf) const
258 {
259 DUNE_THROW(Dune::NotImplemented, "This model does not implement a fluxIntegral method!");
260 }
261
262 // \}
263
268 // \{
269
283 template<class SubControlVolume>
285 const Problem& problem,
286 const Element& element,
287 const ElementDiscretization& elemDisc,
288 const ElementVariables& prevElemVars,
289 const ElementVariables& curElemVars,
290 const SubControlVolume& scv) const
291 {
292 // mass balance within the element. this is the
293 // \f$\frac{\partial S}{\partial t}\f$ term if using implicit or explicit
294 // euler as time discretization.
295 //
296 // TODO: We might need a more explicit way for
297 // doing the time discretization...
298
300 NumEqVector prevStorage = this->asImp().storageIntegral(elemDisc, prevElemVars, scv, /*previous time level?*/true);
301 NumEqVector storage = this->asImp().storageIntegral(elemDisc, curElemVars, scv, /*previous time level?*/false);
302
303 storage -= prevStorage;
304 storage /= timeLoop_->timeStepSize();
305
306 residual[scv.localDofIndex()] += storage;
307 }
308
321 template<class SubControlVolume>
323 const Problem& problem,
324 const Element& element,
325 const ElementDiscretization& elemDisc,
326 const ElementVariables& curElemVars,
327 const SubControlVolume& scv) const
328 {
330 NumEqVector source = this->asImp().sourceIntegral(elemDisc, curElemVars, scv);
332 residual[scv.localDofIndex()] -= source;
333 }
334
345 template<class SubControlVolumeFace>
347 const Problem& problem,
348 const Element& element,
349 const ElementDiscretization& elemDisc,
350 const ElementVariables& elemVars,
351 const SubControlVolumeFace& scvf) const {}
352
363 template<class SubControlVolumeFace>
364 NumEqVector evalFlux(const Problem& problem,
365 const Element& element,
366 const ElementDiscretization& elemDisc,
367 const ElementVariables& elemVars,
368 const SubControlVolumeFace& scvf) const
369 {
370 return asImp().evalFlux(problem, element, elemDisc, elemVars, scvf);
371 }
372
373 //\}
374
378 // \{
379
381 const Problem& problem() const
382 { return *problem_; }
383
386 const TimeLoop& timeLoop() const
387 { return *timeLoop_; }
388
390 bool isStationary() const
391 { return !timeLoop_; }
392
393 // \}
394
395protected:
396 Implementation& asImp()
397 { return *static_cast<Implementation*>(this); }
398
399 const Implementation& asImp() const
400 { return *static_cast<const Implementation*>(this); }
401
402private:
403 const Problem* problem_;
404 const TimeLoop* timeLoop_;
405};
406
407} // end namespace Dumux::Experimental
408
409#endif
GVC GridVariablesCache
export type of the grid variables cache
Definition discretization/gridvariables.hh:33
ElementResidualVector evalStorage(const Element &element, const ElementDiscretization &elemDisc, const ElementVariables &prevElemVars, const ElementVariables &curElemVars) const
Compute the storage local residual, i.e. the deviation of the storage term from zero for instationary...
Definition assembly/localresidual.hh:81
LocalResidual(const Problem *problem, const TimeLoop *timeLoop=nullptr)
the constructor
Definition assembly/localresidual.hh:59
void evalSource(ElementResidualVector &residual, const Problem &problem, const Element &element, const ElementDiscretization &elemDisc, const ElementVariables &curElemVars, const SubControlVolume &scv) const
Compute the source local residual, i.e. the deviation of the source term from zero.
Definition assembly/localresidual.hh:322
NumEqVector sourceIntegral(const ElementDiscretization &elemDisc, const ElementVariables &elemVars, const SubControlVolume &scv) const
Calculate the source term integral of the equation.
Definition assembly/localresidual.hh:220
NumEqVector evalFlux(const Problem &problem, const Element &element, const ElementDiscretization &elemDisc, const ElementVariables &elemVars, const SubControlVolumeFace &scvf) const
Compute the fluxes of the local residual.
Definition assembly/localresidual.hh:364
ElementResidualVector evalFluxAndSource(const Element &element, const ElementDiscretization &elemDisc, const ElementVariables &elemVars) const
Compute the flux and source.
Definition assembly/localresidual.hh:113
const Problem & problem() const
the problem
Definition assembly/localresidual.hh:381
void addToElementStorageResidual(ElementResidualVector &residual, const Problem &problem, const Element &element, const ElementDiscretization &elemDisc, const ElementVariables &prevElemVars, const ElementVariables &curElemVars) const
add additional storage contributions (e.g. hybrid CVFE or FE schemes)
Definition assembly/localresidual.hh:144
const TimeLoop & timeLoop() const
Definition assembly/localresidual.hh:386
bool isStationary() const
returns true if the residual is stationary
Definition assembly/localresidual.hh:390
const Implementation & asImp() const
Definition assembly/localresidual.hh:399
void evalFlux(ElementResidualVector &residual, const Problem &problem, const Element &element, const ElementDiscretization &elemDisc, const ElementVariables &elemVars, const SubControlVolumeFace &scvf) const
Compute the fluxes of the local residual.
Definition assembly/localresidual.hh:346
NumEqVector storageIntegral(const ElementDiscretization &elemDisc, const ElementVariables &elemVars, const SubControlVolume &scv, bool isPreviousTimeLevel) const
Calculate the source term integral of the equation.
Definition assembly/localresidual.hh:201
void addBoundaryFluxIntegral(ElementResidualVector &residual, const Problem &problem, const ElementDiscretization &elemDisc, const ElementVariables &elemVars) const
add boundary flux contributions
Definition assembly/localresidual.hh:161
void addToElementFluxAndSourceResidual(ElementResidualVector &residual, const Problem &problem, const Element &element, const ElementDiscretization &elemDisc, const ElementVariables &curElemVars) const
add additional flux and source contributions (e.g. hybrid CVFE or FE schemes)
Definition assembly/localresidual.hh:153
Implementation & asImp()
Definition assembly/localresidual.hh:396
void evalStorage(ElementResidualVector &residual, const Problem &problem, const Element &element, const ElementDiscretization &elemDisc, const ElementVariables &prevElemVars, const ElementVariables &curElemVars, const SubControlVolume &scv) const
Compute the storage local residual, i.e. the deviation of the storage term from zero for instationary...
Definition assembly/localresidual.hh:284
NumEqVector fluxIntegral(const ElementDiscretization &elemDisc, const ElementVariables &elemVars, const SubControlVolumeFace &scvf) const
Calculate the flux integral of the equation.
Definition assembly/localresidual.hh:255
ReservedBlockVector< NumEqVector, Dumux::Detail::LocalDofs::maxNumLocalDofs< ElementDiscretization >()> ElementResidualVector
the container storing all element residuals
Definition assembly/localresidual.hh:56
Base class for all standard finite volume or finite element problems.
Definition common/problem.hh:39
A arithmetic block vector type based on DUNE's reserved vector.
Definition reservedblockvector.hh:26
Manages the handling of time dependent problems.
Definition common/timeloop.hh:84
The default time loop for instationary simulations.
Definition common/timeloop.hh:139
Defines all properties used in Dumux.
Manages the handling of time dependent problems.
Concept for finite-volume discretizations.
Definition concepts.hh:26
Concepts for discretization types.
Helper classes to compute the integration elements.
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
auto quadratureRule(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolume &scv, QuadratureRules::MidpointQuadrature)
Midpoint quadrature for scv.
Definition quadraturerules.hh:159
Definition assembly/assembler.hh:44
std::ranges::range auto scvs(const FVElementGeometry &fvGeometry, const LocalDof &localDof)
Definition localdof.hh:82
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition extrusion.hh:236
A helper to deduce a vector with the same size as numbers of equations.
Point source types and helpers for handling point sources.
Quadrature rules over sub-control volumes and sub-control volume faces.
A arithmetic block vector type based on DUNE's reserved vector.