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>
35#include "volvardeflectionhelper_.hh"
40namespace Detail::CVFE {
44 template<
class... Args>
45 constexpr void operator()(Args&&...)
const {}
48template<
class X,
class Y>
49using Impl = std::conditional_t<!std::is_same_v<X, void>, X, Y>;
63template<
class TypeTag,
class Assembler,
class Implementation,
bool implicit>
69 using ElementVolumeVariables =
typename GridVariables::GridVolumeVariables::LocalView;
70 using SolutionVector =
typename Assembler::SolutionVector;
77 using ParentType::ParentType;
90 template <
class Res
idualVector,
class PartialReassembler = DefaultPartialReassembler,
class CouplingFunction = Detail::CVFE::NoOperator>
93 const CouplingFunction& maybeAssembleCouplingBlocks = {})
95 this->
asImp_().bindLocalViews();
96 const auto eIdxGlobal = this->
asImp_().problem().gridGeometry().elementMapper().index(this->
element());
97 if (partialReassembler
100 const auto residual = this->
asImp_().evalLocalResidual();
101 for (
const auto& scv : scvs(this->
fvGeometry()))
102 res[scv.dofIndex()] += residual[scv.localDofIndex()];
105 maybeAssembleCouplingBlocks(residual);
109 const auto residual = this->
asImp_().assembleJacobianAndResidualImpl(jac, gridVariables, partialReassembler);
110 for (
const auto& scv : scvs(this->
fvGeometry()))
111 res[scv.dofIndex()] += residual[scv.localDofIndex()];
114 maybeAssembleCouplingBlocks(residual);
122 const auto& gridGeometry = this->
asImp_().problem().gridGeometry();
123 Dune::Hybrid::forEach(std::make_integer_sequence<int, dim>{}, [&](
auto d)
125 constexpr int codim = dim - d;
126 const auto& localCoeffs = gridGeometry.feCache().get(this->
element().type()).localCoefficients();
127 for (
int idx = 0; idx < localCoeffs.size(); ++idx)
129 const auto& localKey = localCoeffs.localKey(idx);
132 if (localKey.codim() != codim)
136 auto entity = this->
element().template subEntity<codim>(localKey.subEntity());
137 if (entity.partitionType() == Dune::InteriorEntity || entity.partitionType() == Dune::BorderEntity)
144 const auto dofIndex = gridGeometry.dofMapper().index(entity);
147 using BlockType =
typename JacobianMatrix::block_type;
148 BlockType &J = jac[dofIndex][dofIndex];
149 for (
int j = 0; j < BlockType::rows; ++j)
158 auto applyDirichlet = [&] (
const auto& scvI,
159 const auto& dirichletValues,
163 res[scvI.dofIndex()][eqIdx] = this->
curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
165 auto& row = jac[scvI.dofIndex()];
166 for (
auto col = row.begin(); col != row.end(); ++col)
167 row[col.index()][eqIdx] = 0.0;
169 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
172 if (this->
asImp_().
problem().gridGeometry().dofOnPeriodicBoundary(scvI.dofIndex()))
174 const auto periodicDof = this->
asImp_().problem().gridGeometry().periodicallyMappedDof(scvI.dofIndex());
175 res[periodicDof][eqIdx] = this->
curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
176 const auto end = jac[periodicDof].end();
177 for (
auto it = jac[periodicDof].begin(); it != end; ++it)
178 (*it) = periodicDof != it.index() ? 0.0 : 1.0;
182 this->
asImp_().enforceDirichletConstraints(applyDirichlet);
191 this->
asImp_().bindLocalViews();
192 this->
asImp_().assembleJacobianAndResidualImpl(jac, gridVariables);
194 auto applyDirichlet = [&] (
const auto& scvI,
195 const auto& dirichletValues,
199 auto& row = jac[scvI.dofIndex()];
200 for (
auto col = row.begin(); col != row.end(); ++col)
201 row[col.index()][eqIdx] = 0.0;
203 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
206 this->
asImp_().enforceDirichletConstraints(applyDirichlet);
212 template <
class Res
idualVector>
215 this->
asImp_().bindLocalViews();
218 for (
const auto& scv : scvs(this->
fvGeometry()))
219 res[scv.dofIndex()] += residual[scv.localDofIndex()];
221 auto applyDirichlet = [&] (
const auto& scvI,
222 const auto& dirichletValues,
226 res[scvI.dofIndex()][eqIdx] = this->
curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
229 this->
asImp_().enforceDirichletConstraints(applyDirichlet);
233 template<
typename ApplyFunction>
237 this->
asImp_().evalDirichletBoundaries(applyDirichlet);
239 this->
asImp_().enforceInternalDirichletConstraints(applyDirichlet);
245 template<
typename ApplyDirichletFunctionType >
252 for (
const auto& scvI : scvs(this->
fvGeometry()))
255 if (bcTypes.hasDirichlet())
257 const auto dirichletValues = this->
asImp_().problem().dirichlet(this->
element(), scvI);
260 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
262 if (bcTypes.isDirichlet(eqIdx))
264 const auto pvIdx = bcTypes.eqToDirichletIndex(eqIdx);
265 assert(0 <= pvIdx && pvIdx < numEq);
266 applyDirichlet(scvI, dirichletValues, eqIdx, pvIdx);
278 template<
class... Args>
285 template<
class... Args>
298template<
class TypeTag,
class Assembler, DiffMethod diffMethod = DiffMethod::numeric,
bool implicit = true,
class Implementation =
void>
306template<
class TypeTag,
class Assembler,
class Implementation>
309 Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true, Implementation>>,
322 static constexpr bool enableGridFluxVarsCache
323 = GridVariables::GridFluxVariablesCache::cachingEnabled;
324 static constexpr bool solutionDependentFluxVarsCache
325 = GridVariables::GridFluxVariablesCache::FluxVariablesCache::isSolDependent;
339 template <
class PartialReassembler = DefaultPartialReassembler>
344 const auto& element = this->element();
345 const auto& fvGeometry = this->fvGeometry();
346 const auto& curSol = this->asImp_().curSol();
348 auto&& curElemVolVars = this->curElemVolVars();
349 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
352 const auto origResiduals = this->evalLocalResidual();
364 static const bool updateAllVolVars = getParamFromGroup<bool>(
365 this->asImp_().problem().paramGroup(),
"Assembly.BoxVolVarsDependOnAllElementDofs",
false
369 auto elemSol =
elementSolution(element, curSol, fvGeometry.gridGeometry());
374 Detail::VolVarsDeflectionHelper deflectionHelper(
375 [&] (
const auto& scv) -> VolumeVariables& {
376 return this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
383 for (
const auto& scv : scvs(fvGeometry))
386 const auto dofIdx = scv.dofIndex();
387 deflectionHelper.setCurrent(scv);
390 for (
int pvIdx = 0; pvIdx < numEq; pvIdx++)
394 auto evalResiduals = [&](Scalar priVar)
397 elemSol[scv.localDofIndex()][pvIdx] = priVar;
398 deflectionHelper.deflect(elemSol, scv, this->asImp_().problem());
399 if constexpr (solutionDependentFluxVarsCache)
401 elemFluxVarsCache.update(element, fvGeometry, curElemVolVars);
402 if constexpr (enableGridFluxVarsCache)
403 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
405 this->asImp_().maybeUpdateCouplingContext(scv, elemSol, pvIdx);
406 return this->evalLocalResidual();
411 static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(),
"Assembly.NumericDifferenceMethod");
413 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
416 for (
const auto& scvJ : scvs(fvGeometry))
419 if (!partialReassembler
422 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
428 A[scvJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[scvJ.localDofIndex()][eqIdx];
434 deflectionHelper.restore(scv);
437 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
438 this->asImp_().maybeUpdateCouplingContext(scv, elemSol, pvIdx);
444 if constexpr (enableGridFluxVarsCache)
445 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
448 this->asImp_().maybeEvalAdditionalDomainDerivatives(origResiduals, A, gridVariables);
450 return origResiduals;
460template<
class TypeTag,
class Assembler,
class Implementation>
463 Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false, Implementation>>,
488 template <
class PartialReassembler = DefaultPartialReassembler>
492 if (partialReassembler)
493 DUNE_THROW(Dune::NotImplemented,
"partial reassembly for explicit time discretization");
496 const auto& element = this->element();
497 const auto& fvGeometry = this->fvGeometry();
498 const auto& curSol = this->asImp_().curSol();
499 auto&& curElemVolVars = this->curElemVolVars();
502 const auto origResiduals = this->evalLocalResidual();
503 const auto origStorageResiduals = this->evalLocalStorageResidual();
514 auto elemSol =
elementSolution(element, curSol, fvGeometry.gridGeometry());
520 for (
const auto& scv : scvs(fvGeometry))
523 const auto dofIdx = scv.dofIndex();
524 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
525 const VolumeVariables origVolVars(curVolVars);
528 for (
int pvIdx = 0; pvIdx < numEq; pvIdx++)
532 auto evalStorage = [&](Scalar priVar)
535 elemSol[scv.localDofIndex()][pvIdx] = priVar;
536 curVolVars.update(elemSol, this->asImp_().problem(), element, scv);
537 return this->evalLocalStorageResidual();
542 static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(),
"Assembly.NumericDifferenceMethod");
544 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
547 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
553 A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
557 curVolVars = origVolVars;
560 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
565 return origResiduals;
574template<
class TypeTag,
class Assembler,
class Implementation>
577 Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true, Implementation>>,
597 template <
class PartialReassembler = DefaultPartialReassembler>
601 if (partialReassembler)
602 DUNE_THROW(Dune::NotImplemented,
"partial reassembly for analytic differentiation");
605 const auto& element = this->element();
606 const auto& fvGeometry = this->fvGeometry();
607 const auto& problem = this->asImp_().problem();
608 const auto& curElemVolVars = this->curElemVolVars();
609 const auto& elemFluxVarsCache = this->elemFluxVarsCache();
612 const auto origResiduals = this->evalLocalResidual();
623 for (
const auto& scv : scvs(fvGeometry))
626 const auto dofIdx = scv.dofIndex();
627 const auto& volVars = curElemVolVars[scv];
632 if (!this->assembler().isStationaryProblem())
633 this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
642 this->localResidual().addSourceDerivatives(A[dofIdx][dofIdx],
651 for (
const auto& scvf : scvfs(fvGeometry))
653 if (!scvf.boundary())
656 this->localResidual().addFluxDerivatives(A,
669 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
670 if (this->elemBcTypes().get(fvGeometry, insideScv).hasNeumann())
673 this->localResidual().addRobinFluxDerivatives(A[insideScv.dofIndex()],
684 return origResiduals;
694template<
class TypeTag,
class Assembler,
class Implementation>
697 Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false, Implementation>>,
717 template <
class PartialReassembler = DefaultPartialReassembler>
721 if (partialReassembler)
722 DUNE_THROW(Dune::NotImplemented,
"partial reassembly for explicit time discretization");
725 const auto& element = this->element();
726 const auto& fvGeometry = this->fvGeometry();
727 const auto& problem = this->asImp_().problem();
728 const auto& curElemVolVars = this->curElemVolVars();
731 const auto origResiduals = this->evalLocalResidual();
742 for (
const auto& scv : scvs(fvGeometry))
745 const auto dofIdx = scv.dofIndex();
746 const auto& volVars = curElemVolVars[scv];
750 this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
758 return origResiduals;
Control volume finite element local assembler using analytic differentiation and explicit time discre...
Definition: assembly/cvfelocalassembler.hh:699
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: assembly/cvfelocalassembler.hh:708
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:718
typename ParentType::LocalResidual LocalResidual
Definition: assembly/cvfelocalassembler.hh:707
Control volume finite element local assembler using analytic differentiation and implicit time discre...
Definition: assembly/cvfelocalassembler.hh:579
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:598
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: assembly/cvfelocalassembler.hh:588
typename ParentType::LocalResidual LocalResidual
Definition: assembly/cvfelocalassembler.hh:587
Control volume finite element local assembler using numeric differentiation and implicit time discret...
Definition: assembly/cvfelocalassembler.hh:311
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: assembly/cvfelocalassembler.hh:330
typename ParentType::LocalResidual LocalResidual
Definition: assembly/cvfelocalassembler.hh:329
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:340
Control volume finite element local assembler using numeric differentiation and explicit time discret...
Definition: assembly/cvfelocalassembler.hh:465
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: assembly/cvfelocalassembler.hh:479
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:489
typename ParentType::LocalResidual LocalResidual
Definition: assembly/cvfelocalassembler.hh:478
A base class for all local CVFE assemblers.
Definition: assembly/cvfelocalassembler.hh:65
void enforceDirichletConstraints(const ApplyFunction &applyDirichlet)
Enforce Dirichlet constraints.
Definition: assembly/cvfelocalassembler.hh:234
void assembleResidual(ResidualVector &res)
Assemble the residual only.
Definition: assembly/cvfelocalassembler.hh:213
void maybeUpdateCouplingContext(Args &&...)
Update the coupling context for coupled models.
Definition: assembly/cvfelocalassembler.hh:279
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:189
void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
Evaluates Dirichlet boundaries.
Definition: assembly/cvfelocalassembler.hh:246
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:91
void bindLocalViews()
Definition: assembly/cvfelocalassembler.hh:79
void maybeEvalAdditionalDomainDerivatives(Args &&...)
Update the additional domain derivatives for coupled models.
Definition: assembly/cvfelocalassembler.hh:286
An assembler for Jacobian and residual contribution per element (CVFE methods)
Definition: assembly/cvfelocalassembler.hh:299
A base class for all local assemblers.
Definition: assembly/fvlocalassemblerbase.hh:36
void bindLocalViews()
Convenience function bind and prepare all relevant variables required for the evaluation of the local...
Definition: assembly/fvlocalassemblerbase.hh:173
ElementVolumeVariables & curElemVolVars()
The current element volume variables.
Definition: assembly/fvlocalassemblerbase.hh:253
ElementBoundaryTypes & elemBcTypes()
The element's boundary types.
Definition: assembly/fvlocalassemblerbase.hh:269
Implementation & asImp_()
Definition: assembly/fvlocalassemblerbase.hh:297
ElementResidualVector evalLocalResidual() const
Convenience function to evaluate the complete local residual for the current element....
Definition: assembly/fvlocalassemblerbase.hh:108
const Problem & problem() const
The problem.
Definition: assembly/fvlocalassemblerbase.hh:229
FVElementGeometry & fvGeometry()
The global finite volume geometry.
Definition: assembly/fvlocalassemblerbase.hh:249
bool elementIsGhost() const
Returns if element is a ghost entity.
Definition: assembly/fvlocalassemblerbase.hh:241
std::decay_t< decltype(std::declval< Assembler >().localResidual())> LocalResidual
Definition: assembly/fvlocalassemblerbase.hh:55
const Element & element() const
The current element.
Definition: assembly/fvlocalassemblerbase.hh:237
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:49
A helper class for local assemblers using numeric differentiation to determine the epsilon.
Definition: numericepsilon.hh:29
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
@ green
does not need to be reassembled
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
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
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.