version 3.8
assembly/fvlocalassemblerbase.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_ASSEMBLER_BASE_HH
13#define DUMUX_FV_LOCAL_ASSEMBLER_BASE_HH
14
15#include <dune/common/reservedvector.hh>
16#include <dune/grid/common/gridenums.hh> // for GhostEntity
17#include <dune/istl/matrixindexset.hh>
18
23
24namespace Dumux {
25
34template<class TypeTag, class Assembler, class Implementation, bool useImplicitAssembly>
36{
42 using SolutionVector = typename Assembler::SolutionVector;
44 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
45 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
46 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
48 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
50 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
51 using Element = typename GridView::template Codim<0>::Entity;
52 static constexpr auto numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
53
54public:
55 using LocalResidual = std::decay_t<decltype(std::declval<Assembler>().localResidual())>;
56 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
57
61 explicit FVLocalAssemblerBase(const Assembler& assembler,
62 const Element& element,
63 const SolutionVector& curSol)
65 element,
66 curSol,
67 localView(assembler.gridGeometry()),
68 localView(assembler.gridVariables().curGridVolVars()),
69 localView(assembler.gridVariables().prevGridVolVars()),
70 localView(assembler.gridVariables().gridFluxVarsCache()),
72 element.partitionType() == Dune::GhostEntity)
73 {}
74
78 explicit FVLocalAssemblerBase(const Assembler& assembler,
79 const Element& element,
80 const SolutionVector& curSol,
81 const FVElementGeometry& fvGeometry,
82 const ElementVolumeVariables& curElemVolVars,
83 const ElementVolumeVariables& prevElemVolVars,
84 const ElementFluxVariablesCache& elemFluxVarsCache,
86 const bool elementIsGhost)
87 : assembler_(assembler)
88 , element_(element)
89 , curSol_(curSol)
90 , fvGeometry_(fvGeometry)
91 , curElemVolVars_(curElemVolVars)
92 , prevElemVolVars_(prevElemVolVars)
93 , elemFluxVarsCache_(elemFluxVarsCache)
94 , localResidual_(localResidual)
95 , elementIsGhost_(elementIsGhost)
96 {}
97
101 static constexpr bool isImplicit()
102 { return useImplicitAssembly; }
103
109 {
110 if (!isImplicit())
111 if (this->assembler().isStationaryProblem())
112 DUNE_THROW(Dune::InvalidStateException, "Using explicit jacobian assembler with stationary local residual");
113
114 if (elementIsGhost())
115 return ElementResidualVector(0.0);
116
119 }
120
125 ElementResidualVector evalLocalResidual(const ElementVolumeVariables& elemVolVars) const
126 {
127 if (!assembler().isStationaryProblem())
128 {
130 residual += evalLocalStorageResidual();
131 return residual;
132 }
133 else
134 return evalLocalFluxAndSourceResidual(elemVolVars);
135 }
136
143 {
146 }
147
154 ElementResidualVector evalLocalFluxAndSourceResidual(const ElementVolumeVariables& elemVolVars) const
155 {
156 return localResidual_.evalFluxAndSource(element_, fvGeometry_, elemVolVars, elemFluxVarsCache_, elemBcTypes_);
157 }
158
165 {
166 return localResidual_.evalStorage(element_, fvGeometry_, prevElemVolVars_, curElemVolVars_);
167 }
168
174 {
175 // get some references for convenience
176 const auto& element = this->element();
177 const auto& curSol = this->curSol();
178 const auto& prevSol = this->assembler().prevSol();
179 auto&& fvGeometry = this->fvGeometry();
180 auto&& curElemVolVars = this->curElemVolVars();
181 auto&& prevElemVolVars = this->prevElemVolVars();
182 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
183
184 // bind the caches
185 fvGeometry.bind(element);
186
187 if (isImplicit())
188 {
191 if (!this->assembler().isStationaryProblem())
192 prevElemVolVars.bindElement(element, fvGeometry, this->assembler().prevSol());
193 }
194 else
195 {
197 prevElemVolVars.bind(element, fvGeometry, prevSol);
199 }
200 }
201
205 template<typename ApplyFunction, class P = Problem, typename std::enable_if_t<P::enableInternalDirichletConstraints(), int> = 0>
206 void enforceInternalDirichletConstraints(const ApplyFunction& applyDirichlet)
207 {
208 // enforce Dirichlet constraints strongly by overwriting partial derivatives with 1 or 0
209 // and set the residual to (privar - dirichletvalue)
210 for (const auto& scvI : scvs(this->fvGeometry()))
211 {
212 const auto internalDirichletConstraints = asImp_().problem().hasInternalDirichletConstraint(this->element(), scvI);
213 if (internalDirichletConstraints.any())
214 {
215 const auto dirichletValues = asImp_().problem().internalDirichlet(this->element(), scvI);
216 // set the Dirichlet conditions in residual and jacobian
217 for (int eqIdx = 0; eqIdx < internalDirichletConstraints.size(); ++eqIdx)
218 if (internalDirichletConstraints[eqIdx])
219 applyDirichlet(scvI, dirichletValues, eqIdx, eqIdx);
220 }
221 }
222 }
223
224 template<typename ApplyFunction, class P = Problem, typename std::enable_if_t<!P::enableInternalDirichletConstraints(), int> = 0>
225 void enforceInternalDirichletConstraints(const ApplyFunction& applyDirichlet)
226 {}
227
229 const Problem& problem() const
230 { return assembler_.problem(); }
231
233 const Assembler& assembler() const
234 { return assembler_; }
235
237 const Element& element() const
238 { return element_; }
239
241 bool elementIsGhost() const
242 { return elementIsGhost_; }
243
245 const SolutionVector& curSol() const
246 { return curSol_; }
247
249 FVElementGeometry& fvGeometry()
250 { return fvGeometry_; }
251
253 ElementVolumeVariables& curElemVolVars()
254 { return curElemVolVars_; }
255
257 ElementVolumeVariables& prevElemVolVars()
258 { return prevElemVolVars_; }
259
261 ElementFluxVariablesCache& elemFluxVarsCache()
262 { return elemFluxVarsCache_; }
263
266 { return localResidual_; }
267
269 ElementBoundaryTypes& elemBcTypes()
270 { return elemBcTypes_; }
271
273 const FVElementGeometry& fvGeometry() const
274 { return fvGeometry_; }
275
277 const ElementVolumeVariables& curElemVolVars() const
278 { return curElemVolVars_; }
279
281 const ElementVolumeVariables& prevElemVolVars() const
282 { return prevElemVolVars_; }
283
285 const ElementFluxVariablesCache& elemFluxVarsCache() const
286 { return elemFluxVarsCache_; }
287
289 const ElementBoundaryTypes& elemBcTypes() const
290 { return elemBcTypes_; }
291
294 { return localResidual_; }
295
296protected:
297 Implementation &asImp_()
298 { return *static_cast<Implementation*>(this); }
299
300 const Implementation &asImp_() const
301 { return *static_cast<const Implementation*>(this); }
302
303 template<class T = TypeTag, typename std::enable_if_t<!GetPropType<T, Properties::GridVariables>::GridVolumeVariables::cachingEnabled, int> = 0>
304 VolumeVariables& getVolVarAccess(GridVolumeVariables& gridVolVars, ElementVolumeVariables& elemVolVars, const SubControlVolume& scv)
305 { return elemVolVars[scv]; }
306
307 template<class T = TypeTag, typename std::enable_if_t<GetPropType<T, Properties::GridVariables>::GridVolumeVariables::cachingEnabled, int> = 0>
308 VolumeVariables& getVolVarAccess(GridVolumeVariables& gridVolVars, ElementVolumeVariables& elemVolVars, const SubControlVolume& scv)
309 { return gridVolVars.volVars(scv); }
310
311private:
312
313 const Assembler& assembler_;
314 const Element& element_;
315 const SolutionVector& curSol_;
316
317 FVElementGeometry fvGeometry_;
318 ElementVolumeVariables curElemVolVars_;
319 ElementVolumeVariables prevElemVolVars_;
320 ElementFluxVariablesCache elemFluxVarsCache_;
321 ElementBoundaryTypes elemBcTypes_;
322
323 LocalResidual localResidual_;
324 bool elementIsGhost_;
325};
326
327
328} // end namespace Dumux
329
330#endif
A base class for all local assemblers.
Definition: assembly/fvlocalassemblerbase.hh:36
const ElementVolumeVariables & prevElemVolVars() const
The element volume variables of the provious time step.
Definition: assembly/fvlocalassemblerbase.hh:281
void bindLocalViews()
Convenience function bind and prepare all relevant variables required for the evaluation of the local...
Definition: assembly/fvlocalassemblerbase.hh:173
ElementResidualVector evalLocalResidual(const ElementVolumeVariables &elemVolVars) const
Evaluates the complete local residual for the current element.
Definition: assembly/fvlocalassemblerbase.hh:125
ElementVolumeVariables & curElemVolVars()
The current element volume variables.
Definition: assembly/fvlocalassemblerbase.hh:253
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: assembly/fvlocalassemblerbase.hh:56
ElementBoundaryTypes & elemBcTypes()
The element's boundary types.
Definition: assembly/fvlocalassemblerbase.hh:269
Implementation & asImp_()
Definition: assembly/fvlocalassemblerbase.hh:297
ElementVolumeVariables & prevElemVolVars()
The element volume variables of the provious time step.
Definition: assembly/fvlocalassemblerbase.hh:257
ElementResidualVector evalLocalResidual() const
Convenience function to evaluate the complete local residual for the current element....
Definition: assembly/fvlocalassemblerbase.hh:108
const Problem & problem() const
The problem.
Definition: assembly/fvlocalassemblerbase.hh:229
const ElementBoundaryTypes & elemBcTypes() const
The element's boundary types.
Definition: assembly/fvlocalassemblerbase.hh:289
const ElementFluxVariablesCache & elemFluxVarsCache() const
The element flux variables cache.
Definition: assembly/fvlocalassemblerbase.hh:285
const Implementation & asImp_() const
Definition: assembly/fvlocalassemblerbase.hh:300
FVLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementVolumeVariables &prevElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const LocalResidual &localResidual, const bool elementIsGhost)
The constructor. General version explicitly expecting each argument.
Definition: assembly/fvlocalassemblerbase.hh:78
const ElementVolumeVariables & curElemVolVars() const
The current element volume variables.
Definition: assembly/fvlocalassemblerbase.hh:277
FVLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol)
The constructor. Delegates to the general constructor.
Definition: assembly/fvlocalassemblerbase.hh:61
const LocalResidual & localResidual() const
The local residual for the current element.
Definition: assembly/fvlocalassemblerbase.hh:293
FVElementGeometry & fvGeometry()
The global finite volume geometry.
Definition: assembly/fvlocalassemblerbase.hh:249
ElementResidualVector evalLocalFluxAndSourceResidual(const ElementVolumeVariables &elemVolVars) const
Evaluates the flux and source terms (i.e, the terms without a time derivative) of the local residual ...
Definition: assembly/fvlocalassemblerbase.hh:154
const Assembler & assembler() const
The assembler.
Definition: assembly/fvlocalassemblerbase.hh:233
ElementFluxVariablesCache & elemFluxVarsCache()
The element flux variables cache.
Definition: assembly/fvlocalassemblerbase.hh:261
bool elementIsGhost() const
Returns if element is a ghost entity.
Definition: assembly/fvlocalassemblerbase.hh:241
std::decay_t< decltype(std::declval< Assembler >().localResidual())> LocalResidual
Definition: assembly/fvlocalassemblerbase.hh:55
LocalResidual & localResidual()
The local residual for the current element.
Definition: assembly/fvlocalassemblerbase.hh:265
static constexpr bool isImplicit()
Returns true if the assembler considers implicit assembly.
Definition: assembly/fvlocalassemblerbase.hh:101
ElementResidualVector evalLocalFluxAndSourceResidual() const
Convenience function to evaluate the flux and source terms (i.e, the terms without a time derivative)...
Definition: assembly/fvlocalassemblerbase.hh:142
void enforceInternalDirichletConstraints(const ApplyFunction &applyDirichlet)
Enforces Dirichlet constraints if enabled in the problem.
Definition: assembly/fvlocalassemblerbase.hh:206
const Element & element() const
The current element.
Definition: assembly/fvlocalassemblerbase.hh:237
VolumeVariables & getVolVarAccess(GridVolumeVariables &gridVolVars, ElementVolumeVariables &elemVolVars, const SubControlVolume &scv)
Definition: assembly/fvlocalassemblerbase.hh:304
const FVElementGeometry & fvGeometry() const
The finite volume geometry.
Definition: assembly/fvlocalassemblerbase.hh:273
ElementResidualVector evalLocalStorageResidual() const
Convenience function to evaluate storage term (i.e, the term with a time derivative) of the local res...
Definition: assembly/fvlocalassemblerbase.hh:164
const SolutionVector & curSol() const
The current solution.
Definition: assembly/fvlocalassemblerbase.hh:245
Defines all properties used in Dumux.
An enum class to define various differentiation methods available in order to compute the derivatives...
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Definition: adapt.hh:17
Definition: common/pdesolver.hh:24
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A arithmetic block vector type based on DUNE's reserved vector.