14#ifndef DUMUX_MULTIDOMAIN_CC_LOCAL_ASSEMBLER_HH
15#define DUMUX_MULTIDOMAIN_CC_LOCAL_ASSEMBLER_HH
17#include <dune/common/indices.hh>
18#include <dune/common/reservedvector.hh>
19#include <dune/common/hybridutilities.hh>
20#include <dune/grid/common/gridenums.hh>
21#include <dune/istl/matrixindexset.hh>
47template<std::
size_t id,
class TypeTag,
class Assembler,
class Implementation,
bool implicit>
55 using SolutionVector =
typename Assembler::SolutionVector;
59 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
60 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
61 using ElementFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache::LocalView;
62 using Scalar =
typename GridVariables::Scalar;
64 using GridGeometry =
typename GridVariables::GridGeometry;
65 using FVElementGeometry =
typename GridGeometry::LocalView;
66 using SubControlVolume =
typename GridGeometry::SubControlVolume;
67 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
69 using GridView =
typename GridGeometry::GridView;
70 using Element =
typename GridView::template Codim<0>::Entity;
72 using CouplingManager =
typename Assembler::CouplingManager;
77 static constexpr auto domainId =
typename Dune::index_constant<id>();
81 using ParentType::ParentType;
86 const SolutionVector&
curSol,
104 template<
class JacobianMatrixRow,
class SubRes
idualVector,
class Gr
idVariablesTuple>
107 this->
asImp_().bindLocalViews();
110 const auto globalI = this->
fvGeometry().gridGeometry().elementMapper().index(this->
element());
111 res[globalI] = this->
asImp_().assembleJacobianAndResidualImplInverse(jacRow[
domainId], *std::get<domainId>(gridVariables));
114 using namespace Dune::Hybrid;
115 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](
auto&& i)
126 template<std::size_t otherId,
class JacRow,
class GridVariables,
127 typename std::enable_if_t<(otherId == id),
int> = 0>
129 const LocalResidualValues& res, GridVariables& gridVariables)
135 template<std::size_t otherId,
class JacRow,
class GridVariables,
136 typename std::enable_if_t<(otherId != id),
int> = 0>
138 const LocalResidualValues& res, GridVariables& gridVariables)
140 this->
asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
146 template<
class SubRes
idualVector>
149 this->
asImp_().bindLocalViews();
150 const auto globalI = this->
fvGeometry().gridGeometry().elementMapper().index(this->
element());
153 if constexpr (Problem::enableInternalDirichletConstraints())
155 const auto applyDirichlet = [&] (
const auto& scvI,
156 const auto& dirichletValues,
160 res[scvI.dofIndex()][eqIdx] = this->
curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
163 this->
asImp_().enforceInternalDirichletConstraints(applyDirichlet);
173 ElementResidualVector residual(this->
fvGeometry().numScv());
179 const auto& curVolVars = elemVolVars[scv];
181 source *= -Extrusion::volume(this->
fvGeometry(), scv)*curVolVars.extrusionFactor();
182 residual[scv.indexInElement()] = std::move(source);
206 const SubControlVolumeFace& scvf)
const
232 if (!this->
assembler().isStationaryProblem())
252 {
return couplingManager_; }
269template<std::
size_t id,
class TypeTag,
class Assembler, DiffMethod DM = DiffMethod::numeric,
bool implicit = true>
278template<std::
size_t id,
class TypeTag,
class Assembler>
291 using FVElementGeometry =
typename GridGeometry::LocalView;
292 using GridView =
typename GridGeometry::GridView;
293 using Element =
typename GridView::template Codim<0>::Entity;
296 enum { dim = GridView::dimension };
300 static constexpr int maxElementStencilSize = GridGeometry::maxElementStencilSize;
301 static constexpr auto domainI = Dune::index_constant<id>();
304 using ParentType::ParentType;
312 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
324 const auto& gridGeometry =
fvGeometry.gridGeometry();
329 const auto globalI = gridGeometry.elementMapper().index(
element);
330 const auto& connectivityMap = gridGeometry.connectivityMap();
331 const auto numNeighbors = connectivityMap[globalI].size();
334 Dune::ReservedVector<Element, maxElementStencilSize> neighborElements;
335 neighborElements.resize(numNeighbors);
339 Residuals origResiduals(numNeighbors + 1); origResiduals = 0.0;
345 const auto evalNeighborFlux = [&] (
const auto& neighbor,
const auto& scvf)
347 if (neighbor.partitionType() == Dune::GhostEntity)
348 return LocalResidualValues(0.0);
356 for (
const auto& dataJ : connectivityMap[globalI])
358 neighborElements[j-1] = gridGeometry.element(dataJ.globalJ);
359 for (
const auto scvfIdx : dataJ.scvfsJ)
360 origResiduals[j] += evalNeighborFlux(neighborElements[j-1],
fvGeometry.scvf(scvfIdx));
372 const auto origPriVars =
curSol[globalI];
373 const auto origVolVars = curVolVars;
380 Residuals partialDerivs(numNeighbors + 1);
382 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
386 auto evalResiduals = [&](Scalar priVar)
388 Residuals partialDerivsTmp(numNeighbors + 1);
389 partialDerivsTmp = 0.0;
391 elemSol[0][pvIdx] = priVar;
392 this->
couplingManager().updateCouplingContext(domainI, *
this, domainI, globalI, elemSol[0], pvIdx);
395 if (enableGridFluxVarsCache)
402 for (std::size_t k = 0; k < connectivityMap[globalI].size(); ++k)
403 for (
auto scvfIdx : connectivityMap[globalI][k].scvfsJ)
404 partialDerivsTmp[k+1] += evalNeighborFlux(neighborElements[k],
fvGeometry.scvf(scvfIdx));
406 return partialDerivsTmp;
413 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
420 partialDerivs[0] = 0.0;
421 partialDerivs[0][pvIdx] = 1.0;
425 curVolVars = origVolVars;
428 elemSol[0][pvIdx] = origPriVars[pvIdx];
431 this->
couplingManager().updateCouplingContext(domainI, *
this, domainI, globalI, elemSol[0], pvIdx);
434 if constexpr (Problem::enableInternalDirichletConstraints())
437 const auto internalDirichletConstraintsOwnElement = this->
problem().hasInternalDirichletConstraint(this->
element(), scv);
438 const auto dirichletValues = this->
problem().internalDirichlet(this->
element(), scv);
440 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
442 if (internalDirichletConstraintsOwnElement[eqIdx])
444 origResiduals[0][eqIdx] = origVolVars.priVars()[eqIdx] - dirichletValues[eqIdx];
445 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
448 A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx];
453 for (
const auto& dataJ : connectivityMap[globalI])
455 const auto& neighborElement = neighborElements[j-1];
456 const auto& neighborScv =
fvGeometry.scv(dataJ.globalJ);
457 const auto internalDirichletConstraintsNeighbor = this->
problem().hasInternalDirichletConstraint(neighborElement, neighborScv);
459 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
461 if (internalDirichletConstraintsNeighbor[eqIdx])
462 A[dataJ.globalJ][globalI][eqIdx][pvIdx] = 0.0;
464 A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j][eqIdx];
472 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
475 A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx];
479 for (
const auto& dataJ : connectivityMap[globalI])
480 A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j++][eqIdx];
491 if (enableGridFluxVarsCache)
495 typename ParentType::LocalResidual::ElementResidualVector origElementResidual({origResiduals[0]});
496 this->
couplingManager().evalAdditionalDomainDerivatives(domainI, *
this, origElementResidual, A, gridVariables);
499 return origResiduals[0];
506 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
508 const LocalResidualValues& res, GridVariables& gridVariables)
517 const auto& gridGeometry =
fvGeometry.gridGeometry();
522 const auto globalI = gridGeometry.elementMapper().index(
element);
524 const auto& curSolJ = this->
curSol()[domainJ];
530 auto updateCoupledVariables = [&] ()
534 if (enableGridFluxVarsCache)
536 if (enableGridVolVarsCache)
538 this->
couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
549 if (enableGridVolVarsCache)
556 for (
const auto& globalJ : stencil)
559 const auto origPriVarsJ = curSolJ[globalJ];
560 auto priVarsJ = origPriVarsJ;
562 const auto origResidual = this->
couplingManager().evalCouplingResidual(domainI, *
this, domainJ, globalJ)[0];
564 for (
int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
566 auto evalCouplingResidual = [&](Scalar priVar)
569 priVarsJ[pvIdx] = priVar;
570 this->
couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
571 updateCoupledVariables();
572 return this->
couplingManager().evalCouplingResidual(domainI, *
this, domainJ, globalJ)[0];
576 LocalResidualValues partialDeriv(0.0);
577 const auto& paramGroup = this->
assembler().problem(domainJ).paramGroup();
579 static const auto epsCoupl = this->
couplingManager().numericEpsilon(domainJ, paramGroup);
581 epsCoupl(priVarsJ[pvIdx], pvIdx), numDiffMethod);
584 if constexpr (Problem::enableInternalDirichletConstraints())
587 const auto internalDirichletConstraints = this->
problem().hasInternalDirichletConstraint(this->
element(), scv);
588 if (internalDirichletConstraints.none())
590 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
591 A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
595 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
597 if (internalDirichletConstraints[eqIdx])
598 A[globalI][globalJ][eqIdx][pvIdx] = 0.0;
600 A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
606 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
607 A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
611 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
614 this->
couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
620 updateCoupledVariables();
631template<std::
size_t id,
class TypeTag,
class Assembler>
634 SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, false>, false>
644 static constexpr auto domainI = Dune::index_constant<id>();
647 using ParentType::ParentType;
659 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
662 if (this->
assembler().isStationaryProblem())
663 DUNE_THROW(Dune::InvalidStateException,
"Using explicit jacobian assembler with stationary local residual");
678 const auto& gridGeometry =
fvGeometry.gridGeometry();
682 const auto globalI = gridGeometry.elementMapper().index(
element);
689 const auto origPriVars =
curSol[globalI];
690 const auto origVolVars = curVolVars;
696 LocalResidualValues partialDeriv;
697 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
702 auto evalStorage = [&](Scalar priVar)
706 elemSol[0][pvIdx] = priVar;
717 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
723 else partialDeriv[pvIdx] = 1.0;
726 curVolVars = origVolVars;
729 elemSol[0][pvIdx] = origPriVars[pvIdx];
733 if constexpr (Problem::enableInternalDirichletConstraints())
736 const auto internalDirichletConstraints = this->
problem().hasInternalDirichletConstraint(this->
element(), scv);
737 const auto dirichletValues = this->
problem().internalDirichlet(this->
element(), scv);
739 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
741 if (internalDirichletConstraints[eqIdx])
743 residual[eqIdx] = origVolVars.priVars()[eqIdx] - dirichletValues[eqIdx];
744 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
747 A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx];
752 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
753 A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx];
767 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
769 const LocalResidualValues& res, GridVariables& gridVariables)
779template<std::
size_t id,
class TypeTag,
class Assembler>
782 SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::analytic, true>, true>
789 using Element =
typename GridView::template Codim<0>::Entity;
793 enum { dim = GridView::dimension };
795 static constexpr auto domainI = Dune::index_constant<id>();
798 using ParentType::ParentType;
806 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
812 const auto globalI = this->
assembler().gridGeometry(domainI).elementMapper().index(this->
element());
813 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
814 A[globalI][globalI][pvIdx][pvIdx] = 1.0;
817 return LocalResidualValues(0.0);
828 const auto globalI = this->
assembler().gridGeometry(domainI).elementMapper().index(
element);
833 if (!this->
assembler().isStationaryProblem())
843 if (!scvf.boundary())
852 if (bcTypes.hasDirichlet() && !bcTypes.hasNeumann())
856 else if (bcTypes.hasNeumann() && !bcTypes.hasDirichlet())
860 DUNE_THROW(Dune::NotImplemented,
"Mixed boundary conditions. Use pure boundary conditions by converting Dirichlet BCs to Robin BCs");
864 if constexpr (Problem::enableInternalDirichletConstraints())
867 const auto internalDirichletConstraints = this->
problem().hasInternalDirichletConstraint(this->
element(), scv);
868 const auto dirichletValues = this->
problem().internalDirichlet(this->
element(), scv);
872 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
874 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
876 if (internalDirichletConstraints[eqIdx])
878 residual[eqIdx] = volVars.priVars()[eqIdx] - dirichletValues[eqIdx];
879 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
883 if (!scvf.boundary())
884 A[globalI][
fvGeometry.scv(scvf.outsideScvIdx()).dofIndex()][eqIdx][pvIdx] = 0.0;
899 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
901 const LocalResidualValues& res, GridVariables& gridVariables)
910 const auto& gridGeometry =
fvGeometry.gridGeometry();
915 const auto globalI = gridGeometry.elementMapper().index(
element);
918 for (
const auto globalJ : stencil)
920 const auto& elementJ = this->
assembler().gridGeometry(domainJ).element(globalJ);
924 if constexpr (Problem::enableInternalDirichletConstraints())
927 const auto internalDirichletConstraints = this->
problem().hasInternalDirichletConstraint(this->
element(), scv);
928 if (internalDirichletConstraints.any())
930 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
931 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
932 if (internalDirichletConstraints[eqIdx])
933 A[globalI][globalJ][eqIdx][pvIdx] = 0.0;
A base class for all local assemblers.
The interface of the coupling manager for multi domain problems.
Definition multidomain/couplingmanager.hh:37
ElementVolumeVariables & curElemVolVars()
Definition assembly/fvlocalassemblerbase.hh:253
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition assembly/fvlocalassemblerbase.hh:56
Implementation & asImp_()
Definition assembly/fvlocalassemblerbase.hh:297
ElementVolumeVariables & prevElemVolVars()
Definition assembly/fvlocalassemblerbase.hh:257
ElementResidualVector evalLocalResidual() const
Definition assembly/fvlocalassemblerbase.hh:108
FVLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol)
Definition assembly/fvlocalassemblerbase.hh:61
FVElementGeometry & fvGeometry()
Definition assembly/fvlocalassemblerbase.hh:249
const Assembler & assembler() const
Definition assembly/fvlocalassemblerbase.hh:233
ElementFluxVariablesCache & elemFluxVarsCache()
Definition assembly/fvlocalassemblerbase.hh:261
bool elementIsGhost() const
Definition assembly/fvlocalassemblerbase.hh:241
std::decay_t< decltype(std::declval< Assembler >().localResidual())> LocalResidual
Definition assembly/fvlocalassemblerbase.hh:55
LocalResidual & localResidual()
Definition assembly/fvlocalassemblerbase.hh:265
const Element & element() const
Definition assembly/fvlocalassemblerbase.hh:237
VolumeVariables & getVolVarAccess(GridVolumeVariables &gridVolVars, ElementVolumeVariables &elemVolVars, const SubControlVolume &scv)
Definition assembly/fvlocalassemblerbase.hh:304
const SolutionVector & curSol() const
Definition assembly/fvlocalassemblerbase.hh:245
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
A helper class for local assemblers using numeric differentiation to determine the epsilon.
Definition numericepsilon.hh:29
A arithmetic block vector type based on DUNE's reserved vector.
Definition reservedblockvector.hh:26
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const LocalResidualValues &res, GridVariables &gridVariables)
Computes the derivatives of the residual of the given element with respect to primary variables of do...
Definition multidomain/subdomaincclocalassembler.hh:507
LocalResidualValues assembleJacobianAndResidualImplInverse(JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition multidomain/subdomaincclocalassembler.hh:313
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const LocalResidualValues &res, GridVariables &gridVariables)
Computes the derivatives of the residual of the given element with respect to primary variables of do...
Definition multidomain/subdomaincclocalassembler.hh:900
LocalResidualValues assembleJacobianAndResidualImpl(JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition multidomain/subdomaincclocalassembler.hh:807
LocalResidualValues assembleJacobianAndResidualImplInverse(JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition multidomain/subdomaincclocalassembler.hh:660
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const LocalResidualValues &res, GridVariables &gridVariables)
Computes the derivatives of the residual of the given element with respect to primary variables of do...
Definition multidomain/subdomaincclocalassembler.hh:768
typename ParentType::LocalResidual LocalResidual
the local residual type of this domain
Definition multidomain/subdomaincclocalassembler.hh:79
void bindLocalViews()
Prepares all local views necessary for local assembly.
Definition multidomain/subdomaincclocalassembler.hh:215
void assembleJacobianAndResidual(JacobianMatrixRow &jacRow, SubResidualVector &res, GridVariablesTuple &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition multidomain/subdomaincclocalassembler.hh:105
SubDomainCCLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, CouplingManager &couplingManager)
the constructor
Definition multidomain/subdomaincclocalassembler.hh:84
ElementResidualVector evalLocalSourceResidual(const Element &neighbor) const
Evaluates the local source term depending on time discretization scheme.
Definition multidomain/subdomaincclocalassembler.hh:191
LocalResidualValues evalFluxResidual(const Element &neighbor, const SubControlVolumeFace &scvf) const
Evaluates the fluxes depending on the chose time discretization scheme.
Definition multidomain/subdomaincclocalassembler.hh:205
void assembleResidual(SubResidualVector &res)
Assemble the residual only.
Definition multidomain/subdomaincclocalassembler.hh:147
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacRow &jacRow, const LocalResidualValues &res, GridVariables &gridVariables)
Assemble the entries in a coupling block of the jacobian. There is no coupling block between a domain...
Definition multidomain/subdomaincclocalassembler.hh:128
const Problem & problem() const
return reference to the underlying problem
Definition multidomain/subdomaincclocalassembler.hh:247
CouplingManager & couplingManager()
Definition multidomain/subdomaincclocalassembler.hh:251
static constexpr auto domainId
Definition multidomain/subdomaincclocalassembler.hh:77
ElementResidualVector evalLocalSourceResidual(const Element &element, const ElementVolumeVariables &elemVolVars) const
Evaluates the local source term for an element and given element volume variables.
Definition multidomain/subdomaincclocalassembler.hh:170
LocalResidualValues evalLocalStorageResidual() const
Evaluates the storage terms within the element.
Definition multidomain/subdomaincclocalassembler.hh:197
The cell-centered scheme multidomain local assembler.
Definition multidomain/subdomaincclocalassembler.hh:270
Defines all properties used in Dumux.
An enum class to define various differentiation methods available in order to compute the derivatives...
Element solution classes and factory functions.
Helper classes to compute the integration elements.
DiffMethod
Differentiation methods in order to compute the derivatives of the residual i.e. the entries in the j...
Definition diffmethod.hh:25
@ analytic
Definition diffmethod.hh:26
@ numeric
Definition diffmethod.hh:26
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
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition localview.hh:26
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
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
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition extrusion.hh:166
Definition common/pdesolver.hh:24
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.
A arithmetic block vector type based on DUNE's reserved vector.