version 3.11-dev
Loading...
Searching...
No Matches
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
29#include <dumux/common/typetraits/localdofs_.hh>
30
36
39
40#include <dumux/assembly/cvfevolvarsdeflectionpolicy_.hh>
41
42namespace Dumux::Experimental {
43
53template<class TypeTag, class Assembler, class Implementation>
54class CVFELocalAssemblerBase : public Experimental::FVLocalAssemblerBase<TypeTag, Assembler, Implementation>
55{
60 using ElementVolumeVariables = typename GridVariables::GridVolumeVariables::LocalView;
61 using SolutionVector = typename Assembler::SolutionVector;
62
63 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
65
66public:
67
68 using ParentType::ParentType;
69
72
73
75 {
77 this->elemBcTypes().update(this->asImp_().problem(), this->element(), this->fvGeometry());
78 }
79
84 template <class ResidualVector, class StageParams, class PartialReassembler = DefaultPartialReassembler, class CouplingFunction = Noop>
85 void assembleJacobianAndResidual(JacobianMatrix& jac, ResidualVector& res, GridVariables& gridVariables,
86 const StageParams& stageParams, ResidualVector& temporal, ResidualVector& spatial,
87 ResidualVector& constrainedDofs,
88 const PartialReassembler* partialReassembler = nullptr,
89 const CouplingFunction& maybeAssembleCouplingBlocks = noop)
90 {
91 this->asImp_().bindLocalViews();
92 const auto eIdxGlobal = this->asImp_().problem().gridGeometry().elementMapper().index(this->element());
93
94 this->localResidual().spatialWeight(1.0);
95 this->localResidual().temporalWeight(1.0);
96
97 const auto sWeight = stageParams.spatialWeight(stageParams.size()-1);
98 const auto tWeight = stageParams.temporalWeight(stageParams.size()-1);
99
100 const auto flux = this->evalLocalFluxAndSourceResidual(this->curElemVolVars());
101 const auto storage = this->localResidual().evalStorage(this->fvGeometry(), this->curElemVolVars());
102 ElementResidualVector origResidual(flux.size());
103 origResidual = 0.0;
104 for (const auto& scv : scvs(this->fvGeometry()))
105 {
106 spatial[scv.dofIndex()] += flux[scv.localDofIndex()];
107 temporal[scv.dofIndex()] += storage[scv.localDofIndex()];
108 origResidual[scv.localDofIndex()] += flux[scv.localDofIndex()]*sWeight + storage[scv.localDofIndex()]*tWeight;
109 res[scv.dofIndex()] += origResidual[scv.localDofIndex()];
110 }
111
112 this->localResidual().spatialWeight(sWeight);
113 this->localResidual().temporalWeight(tWeight);
114
115 if (partialReassembler && partialReassembler->elementColor(eIdxGlobal) == EntityColor::green)
116 {
117 // assemble the coupling blocks for coupled models (does nothing if not coupled)
118 maybeAssembleCouplingBlocks(origResidual);
119 }
120 else if (!this->elementIsGhost())
121 {
122 this->asImp_().assembleJacobian(jac, gridVariables, origResidual, partialReassembler); // forward to the internal implementation
123
124 // assemble the coupling blocks for coupled models (does nothing if not coupled)
125 maybeAssembleCouplingBlocks(origResidual);
126 }
127 else
128 {
129 // Treatment of ghost elements
130 assert(this->elementIsGhost());
131
132 // handle dofs per codimension
133 const auto& gridGeometry = this->asImp_().problem().gridGeometry();
134 Dune::Hybrid::forEach(std::make_integer_sequence<int, dim+1>{}, [&](auto d)
135 {
136 constexpr int codim = dim - d;
137 const auto& localCoeffs = gridGeometry.feCache().get(this->element().type()).localCoefficients();
138 for (int idx = 0; idx < localCoeffs.size(); ++idx)
139 {
140 const auto& localKey = localCoeffs.localKey(idx);
141
142 // skip if we are not handling this codim right now
143 if (localKey.codim() != codim)
144 continue;
145
146 // do not change the non-ghost entities
147 auto entity = this->element().template subEntity<codim>(localKey.subEntity());
148 if (entity.partitionType() == Dune::InteriorEntity || entity.partitionType() == Dune::BorderEntity)
149 continue;
150
151 // Set identity rows for ALL DOFs of this ghost entity.
152 // Entities with multiple DOFs (e.g. PQ3 edge interior DOFs with 2 per edge)
153 // require iterating over all DOF indices via asMultiMapper(dofMapper()).indices(entity).
154 using BlockType = typename JacobianMatrix::block_type;
155 for (const auto dofIndex : asMultiMapper(gridGeometry.dofMapper()).indices(entity))
156 {
157 BlockType &J = jac[dofIndex][dofIndex];
158 for (int j = 0; j < BlockType::rows; ++j)
159 J[j][j] = 1.0;
160 res[dofIndex] = 0;
161 constrainedDofs[dofIndex] = 1;
162 }
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& localDof : localDofs(this->fvGeometry()))
230 res[localDof.dofIndex()] += residual[localDof.index()];
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{
343 using ElementVolumeVariables = typename GridVariables::GridVolumeVariables::LocalView;
346
347 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
349
350 static constexpr bool enableGridFluxVarsCache
351 = GridVariables::GridFluxVariablesCache::cachingEnabled;
352 static constexpr bool solutionDependentFluxVarsCache
353 = GridVariables::GridFluxVariablesCache::FluxVariablesCache::isSolDependent;
354
355public:
356
359 using ParentType::ParentType;
360
361 template <class PartialReassembler = DefaultPartialReassembler>
362 void assembleJacobian(JacobianMatrix& A, GridVariables& gridVariables,
363 const ElementResidualVector& origResiduals,
364 const PartialReassembler* partialReassembler = nullptr)
365 {
366 if (this->isImplicit())
367 assembleJacobianImplicit_(A, gridVariables, origResiduals, partialReassembler);
368 else
369 assembleJacobianExplicit_(A, gridVariables, origResiduals, partialReassembler);
370 }
371
372private:
379 template <class PartialReassembler = DefaultPartialReassembler>
380 void assembleJacobianImplicit_(JacobianMatrix& A, GridVariables& gridVariables,
381 const ElementResidualVector& origResiduals,
382 const PartialReassembler* partialReassembler = nullptr)
383 {
384 // get some aliases for convenience
385 const auto& element = this->element();
386 const auto& fvGeometry = this->fvGeometry();
387 const auto& curSol = this->asImp_().curSol();
388
389 auto&& curElemVolVars = this->curElemVolVars();
390 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
391
393 // //
394 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
395 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
396 // actual element. In the actual element we evaluate the derivative of the entire residual. //
397 // //
399
400 // if all volvars in the stencil have to be updated or if it's enough to only update the
401 // volVars for the scv whose associated dof has been deflected
402 static const bool updateAllVolVars = getParamFromGroup<bool>(
403 this->asImp_().problem().paramGroup(), "Assembly.BoxVolVarsDependOnAllElementDofs", false
404 );
405
406 // create the element solution
407 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
408
409 // create the vector storing the partial derivatives
410 ElementResidualVector partialDerivs(Dumux::Detail::LocalDofs::numLocalDofs(fvGeometry));
411
412 auto deflectionPolicy = Dumux::Detail::CVFE::makeVariablesDeflectionPolicy(
413 gridVariables.curGridVolVars(),
414 curElemVolVars,
415 fvGeometry,
416 updateAllVolVars
417 );
418
419 auto assembleDerivative = [&, this](const auto& localDof)
420 {
421 // dof index and corresponding actual pri vars
422 const auto dofIdx = localDof.dofIndex();
423 const auto localIdx = localDof.index();
424 deflectionPolicy.store(localDof);
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[localIdx][pvIdx] = priVar;
435 deflectionPolicy.update(elemSol, localDof, 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(localDof, 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[localIdx][pvIdx], partialDerivs, origResiduals,
450 eps_(elemSol[localIdx][pvIdx], pvIdx), numDiffMethod);
451
452 // update the global stiffness matrix with the current partial derivatives
453 for (const auto& localDofJ : localDofs(fvGeometry))
454 {
455 // don't add derivatives for green dofs
456 if (!partialReassembler
457 || partialReassembler->dofColor(localDofJ.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[localDofJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[localDofJ.index()][eqIdx];
466 }
467 }
468 }
469
470 // restore the original state of the scv's volume variables
471 deflectionPolicy.restore(localDof);
472
473 // restore the original element solution
474 elemSol[localIdx][pvIdx] = curSol[localDof.dofIndex()][pvIdx];
475 this->asImp_().maybeUpdateCouplingContext(localDof, elemSol, pvIdx);
476 }
477 };
478
479 // calculation of the derivatives
480 for (const auto& localDof : localDofs(fvGeometry))
481 assembleDerivative(localDof);
482
483 // restore original state of the flux vars cache in case of global caching.
484 // In the case of local caching this is obsolete because the elemFluxVarsCache used here goes out of scope after this.
485 if constexpr (enableGridFluxVarsCache)
486 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
487
488 // evaluate additional derivatives that might arise from the coupling (no-op if not coupled)
489 this->asImp_().maybeEvalAdditionalDomainDerivatives(origResiduals, A, gridVariables);
490 }
491
497 template <class PartialReassembler = DefaultPartialReassembler>
498 void assembleJacobianExplicit_(JacobianMatrix& A, GridVariables& gridVariables,
499 const ElementResidualVector& origResiduals,
500 const PartialReassembler* partialReassembler = nullptr)
501 {
502 if (partialReassembler)
503 DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization");
504
505 // get some aliases for convenience
506 const auto& element = this->element();
507 const auto& fvGeometry = this->fvGeometry();
508 const auto& curSol = this->asImp_().curSol();
509 auto&& curElemVolVars = this->curElemVolVars();
510
511 // create the element solution
512 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
513
514 // create the vector storing the partial derivatives
515 ElementResidualVector partialDerivs(Dumux::Detail::LocalDofs::numLocalDofs(fvGeometry));
516
517 // calculation of the derivatives
518 for (const auto& scv : scvs(fvGeometry))
519 {
520 // dof index and corresponding actual pri vars
521 const auto dofIdx = scv.dofIndex();
522 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
523 const VolumeVariables origVolVars(curVolVars);
524
525 // calculate derivatives w.r.t to the privars at the dof at hand
526 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
527 {
528 partialDerivs = 0.0;
529
530 auto evalStorage = [&](Scalar priVar)
531 {
532 elemSol[scv.localDofIndex()][pvIdx] = priVar;
533 curVolVars.update(elemSol, this->asImp_().problem(), element, scv);
534 return this->evalStorage();
535 };
536
537 // derive the residuals numerically
538 static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()};
539 static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(), "Assembly.NumericDifferenceMethod");
540 NumericDifferentiation::partialDerivative(evalStorage, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origResiduals,
541 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
542
543 // update the global stiffness matrix with the current partial derivatives
544 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
545 {
546 // A[i][col][eqIdx][pvIdx] is the rate of change of
547 // the residual of equation 'eqIdx' at dof 'i'
548 // depending on the primary variable 'pvIdx' at dof
549 // 'col'.
550 A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
551 }
552
553 // restore the original state of the scv's volume variables
554 curVolVars = origVolVars;
555
556 // restore the original element solution
557 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
558 }
559 }
560 }
561};
562
563} // end namespace Dumux
564
565#endif
A linear system assembler (residual and Jacobian) for general discretization schemes.
Definition assembly/assembler.hh:83
GetPropType< TypeTag, Properties::SolutionVector > SolutionVector
Definition assembly/assembler.hh:96
void assembleJacobian(JacobianMatrix &A, GridVariables &gridVariables, const ElementResidualVector &origResiduals, const PartialReassembler *partialReassembler=nullptr)
Definition experimental/assembly/cvfelocalassembler.hh:362
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition experimental/assembly/cvfelocalassembler.hh:358
typename ParentType::LocalResidual LocalResidual
Definition experimental/assembly/cvfelocalassembler.hh:357
A base class for all local CVFE assemblers.
Definition experimental/assembly/cvfelocalassembler.hh:55
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:74
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:85
typename ParentType::LocalResidual LocalResidual
Definition experimental/assembly/cvfelocalassembler.hh:70
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition experimental/assembly/cvfelocalassembler.hh:71
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
The grid variables class for general schemes, storing variables and data.
Definition discretization/gridvariables.hh:27
ReservedBlockVector< NumEqVector, Dumux::Detail::LocalDofs::maxNumLocalDofs< ElementDiscretization >()> ElementResidualVector
the container storing all element residuals
Definition assembly/localresidual.hh:56
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.
A base class for all local assemblers.
DiffMethod
Differentiation methods in order to compute the derivatives of the residual i.e. the entries in the j...
Definition diffmethod.hh:25
@ numeric
Definition diffmethod.hh:26
@ green
does not need to be reassembled
Definition entitycolor.hh:40
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
@ element
Definition fieldtype.hh:23
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition parameters.hh:149
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.
Adapter to expose a multi-DOF mapper interface for single- and multi-DOF mappers.
Definition assembly/assembler.hh:44
constexpr auto asMultiMapper(const Mapper &mapper)
Definition multimapperview.hh:48
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:82
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.