13#ifndef DUMUX_CVFE_LOCAL_ASSEMBLER_HH
14#define DUMUX_CVFE_LOCAL_ASSEMBLER_HH
16#include <dune/common/exceptions.hh>
17#include <dune/common/hybridutilities.hh>
18#include <dune/common/reservedvector.hh>
19#include <dune/grid/common/gridenums.hh>
20#include <dune/istl/matrixindexset.hh>
21#include <dune/istl/bvector.hh>
27#include <dumux/common/typetraits/localdofs_.hh>
28#include <dumux/common/typetraits/boundary_.hh>
39#include "cvfevolvarsdeflectionpolicy_.hh"
48 template<
class... Args>
49 constexpr void operator()(Args&&...)
const {}
52template<
class X,
class Y>
53using Impl = std::conditional_t<!std::is_same_v<X, void>, X, Y>;
67template<
class TypeTag,
class Assembler,
class Implementation,
bool implicit>
74 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
75 using ElementVolumeVariables =
typename GridVariables::GridVolumeVariables::LocalView;
77 using SolutionVector =
typename Assembler::SolutionVector;
79 using FVElementGeometry =
typename GridGeometry::LocalView;
82 static constexpr int dim = GridGeometry::GridView::dimension;
86 using ParentType::ParentType;
99 template <
class Res
idualVector,
class PartialReassembler = DefaultPartialReassembler,
class CouplingFunction = Detail::CVFE::NoOperator>
102 const CouplingFunction& maybeAssembleCouplingBlocks = {})
104 this->
asImp_().bindLocalViews();
105 const auto eIdxGlobal = this->
asImp_().problem().gridGeometry().elementMapper().index(this->
element());
106 if (partialReassembler
109 const auto residual = this->
asImp_().evalLocalResidual();
112 res[localDof.dofIndex()] += residual[localDof.index()];
115 maybeAssembleCouplingBlocks(residual);
119 const auto residual = this->
asImp_().assembleJacobianAndResidualImpl(jac, gridVariables, partialReassembler);
122 res[localDof.dofIndex()] += residual[localDof.index()];
125 maybeAssembleCouplingBlocks(residual);
133 const auto& gridGeometry = this->
asImp_().problem().gridGeometry();
134 Dune::Hybrid::forEach(std::make_integer_sequence<int, dim+1>{}, [&](
auto d)
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)
140 const auto& localKey = localCoeffs.localKey(idx);
143 if (localKey.codim() != codim)
147 auto entity = this->
element().template subEntity<codim>(localKey.subEntity());
148 if (entity.partitionType() == Dune::InteriorEntity || entity.partitionType() == Dune::BorderEntity)
154 using BlockType =
typename JacobianMatrix::block_type;
155 for (
const auto dofIndex :
asMultiMapper(gridGeometry.dofMapper()).indices(entity))
157 BlockType &J = jac[dofIndex][dofIndex];
158 for (
int j = 0; j < BlockType::rows; ++j)
166 auto applyDirichlet = [&] (
const auto& scvOrLocalDofI,
167 const auto& dirichletValues,
171 res[scvOrLocalDofI.dofIndex()][eqIdx] = this->
curElemVolVars()[scvOrLocalDofI].priVars()[pvIdx] - dirichletValues[pvIdx];
173 auto& row = jac[scvOrLocalDofI.dofIndex()];
174 for (
auto col = row.begin(); col != row.end(); ++col)
175 row[col.index()][eqIdx] = 0.0;
177 jac[scvOrLocalDofI.dofIndex()][scvOrLocalDofI.dofIndex()][eqIdx][pvIdx] = 1.0;
180 if (this->
asImp_().
problem().gridGeometry().dofOnPeriodicBoundary(scvOrLocalDofI.dofIndex()))
182 const auto periodicDof = this->
asImp_().problem().gridGeometry().periodicallyMappedDof(scvOrLocalDofI.dofIndex());
183 res[periodicDof][eqIdx] = this->
asImp_().curSol()[periodicDof][pvIdx] - dirichletValues[pvIdx];
185 auto& rowP = jac[periodicDof];
186 for (
auto col = rowP.begin(); col != rowP.end(); ++col)
187 rowP[col.index()][eqIdx] = 0.0;
189 rowP[periodicDof][eqIdx][pvIdx] = 1.0;
193 this->
asImp_().enforceDirichletConstraints(applyDirichlet);
202 this->
asImp_().bindLocalViews();
203 this->
asImp_().assembleJacobianAndResidualImpl(jac, gridVariables);
205 auto applyDirichlet = [&] (
const auto& scvOrLocalDofI,
206 const auto& dirichletValues,
210 auto& row = jac[scvOrLocalDofI.dofIndex()];
211 for (
auto col = row.begin(); col != row.end(); ++col)
212 row[col.index()][eqIdx] = 0.0;
214 jac[scvOrLocalDofI.dofIndex()][scvOrLocalDofI.dofIndex()][eqIdx][pvIdx] = 1.0;
217 this->
asImp_().enforceDirichletConstraints(applyDirichlet);
223 template <
class Res
idualVector>
226 this->
asImp_().bindLocalViews();
230 res[localDof.dofIndex()] += residual[localDof.index()];
232 auto applyDirichlet = [&] (
const auto& scvOrLocalDofI,
233 const auto& dirichletValues,
237 res[scvOrLocalDofI.dofIndex()][eqIdx] = this->
curElemVolVars()[scvOrLocalDofI].priVars()[pvIdx] - dirichletValues[pvIdx];
240 this->
asImp_().enforceDirichletConstraints(applyDirichlet);
244 template<
typename ApplyFunction>
248 this->
asImp_().evalDirichletBoundaries(applyDirichlet);
250 this->
asImp_().enforceInternalDirichletConstraints(applyDirichlet);
256 template<
typename ApplyDirichletFunctionType >
262 if constexpr (!Detail::hasProblemBoundaryTypesForFaceFunction<Problem, FVElementGeometry>())
269 if (bcTypes.hasDirichlet())
271 const auto dirichletValues = this->
asImp_().problem().dirichlet(this->
element(), scvI);
274 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
276 if (bcTypes.isDirichlet(eqIdx))
278 const auto pvIdx = bcTypes.eqToDirichletIndex(eqIdx);
279 assert(0 <= pvIdx && pvIdx < numEq);
280 applyDirichlet(scvI, dirichletValues, eqIdx, pvIdx);
293 template<
class... Args>
300 template<
class... Args>
314template<
class TypeTag,
class Assembler, DiffMethod diffMethod = DiffMethod::numeric,
bool implicit = true,
class Implementation =
void>
322template<
class TypeTag,
class Assembler,
class Implementation>
325 Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true, Implementation>>,
333 using ElementVolumeVariables =
typename GridVariables::GridVolumeVariables::LocalView;
340 static constexpr bool enableGridFluxVarsCache
341 = GridVariables::GridFluxVariablesCache::cachingEnabled;
342 static constexpr bool solutionDependentFluxVarsCache
343 = GridVariables::GridFluxVariablesCache::FluxVariablesCache::isSolDependent;
349 using ParentType::ParentType;
357 template <
class PartialReassembler = DefaultPartialReassembler>
383 this->
asImp_().
problem().paramGroup(),
"Assembly.BoxVolVarsDependOnAllElementDofs",
false
392 auto deflectionPolicy = Detail::CVFE::makeVariablesDeflectionPolicy(
393 gridVariables.curGridVolVars(),
399 auto assembleDerivative = [&,
this](
const auto& localDof)
402 const auto dofIdx = localDof.dofIndex();
403 const auto localIdx = localDof.index();
404 deflectionPolicy.store(localDof);
407 for (
int pvIdx = 0; pvIdx < numEq; pvIdx++)
411 auto evalResiduals = [&](PrimaryVariable priVar)
414 elemSol[localIdx][pvIdx] = priVar;
415 deflectionPolicy.update(elemSol, localDof, this->
asImp_().
problem());
416 if constexpr (solutionDependentFluxVarsCache)
419 if constexpr (enableGridFluxVarsCache)
422 this->
asImp_().maybeUpdateCouplingContext(localDof, elemSol, pvIdx);
430 eps_(elemSol[localIdx][pvIdx], pvIdx), numDiffMethod);
436 if (!partialReassembler
439 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
445 A[localDofJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[localDofJ.index()][eqIdx];
451 deflectionPolicy.restore(localDof);
454 elemSol[localIdx][pvIdx] =
curSol[localDof.dofIndex()][pvIdx];
455 this->
asImp_().maybeUpdateCouplingContext(localDof, elemSol, pvIdx);
461 assembleDerivative(localDof);
465 if constexpr (enableGridFluxVarsCache)
469 this->
asImp_().maybeEvalAdditionalDomainDerivatives(origResiduals, A, gridVariables);
471 return origResiduals;
481template<
class TypeTag,
class Assembler,
class Implementation>
484 Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false, Implementation>>,
501 using ParentType::ParentType;
509 template <
class PartialReassembler = DefaultPartialReassembler>
513 if (partialReassembler)
514 DUNE_THROW(Dune::NotImplemented,
"partial reassembly for explicit time discretization");
544 const auto localIdx = scv.localDofIndex();
545 const auto dofIdx = scv.dofIndex();
547 const VolumeVariables origVolVars(curVolVars);
550 for (
int pvIdx = 0; pvIdx < numEq; pvIdx++)
554 auto evalStorage = [&](Scalar priVar)
557 elemSol[localIdx][pvIdx] = priVar;
566 eps_(elemSol[localIdx][pvIdx], pvIdx), numDiffMethod);
569 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
575 A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[localIdx][eqIdx];
579 curVolVars = origVolVars;
582 elemSol[localIdx][pvIdx] =
curSol[dofIdx][pvIdx];
587 return origResiduals;
596template<
class TypeTag,
class Assembler,
class Implementation>
599 Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true, Implementation>>,
611 using ParentType::ParentType;
619 template <
class PartialReassembler = DefaultPartialReassembler>
623 if (partialReassembler)
624 DUNE_THROW(Dune::NotImplemented,
"partial reassembly for analytic differentiation");
648 const auto dofIdx = scv.dofIndex();
654 if (!this->
assembler().isStationaryProblem())
655 this->
localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
664 this->
localResidual().addSourceDerivatives(A[dofIdx][dofIdx],
675 if (!scvf.boundary())
691 const auto& insideScv =
fvGeometry.scv(scvf.insideScvIdx());
695 this->
localResidual().addRobinFluxDerivatives(A[insideScv.dofIndex()],
706 return origResiduals;
716template<
class TypeTag,
class Assembler,
class Implementation>
719 Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false, Implementation>>,
731 using ParentType::ParentType;
739 template <
class PartialReassembler = DefaultPartialReassembler>
743 if (partialReassembler)
744 DUNE_THROW(Dune::NotImplemented,
"partial reassembly for explicit time discretization");
767 const auto dofIdx = scv.dofIndex();
772 this->
localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
780 return origResiduals;
A base class for all local assemblers.
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition assembly/cvfelocalassembler.hh:730
ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix &A, GridVariables &gridVariables, const PartialReassembler *partialReassembler=nullptr)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition assembly/cvfelocalassembler.hh:740
typename ParentType::LocalResidual LocalResidual
Definition assembly/cvfelocalassembler.hh:729
ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix &A, GridVariables &gridVariables, const PartialReassembler *partialReassembler=nullptr)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition assembly/cvfelocalassembler.hh:620
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition assembly/cvfelocalassembler.hh:610
typename ParentType::LocalResidual LocalResidual
Definition assembly/cvfelocalassembler.hh:609
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition assembly/cvfelocalassembler.hh:348
typename ParentType::LocalResidual LocalResidual
Definition assembly/cvfelocalassembler.hh:347
ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix &A, GridVariables &gridVariables, const PartialReassembler *partialReassembler=nullptr)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition assembly/cvfelocalassembler.hh:358
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition assembly/cvfelocalassembler.hh:500
ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix &A, GridVariables &gridVariables, const PartialReassembler *partialReassembler=nullptr)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition assembly/cvfelocalassembler.hh:510
typename ParentType::LocalResidual LocalResidual
Definition assembly/cvfelocalassembler.hh:499
A base class for all local CVFE assemblers.
Definition assembly/cvfelocalassembler.hh:69
void enforceDirichletConstraints(const ApplyFunction &applyDirichlet)
Enforce Dirichlet constraints.
Definition assembly/cvfelocalassembler.hh:245
void assembleResidual(ResidualVector &res)
Assemble the residual only.
Definition assembly/cvfelocalassembler.hh:224
void maybeUpdateCouplingContext(Args &&...)
Update the coupling context for coupled models.
Definition assembly/cvfelocalassembler.hh:294
void assembleJacobian(JacobianMatrix &jac, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition assembly/cvfelocalassembler.hh:200
void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
Evaluates Dirichlet boundaries.
Definition assembly/cvfelocalassembler.hh:257
void assembleJacobianAndResidual(JacobianMatrix &jac, ResidualVector &res, GridVariables &gridVariables, const PartialReassembler *partialReassembler=nullptr, const CouplingFunction &maybeAssembleCouplingBlocks={})
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition assembly/cvfelocalassembler.hh:100
void bindLocalViews()
Definition assembly/cvfelocalassembler.hh:88
void maybeEvalAdditionalDomainDerivatives(Args &&...)
Update the additional domain derivatives for coupled models.
Definition assembly/cvfelocalassembler.hh:301
An assembler for Jacobian and residual contribution per element (CVFE methods).
Definition assembly/cvfelocalassembler.hh:315
void bindLocalViews()
Definition assembly/fvlocalassemblerbase.hh:173
ElementVolumeVariables & curElemVolVars()
Definition assembly/fvlocalassemblerbase.hh:253
ElementBoundaryTypes & elemBcTypes()
Definition assembly/fvlocalassemblerbase.hh:269
Implementation & asImp_()
Definition assembly/fvlocalassemblerbase.hh:297
ElementResidualVector evalLocalResidual() const
Definition assembly/fvlocalassemblerbase.hh:108
const Problem & problem() const
Definition assembly/fvlocalassemblerbase.hh:229
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
ElementResidualVector evalLocalStorageResidual() const
Definition assembly/fvlocalassemblerbase.hh:164
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:32
detects which entries in the Jacobian have to be recomputed
Definition partialreassembler.hh:420
Defines all properties used in Dumux.
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
@ analytic
Definition diffmethod.hh:26
@ 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
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 elementvariables.hh:24
constexpr auto asMultiMapper(const Mapper &mapper)
Definition multimapperview.hh:48
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.