14#ifndef DUMUX_EXPERIMENTAL_CC_LOCAL_ASSEMBLER_HH
15#define DUMUX_EXPERIMENTAL_CC_LOCAL_ASSEMBLER_HH
17#include <dune/common/reservedvector.hh>
18#include <dune/grid/common/gridenums.hh>
19#include <dune/istl/matrixindexset.hh>
47template<
class TypeTag,
class Assembler,
class Implementation>
52 using FVElementGeometry =
typename GridGeometry::LocalView;
53 using Element =
typename FVElementGeometry::Element;
54 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
55 using GridView =
typename GridGeometry::GridView;
64 using ParentType::ParentType;
66 template <
class Res
idualVector,
class StageParams,
class PartialReassembler = DefaultPartialReassembler,
class CouplingFunction = Noop>
68 const StageParams& stageParams, ResidualVector& temporal, ResidualVector& spatial,
69 ResidualVector& constrainedDofs,
71 const CouplingFunction& maybeAssembleCouplingBlocks =
noop)
73 this->
asImp_().bindLocalViews();
74 const auto globalI = this->
fvGeometry().gridGeometry().elementMapper().index(this->
element());
81 res[globalI] = spatial[globalI]*stageParams.spatialWeight(stageParams.size()-1)
82 + temporal[globalI]*stageParams.temporalWeight(stageParams.size()-1);
84 this->
localResidual().spatialWeight(stageParams.spatialWeight(stageParams.size()-1));
85 this->
localResidual().temporalWeight(stageParams.temporalWeight(stageParams.size()-1));
88 if (partialReassembler
92 maybeAssembleCouplingBlocks(res[globalI]);
96 this->
asImp_().assembleJacobian(jac, gridVariables, res[globalI]);
99 maybeAssembleCouplingBlocks(res[globalI]);
106 template <
class Res
idualVector>
109 this->
asImp_().bindLocalViews();
110 const auto globalI = this->
assembler().gridGeometry().elementMapper().index(this->
element());
111 res[globalI] = this->
asImp_().evalLocalResidual()[0];
114 if constexpr (Problem::enableInternalDirichletConstraints())
116 const auto applyDirichlet = [&] (
const auto& scvI,
117 const auto& dirichletValues,
121 res[scvI.dofIndex()][eqIdx] = this->
curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
124 this->
asImp_().enforceInternalDirichletConstraints(applyDirichlet);
132 const SubControlVolumeFace& scvf)
const
144 template<
class SubRes
idualVector>
147 this->
asImp_().bindLocalViews();
148 const auto globalI = this->
fvGeometry().gridGeometry().elementMapper().index(this->
element());
152 if constexpr (Problem::enableInternalDirichletConstraints())
153 DUNE_THROW(Dune::NotImplemented,
"Not implemented");
160 template<
class... Args>
167 template<
class... Args>
180template<
class TypeTag,
class Assembler, DiffMethod diffMethod = DiffMethod::numeric,
class Implementation =
void>
189template<
class TypeTag,
class Assembler,
class Implementation>
192 NonVoidOr<CCLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, Implementation>, Implementation>>
199 using FVElementGeometry =
typename GridGeometry::LocalView;
200 using Element =
typename FVElementGeometry::Element;
203 using Problem =
typename GridVariables::GridVolumeVariables::Problem;
209 static constexpr int maxElementStencilSize = GridGeometry::maxElementStencilSize;
222 void assembleJacobian(JacobianMatrix& A, GridVariables& gridVariables,
const NumEqVector& origResidual)
224 if (this->isImplicit())
225 assembleJacobianImplicit_(A, gridVariables, origResidual);
227 assembleJacobianExplicit_(A, gridVariables, origResidual);
245 const auto& element = this->element();
246 const auto& fvGeometry = this->fvGeometry();
247 const auto& gridGeometry = this->fvGeometry().gridGeometry();
248 auto&& curElemVolVars = this->curElemVolVars();
249 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
252 const auto globalI = gridGeometry.elementMapper().index(element);
253 const auto& connectivityMap = gridGeometry.connectivityMap();
254 const auto numNeighbors = connectivityMap[globalI].size();
257 Dune::ReservedVector<Element, maxElementStencilSize> neighborElements;
258 neighborElements.resize(numNeighbors);
262 Residuals origResiduals(numNeighbors + 1); origResiduals = 0.0;
263 origResiduals[0] = origResidual;
268 auto evalNeighborFlux = [&] (
const auto& neighbor,
const auto& scvf)
270 if (neighbor.partitionType() == Dune::GhostEntity)
273 return this->evalFlux(neighbor, scvf);
279 for (
const auto& dataJ : connectivityMap[globalI])
281 neighborElements[j-1] = gridGeometry.element(dataJ.globalJ);
282 for (
const auto scvfIdx : dataJ.scvfsJ)
283 origResiduals[j] += evalNeighborFlux(neighborElements[j-1], fvGeometry.scvf(scvfIdx));
289 const auto& scv = fvGeometry.scv(globalI);
290 auto& curVolVars = ParentType::getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
294 const auto& curSol = this->asImp_().curSol();
295 const auto origPriVars = curSol[globalI];
296 const auto origVolVars = curVolVars;
299 auto elemSol = elementSolution<FVElementGeometry>(origPriVars);
303 Residuals partialDerivs(numNeighbors + 1);
305 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
309 auto evalResiduals = [&](Scalar priVar)
311 Residuals partialDerivsTmp(numNeighbors + 1);
312 partialDerivsTmp = 0.0;
314 elemSol[0][pvIdx] = priVar;
315 this->asImp_().maybeUpdateCouplingContext(scv, elemSol, pvIdx);
316 curVolVars.update(elemSol, this->asImp_().problem(), element, scv);
317 elemFluxVarsCache.update(element, fvGeometry, curElemVolVars);
318 if (enableGridFluxVarsCache)
319 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
322 partialDerivsTmp[0] = this->evalLocalResidual()[0];
325 for (std::size_t k = 0; k < numNeighbors; ++k)
326 for (
auto scvfIdx : connectivityMap[globalI][k].scvfsJ)
327 partialDerivsTmp[k+1] += evalNeighborFlux(neighborElements[k], fvGeometry.scvf(scvfIdx));
329 return partialDerivsTmp;
333 static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()};
334 static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(),
"Assembly.NumericDifferenceMethod");
336 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
341 if (this->elementIsGhost())
343 partialDerivs[0] = 0.0;
344 partialDerivs[0][pvIdx] = 1.0;
348 curVolVars = origVolVars;
351 elemSol[0][pvIdx] = origPriVars[pvIdx];
354 this->asImp_().maybeUpdateCouplingContext(scv, elemSol, pvIdx);
358 if constexpr (Problem::enableInternalDirichletConstraints())
361 const auto internalDirichletConstraintsOwnElement = this->asImp_().problem().hasInternalDirichletConstraint(this->
element(), scv);
362 const auto dirichletValues = this->asImp_().problem().internalDirichlet(this->
element(), scv);
364 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
366 if (internalDirichletConstraintsOwnElement[eqIdx])
368 origResiduals[0][eqIdx] = origVolVars.priVars()[eqIdx] - dirichletValues[eqIdx];
369 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
372 A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx];
377 for (
const auto& dataJ : connectivityMap[globalI])
379 const auto& neighborElement = neighborElements[j-1];
380 const auto& neighborScv = fvGeometry.scv(dataJ.globalJ);
381 const auto internalDirichletConstraintsNeighbor = this->asImp_().problem().hasInternalDirichletConstraint(neighborElement, neighborScv);
383 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
385 if (internalDirichletConstraintsNeighbor[eqIdx])
386 A[dataJ.globalJ][globalI][eqIdx][pvIdx] = 0.0;
388 A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j][eqIdx];
396 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
399 A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx];
403 for (
const auto& dataJ : connectivityMap[globalI])
404 A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j++][eqIdx];
415 if (enableGridFluxVarsCache)
416 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
419 ElementResidualVector orig{origResidual};
420 this->asImp_().maybeEvalAdditionalDomainDerivatives(orig, A, gridVariables);
427 void assembleJacobianExplicit_(JacobianMatrix& A, GridVariables& gridVariables,
const NumEqVector& origResidual)
429 if (this->assembler().isStationaryProblem())
430 DUNE_THROW(Dune::InvalidStateException,
"Using explicit jacobian assembler with stationary local residual");
433 auto storageResidual = origResidual;
443 const auto& fvGeometry = this->fvGeometry();
444 const auto& gridGeometry = this->fvGeometry().gridGeometry();
445 auto&& curElemVolVars = this->curElemVolVars();
448 const auto globalI = gridGeometry.elementMapper().index(element);
449 const auto& scv = fvGeometry.scv(globalI);
450 auto& curVolVars = ParentType::getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
454 const auto& curSol = this->asImp_().curSol();
455 const auto origPriVars = curSol[globalI];
456 const auto origVolVars = curVolVars;
459 auto elemSol = elementSolution<FVElementGeometry>(origPriVars);
463 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
468 auto evalStorage = [&](Scalar priVar)
472 elemSol[0][pvIdx] = priVar;
473 curVolVars.update(elemSol, this->asImp_().problem(),
element, scv);
474 return this->evalStorage()[0];
478 if (!this->elementIsGhost())
480 static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()};
481 static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(),
"Assembly.NumericDifferenceMethod");
483 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
489 else partialDeriv[pvIdx] = 1.0;
492 curVolVars = origVolVars;
495 elemSol[0][pvIdx] = origPriVars[pvIdx];
499 if constexpr (Problem::enableInternalDirichletConstraints())
502 const auto internalDirichletConstraints = this->asImp_().problem().hasInternalDirichletConstraint(this->
element(), scv);
503 const auto dirichletValues = this->asImp_().problem().internalDirichlet(this->
element(), scv);
505 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
509 if (internalDirichletConstraints[eqIdx])
511 storageResidual[eqIdx] = origVolVars.priVars()[eqIdx] - dirichletValues[eqIdx];
512 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
515 A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx];
520 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
521 A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx];
The local element solution class for cell-centered methods.
Cell-centered scheme local assembler using numeric differentiation.
Definition: experimental/assembly/cclocalassembler.hh:193
void assembleJacobian(JacobianMatrix &A, GridVariables &gridVariables, const NumEqVector &origResidual)
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition: experimental/assembly/cclocalassembler.hh:222
typename ParentType::LocalResidual LocalResidual
Definition: experimental/assembly/cclocalassembler.hh:215
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: experimental/assembly/cclocalassembler.hh:216
A base class for all local cell-centered assemblers.
Definition: experimental/assembly/cclocalassembler.hh:49
void assembleResidual(ResidualVector &res)
Assemble the residual only.
Definition: experimental/assembly/cclocalassembler.hh:107
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)
Definition: experimental/assembly/cclocalassembler.hh:67
void maybeUpdateCouplingContext(Args &&...)
Update the coupling context for coupled models.
Definition: experimental/assembly/cclocalassembler.hh:161
void maybeEvalAdditionalDomainDerivatives(Args &&...)
Update the additional domain derivatives for coupled models.
Definition: experimental/assembly/cclocalassembler.hh:168
void assembleCurrentResidual(SubResidualVector &spatialRes, SubResidualVector &temporalRes)
Assemble the residual only.
Definition: experimental/assembly/cclocalassembler.hh:145
NumEqVector evalFlux(const Element &neighbor, const SubControlVolumeFace &scvf) const
Evaluates the fluxes (element can potentially be a neighbor)
Definition: experimental/assembly/cclocalassembler.hh:131
An assembler for Jacobian and residual contribution per element (cell-centered methods)
Definition: experimental/assembly/cclocalassembler.hh:181
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
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
const Assembler & assembler() const
The assembler.
Definition: experimental/assembly/fvlocalassemblerbase.hh:214
FVElementGeometry & fvGeometry()
The global finite volume geometry.
Definition: experimental/assembly/fvlocalassemblerbase.hh:230
Implementation & asImp_()
Definition: experimental/assembly/fvlocalassemblerbase.hh:274
ElementFluxVariablesCache & elemFluxVarsCache()
The element flux variables cache.
Definition: experimental/assembly/fvlocalassemblerbase.hh:238
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
The flux stencil specialized for different discretization schemes.
Definition: fluxstencil.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:49
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:420
A arithmetic block vector type based on DUNE's reserved vector.
Definition: reservedblockvector.hh:26
Defines all properties used in Dumux.
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.
The flux stencil specialized for different discretization schemes.
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
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition: numeqvector.hh:34
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 helper to deduce a vector with the same size as numbers of equations.
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.
A arithmetic block vector type based on DUNE's reserved vector.