14#ifndef DUMUX_EXPERIMENTAL_CVFE_LOCAL_ASSEMBLER_HH
15#define DUMUX_EXPERIMENTAL_CVFE_LOCAL_ASSEMBLER_HH
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>
29#include <dumux/common/typetraits/localdofs_.hh>
40#include <dumux/assembly/cvfevolvarsdeflectionpolicy_.hh>
53template<
class TypeTag,
class Assembler,
class Implementation>
60 using ElementVolumeVariables =
typename GridVariables::GridVolumeVariables::LocalView;
68 using ParentType::ParentType;
84 template <
class Res
idualVector,
class StageParams,
class PartialReassembler = DefaultPartialReassembler,
class CouplingFunction = Noop>
86 const StageParams& stageParams, ResidualVector& temporal, ResidualVector& spatial,
87 ResidualVector& constrainedDofs,
89 const CouplingFunction& maybeAssembleCouplingBlocks =
noop)
91 this->
asImp_().bindLocalViews();
92 const auto eIdxGlobal = this->
asImp_().problem().gridGeometry().elementMapper().index(this->
element());
97 const auto sWeight = stageParams.spatialWeight(stageParams.size()-1);
98 const auto tWeight = stageParams.temporalWeight(stageParams.size()-1);
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()];
115 if (partialReassembler && partialReassembler->elementColor(eIdxGlobal) ==
EntityColor::green)
118 maybeAssembleCouplingBlocks(origResidual);
122 this->
asImp_().assembleJacobian(jac, gridVariables, origResidual, partialReassembler);
125 maybeAssembleCouplingBlocks(origResidual);
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)
161 constrainedDofs[dofIndex] = 1;
167 auto applyDirichlet = [&] (
const auto& scvI,
168 const auto& dirichletValues,
172 res[scvI.dofIndex()][eqIdx] = this->
curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
173 constrainedDofs[scvI.dofIndex()][eqIdx] = 1;
175 auto& row = jac[scvI.dofIndex()];
176 for (
auto col = row.begin(); col != row.end(); ++col)
177 row[col.index()][eqIdx] = 0.0;
179 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
182 if (this->
asImp_().
problem().gridGeometry().dofOnPeriodicBoundary(scvI.dofIndex()))
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;
193 this->
asImp_().enforceDirichletConstraints(applyDirichlet);
202 this->
asImp_().bindLocalViews();
203 this->
asImp_().assembleJacobian(jac, gridVariables);
205 auto applyDirichlet = [&] (
const auto& scvI,
206 const auto& dirichletValues,
210 auto& row = jac[scvI.dofIndex()];
211 for (
auto col = row.begin(); col != row.end(); ++col)
212 row[col.index()][eqIdx] = 0.0;
214 jac[scvI.dofIndex()][scvI.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& scvI,
233 const auto& dirichletValues,
237 res[scvI.dofIndex()][eqIdx] = this->
curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
240 this->
asImp_().enforceDirichletConstraints(applyDirichlet);
246 template<
class Res
idualVector>
249 this->
asImp_().bindLocalViews();
254 spatialRes[scv.dofIndex()] += flux[scv.localDofIndex()];
255 temporalRes[scv.dofIndex()] += storage[scv.localDofIndex()];
260 template<
typename ApplyFunction>
264 this->
asImp_().evalDirichletBoundaries(applyDirichlet);
266 this->
asImp_().enforceInternalDirichletConstraints(applyDirichlet);
272 template<
typename ApplyDirichletFunctionType >
282 if (bcTypes.hasDirichlet())
284 const auto dirichletValues = this->
asImp_().problem().dirichlet(this->
element(), scvI);
287 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
289 if (bcTypes.isDirichlet(eqIdx))
291 const auto pvIdx = bcTypes.eqToDirichletIndex(eqIdx);
292 assert(0 <= pvIdx && pvIdx < numEq);
293 applyDirichlet(scvI, dirichletValues, eqIdx, pvIdx);
305 template<
class... Args>
312 template<
class... Args>
325template<
class TypeTag,
class Assembler, DiffMethod diffMethod = DiffMethod::numeric,
class Implementation =
void>
334template<
class TypeTag,
class Assembler,
class Implementation>
337 NonVoidOr<CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, Implementation>, Implementation>>
343 using ElementVolumeVariables =
typename GridVariables::GridVolumeVariables::LocalView;
350 static constexpr bool enableGridFluxVarsCache
351 = GridVariables::GridFluxVariablesCache::cachingEnabled;
352 static constexpr bool solutionDependentFluxVarsCache
353 = GridVariables::GridFluxVariablesCache::FluxVariablesCache::isSolDependent;
359 using ParentType::ParentType;
361 template <
class PartialReassembler = DefaultPartialReassembler>
367 assembleJacobianImplicit_(A, gridVariables, origResiduals, partialReassembler);
369 assembleJacobianExplicit_(A, gridVariables, origResiduals, partialReassembler);
379 template <
class PartialReassembler = DefaultPartialReassembler>
380 void assembleJacobianImplicit_(JacobianMatrix& A,
GridVariables& gridVariables,
381 const ElementResidualVector& origResiduals,
385 const auto& element = this->element();
386 const auto& fvGeometry = this->fvGeometry();
387 const auto& curSol = this->asImp_().curSol();
389 auto&& curElemVolVars = this->curElemVolVars();
390 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
403 this->asImp_().problem().paramGroup(),
"Assembly.BoxVolVarsDependOnAllElementDofs",
false
407 auto elemSol =
elementSolution(element, curSol, fvGeometry.gridGeometry());
410 ElementResidualVector partialDerivs(Dumux::Detail::LocalDofs::numLocalDofs(fvGeometry));
412 auto deflectionPolicy = Dumux::Detail::CVFE::makeVariablesDeflectionPolicy(
413 gridVariables.curGridVolVars(),
419 auto assembleDerivative = [&,
this](
const auto& localDof)
422 const auto dofIdx = localDof.dofIndex();
423 const auto localIdx = localDof.index();
424 deflectionPolicy.store(localDof);
427 for (
int pvIdx = 0; pvIdx < numEq; pvIdx++)
431 auto evalResiduals = [&](Scalar priVar)
434 elemSol[localIdx][pvIdx] = priVar;
435 deflectionPolicy.update(elemSol, localDof, this->asImp_().problem());
436 if constexpr (solutionDependentFluxVarsCache)
438 elemFluxVarsCache.update(element, fvGeometry, curElemVolVars);
439 if constexpr (enableGridFluxVarsCache)
440 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
442 this->asImp_().maybeUpdateCouplingContext(localDof, elemSol, pvIdx);
443 return this->evalLocalResidual();
447 static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()};
448 static const int numDiffMethod =
getParamFromGroup<int>(this->asImp_().problem().paramGroup(),
"Assembly.NumericDifferenceMethod");
450 eps_(elemSol[localIdx][pvIdx], pvIdx), numDiffMethod);
453 for (
const auto& localDofJ :
localDofs(fvGeometry))
456 if (!partialReassembler
459 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
465 A[localDofJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[localDofJ.index()][eqIdx];
471 deflectionPolicy.restore(localDof);
474 elemSol[localIdx][pvIdx] = curSol[localDof.dofIndex()][pvIdx];
475 this->asImp_().maybeUpdateCouplingContext(localDof, elemSol, pvIdx);
480 for (
const auto& localDof :
localDofs(fvGeometry))
481 assembleDerivative(localDof);
485 if constexpr (enableGridFluxVarsCache)
486 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
489 this->asImp_().maybeEvalAdditionalDomainDerivatives(origResiduals, A, gridVariables);
497 template <
class PartialReassembler = DefaultPartialReassembler>
498 void assembleJacobianExplicit_(JacobianMatrix& A, GridVariables& gridVariables,
499 const ElementResidualVector& origResiduals,
500 const PartialReassembler* partialReassembler =
nullptr)
502 if (partialReassembler)
503 DUNE_THROW(Dune::NotImplemented,
"partial reassembly for explicit time discretization");
507 const auto& fvGeometry = this->fvGeometry();
508 const auto& curSol = this->asImp_().curSol();
509 auto&& curElemVolVars = this->curElemVolVars();
512 auto elemSol =
elementSolution(element, curSol, fvGeometry.gridGeometry());
515 ElementResidualVector partialDerivs(Dumux::Detail::LocalDofs::numLocalDofs(fvGeometry));
518 for (
const auto& scv :
scvs(fvGeometry))
521 const auto dofIdx = scv.dofIndex();
522 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
523 const VolumeVariables origVolVars(curVolVars);
526 for (
int pvIdx = 0; pvIdx < numEq; pvIdx++)
530 auto evalStorage = [&](Scalar priVar)
532 elemSol[scv.localDofIndex()][pvIdx] = priVar;
533 curVolVars.update(elemSol, this->asImp_().problem(), element, scv);
534 return this->evalStorage();
538 static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()};
539 static const int numDiffMethod =
getParamFromGroup<int>(this->asImp_().problem().paramGroup(),
"Assembly.NumericDifferenceMethod");
541 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
544 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
550 A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
554 curVolVars = origVolVars;
557 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
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
bool isImplicit() const
Definition experimental/assembly/fvlocalassemblerbase.hh:270
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.
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
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.