version 3.11-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-FileCopyrightText: 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#include <dumux/common/typetraits/localdofs_.hh>
29
35
38
39#include <dumux/assembly/cvfevolvarsdeflectionpolicy_.hh>
40
41namespace Dumux::Experimental {
42
52template<class TypeTag, class Assembler, class Implementation>
53class CVFELocalAssemblerBase : public Experimental::FVLocalAssemblerBase<TypeTag, Assembler, Implementation>
54{
59 using ElementVolumeVariables = typename GridVariables::GridVolumeVariables::LocalView;
60 using SolutionVector = typename Assembler::SolutionVector;
61
62 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
64
65public:
66
67 using ParentType::ParentType;
68
70 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
71
72
74 {
76 this->elemBcTypes().update(this->asImp_().problem(), this->element(), this->fvGeometry());
77 }
78
83 template <class ResidualVector, class StageParams, class PartialReassembler = DefaultPartialReassembler, class CouplingFunction = Noop>
84 void assembleJacobianAndResidual(JacobianMatrix& jac, ResidualVector& res, GridVariables& gridVariables,
85 const StageParams& stageParams, ResidualVector& temporal, ResidualVector& spatial,
86 ResidualVector& constrainedDofs,
87 const PartialReassembler* partialReassembler = nullptr,
88 const CouplingFunction& maybeAssembleCouplingBlocks = noop)
89 {
90 this->asImp_().bindLocalViews();
91 const auto eIdxGlobal = this->asImp_().problem().gridGeometry().elementMapper().index(this->element());
92
93 this->localResidual().spatialWeight(1.0);
94 this->localResidual().temporalWeight(1.0);
95
96 const auto sWeight = stageParams.spatialWeight(stageParams.size()-1);
97 const auto tWeight = stageParams.temporalWeight(stageParams.size()-1);
98
99 const auto flux = this->evalLocalFluxAndSourceResidual(this->curElemVolVars());
100 const auto storage = this->localResidual().evalStorage(this->fvGeometry(), this->curElemVolVars());
101 ElementResidualVector origResidual(flux.size());
102 origResidual = 0.0;
103 for (const auto& scv : scvs(this->fvGeometry()))
104 {
105 spatial[scv.dofIndex()] += flux[scv.localDofIndex()];
106 temporal[scv.dofIndex()] += storage[scv.localDofIndex()];
107 origResidual[scv.localDofIndex()] += flux[scv.localDofIndex()]*sWeight + storage[scv.localDofIndex()]*tWeight;
108 res[scv.dofIndex()] += origResidual[scv.localDofIndex()];
109 }
110
111 this->localResidual().spatialWeight(sWeight);
112 this->localResidual().temporalWeight(tWeight);
113
114 if (partialReassembler && partialReassembler->elementColor(eIdxGlobal) == EntityColor::green)
115 {
116 // assemble the coupling blocks for coupled models (does nothing if not coupled)
117 maybeAssembleCouplingBlocks(origResidual);
118 }
119 else if (!this->elementIsGhost())
120 {
121 this->asImp_().assembleJacobian(jac, gridVariables, origResidual, partialReassembler); // forward to the internal implementation
122
123 // assemble the coupling blocks for coupled models (does nothing if not coupled)
124 maybeAssembleCouplingBlocks(origResidual);
125 }
126 else
127 {
128 // Treatment of ghost elements
129 assert(this->elementIsGhost());
130
131 // handle dofs per codimension
132 const auto& gridGeometry = this->asImp_().problem().gridGeometry();
133 Dune::Hybrid::forEach(std::make_integer_sequence<int, dim>{}, [&](auto d)
134 {
135 constexpr int codim = dim - d;
136 const auto& localCoeffs = gridGeometry.feCache().get(this->element().type()).localCoefficients();
137 for (int idx = 0; idx < localCoeffs.size(); ++idx)
138 {
139 const auto& localKey = localCoeffs.localKey(idx);
140
141 // skip if we are not handling this codim right now
142 if (localKey.codim() != codim)
143 continue;
144
145 // do not change the non-ghost entities
146 auto entity = this->element().template subEntity<codim>(localKey.subEntity());
147 if (entity.partitionType() == Dune::InteriorEntity || entity.partitionType() == Dune::BorderEntity)
148 continue;
149
150 // WARNING: this only works if the mapping from codim+subEntity to
151 // global dofIndex is unique (on dof per entity of this codim).
152 // For more general mappings, we should use a proper local-global mapping here.
153 // For example through dune-functions.
154 const auto dofIndex = gridGeometry.dofMapper().index(entity);
155
156 // this might be a vector-valued dof
157 using BlockType = typename JacobianMatrix::block_type;
158 BlockType &J = jac[dofIndex][dofIndex];
159 for (int j = 0; j < BlockType::rows; ++j)
160 J[j][j] = 1.0;
161
162 // set residual for the ghost dof
163 res[dofIndex] = 0;
164 constrainedDofs[dofIndex] = 1;
165 }
166 });
167 }
168
169 auto applyDirichlet = [&] (const auto& scvI,
170 const auto& dirichletValues,
171 const auto eqIdx,
172 const auto pvIdx)
173 {
174 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
175 constrainedDofs[scvI.dofIndex()][eqIdx] = 1;
176
177 auto& row = jac[scvI.dofIndex()];
178 for (auto col = row.begin(); col != row.end(); ++col)
179 row[col.index()][eqIdx] = 0.0;
180
181 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
182
183 // if a periodic dof has Dirichlet values also apply the same Dirichlet values to the other dof
184 if (this->asImp_().problem().gridGeometry().dofOnPeriodicBoundary(scvI.dofIndex()))
185 {
186 const auto periodicDof = this->asImp_().problem().gridGeometry().periodicallyMappedDof(scvI.dofIndex());
187 res[periodicDof][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
188 constrainedDofs[periodicDof][eqIdx] = 1;
189 const auto end = jac[periodicDof].end();
190 for (auto it = jac[periodicDof].begin(); it != end; ++it)
191 (*it) = periodicDof != it.index() ? 0.0 : 1.0;
192 }
193 };
194
195 this->asImp_().enforceDirichletConstraints(applyDirichlet);
196 }
197
202 void assembleJacobian(JacobianMatrix& jac, GridVariables& gridVariables)
203 {
204 this->asImp_().bindLocalViews();
205 this->asImp_().assembleJacobian(jac, gridVariables); // forward to the internal implementation
206
207 auto applyDirichlet = [&] (const auto& scvI,
208 const auto& dirichletValues,
209 const auto eqIdx,
210 const auto pvIdx)
211 {
212 auto& row = jac[scvI.dofIndex()];
213 for (auto col = row.begin(); col != row.end(); ++col)
214 row[col.index()][eqIdx] = 0.0;
215
216 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
217 };
218
219 this->asImp_().enforceDirichletConstraints(applyDirichlet);
220 }
221
225 template <class ResidualVector>
226 void assembleResidual(ResidualVector& res)
227 {
228 this->asImp_().bindLocalViews();
229 const auto residual = this->evalLocalResidual();
230
231 for (const auto& localDof : localDofs(this->fvGeometry()))
232 res[localDof.dofIndex()] += residual[localDof.index()];
233
234 auto applyDirichlet = [&] (const auto& scvI,
235 const auto& dirichletValues,
236 const auto eqIdx,
237 const auto pvIdx)
238 {
239 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
240 };
241
242 this->asImp_().enforceDirichletConstraints(applyDirichlet);
243 }
244
248 template<class ResidualVector>
249 void assembleCurrentResidual(ResidualVector& spatialRes, ResidualVector& temporalRes)
250 {
251 this->asImp_().bindLocalViews();
252 const auto flux = this->evalLocalFluxAndSourceResidual(this->curElemVolVars());
253 const auto storage = this->localResidual().evalStorage(this->fvGeometry(), this->curElemVolVars());
254 for (const auto& scv : scvs(this->fvGeometry()))
255 {
256 spatialRes[scv.dofIndex()] += flux[scv.localDofIndex()];
257 temporalRes[scv.dofIndex()] += storage[scv.localDofIndex()];
258 }
259 }
260
262 template<typename ApplyFunction>
263 void enforceDirichletConstraints(const ApplyFunction& applyDirichlet)
264 {
265 // enforce Dirichlet boundary conditions
266 this->asImp_().evalDirichletBoundaries(applyDirichlet);
267 // take care of internal Dirichlet constraints (if enabled)
268 this->asImp_().enforceInternalDirichletConstraints(applyDirichlet);
269 }
270
274 template< typename ApplyDirichletFunctionType >
275 void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
276 {
277 // enforce Dirichlet boundaries by overwriting partial derivatives with 1 or 0
278 // and set the residual to (privar - dirichletvalue)
279 if (this->elemBcTypes().hasDirichlet())
280 {
281 for (const auto& scvI : scvs(this->fvGeometry()))
282 {
283 const auto bcTypes = this->elemBcTypes().get(this->fvGeometry(), scvI);
284 if (bcTypes.hasDirichlet())
285 {
286 const auto dirichletValues = this->asImp_().problem().dirichlet(this->element(), scvI);
287
288 // set the Dirichlet conditions in residual and jacobian
289 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
290 {
291 if (bcTypes.isDirichlet(eqIdx))
292 {
293 const auto pvIdx = bcTypes.eqToDirichletIndex(eqIdx);
294 assert(0 <= pvIdx && pvIdx < numEq);
295 applyDirichlet(scvI, dirichletValues, eqIdx, pvIdx);
296 }
297 }
298 }
299 }
300 }
301 }
302
307 template<class... Args>
308 void maybeUpdateCouplingContext(Args&&...) {}
309
314 template<class... Args>
316};
317
327template<class TypeTag, class Assembler, DiffMethod diffMethod = DiffMethod::numeric, class Implementation = void>
329
336template<class TypeTag, class Assembler, class Implementation>
337class CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, Implementation>
338: public CVFELocalAssemblerBase<TypeTag, Assembler,
339 NonVoidOr<CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, Implementation>, Implementation>>
340{
345 using ElementVolumeVariables = typename GridVariables::GridVolumeVariables::LocalView;
348
349 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
351
352 static constexpr bool enableGridFluxVarsCache
353 = GridVariables::GridFluxVariablesCache::cachingEnabled;
354 static constexpr bool solutionDependentFluxVarsCache
355 = GridVariables::GridFluxVariablesCache::FluxVariablesCache::isSolDependent;
356
357public:
358
360 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
362
363 template <class PartialReassembler = DefaultPartialReassembler>
364 void assembleJacobian(JacobianMatrix& A, GridVariables& gridVariables,
365 const ElementResidualVector& origResiduals,
366 const PartialReassembler* partialReassembler = nullptr)
367 {
368 if (this->isImplicit())
369 assembleJacobianImplicit_(A, gridVariables, origResiduals, partialReassembler);
370 else
371 assembleJacobianExplicit_(A, gridVariables, origResiduals, partialReassembler);
372 }
373
374private:
381 template <class PartialReassembler = DefaultPartialReassembler>
382 void assembleJacobianImplicit_(JacobianMatrix& A, GridVariables& gridVariables,
383 const ElementResidualVector& origResiduals,
384 const PartialReassembler* partialReassembler = nullptr)
385 {
386 // get some aliases for convenience
387 const auto& element = this->element();
388 const auto& fvGeometry = this->fvGeometry();
389 const auto& curSol = this->asImp_().curSol();
390
391 auto&& curElemVolVars = this->curElemVolVars();
392 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
393
395 // //
396 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
397 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
398 // actual element. In the actual element we evaluate the derivative of the entire residual. //
399 // //
401
402 // if all volvars in the stencil have to be updated or if it's enough to only update the
403 // volVars for the scv whose associated dof has been deflected
404 static const bool updateAllVolVars = getParamFromGroup<bool>(
405 this->asImp_().problem().paramGroup(), "Assembly.BoxVolVarsDependOnAllElementDofs", false
406 );
407
408 // create the element solution
409 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
410
411 // create the vector storing the partial derivatives
412 ElementResidualVector partialDerivs(Dumux::Detail::LocalDofs::numLocalDofs(fvGeometry));
413
414 auto deflectionPolicy = Detail::CVFE::makeVariablesDeflectionPolicy(
415 gridVariables.curGridVolVars(),
416 curElemVolVars,
417 fvGeometry,
418 updateAllVolVars
419 );
420
421 auto assembleDerivative = [&, this](const auto& localDof)
422 {
423 // dof index and corresponding actual pri vars
424 const auto dofIdx = localDof.dofIndex();
425 const auto localIdx = localDof.index();
426 deflectionPolicy.store(localDof);
427
428 // calculate derivatives w.r.t to the privars at the dof at hand
429 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
430 {
431 partialDerivs = 0.0;
432
433 auto evalResiduals = [&](Scalar priVar)
434 {
435 // update the volume variables and compute element residual
436 elemSol[localIdx][pvIdx] = priVar;
437 deflectionPolicy.update(elemSol, localDof, this->asImp_().problem());
438 if constexpr (solutionDependentFluxVarsCache)
439 {
440 elemFluxVarsCache.update(element, fvGeometry, curElemVolVars);
441 if constexpr (enableGridFluxVarsCache)
442 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
443 }
444 this->asImp_().maybeUpdateCouplingContext(localDof, elemSol, pvIdx);
445 return this->evalLocalResidual();
446 };
447
448 // derive the residuals numerically
449 static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()};
450 static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(), "Assembly.NumericDifferenceMethod");
451 NumericDifferentiation::partialDerivative(evalResiduals, elemSol[localIdx][pvIdx], partialDerivs, origResiduals,
452 eps_(elemSol[localIdx][pvIdx], pvIdx), numDiffMethod);
453
454 // update the global stiffness matrix with the current partial derivatives
455 for (const auto& localDofJ : localDofs(fvGeometry))
456 {
457 // don't add derivatives for green dofs
458 if (!partialReassembler
459 || partialReassembler->dofColor(localDofJ.dofIndex()) != EntityColor::green)
460 {
461 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
462 {
463 // A[i][col][eqIdx][pvIdx] is the rate of change of
464 // the residual of equation 'eqIdx' at dof 'i'
465 // depending on the primary variable 'pvIdx' at dof
466 // 'col'.
467 A[localDofJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[localDofJ.index()][eqIdx];
468 }
469 }
470 }
471
472 // restore the original state of the scv's volume variables
473 deflectionPolicy.restore(localDof);
474
475 // restore the original element solution
476 elemSol[localIdx][pvIdx] = curSol[localDof.dofIndex()][pvIdx];
477 this->asImp_().maybeUpdateCouplingContext(localDof, elemSol, pvIdx);
478 }
479 };
480
481 // calculation of the derivatives
482 for (const auto& localDof : localDofs(fvGeometry))
483 assembleDerivative(localDof);
484
485 // restore original state of the flux vars cache in case of global caching.
486 // In the case of local caching this is obsolete because the elemFluxVarsCache used here goes out of scope after this.
487 if constexpr (enableGridFluxVarsCache)
488 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
489
490 // evaluate additional derivatives that might arise from the coupling (no-op if not coupled)
491 this->asImp_().maybeEvalAdditionalDomainDerivatives(origResiduals, A, gridVariables);
492 }
493
499 template <class PartialReassembler = DefaultPartialReassembler>
500 void assembleJacobianExplicit_(JacobianMatrix& A, GridVariables& gridVariables,
501 const ElementResidualVector& origResiduals,
502 const PartialReassembler* partialReassembler = nullptr)
503 {
504 if (partialReassembler)
505 DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization");
506
507 // get some aliases for convenience
508 const auto& element = this->element();
509 const auto& fvGeometry = this->fvGeometry();
510 const auto& curSol = this->asImp_().curSol();
511 auto&& curElemVolVars = this->curElemVolVars();
512
513 // create the element solution
514 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
515
516 // create the vector storing the partial derivatives
517 ElementResidualVector partialDerivs(Dumux::Detail::LocalDofs::numLocalDofs(fvGeometry));
518
519 // calculation of the derivatives
520 for (const auto& scv : scvs(fvGeometry))
521 {
522 // dof index and corresponding actual pri vars
523 const auto dofIdx = scv.dofIndex();
524 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
525 const VolumeVariables origVolVars(curVolVars);
526
527 // calculate derivatives w.r.t to the privars at the dof at hand
528 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
529 {
530 partialDerivs = 0.0;
531
532 auto evalStorage = [&](Scalar priVar)
533 {
534 elemSol[scv.localDofIndex()][pvIdx] = priVar;
535 curVolVars.update(elemSol, this->asImp_().problem(), element, scv);
536 return this->evalStorage();
537 };
538
539 // derive the residuals numerically
540 static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()};
541 static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(), "Assembly.NumericDifferenceMethod");
542 NumericDifferentiation::partialDerivative(evalStorage, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origResiduals,
543 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
544
545 // update the global stiffness matrix with the current partial derivatives
546 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
547 {
548 // A[i][col][eqIdx][pvIdx] is the rate of change of
549 // the residual of equation 'eqIdx' at dof 'i'
550 // depending on the primary variable 'pvIdx' at dof
551 // 'col'.
552 A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
553 }
554
555 // restore the original state of the scv's volume variables
556 curVolVars = origVolVars;
557
558 // restore the original element solution
559 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
560 }
561 }
562 }
563};
564
565} // end namespace Dumux
566
567#endif
Control volume finite element local assembler using numeric differentiation.
Definition: experimental/assembly/cvfelocalassembler.hh:340
void assembleJacobian(JacobianMatrix &A, GridVariables &gridVariables, const ElementResidualVector &origResiduals, const PartialReassembler *partialReassembler=nullptr)
Definition: experimental/assembly/cvfelocalassembler.hh:364
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: experimental/assembly/cvfelocalassembler.hh:360
typename ParentType::LocalResidual LocalResidual
Definition: experimental/assembly/cvfelocalassembler.hh:359
A base class for all local CVFE assemblers.
Definition: experimental/assembly/cvfelocalassembler.hh:54
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:202
void bindLocalViews()
Definition: experimental/assembly/cvfelocalassembler.hh:73
void assembleResidual(ResidualVector &res)
Assemble the residual only.
Definition: experimental/assembly/cvfelocalassembler.hh:226
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:84
typename ParentType::LocalResidual LocalResidual
Definition: experimental/assembly/cvfelocalassembler.hh:69
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: experimental/assembly/cvfelocalassembler.hh:70
void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
Evaluates Dirichlet boundaries.
Definition: experimental/assembly/cvfelocalassembler.hh:275
void assembleCurrentResidual(ResidualVector &spatialRes, ResidualVector &temporalRes)
Assemble the residual only.
Definition: experimental/assembly/cvfelocalassembler.hh:249
void maybeUpdateCouplingContext(Args &&...)
Update the coupling context for coupled models.
Definition: experimental/assembly/cvfelocalassembler.hh:308
void maybeEvalAdditionalDomainDerivatives(Args &&...)
Update the additional domain derivatives for coupled models.
Definition: experimental/assembly/cvfelocalassembler.hh:315
void enforceDirichletConstraints(const ApplyFunction &applyDirichlet)
Enforce Dirichlet constraints.
Definition: experimental/assembly/cvfelocalassembler.hh:263
An assembler for Jacobian and residual contribution per element (CVFE methods)
Definition: experimental/assembly/cvfelocalassembler.hh:328
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
Class representing dofs on elements for control-volume finite element schemes.
Definition: experimental/assembly/cclocalassembler.hh:36
constexpr auto noop
Function that performs no operation.
Definition: common/typetraits/typetraits.hh:29
std::ranges::range auto scvs(const FVElementGeometry &fvGeometry, const LocalDof &localDof)
Definition: localdof.hh:79
auto localDofs(const FVElementGeometry &fvGeometry)
range over local dofs
Definition: localdof.hh:50
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.