version 3.10-dev
experimental/assembly/cclocalassembler.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_CC_LOCAL_ASSEMBLER_HH
15#define DUMUX_EXPERIMENTAL_CC_LOCAL_ASSEMBLER_HH
16
17#include <dune/common/reservedvector.hh>
18#include <dune/grid/common/gridenums.hh> // for GhostEntity
19#include <dune/istl/matrixindexset.hh>
20
33
35
37
47template<class TypeTag, class Assembler, class Implementation>
48class CCLocalAssemblerBase : public Experimental::FVLocalAssemblerBase<TypeTag, Assembler, Implementation>
49{
52 using FVElementGeometry = typename GridGeometry::LocalView;
53 using Element = typename FVElementGeometry::Element;
54 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
55 using GridView = typename GridGeometry::GridView;
59 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
61
62public:
63
64 using ParentType::ParentType;
65
66 template <class ResidualVector, class StageParams, class PartialReassembler = DefaultPartialReassembler, class CouplingFunction = Noop>
67 void assembleJacobianAndResidual(JacobianMatrix& jac, ResidualVector& res, GridVariables& gridVariables,
68 const StageParams& stageParams, ResidualVector& temporal, ResidualVector& spatial,
69 ResidualVector& constrainedDofs,
70 const PartialReassembler* partialReassembler = nullptr,
71 const CouplingFunction& maybeAssembleCouplingBlocks = noop)
72 {
73 this->asImp_().bindLocalViews();
74 const auto globalI = this->fvGeometry().gridGeometry().elementMapper().index(this->element());
75
76 this->localResidual().spatialWeight(1.0);
77 this->localResidual().temporalWeight(1.0);
78
79 spatial[globalI] = this->evalLocalFluxAndSourceResidual(this->curElemVolVars())[0];
80 temporal[globalI] = this->localResidual().evalStorage(this->fvGeometry(), this->curElemVolVars())[0];
81 res[globalI] = spatial[globalI]*stageParams.spatialWeight(stageParams.size()-1)
82 + temporal[globalI]*stageParams.temporalWeight(stageParams.size()-1);
83
84 this->localResidual().spatialWeight(stageParams.spatialWeight(stageParams.size()-1));
85 this->localResidual().temporalWeight(stageParams.temporalWeight(stageParams.size()-1));
86
87
88 if (partialReassembler
89 && partialReassembler->elementColor(globalI) == EntityColor::green)
90 {
91 // assemble the coupling blocks for coupled models (does nothing if not coupled)
92 maybeAssembleCouplingBlocks(res[globalI]);
93 }
94 else
95 {
96 this->asImp_().assembleJacobian(jac, gridVariables, res[globalI]); // forward to the internal implementation
97
98 // assemble the coupling blocks for coupled models (does nothing if not coupled)
99 maybeAssembleCouplingBlocks(res[globalI]);
100 }
101 }
102
106 template <class ResidualVector>
107 void assembleResidual(ResidualVector& res)
108 {
109 this->asImp_().bindLocalViews();
110 const auto globalI = this->assembler().gridGeometry().elementMapper().index(this->element());
111 res[globalI] = this->asImp_().evalLocalResidual()[0]; // forward to the internal implementation
112
114 if constexpr (Problem::enableInternalDirichletConstraints())
115 {
116 const auto applyDirichlet = [&] (const auto& scvI,
117 const auto& dirichletValues,
118 const auto eqIdx,
119 const auto pvIdx)
120 {
121 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
122 };
123
124 this->asImp_().enforceInternalDirichletConstraints(applyDirichlet);
125 }
126 }
127
131 NumEqVector evalFlux(const Element& neighbor,
132 const SubControlVolumeFace& scvf) const
133 {
134 return this->localResidual().evalFlux(
135 this->asImp_().problem(), neighbor,
136 this->fvGeometry(), this->curElemVolVars(),
137 this->elemFluxVarsCache(), scvf
138 );
139 }
140
144 template<class SubResidualVector>
145 void assembleCurrentResidual(SubResidualVector& spatialRes, SubResidualVector& temporalRes)
146 {
147 this->asImp_().bindLocalViews();
148 const auto globalI = this->fvGeometry().gridGeometry().elementMapper().index(this->element());
149 spatialRes[globalI] = this->evalLocalFluxAndSourceResidual(this->curElemVolVars())[0];
150 temporalRes[globalI] = this->localResidual().evalStorage(this->fvGeometry(), this->curElemVolVars())[0];
151
152 if constexpr (Problem::enableInternalDirichletConstraints())
153 DUNE_THROW(Dune::NotImplemented, "Not implemented");
154 }
155
160 template<class... Args>
161 void maybeUpdateCouplingContext(Args&&...) {}
162
167 template<class... Args>
169};
170
180template<class TypeTag, class Assembler, DiffMethod diffMethod = DiffMethod::numeric, class Implementation = void>
182
189template<class TypeTag, class Assembler, class Implementation>
190class CCLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, Implementation>
191: public CCLocalAssemblerBase<TypeTag, Assembler,
192 NonVoidOr<CCLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, Implementation>, Implementation>>
193{
199 using FVElementGeometry = typename GridGeometry::LocalView;
200 using Element = typename FVElementGeometry::Element;
203 using Problem = typename GridVariables::GridVolumeVariables::Problem;
204
205 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
207
209 static constexpr int maxElementStencilSize = GridGeometry::maxElementStencilSize;
211
212public:
214
216 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
217
222 void assembleJacobian(JacobianMatrix& A, GridVariables& gridVariables, const NumEqVector& origResidual)
223 {
224 if (this->isImplicit())
225 assembleJacobianImplicit_(A, gridVariables, origResidual);
226 else
227 assembleJacobianExplicit_(A, gridVariables, origResidual);
228 }
229
230private:
231
236 void assembleJacobianImplicit_(JacobianMatrix& A, GridVariables& gridVariables, const NumEqVector& origResidual)
237 {
239 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
240 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
241 // actual element. In the actual element we evaluate the derivative of the entire residual. //
243
244 // get some aliases for convenience
245 const auto& element = this->element();
246 const auto& fvGeometry = this->fvGeometry();
247 const auto& gridGeometry = this->fvGeometry().gridGeometry();
248 auto&& curElemVolVars = this->curElemVolVars();
249 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
250
251 // get stencil information
252 const auto globalI = gridGeometry.elementMapper().index(element);
253 const auto& connectivityMap = gridGeometry.connectivityMap();
254 const auto numNeighbors = connectivityMap[globalI].size();
255
256 // container to store the neighboring elements
257 Dune::ReservedVector<Element, maxElementStencilSize> neighborElements;
258 neighborElements.resize(numNeighbors);
259
260 // assemble the undeflected residual
262 Residuals origResiduals(numNeighbors + 1); origResiduals = 0.0;
263 origResiduals[0] = origResidual;
264
265 // lambda for convenient evaluation of the fluxes across scvfs in the neighbors
266 // if the neighbor is a ghost we don't want to add anything to their residual
267 // so we return 0 and omit computing the flux
268 auto evalNeighborFlux = [&] (const auto& neighbor, const auto& scvf)
269 {
270 if (neighbor.partitionType() == Dune::GhostEntity)
271 return NumEqVector(0.0);
272 else
273 return this->evalFlux(neighbor, scvf);
274 };
275
276 // get the elements in which we need to evaluate the fluxes
277 // and calculate these in the undeflected state
278 unsigned int j = 1;
279 for (const auto& dataJ : connectivityMap[globalI])
280 {
281 neighborElements[j-1] = gridGeometry.element(dataJ.globalJ);
282 for (const auto scvfIdx : dataJ.scvfsJ)
283 origResiduals[j] += evalNeighborFlux(neighborElements[j-1], fvGeometry.scvf(scvfIdx));
284
285 ++j;
286 }
287
288 // reference to the element's scv (needed later) and corresponding vol vars
289 const auto& scv = fvGeometry.scv(globalI);
290 auto& curVolVars = ParentType::getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
291
292 // save a copy of the original privars and vol vars in order
293 // to restore the original solution after deflection
294 const auto& curSol = this->asImp_().curSol();
295 const auto origPriVars = curSol[globalI];
296 const auto origVolVars = curVolVars;
297
298 // element solution container to be deflected
299 auto elemSol = elementSolution<FVElementGeometry>(origPriVars);
300
301 // derivatives in the neighbors with respect to the current elements
302 // in index 0 we save the derivative of the element residual with respect to it's own dofs
303 Residuals partialDerivs(numNeighbors + 1);
304
305 for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
306 {
307 partialDerivs = 0.0;
308
309 auto evalResiduals = [&](Scalar priVar)
310 {
311 Residuals partialDerivsTmp(numNeighbors + 1);
312 partialDerivsTmp = 0.0;
313 // update the volume variables and the flux var cache
314 elemSol[0][pvIdx] = priVar;
315 this->asImp_().maybeUpdateCouplingContext(scv, elemSol, pvIdx);
316 curVolVars.update(elemSol, this->asImp_().problem(), element, scv);
317 elemFluxVarsCache.update(element, fvGeometry, curElemVolVars);
318 if (enableGridFluxVarsCache)
319 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
320
321 // calculate the residual with the deflected primary variables
322 partialDerivsTmp[0] = this->evalLocalResidual()[0];
323
324 // calculate the fluxes in the neighbors with the deflected primary variables
325 for (std::size_t k = 0; k < numNeighbors; ++k)
326 for (auto scvfIdx : connectivityMap[globalI][k].scvfsJ)
327 partialDerivsTmp[k+1] += evalNeighborFlux(neighborElements[k], fvGeometry.scvf(scvfIdx));
328
329 return partialDerivsTmp;
330 };
331
332 // derive the residuals numerically
333 static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()};
334 static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(), "Assembly.NumericDifferenceMethod");
335 NumericDifferentiation::partialDerivative(evalResiduals, elemSol[0][pvIdx], partialDerivs, origResiduals,
336 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
337
338 // Correct derivative for ghost elements, i.e. set a 1 for the derivative w.r.t. the
339 // current primary variable and a 0 elsewhere. As we always solve for a delta of the
340 // solution with respect to the initial one, this results in a delta of 0 for ghosts.
341 if (this->elementIsGhost())
342 {
343 partialDerivs[0] = 0.0;
344 partialDerivs[0][pvIdx] = 1.0;
345 }
346
347 // restore the original state of the scv's volume variables
348 curVolVars = origVolVars;
349
350 // restore the current element solution
351 elemSol[0][pvIdx] = origPriVars[pvIdx];
352
353 // restore the undeflected state of the coupling context
354 this->asImp_().maybeUpdateCouplingContext(scv, elemSol, pvIdx);
355
356 // add the current partial derivatives to the global jacobian matrix
357 // no special treatment is needed if globalJ is a ghost because then derivatives have been assembled to 0 above
358 if constexpr (Problem::enableInternalDirichletConstraints())
359 {
360 // check if own element has internal Dirichlet constraint
361 const auto internalDirichletConstraintsOwnElement = this->asImp_().problem().hasInternalDirichletConstraint(this->element(), scv);
362 const auto dirichletValues = this->asImp_().problem().internalDirichlet(this->element(), scv);
363
364 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
365 {
366 if (internalDirichletConstraintsOwnElement[eqIdx])
367 {
368 origResiduals[0][eqIdx] = origVolVars.priVars()[eqIdx] - dirichletValues[eqIdx];
369 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
370 }
371 else
372 A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx];
373 }
374
375 // off-diagonal entries
376 j = 1;
377 for (const auto& dataJ : connectivityMap[globalI])
378 {
379 const auto& neighborElement = neighborElements[j-1];
380 const auto& neighborScv = fvGeometry.scv(dataJ.globalJ);
381 const auto internalDirichletConstraintsNeighbor = this->asImp_().problem().hasInternalDirichletConstraint(neighborElement, neighborScv);
382
383 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
384 {
385 if (internalDirichletConstraintsNeighbor[eqIdx])
386 A[dataJ.globalJ][globalI][eqIdx][pvIdx] = 0.0;
387 else
388 A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j][eqIdx];
389 }
390
391 ++j;
392 }
393 }
394 else // no internal Dirichlet constraints specified
395 {
396 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
397 {
398 // the diagonal entries
399 A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx];
400
401 // off-diagonal entries
402 j = 1;
403 for (const auto& dataJ : connectivityMap[globalI])
404 A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j++][eqIdx];
405 }
406 }
407 }
408
409 // restore original state of the flux vars cache in case of global caching.
410 // This has to be done in order to guarantee that everything is in an undeflected
411 // state before the assembly of another element is called. In the case of local caching
412 // this is obsolete because the elemFluxVarsCache used here goes out of scope after this.
413 // We only have to do this for the last primary variable, for all others the flux var cache
414 // is updated with the correct element volume variables before residual evaluations
415 if (enableGridFluxVarsCache)
416 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
417
418 // evaluate additional derivatives that might arise from the coupling (no-op if not coupled)
419 ElementResidualVector orig{origResidual};
420 this->asImp_().maybeEvalAdditionalDomainDerivatives(orig, A, gridVariables);
421 }
422
427 void assembleJacobianExplicit_(JacobianMatrix& A, GridVariables& gridVariables, const NumEqVector& origResidual)
428 {
429 if (this->assembler().isStationaryProblem())
430 DUNE_THROW(Dune::InvalidStateException, "Using explicit jacobian assembler with stationary local residual");
431
432 // assemble the undeflected residual
433 auto storageResidual = origResidual;
434
436 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
437 // neighboring elements all derivatives are zero. For the assembled element only the storage //
438 // derivatives are non-zero. //
440
441 // get some aliases for convenience
442 const auto& element = this->element();
443 const auto& fvGeometry = this->fvGeometry();
444 const auto& gridGeometry = this->fvGeometry().gridGeometry();
445 auto&& curElemVolVars = this->curElemVolVars();
446
447 // reference to the element's scv (needed later) and corresponding vol vars
448 const auto globalI = gridGeometry.elementMapper().index(element);
449 const auto& scv = fvGeometry.scv(globalI);
450 auto& curVolVars = ParentType::getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
451
452 // save a copy of the original privars and vol vars in order
453 // to restore the original solution after deflection
454 const auto& curSol = this->asImp_().curSol();
455 const auto origPriVars = curSol[globalI];
456 const auto origVolVars = curVolVars;
457
458 // element solution container to be deflected
459 auto elemSol = elementSolution<FVElementGeometry>(origPriVars);
460
461 // derivatives in the neighbors with respect to the current elements
462 NumEqVector partialDeriv;
463 for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
464 {
465 // reset derivatives of element dof with respect to itself
466 partialDeriv = 0.0;
467
468 auto evalStorage = [&](Scalar priVar)
469 {
470 // update the volume variables and calculate
471 // the residual with the deflected primary variables
472 elemSol[0][pvIdx] = priVar;
473 curVolVars.update(elemSol, this->asImp_().problem(), element, scv);
474 return this->evalStorage()[0];
475 };
476
477 // for non-ghosts compute the derivative numerically
478 if (!this->elementIsGhost())
479 {
480 static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()};
481 static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(), "Assembly.NumericDifferenceMethod");
482 NumericDifferentiation::partialDerivative(evalStorage, elemSol[0][pvIdx], partialDeriv, storageResidual,
483 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
484 }
485
486 // for ghost elements we assemble a 1.0 where the primary variable and zero everywhere else
487 // as we always solve for a delta of the solution with respect to the initial solution this
488 // results in a delta of zero for ghosts
489 else partialDeriv[pvIdx] = 1.0;
490
491 // restore the original state of the scv's volume variables
492 curVolVars = origVolVars;
493
494 // restore the current element solution
495 elemSol[0][pvIdx] = origPriVars[pvIdx];
496
497 // add the current partial derivatives to the global jacobian matrix
498 // only diagonal entries for explicit jacobians
499 if constexpr (Problem::enableInternalDirichletConstraints())
500 {
501 // check if own element has internal Dirichlet constraint
502 const auto internalDirichletConstraints = this->asImp_().problem().hasInternalDirichletConstraint(this->element(), scv);
503 const auto dirichletValues = this->asImp_().problem().internalDirichlet(this->element(), scv);
504
505 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
506 {
507 // TODO: we need to set constrained dof flags for these dofs and not modify the residual here
508 // this function only takes care of the Jacobian
509 if (internalDirichletConstraints[eqIdx])
510 {
511 storageResidual[eqIdx] = origVolVars.priVars()[eqIdx] - dirichletValues[eqIdx];
512 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
513 }
514 else
515 A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx];
516 }
517 }
518 else
519 {
520 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
521 A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx];
522 }
523 }
524 }
525};
526
527} // end namespace Dumux
528
529#endif
The local element solution class for cell-centered methods.
Cell-centered scheme local assembler using numeric differentiation.
Definition: experimental/assembly/cclocalassembler.hh:193
void assembleJacobian(JacobianMatrix &A, GridVariables &gridVariables, const NumEqVector &origResidual)
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition: experimental/assembly/cclocalassembler.hh:222
typename ParentType::LocalResidual LocalResidual
Definition: experimental/assembly/cclocalassembler.hh:215
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: experimental/assembly/cclocalassembler.hh:216
A base class for all local cell-centered assemblers.
Definition: experimental/assembly/cclocalassembler.hh:49
void assembleResidual(ResidualVector &res)
Assemble the residual only.
Definition: experimental/assembly/cclocalassembler.hh:107
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)
Definition: experimental/assembly/cclocalassembler.hh:67
void maybeUpdateCouplingContext(Args &&...)
Update the coupling context for coupled models.
Definition: experimental/assembly/cclocalassembler.hh:161
void maybeEvalAdditionalDomainDerivatives(Args &&...)
Update the additional domain derivatives for coupled models.
Definition: experimental/assembly/cclocalassembler.hh:168
void assembleCurrentResidual(SubResidualVector &spatialRes, SubResidualVector &temporalRes)
Assemble the residual only.
Definition: experimental/assembly/cclocalassembler.hh:145
NumEqVector evalFlux(const Element &neighbor, const SubControlVolumeFace &scvf) const
Evaluates the fluxes (element can potentially be a neighbor)
Definition: experimental/assembly/cclocalassembler.hh:131
An assembler for Jacobian and residual contribution per element (cell-centered methods)
Definition: experimental/assembly/cclocalassembler.hh:181
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
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
const Assembler & assembler() const
The assembler.
Definition: experimental/assembly/fvlocalassemblerbase.hh:214
FVElementGeometry & fvGeometry()
The global finite volume geometry.
Definition: experimental/assembly/fvlocalassemblerbase.hh:230
Implementation & asImp_()
Definition: experimental/assembly/fvlocalassemblerbase.hh:274
ElementFluxVariablesCache & elemFluxVarsCache()
The element flux variables cache.
Definition: experimental/assembly/fvlocalassemblerbase.hh:238
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
The flux stencil specialized for different discretization schemes.
Definition: fluxstencil.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
A arithmetic block vector type based on DUNE's reserved vector.
Definition: reservedblockvector.hh:26
Defines all properties used in Dumux.
Type traits.
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.
The flux stencil specialized for different discretization schemes.
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
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
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 helper to deduce a vector with the same size as numbers of equations.
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.
A arithmetic block vector type based on DUNE's reserved vector.