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->
asImp_().curSol()[periodicDof][pvIdx] - dirichletValues[pvIdx];
177 auto& rowP = jac[periodicDof];
178 for (
auto col = rowP.begin(); col != rowP.end(); ++col)
179 rowP[col.index()][eqIdx] = 0.0;
181 rowP[periodicDof][eqIdx][pvIdx] = 1.0;
185 this->
asImp_().enforceDirichletConstraints(applyDirichlet);
194 this->
asImp_().bindLocalViews();
195 this->
asImp_().assembleJacobianAndResidualImpl(jac, gridVariables);
197 auto applyDirichlet = [&] (
const auto& scvI,
198 const auto& dirichletValues,
202 auto& row = jac[scvI.dofIndex()];
203 for (
auto col = row.begin(); col != row.end(); ++col)
204 row[col.index()][eqIdx] = 0.0;
206 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
209 this->
asImp_().enforceDirichletConstraints(applyDirichlet);
215 template <
class Res
idualVector>
218 this->
asImp_().bindLocalViews();
221 for (
const auto& scv : scvs(this->
fvGeometry()))
222 res[scv.dofIndex()] += residual[scv.localDofIndex()];
224 auto applyDirichlet = [&] (
const auto& scvI,
225 const auto& dirichletValues,
229 res[scvI.dofIndex()][eqIdx] = this->
curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
232 this->
asImp_().enforceDirichletConstraints(applyDirichlet);
236 template<
typename ApplyFunction>
240 this->
asImp_().evalDirichletBoundaries(applyDirichlet);
242 this->
asImp_().enforceInternalDirichletConstraints(applyDirichlet);
248 template<
typename ApplyDirichletFunctionType >
255 for (
const auto& scvI : scvs(this->
fvGeometry()))
258 if (bcTypes.hasDirichlet())
260 const auto dirichletValues = this->
asImp_().problem().dirichlet(this->
element(), scvI);
263 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
265 if (bcTypes.isDirichlet(eqIdx))
267 const auto pvIdx = bcTypes.eqToDirichletIndex(eqIdx);
268 assert(0 <= pvIdx && pvIdx < numEq);
269 applyDirichlet(scvI, dirichletValues, eqIdx, pvIdx);
281 template<
class... Args>
288 template<
class... Args>
301template<
class TypeTag,
class Assembler, DiffMethod diffMethod = DiffMethod::numeric,
bool implicit = true,
class Implementation =
void>
309template<
class TypeTag,
class Assembler,
class Implementation>
312 Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true, Implementation>>,
325 static constexpr bool enableGridFluxVarsCache
326 = GridVariables::GridFluxVariablesCache::cachingEnabled;
327 static constexpr bool solutionDependentFluxVarsCache
328 = GridVariables::GridFluxVariablesCache::FluxVariablesCache::isSolDependent;
334 using ParentType::ParentType;
342 template <
class PartialReassembler = DefaultPartialReassembler>
368 this->
asImp_().
problem().paramGroup(),
"Assembly.BoxVolVarsDependOnAllElementDofs",
false
377 Detail::VolVarsDeflectionHelper deflectionHelper(
378 [&] (
const auto& scv) -> VolumeVariables& {
389 const auto dofIdx = scv.dofIndex();
390 deflectionHelper.setCurrent(scv);
393 for (
int pvIdx = 0; pvIdx < numEq; pvIdx++)
397 auto evalResiduals = [&](Scalar priVar)
400 elemSol[scv.localDofIndex()][pvIdx] = priVar;
401 deflectionHelper.deflect(elemSol, scv, this->
asImp_().
problem());
402 if constexpr (solutionDependentFluxVarsCache)
405 if constexpr (enableGridFluxVarsCache)
408 this->
asImp_().maybeUpdateCouplingContext(scv, elemSol, pvIdx);
416 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
422 if (!partialReassembler
425 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
431 A[scvJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[scvJ.localDofIndex()][eqIdx];
437 deflectionHelper.restore(scv);
440 elemSol[scv.localDofIndex()][pvIdx] =
curSol[scv.dofIndex()][pvIdx];
441 this->
asImp_().maybeUpdateCouplingContext(scv, elemSol, pvIdx);
447 if constexpr (enableGridFluxVarsCache)
451 this->
asImp_().maybeEvalAdditionalDomainDerivatives(origResiduals, A, gridVariables);
453 return origResiduals;
463template<
class TypeTag,
class Assembler,
class Implementation>
466 Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false, Implementation>>,
483 using ParentType::ParentType;
491 template <
class PartialReassembler = DefaultPartialReassembler>
495 if (partialReassembler)
496 DUNE_THROW(Dune::NotImplemented,
"partial reassembly for explicit time discretization");
526 const auto dofIdx = scv.dofIndex();
528 const VolumeVariables origVolVars(curVolVars);
531 for (
int pvIdx = 0; pvIdx < numEq; pvIdx++)
535 auto evalStorage = [&](Scalar priVar)
538 elemSol[scv.localDofIndex()][pvIdx] = priVar;
547 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
550 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
556 A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
560 curVolVars = origVolVars;
563 elemSol[scv.localDofIndex()][pvIdx] =
curSol[scv.dofIndex()][pvIdx];
568 return origResiduals;
577template<
class TypeTag,
class Assembler,
class Implementation>
580 Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true, Implementation>>,
592 using ParentType::ParentType;
600 template <
class PartialReassembler = DefaultPartialReassembler>
604 if (partialReassembler)
605 DUNE_THROW(Dune::NotImplemented,
"partial reassembly for analytic differentiation");
629 const auto dofIdx = scv.dofIndex();
635 if (!this->
assembler().isStationaryProblem())
636 this->
localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
645 this->
localResidual().addSourceDerivatives(A[dofIdx][dofIdx],
656 if (!scvf.boundary())
672 const auto& insideScv =
fvGeometry.scv(scvf.insideScvIdx());
676 this->
localResidual().addRobinFluxDerivatives(A[insideScv.dofIndex()],
687 return origResiduals;
697template<
class TypeTag,
class Assembler,
class Implementation>
700 Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false, Implementation>>,
712 using ParentType::ParentType;
720 template <
class PartialReassembler = DefaultPartialReassembler>
724 if (partialReassembler)
725 DUNE_THROW(Dune::NotImplemented,
"partial reassembly for explicit time discretization");
748 const auto dofIdx = scv.dofIndex();
753 this->
localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
761 return origResiduals;
A base class for all local assemblers.
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition assembly/cvfelocalassembler.hh:711
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:721
typename ParentType::LocalResidual LocalResidual
Definition assembly/cvfelocalassembler.hh:710
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:601
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition assembly/cvfelocalassembler.hh:591
typename ParentType::LocalResidual LocalResidual
Definition assembly/cvfelocalassembler.hh:590
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition assembly/cvfelocalassembler.hh:333
typename ParentType::LocalResidual LocalResidual
Definition assembly/cvfelocalassembler.hh:332
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:343
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition assembly/cvfelocalassembler.hh:482
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:492
typename ParentType::LocalResidual LocalResidual
Definition assembly/cvfelocalassembler.hh:481
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:237
void assembleResidual(ResidualVector &res)
Assemble the residual only.
Definition assembly/cvfelocalassembler.hh:216
void maybeUpdateCouplingContext(Args &&...)
Update the coupling context for coupled models.
Definition assembly/cvfelocalassembler.hh:282
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:192
void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
Evaluates Dirichlet boundaries.
Definition assembly/cvfelocalassembler.hh:249
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:289
An assembler for Jacobian and residual contribution per element (CVFE methods).
Definition assembly/cvfelocalassembler.hh:302
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: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
@ 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
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.