version 3.10-dev
experimental/assembly/cvfelocalassembler.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//
14#ifndef DUMUX_EXPERIMENTAL_CVFE_LOCAL_ASSEMBLER_HH
15#define DUMUX_EXPERIMENTAL_CVFE_LOCAL_ASSEMBLER_HH
16
17#include <dune/common/exceptions.hh>
18#include <dune/common/hybridutilities.hh>
19#include <dune/common/reservedvector.hh>
20#include <dune/grid/common/gridenums.hh>
21#include <dune/istl/matrixindexset.hh>
22#include <dune/istl/bvector.hh>
23
28
34
36
37#include <dumux/assembly/volvardeflectionhelper_.hh>
38
39namespace Dumux::Experimental {
40
50template<class TypeTag, class Assembler, class Implementation>
51class CVFELocalAssemblerBase : public Experimental::FVLocalAssemblerBase<TypeTag, Assembler, Implementation>
52{
57 using ElementVolumeVariables = typename GridVariables::GridVolumeVariables::LocalView;
58 using SolutionVector = typename Assembler::SolutionVector;
59
60 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
62
63public:
64
65 using ParentType::ParentType;
66
68 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
69
70
72 {
74 this->elemBcTypes().update(this->asImp_().problem(), this->element(), this->fvGeometry());
75 }
76
81 template <class ResidualVector, class StageParams, class PartialReassembler = DefaultPartialReassembler, class CouplingFunction = Noop>
82 void assembleJacobianAndResidual(JacobianMatrix& jac, ResidualVector& res, GridVariables& gridVariables,
83 const StageParams& stageParams, ResidualVector& temporal, ResidualVector& spatial,
84 ResidualVector& constrainedDofs,
85 const PartialReassembler* partialReassembler = nullptr,
86 const CouplingFunction& maybeAssembleCouplingBlocks = noop)
87 {
88 this->asImp_().bindLocalViews();
89 const auto eIdxGlobal = this->asImp_().problem().gridGeometry().elementMapper().index(this->element());
90
91 this->localResidual().spatialWeight(1.0);
92 this->localResidual().temporalWeight(1.0);
93
94 const auto sWeight = stageParams.spatialWeight(stageParams.size()-1);
95 const auto tWeight = stageParams.temporalWeight(stageParams.size()-1);
96
97 const auto flux = this->evalLocalFluxAndSourceResidual(this->curElemVolVars());
98 const auto storage = this->localResidual().evalStorage(this->fvGeometry(), this->curElemVolVars());
99 ElementResidualVector origResidual(flux.size());
100 origResidual = 0.0;
101 for (const auto& scv : scvs(this->fvGeometry()))
102 {
103 spatial[scv.dofIndex()] += flux[scv.localDofIndex()];
104 temporal[scv.dofIndex()] += storage[scv.localDofIndex()];
105 origResidual[scv.localDofIndex()] += flux[scv.localDofIndex()]*sWeight + storage[scv.localDofIndex()]*tWeight;
106 res[scv.dofIndex()] += origResidual[scv.localDofIndex()];
107 }
108
109 this->localResidual().spatialWeight(sWeight);
110 this->localResidual().temporalWeight(tWeight);
111
112 if (partialReassembler && partialReassembler->elementColor(eIdxGlobal) == EntityColor::green)
113 {
114 // assemble the coupling blocks for coupled models (does nothing if not coupled)
115 maybeAssembleCouplingBlocks(origResidual);
116 }
117 else if (!this->elementIsGhost())
118 {
119 this->asImp_().assembleJacobian(jac, gridVariables, origResidual, partialReassembler); // forward to the internal implementation
120
121 // assemble the coupling blocks for coupled models (does nothing if not coupled)
122 maybeAssembleCouplingBlocks(origResidual);
123 }
124 else
125 {
126 // Treatment of ghost elements
127 assert(this->elementIsGhost());
128
129 // handle dofs per codimension
130 const auto& gridGeometry = this->asImp_().problem().gridGeometry();
131 Dune::Hybrid::forEach(std::make_integer_sequence<int, dim>{}, [&](auto d)
132 {
133 constexpr int codim = dim - d;
134 const auto& localCoeffs = gridGeometry.feCache().get(this->element().type()).localCoefficients();
135 for (int idx = 0; idx < localCoeffs.size(); ++idx)
136 {
137 const auto& localKey = localCoeffs.localKey(idx);
138
139 // skip if we are not handling this codim right now
140 if (localKey.codim() != codim)
141 continue;
142
143 // do not change the non-ghost entities
144 auto entity = this->element().template subEntity<codim>(localKey.subEntity());
145 if (entity.partitionType() == Dune::InteriorEntity || entity.partitionType() == Dune::BorderEntity)
146 continue;
147
148 // WARNING: this only works if the mapping from codim+subEntity to
149 // global dofIndex is unique (on dof per entity of this codim).
150 // For more general mappings, we should use a proper local-global mapping here.
151 // For example through dune-functions.
152 const auto dofIndex = gridGeometry.dofMapper().index(entity);
153
154 // this might be a vector-valued dof
155 using BlockType = typename JacobianMatrix::block_type;
156 BlockType &J = jac[dofIndex][dofIndex];
157 for (int j = 0; j < BlockType::rows; ++j)
158 J[j][j] = 1.0;
159
160 // set residual for the ghost dof
161 res[dofIndex] = 0;
162 constrainedDofs[dofIndex] = 1;
163 }
164 });
165 }
166
167 auto applyDirichlet = [&] (const auto& scvI,
168 const auto& dirichletValues,
169 const auto eqIdx,
170 const auto pvIdx)
171 {
172 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
173 constrainedDofs[scvI.dofIndex()][eqIdx] = 1;
174
175 auto& row = jac[scvI.dofIndex()];
176 for (auto col = row.begin(); col != row.end(); ++col)
177 row[col.index()][eqIdx] = 0.0;
178
179 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
180
181 // if a periodic dof has Dirichlet values also apply the same Dirichlet values to the other dof
182 if (this->asImp_().problem().gridGeometry().dofOnPeriodicBoundary(scvI.dofIndex()))
183 {
184 const auto periodicDof = this->asImp_().problem().gridGeometry().periodicallyMappedDof(scvI.dofIndex());
185 res[periodicDof][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
186 constrainedDofs[periodicDof][eqIdx] = 1;
187 const auto end = jac[periodicDof].end();
188 for (auto it = jac[periodicDof].begin(); it != end; ++it)
189 (*it) = periodicDof != it.index() ? 0.0 : 1.0;
190 }
191 };
192
193 this->asImp_().enforceDirichletConstraints(applyDirichlet);
194 }
195
200 void assembleJacobian(JacobianMatrix& jac, GridVariables& gridVariables)
201 {
202 this->asImp_().bindLocalViews();
203 this->asImp_().assembleJacobian(jac, gridVariables); // forward to the internal implementation
204
205 auto applyDirichlet = [&] (const auto& scvI,
206 const auto& dirichletValues,
207 const auto eqIdx,
208 const auto pvIdx)
209 {
210 auto& row = jac[scvI.dofIndex()];
211 for (auto col = row.begin(); col != row.end(); ++col)
212 row[col.index()][eqIdx] = 0.0;
213
214 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
215 };
216
217 this->asImp_().enforceDirichletConstraints(applyDirichlet);
218 }
219
223 template <class ResidualVector>
224 void assembleResidual(ResidualVector& res)
225 {
226 this->asImp_().bindLocalViews();
227 const auto residual = this->evalLocalResidual();
228
229 for (const auto& scv : scvs(this->fvGeometry()))
230 res[scv.dofIndex()] += residual[scv.localDofIndex()];
231
232 auto applyDirichlet = [&] (const auto& scvI,
233 const auto& dirichletValues,
234 const auto eqIdx,
235 const auto pvIdx)
236 {
237 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
238 };
239
240 this->asImp_().enforceDirichletConstraints(applyDirichlet);
241 }
242
246 template<class ResidualVector>
247 void assembleCurrentResidual(ResidualVector& spatialRes, ResidualVector& temporalRes)
248 {
249 this->asImp_().bindLocalViews();
250 const auto flux = this->evalLocalFluxAndSourceResidual(this->curElemVolVars());
251 const auto storage = this->localResidual().evalStorage(this->fvGeometry(), this->curElemVolVars());
252 for (const auto& scv : scvs(this->fvGeometry()))
253 {
254 spatialRes[scv.dofIndex()] += flux[scv.localDofIndex()];
255 temporalRes[scv.dofIndex()] += storage[scv.localDofIndex()];
256 }
257 }
258
260 template<typename ApplyFunction>
261 void enforceDirichletConstraints(const ApplyFunction& applyDirichlet)
262 {
263 // enforce Dirichlet boundary conditions
264 this->asImp_().evalDirichletBoundaries(applyDirichlet);
265 // take care of internal Dirichlet constraints (if enabled)
266 this->asImp_().enforceInternalDirichletConstraints(applyDirichlet);
267 }
268
272 template< typename ApplyDirichletFunctionType >
273 void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
274 {
275 // enforce Dirichlet boundaries by overwriting partial derivatives with 1 or 0
276 // and set the residual to (privar - dirichletvalue)
277 if (this->elemBcTypes().hasDirichlet())
278 {
279 for (const auto& scvI : scvs(this->fvGeometry()))
280 {
281 const auto bcTypes = this->elemBcTypes().get(this->fvGeometry(), scvI);
282 if (bcTypes.hasDirichlet())
283 {
284 const auto dirichletValues = this->asImp_().problem().dirichlet(this->element(), scvI);
285
286 // set the Dirichlet conditions in residual and jacobian
287 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
288 {
289 if (bcTypes.isDirichlet(eqIdx))
290 {
291 const auto pvIdx = bcTypes.eqToDirichletIndex(eqIdx);
292 assert(0 <= pvIdx && pvIdx < numEq);
293 applyDirichlet(scvI, dirichletValues, eqIdx, pvIdx);
294 }
295 }
296 }
297 }
298 }
299 }
300
305 template<class... Args>
306 void maybeUpdateCouplingContext(Args&&...) {}
307
312 template<class... Args>
314};
315
325template<class TypeTag, class Assembler, DiffMethod diffMethod = DiffMethod::numeric, class Implementation = void>
327
334template<class TypeTag, class Assembler, class Implementation>
335class CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, Implementation>
336: public CVFELocalAssemblerBase<TypeTag, Assembler,
337 NonVoidOr<CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, Implementation>, Implementation>>
338{
345
346 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
348
349 static constexpr bool enableGridFluxVarsCache
350 = GridVariables::GridFluxVariablesCache::cachingEnabled;
351 static constexpr bool solutionDependentFluxVarsCache
352 = GridVariables::GridFluxVariablesCache::FluxVariablesCache::isSolDependent;
353
354public:
355
357 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
359
360 template <class PartialReassembler = DefaultPartialReassembler>
361 void assembleJacobian(JacobianMatrix& A, GridVariables& gridVariables,
362 const ElementResidualVector& origResiduals,
363 const PartialReassembler* partialReassembler = nullptr)
364 {
365 if (this->isImplicit())
366 assembleJacobianImplicit_(A, gridVariables, origResiduals, partialReassembler);
367 else
368 assembleJacobianExplicit_(A, gridVariables, origResiduals, partialReassembler);
369 }
370
371private:
378 template <class PartialReassembler = DefaultPartialReassembler>
379 void assembleJacobianImplicit_(JacobianMatrix& A, GridVariables& gridVariables,
380 const ElementResidualVector& origResiduals,
381 const PartialReassembler* partialReassembler = nullptr)
382 {
383 // get some aliases for convenience
384 const auto& element = this->element();
385 const auto& fvGeometry = this->fvGeometry();
386 const auto& curSol = this->asImp_().curSol();
387
388 auto&& curElemVolVars = this->curElemVolVars();
389 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
390
392 // //
393 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
394 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
395 // actual element. In the actual element we evaluate the derivative of the entire residual. //
396 // //
398
399 // if all volvars in the stencil have to be updated or if it's enough to only update the
400 // volVars for the scv whose associated dof has been deflected
401 static const bool updateAllVolVars = getParamFromGroup<bool>(
402 this->asImp_().problem().paramGroup(), "Assembly.BoxVolVarsDependOnAllElementDofs", false
403 );
404
405 // create the element solution
406 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
407
408 // create the vector storing the partial derivatives
409 ElementResidualVector partialDerivs(fvGeometry.numScv());
410
411 Dumux::Detail::VolVarsDeflectionHelper deflectionHelper(
412 [&] (const auto& scv) -> VolumeVariables& {
413 return this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
414 },
415 fvGeometry,
416 updateAllVolVars
417 );
418
419 // calculation of the derivatives
420 for (const auto& scv : scvs(fvGeometry))
421 {
422 // dof index and corresponding actual pri vars
423 const auto dofIdx = scv.dofIndex();
424 deflectionHelper.setCurrent(scv);
425
426 // calculate derivatives w.r.t to the privars at the dof at hand
427 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
428 {
429 partialDerivs = 0.0;
430
431 auto evalResiduals = [&](Scalar priVar)
432 {
433 // update the volume variables and compute element residual
434 elemSol[scv.localDofIndex()][pvIdx] = priVar;
435 deflectionHelper.deflect(elemSol, scv, this->asImp_().problem());
436 if constexpr (solutionDependentFluxVarsCache)
437 {
438 elemFluxVarsCache.update(element, fvGeometry, curElemVolVars);
439 if constexpr (enableGridFluxVarsCache)
440 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
441 }
442 this->asImp_().maybeUpdateCouplingContext(scv, elemSol, pvIdx);
443 return this->evalLocalResidual();
444 };
445
446 // derive the residuals numerically
447 static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()};
448 static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(), "Assembly.NumericDifferenceMethod");
449 NumericDifferentiation::partialDerivative(evalResiduals, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origResiduals,
450 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
451
452 // update the global stiffness matrix with the current partial derivatives
453 for (const auto& scvJ : scvs(fvGeometry))
454 {
455 // don't add derivatives for green dofs
456 if (!partialReassembler
457 || partialReassembler->dofColor(scvJ.dofIndex()) != EntityColor::green)
458 {
459 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
460 {
461 // A[i][col][eqIdx][pvIdx] is the rate of change of
462 // the residual of equation 'eqIdx' at dof 'i'
463 // depending on the primary variable 'pvIdx' at dof
464 // 'col'.
465 A[scvJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[scvJ.localDofIndex()][eqIdx];
466 }
467 }
468 }
469
470 // restore the original state of the scv's volume variables
471 deflectionHelper.restore(scv);
472
473 // restore the original element solution
474 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
475 this->asImp_().maybeUpdateCouplingContext(scv, elemSol, pvIdx);
476 }
477 }
478
479 // restore original state of the flux vars cache in case of global caching.
480 // In the case of local caching this is obsolete because the elemFluxVarsCache used here goes out of scope after this.
481 if constexpr (enableGridFluxVarsCache)
482 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
483
484 // evaluate additional derivatives that might arise from the coupling (no-op if not coupled)
485 this->asImp_().maybeEvalAdditionalDomainDerivatives(origResiduals, A, gridVariables);
486 }
487
493 template <class PartialReassembler = DefaultPartialReassembler>
494 void assembleJacobianExplicit_(JacobianMatrix& A, GridVariables& gridVariables,
495 const ElementResidualVector& origResiduals,
496 const PartialReassembler* partialReassembler = nullptr)
497 {
498 if (partialReassembler)
499 DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization");
500
501 // get some aliases for convenience
502 const auto& element = this->element();
503 const auto& fvGeometry = this->fvGeometry();
504 const auto& curSol = this->asImp_().curSol();
505 auto&& curElemVolVars = this->curElemVolVars();
506
507 // create the element solution
508 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
509
510 // create the vector storing the partial derivatives
511 ElementResidualVector partialDerivs(fvGeometry.numScv());
512
513 // calculation of the derivatives
514 for (const auto& scv : scvs(fvGeometry))
515 {
516 // dof index and corresponding actual pri vars
517 const auto dofIdx = scv.dofIndex();
518 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
519 const VolumeVariables origVolVars(curVolVars);
520
521 // calculate derivatives w.r.t to the privars at the dof at hand
522 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
523 {
524 partialDerivs = 0.0;
525
526 auto evalStorage = [&](Scalar priVar)
527 {
528 elemSol[scv.localDofIndex()][pvIdx] = priVar;
529 curVolVars.update(elemSol, this->asImp_().problem(), element, scv);
530 return this->evalStorage();
531 };
532
533 // derive the residuals numerically
534 static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()};
535 static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(), "Assembly.NumericDifferenceMethod");
536 NumericDifferentiation::partialDerivative(evalStorage, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origResiduals,
537 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
538
539 // update the global stiffness matrix with the current partial derivatives
540 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
541 {
542 // A[i][col][eqIdx][pvIdx] is the rate of change of
543 // the residual of equation 'eqIdx' at dof 'i'
544 // depending on the primary variable 'pvIdx' at dof
545 // 'col'.
546 A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
547 }
548
549 // restore the original state of the scv's volume variables
550 curVolVars = origVolVars;
551
552 // restore the original element solution
553 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
554 }
555 }
556 }
557};
558
559} // end namespace Dumux
560
561#endif
Control volume finite element local assembler using numeric differentiation.
Definition: experimental/assembly/cvfelocalassembler.hh:338
void assembleJacobian(JacobianMatrix &A, GridVariables &gridVariables, const ElementResidualVector &origResiduals, const PartialReassembler *partialReassembler=nullptr)
Definition: experimental/assembly/cvfelocalassembler.hh:361
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: experimental/assembly/cvfelocalassembler.hh:357
typename ParentType::LocalResidual LocalResidual
Definition: experimental/assembly/cvfelocalassembler.hh:356
A base class for all local CVFE assemblers.
Definition: experimental/assembly/cvfelocalassembler.hh:52
void assembleJacobian(JacobianMatrix &jac, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: experimental/assembly/cvfelocalassembler.hh:200
void bindLocalViews()
Definition: experimental/assembly/cvfelocalassembler.hh:71
void assembleResidual(ResidualVector &res)
Assemble the residual only.
Definition: experimental/assembly/cvfelocalassembler.hh:224
void assembleJacobianAndResidual(JacobianMatrix &jac, ResidualVector &res, GridVariables &gridVariables, const StageParams &stageParams, ResidualVector &temporal, ResidualVector &spatial, ResidualVector &constrainedDofs, const PartialReassembler *partialReassembler=nullptr, const CouplingFunction &maybeAssembleCouplingBlocks=noop)
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition: experimental/assembly/cvfelocalassembler.hh:82
typename ParentType::LocalResidual LocalResidual
Definition: experimental/assembly/cvfelocalassembler.hh:67
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: experimental/assembly/cvfelocalassembler.hh:68
void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
Evaluates Dirichlet boundaries.
Definition: experimental/assembly/cvfelocalassembler.hh:273
void assembleCurrentResidual(ResidualVector &spatialRes, ResidualVector &temporalRes)
Assemble the residual only.
Definition: experimental/assembly/cvfelocalassembler.hh:247
void maybeUpdateCouplingContext(Args &&...)
Update the coupling context for coupled models.
Definition: experimental/assembly/cvfelocalassembler.hh:306
void maybeEvalAdditionalDomainDerivatives(Args &&...)
Update the additional domain derivatives for coupled models.
Definition: experimental/assembly/cvfelocalassembler.hh:313
void enforceDirichletConstraints(const ApplyFunction &applyDirichlet)
Enforce Dirichlet constraints.
Definition: experimental/assembly/cvfelocalassembler.hh:261
An assembler for Jacobian and residual contribution per element (CVFE methods)
Definition: experimental/assembly/cvfelocalassembler.hh:326
A base class for all local assemblers.
Definition: experimental/assembly/fvlocalassemblerbase.hh:38
ElementVolumeVariables & curElemVolVars()
The current element volume variables.
Definition: experimental/assembly/fvlocalassemblerbase.hh:234
void bindLocalViews()
Convenience function bind and prepare all relevant variables required for the evaluation of the local...
Definition: experimental/assembly/fvlocalassemblerbase.hh:163
ElementResidualVector evalLocalResidual() const
Convenience function to evaluate the complete local residual for the current element....
Definition: experimental/assembly/fvlocalassemblerbase.hh:104
bool elementIsGhost() const
Returns if element is a ghost entity.
Definition: experimental/assembly/fvlocalassemblerbase.hh:222
ElementResidualVector evalLocalFluxAndSourceResidual() const
Convenience function to evaluate the flux and source terms (i.e, the terms without a time derivative)...
Definition: experimental/assembly/fvlocalassemblerbase.hh:133
const Problem & problem() const
The problem.
Definition: experimental/assembly/fvlocalassemblerbase.hh:210
LocalResidual & localResidual()
The local residual for the current element.
Definition: experimental/assembly/fvlocalassemblerbase.hh:242
const Element & element() const
The current element.
Definition: experimental/assembly/fvlocalassemblerbase.hh:218
FVElementGeometry & fvGeometry()
The global finite volume geometry.
Definition: experimental/assembly/fvlocalassemblerbase.hh:230
Implementation & asImp_()
Definition: experimental/assembly/fvlocalassemblerbase.hh:274
ElementBoundaryTypes & elemBcTypes()
The element's boundary types.
Definition: experimental/assembly/fvlocalassemblerbase.hh:246
std::decay_t< decltype(std::declval< Assembler >().localResidual())> LocalResidual
Definition: experimental/assembly/fvlocalassemblerbase.hh:57
Base class for grid variables.
Definition: experimental/discretization/gridvariables.hh:33
static void partialDerivative(const Function &function, Scalar x0, FunctionEvalType &derivative, const FunctionEvalType &fx0, const int numericDifferenceMethod=1)
Computes the derivative of a function with respect to a function parameter.
Definition: numericdifferentiation.hh:50
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:420
Defines all properties used in Dumux.
Type traits.
The local element solution class for control-volume finite element methods.
An enum class to define various differentiation methods available in order to compute the derivatives...
An enum class to define the colors of elements and vertices required for partial Jacobian reassembly.
DiffMethod
Differentiation methods in order to compute the derivatives of the residual i.e. the entries in the j...
Definition: diffmethod.hh:25
@ green
does not need to be reassembled
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::cctpfa||GridGeometry::discMethod==DiscretizationMethods::ccmpfa, CCElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for cell-centered schemes.
Definition: cellcentered/elementsolution.hh:101
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Definition: experimental/assembly/cclocalassembler.hh:36
constexpr auto noop
Function that performs no operation.
Definition: common/typetraits/typetraits.hh:29
A class for numeric differentiation.
An adapter class for local assemblers using numeric differentiation.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Detects which entries in the Jacobian have to be recomputed.