version 3.8
fvlocalresidual.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//
12#ifndef DUMUX_FV_LOCAL_RESIDUAL_HH
13#define DUMUX_FV_LOCAL_RESIDUAL_HH
14
15#include <dune/common/exceptions.hh>
16#include <dune/istl/bvector.hh>
17
24
25namespace Dumux {
26
33template<class TypeTag>
35{
40 using Element = typename GridView::template Codim<0>::Entity;
41 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
44 using SubControlVolume = typename GridGeometry::SubControlVolume;
45 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
46 using Extrusion = Extrusion_t<GridGeometry>;
49 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
51 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
54
55public:
58
60 FVLocalResidual(const Problem* problem,
61 const TimeLoop* timeLoop = nullptr)
62 : problem_(problem)
63 , timeLoop_(timeLoop)
64 {}
65
71 // \{
72
87 const Element &element,
88 const GridGeometry& gridGeometry,
89 const GridVariables& gridVariables,
90 const SolutionVector& sol) const
91 {
92 // make sure FVElementGeometry and volume variables are bound to the element
93 const auto fvGeometry = localView(gridGeometry).bind(element);
94 const auto elemVolVars = localView(gridVariables.curGridVolVars()).bind(element, fvGeometry, sol);
95
96 ElementResidualVector storage(fvGeometry.numScv());
97
98 // calculate the amount of conservation each quantity inside
99 // all sub control volumes
100 for (auto&& scv : scvs(fvGeometry))
101 {
102 auto localScvIdx = scv.localDofIndex();
103 const auto& volVars = elemVolVars[scv];
104 storage[localScvIdx] = asImp().computeStorage(problem, scv, volVars);
105 storage[localScvIdx] *= Extrusion::volume(fvGeometry, scv) * volVars.extrusionFactor();
106 }
107
108 return storage;
109 }
110
111 // \}
112
117 // \{
118
131 ElementResidualVector evalStorage(const Element& element,
132 const FVElementGeometry& fvGeometry,
133 const ElementVolumeVariables& prevElemVolVars,
134 const ElementVolumeVariables& curElemVolVars) const
135 {
136 assert(timeLoop_ && "no time loop set for storage term evaluation");
137
138 // initialize the residual vector for all scvs in this element
139 ElementResidualVector residual(fvGeometry.numScv());
140
141 // evaluate the volume terms (storage + source terms)
142 // forward to the local residual specialized for the discretization methods
143 for (auto&& scv : scvs(fvGeometry))
144 asImp().evalStorage(residual, this->problem(), element, fvGeometry, prevElemVolVars, curElemVolVars, scv);
145
146 return residual;
147 }
148
150 const FVElementGeometry& fvGeometry,
151 const ElementVolumeVariables& elemVolVars,
152 const ElementFluxVariablesCache& elemFluxVarsCache,
153 const ElementBoundaryTypes &bcTypes) const
154 {
155 // initialize the residual vector for all scvs in this element
156 ElementResidualVector residual(fvGeometry.numScv());
157
158 // evaluate the volume terms (storage + source terms)
159 // forward to the local residual specialized for the discretization methods
160 for (auto&& scv : scvs(fvGeometry))
161 asImp().evalSource(residual, this->problem(), element, fvGeometry, elemVolVars, scv);
162
163 // forward to the local residual specialized for the discretization methods
164 for (auto&& scvf : scvfs(fvGeometry))
165 asImp().evalFlux(residual, this->problem(), element, fvGeometry, elemVolVars, bcTypes, elemFluxVarsCache, scvf);
166
167 return residual;
168 }
169
170 // \}
171
172
177 // \{
178
188 NumEqVector computeStorage(const Problem& problem,
189 const SubControlVolume& scv,
190 const VolumeVariables& volVars) const
191 {
192 DUNE_THROW(Dune::NotImplemented, "This model does not implement a storage method!");
193 }
194
208 NumEqVector computeSource(const Problem& problem,
209 const Element& element,
210 const FVElementGeometry& fvGeometry,
211 const ElementVolumeVariables& elemVolVars,
212 const SubControlVolume &scv) const
213 {
214 NumEqVector source(0.0);
215
216 // add contributions from volume flux sources
217 source += problem.source(element, fvGeometry, elemVolVars, scv);
218
219 // add contribution from possible point sources
220 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
221
222 return source;
223 }
224
239 NumEqVector computeFlux(const Problem& problem,
240 const Element& element,
241 const FVElementGeometry& fvGeometry,
242 const ElementVolumeVariables& elemVolVars,
243 const SubControlVolumeFace& scvf,
244 const ElementFluxVariablesCache& elemFluxVarsCache) const
245 {
246 DUNE_THROW(Dune::NotImplemented, "This model does not implement a flux method!");
247 }
248
249 // \}
250
255 // \{
256
273 const Problem& problem,
274 const Element& element,
275 const FVElementGeometry& fvGeometry,
276 const ElementVolumeVariables& prevElemVolVars,
277 const ElementVolumeVariables& curElemVolVars,
278 const SubControlVolume& scv) const
279 {
280 const auto& curVolVars = curElemVolVars[scv];
281 const auto& prevVolVars = prevElemVolVars[scv];
282
283 // mass balance within the element. this is the
284 // \f$\frac{m}{\partial t}\f$ term if using implicit or explicit
285 // euler as time discretization.
286 //
287 // TODO: We might need a more explicit way for
288 // doing the time discretization...
289
291 NumEqVector prevStorage = asImp().computeStorage(problem, scv, prevVolVars);
292 NumEqVector storage = asImp().computeStorage(problem, scv, curVolVars);
293
294 prevStorage *= prevVolVars.extrusionFactor();
295 storage *= curVolVars.extrusionFactor();
296
297 storage -= prevStorage;
298 storage *= Extrusion::volume(fvGeometry, scv);
299 storage /= timeLoop_->timeStepSize();
300
301 residual[scv.localDofIndex()] += storage;
302 }
303
318 const Problem& problem,
319 const Element& element,
320 const FVElementGeometry& fvGeometry,
321 const ElementVolumeVariables& curElemVolVars,
322 const SubControlVolume& scv) const
323 {
325 const auto& curVolVars = curElemVolVars[scv];
326 NumEqVector source = asImp().computeSource(problem, element, fvGeometry, curElemVolVars, scv);
327 source *= Extrusion::volume(fvGeometry, scv)*curVolVars.extrusionFactor();
328
330 residual[scv.localDofIndex()] -= source;
331 }
332
348 const Problem& problem,
349 const Element& element,
350 const FVElementGeometry& fvGeometry,
351 const ElementVolumeVariables& elemVolVars,
352 const ElementBoundaryTypes& elemBcTypes,
353 const ElementFluxVariablesCache& elemFluxVarsCache,
354 const SubControlVolumeFace& scvf) const {}
355
369 NumEqVector evalFlux(const Problem& problem,
370 const Element& element,
371 const FVElementGeometry& fvGeometry,
372 const ElementVolumeVariables& elemVolVars,
373 const ElementFluxVariablesCache& elemFluxVarsCache,
374 const SubControlVolumeFace& scvf) const
375 {
376 return asImp().evalFlux(problem, element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
377 }
378
379 //\}
380
384 // \{
385
387 template<class PartialDerivativeMatrix>
388 void addStorageDerivatives(PartialDerivativeMatrix& partialDerivatives,
389 const Problem& problem,
390 const Element& element,
391 const FVElementGeometry& fvGeometry,
392 const VolumeVariables& curVolVars,
393 const SubControlVolume& scv) const
394 {
395 DUNE_THROW(Dune::NotImplemented, "analytic storage derivative");
396 }
397
399 template<class PartialDerivativeMatrix>
400 void addSourceDerivatives(PartialDerivativeMatrix& partialDerivatives,
401 const Problem& problem,
402 const Element& element,
403 const FVElementGeometry& fvGeometry,
404 const VolumeVariables& curVolVars,
405 const SubControlVolume& scv) const
406 {
407 DUNE_THROW(Dune::NotImplemented, "analytic source derivative");
408 }
409
411 template<class PartialDerivativeMatrices, class T = TypeTag>
412 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod != DiscretizationMethods::box, void>
413 addFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
414 const Problem& problem,
415 const Element& element,
416 const FVElementGeometry& fvGeometry,
417 const ElementVolumeVariables& curElemVolVars,
418 const ElementFluxVariablesCache& elemFluxVarsCache,
419 const SubControlVolumeFace& scvf) const
420 {
421 DUNE_THROW(Dune::NotImplemented, "analytic flux derivative for cell-centered models");
422 }
423
425 template<class JacobianMatrix, class T = TypeTag>
426 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethods::box, void>
427 addFluxDerivatives(JacobianMatrix& A,
428 const Problem& problem,
429 const Element& element,
430 const FVElementGeometry& fvGeometry,
431 const ElementVolumeVariables& curElemVolVars,
432 const ElementFluxVariablesCache& elemFluxVarsCache,
433 const SubControlVolumeFace& scvf) const
434 {
435 DUNE_THROW(Dune::NotImplemented, "analytic flux derivative for box models");
436 }
437
439 template<class PartialDerivativeMatrices>
440 void addCCDirichletFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
441 const Problem& problem,
442 const Element& element,
443 const FVElementGeometry& fvGeometry,
444 const ElementVolumeVariables& curElemVolVars,
445 const ElementFluxVariablesCache& elemFluxVarsCache,
446 const SubControlVolumeFace& scvf) const
447 {
448 DUNE_THROW(Dune::NotImplemented, "analytic Dirichlet flux derivative");
449 }
450
452 template<class PartialDerivativeMatrices>
453 void addRobinFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
454 const Problem& problem,
455 const Element& element,
456 const FVElementGeometry& fvGeometry,
457 const ElementVolumeVariables& curElemVolVars,
458 const ElementFluxVariablesCache& elemFluxVarsCache,
459 const SubControlVolumeFace& scvf) const
460 {
461 DUNE_THROW(Dune::NotImplemented, "analytic Robin flux derivative");
462 }
463
464 //\}
465
469 // \{
470
472 const Problem& problem() const
473 { return *problem_; }
474
477 const TimeLoop& timeLoop() const
478 { return *timeLoop_; }
479
481 bool isStationary() const
482 { return !timeLoop_; }
483
484 // \}
485protected:
486
487 Implementation &asImp()
488 { return *static_cast<Implementation*>(this); }
489
490 const Implementation &asImp() const
491 { return *static_cast<const Implementation*>(this); }
492
493private:
494 const Problem* problem_;
495 const TimeLoop* timeLoop_;
496};
497
498} // end namespace Dumux
499
500#endif
The element-wise residual for finite volume schemes.
Definition: fvlocalresidual.hh:35
Implementation & asImp()
Definition: fvlocalresidual.hh:487
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: fvlocalresidual.hh:208
void addRobinFluxDerivatives(PartialDerivativeMatrices &derivativeMatrices, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Compute the derivative of Robin type boundary conditions ("solution dependent Neumann")
Definition: fvlocalresidual.hh:453
NumEqVector computeStorage(const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars) const
Calculate the source term of the equation.
Definition: fvlocalresidual.hh:188
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
Compute the flux local residual, i.e. the deviation of the flux term from zero.
Definition: fvlocalresidual.hh:347
const Problem & problem() const
the problem
Definition: fvlocalresidual.hh:472
void addSourceDerivatives(PartialDerivativeMatrix &partialDerivatives, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const VolumeVariables &curVolVars, const SubControlVolume &scv) const
Compute the derivative of the source residual.
Definition: fvlocalresidual.hh:400
const TimeLoop & timeLoop() const
Definition: fvlocalresidual.hh:477
void addStorageDerivatives(PartialDerivativeMatrix &partialDerivatives, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const VolumeVariables &curVolVars, const SubControlVolume &scv) const
Compute the derivative of the storage residual.
Definition: fvlocalresidual.hh:388
bool isStationary() const
returns true if the residual is stationary
Definition: fvlocalresidual.hh:481
NumEqVector evalFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Compute the flux local residual, i.e. the deviation of the flux term from zero.
Definition: fvlocalresidual.hh:369
const Implementation & asImp() const
Definition: fvlocalresidual.hh:490
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
FVLocalResidual(const Problem *problem, const TimeLoop *timeLoop=nullptr)
the constructor
Definition: fvlocalresidual.hh:60
ElementResidualVector evalStorage(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &prevElemVolVars, const ElementVolumeVariables &curElemVolVars) const
Compute the storage local residual, i.e. the deviation of the storage term from zero for instationary...
Definition: fvlocalresidual.hh:131
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod !=DiscretizationMethods::box, void > addFluxDerivatives(PartialDerivativeMatrices &derivativeMatrices, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Compute the derivative of the flux residual.
Definition: fvlocalresidual.hh:413
void evalSource(ElementResidualVector &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const SubControlVolume &scv) const
Compute the source local residual, i.e. the deviation of the source term from zero.
Definition: fvlocalresidual.hh:317
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod==DiscretizationMethods::box, void > addFluxDerivatives(JacobianMatrix &A, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Compute the derivative of the flux residual for the box method.
Definition: fvlocalresidual.hh:427
ElementResidualVector evalFluxAndSource(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const ElementBoundaryTypes &bcTypes) const
Definition: fvlocalresidual.hh:149
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: fvlocalresidual.hh:272
void addCCDirichletFluxDerivatives(PartialDerivativeMatrices &derivativeMatrices, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Compute the derivative of the Dirichlet flux residual for cell-centered schemes.
Definition: fvlocalresidual.hh:440
NumEqVector computeFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache) const
Calculate the flux term of the equation.
Definition: fvlocalresidual.hh:239
A arithmetic block vector type based on DUNE's reserved vector.
Definition: reservedblockvector.hh:26
virtual Scalar timeStepSize() const =0
Returns the suggested time step length .
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.
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
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
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
The available discretization methods in Dumux.
constexpr Box box
Definition: method.hh:147
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.
A arithmetic block vector type based on DUNE's reserved vector.