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>
37#include <dumux/assembly/volvardeflectionhelper_.hh>
50template<
class TypeTag,
class Assembler,
class Implementation>
57 using ElementVolumeVariables =
typename GridVariables::GridVolumeVariables::LocalView;
58 using SolutionVector =
typename Assembler::SolutionVector;
65 using ParentType::ParentType;
81 template <
class Res
idualVector,
class StageParams,
class PartialReassembler = DefaultPartialReassembler,
class CouplingFunction = Noop>
83 const StageParams& stageParams, ResidualVector& temporal, ResidualVector& spatial,
84 ResidualVector& constrainedDofs,
86 const CouplingFunction& maybeAssembleCouplingBlocks =
noop)
88 this->
asImp_().bindLocalViews();
89 const auto eIdxGlobal = this->
asImp_().problem().gridGeometry().elementMapper().index(this->
element());
94 const auto sWeight = stageParams.spatialWeight(stageParams.size()-1);
95 const auto tWeight = stageParams.temporalWeight(stageParams.size()-1);
101 for (
const auto& scv : scvs(this->
fvGeometry()))
103 spatial[scv.dofIndex()] += flux[scv.localDofIndex()];
104 temporal[scv.dofIndex()] += storage[scv.localDofIndex()];
105 origResidual[scv.localDofIndex()] += flux[scv.localDofIndex()]*sWeight + storage[scv.localDofIndex()]*tWeight;
106 res[scv.dofIndex()] += origResidual[scv.localDofIndex()];
112 if (partialReassembler && partialReassembler->elementColor(eIdxGlobal) ==
EntityColor::green)
115 maybeAssembleCouplingBlocks(origResidual);
119 this->
asImp_().assembleJacobian(jac, gridVariables, origResidual, partialReassembler);
122 maybeAssembleCouplingBlocks(origResidual);
130 const auto& gridGeometry = this->
asImp_().problem().gridGeometry();
131 Dune::Hybrid::forEach(std::make_integer_sequence<int, dim>{}, [&](
auto d)
133 constexpr int codim = dim - d;
134 const auto& localCoeffs = gridGeometry.feCache().get(this->
element().type()).localCoefficients();
135 for (
int idx = 0; idx < localCoeffs.size(); ++idx)
137 const auto& localKey = localCoeffs.localKey(idx);
140 if (localKey.codim() != codim)
144 auto entity = this->
element().template subEntity<codim>(localKey.subEntity());
145 if (entity.partitionType() == Dune::InteriorEntity || entity.partitionType() == Dune::BorderEntity)
152 const auto dofIndex = gridGeometry.dofMapper().index(entity);
155 using BlockType =
typename JacobianMatrix::block_type;
156 BlockType &J = jac[dofIndex][dofIndex];
157 for (
int j = 0; j < BlockType::rows; ++j)
162 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();
229 for (
const auto& scv : scvs(this->
fvGeometry()))
230 res[scv.dofIndex()] += residual[scv.localDofIndex()];
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();
252 for (
const auto& scv : scvs(this->
fvGeometry()))
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 >
279 for (
const auto& scvI : scvs(this->
fvGeometry()))
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>>
349 static constexpr bool enableGridFluxVarsCache
350 = GridVariables::GridFluxVariablesCache::cachingEnabled;
351 static constexpr bool solutionDependentFluxVarsCache
352 = GridVariables::GridFluxVariablesCache::FluxVariablesCache::isSolDependent;
360 template <
class PartialReassembler = DefaultPartialReassembler>
365 if (this->isImplicit())
366 assembleJacobianImplicit_(A, gridVariables, origResiduals, partialReassembler);
368 assembleJacobianExplicit_(A, gridVariables, origResiduals, partialReassembler);
378 template <
class PartialReassembler = DefaultPartialReassembler>
379 void assembleJacobianImplicit_(JacobianMatrix& A,
GridVariables& gridVariables,
380 const ElementResidualVector& origResiduals,
384 const auto& element = this->element();
385 const auto& fvGeometry = this->fvGeometry();
386 const auto& curSol = this->asImp_().curSol();
388 auto&& curElemVolVars = this->curElemVolVars();
389 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
401 static const bool updateAllVolVars = getParamFromGroup<bool>(
402 this->asImp_().problem().paramGroup(),
"Assembly.BoxVolVarsDependOnAllElementDofs",
false
406 auto elemSol =
elementSolution(element, curSol, fvGeometry.gridGeometry());
409 ElementResidualVector partialDerivs(fvGeometry.numScv());
411 Dumux::Detail::VolVarsDeflectionHelper deflectionHelper(
412 [&] (
const auto& scv) -> VolumeVariables& {
413 return this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
420 for (
const auto& scv : scvs(fvGeometry))
423 const auto dofIdx = scv.dofIndex();
424 deflectionHelper.setCurrent(scv);
427 for (
int pvIdx = 0; pvIdx < numEq; pvIdx++)
431 auto evalResiduals = [&](Scalar priVar)
434 elemSol[scv.localDofIndex()][pvIdx] = priVar;
435 deflectionHelper.deflect(elemSol, scv, 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(scv, 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[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
453 for (
const auto& scvJ : scvs(fvGeometry))
456 if (!partialReassembler
459 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
465 A[scvJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[scvJ.localDofIndex()][eqIdx];
471 deflectionHelper.restore(scv);
474 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
475 this->asImp_().maybeUpdateCouplingContext(scv, elemSol, pvIdx);
481 if constexpr (enableGridFluxVarsCache)
482 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
485 this->asImp_().maybeEvalAdditionalDomainDerivatives(origResiduals, A, gridVariables);
493 template <
class PartialReassembler = DefaultPartialReassembler>
494 void assembleJacobianExplicit_(JacobianMatrix& A, GridVariables& gridVariables,
495 const ElementResidualVector& origResiduals,
496 const PartialReassembler* partialReassembler =
nullptr)
498 if (partialReassembler)
499 DUNE_THROW(Dune::NotImplemented,
"partial reassembly for explicit time discretization");
503 const auto& fvGeometry = this->fvGeometry();
504 const auto& curSol = this->asImp_().curSol();
505 auto&& curElemVolVars = this->curElemVolVars();
508 auto elemSol =
elementSolution(element, curSol, fvGeometry.gridGeometry());
511 ElementResidualVector partialDerivs(fvGeometry.numScv());
514 for (
const auto& scv : scvs(fvGeometry))
517 const auto dofIdx = scv.dofIndex();
518 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
519 const VolumeVariables origVolVars(curVolVars);
522 for (
int pvIdx = 0; pvIdx < numEq; pvIdx++)
526 auto evalStorage = [&](Scalar priVar)
528 elemSol[scv.localDofIndex()][pvIdx] = priVar;
529 curVolVars.update(elemSol, this->asImp_().problem(),
element, scv);
530 return this->evalStorage();
534 static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()};
535 static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(),
"Assembly.NumericDifferenceMethod");
537 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
540 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
546 A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
550 curVolVars = origVolVars;
553 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
Control volume finite element local assembler using numeric differentiation.
Definition: experimental/assembly/cvfelocalassembler.hh:338
void assembleJacobian(JacobianMatrix &A, GridVariables &gridVariables, const ElementResidualVector &origResiduals, const PartialReassembler *partialReassembler=nullptr)
Definition: experimental/assembly/cvfelocalassembler.hh:361
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: experimental/assembly/cvfelocalassembler.hh:357
typename ParentType::LocalResidual LocalResidual
Definition: experimental/assembly/cvfelocalassembler.hh:356
A base class for all local CVFE assemblers.
Definition: experimental/assembly/cvfelocalassembler.hh:52
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:71
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:82
typename ParentType::LocalResidual LocalResidual
Definition: experimental/assembly/cvfelocalassembler.hh:67
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: experimental/assembly/cvfelocalassembler.hh:68
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
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
Base class for grid variables.
Definition: experimental/discretization/gridvariables.hh:33
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.
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
Definition: experimental/assembly/cclocalassembler.hh:36
constexpr auto noop
Function that performs no operation.
Definition: common/typetraits/typetraits.hh:29
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.