13#ifndef DUMUX_CC_LOCAL_ASSEMBLER_HH
14#define DUMUX_CC_LOCAL_ASSEMBLER_HH
16#include <dune/common/reservedvector.hh>
17#include <dune/grid/common/gridenums.hh>
18#include <dune/istl/matrixindexset.hh>
44template<
class TypeTag,
class Assembler,
class Implementation,
bool implicit>
56 using ParentType::ParentType;
62 template <
class Res
idualVector,
class PartialReassembler = DefaultPartialReassembler>
66 this->
asImp_().bindLocalViews();
67 const auto globalI = this->
assembler().gridGeometry().elementMapper().index(this->
element());
68 if (partialReassembler
71 res[globalI] = this->
asImp_().evalLocalResidual()[0];
75 res[globalI] = this->
asImp_().assembleJacobianAndResidualImpl(jac, gridVariables);
85 this->
asImp_().bindLocalViews();
86 this->
asImp_().assembleJacobianAndResidualImpl(jac, gridVariables);
92 template <
class Res
idualVector>
95 this->
asImp_().bindLocalViews();
96 const auto globalI = this->
assembler().gridGeometry().elementMapper().index(this->
element());
97 res[globalI] = this->
asImp_().evalLocalResidual()[0];
100 if constexpr (Problem::enableInternalDirichletConstraints())
102 const auto applyDirichlet = [&] (
const auto& scvI,
103 const auto& dirichletValues,
107 res[scvI.dofIndex()][eqIdx] = this->
curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
110 this->
asImp_().enforceInternalDirichletConstraints(applyDirichlet);
123template<
class TypeTag,
class Assembler, DiffMethod diffMethod = DiffMethod::numeric,
bool implicit = true>
131template<
class TypeTag,
class Assembler>
134 CCLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true>, true >
142 using FVElementGeometry =
typename GridGeometry::LocalView;
145 using Problem =
typename GridVariables::GridVolumeVariables::Problem;
151 static constexpr int maxElementStencilSize = GridGeometry::maxElementStencilSize;
173 const auto& element = this->element();
174 const auto& fvGeometry = this->fvGeometry();
175 const auto& gridGeometry = this->assembler().gridGeometry();
176 auto&& curElemVolVars = this->curElemVolVars();
177 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
180 const auto globalI = gridGeometry.elementMapper().index(element);
181 const auto& connectivityMap = gridGeometry.connectivityMap();
182 const auto numNeighbors = connectivityMap[globalI].size();
185 Dune::ReservedVector<Element, maxElementStencilSize> neighborElements;
186 neighborElements.resize(numNeighbors);
190 Residuals origResiduals(numNeighbors + 1); origResiduals = 0.0;
191 origResiduals[0] = this->evalLocalResidual()[0];
196 auto evalNeighborFlux = [&] (
const auto& neighbor,
const auto& scvf)
198 if (neighbor.partitionType() == Dune::GhostEntity)
201 return this->localResidual().evalFlux(this->problem(),
204 this->curElemVolVars(),
205 this->elemFluxVarsCache(), scvf);
211 for (
const auto& dataJ : connectivityMap[globalI])
213 neighborElements[j-1] = gridGeometry.element(dataJ.globalJ);
214 for (
const auto scvfIdx : dataJ.scvfsJ)
215 origResiduals[j] += evalNeighborFlux(neighborElements[j-1], fvGeometry.scvf(scvfIdx));
221 const auto& scv = fvGeometry.scv(globalI);
222 auto& curVolVars = ParentType::getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
226 const auto& curSol = this->curSol();
227 const auto origPriVars = curSol[globalI];
228 const auto origVolVars = curVolVars;
235 Residuals partialDerivs(numNeighbors + 1);
237 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
241 auto evalResiduals = [&](Scalar priVar)
243 Residuals partialDerivsTmp(numNeighbors + 1);
244 partialDerivsTmp = 0.0;
246 elemSol[0][pvIdx] = priVar;
247 curVolVars.update(elemSol, this->problem(), element, scv);
248 elemFluxVarsCache.update(element, fvGeometry, curElemVolVars);
249 if (enableGridFluxVarsCache)
250 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
253 partialDerivsTmp[0] = this->evalLocalResidual()[0];
256 for (std::size_t k = 0; k < numNeighbors; ++k)
257 for (
auto scvfIdx : connectivityMap[globalI][k].scvfsJ)
258 partialDerivsTmp[k+1] += evalNeighborFlux(neighborElements[k], fvGeometry.scvf(scvfIdx));
260 return partialDerivsTmp;
265 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(),
"Assembly.NumericDifferenceMethod");
267 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
272 if (this->elementIsGhost())
274 partialDerivs[0] = 0.0;
275 partialDerivs[0][pvIdx] = 1.0;
279 curVolVars = origVolVars;
282 elemSol[0][pvIdx] = origPriVars[pvIdx];
286 if constexpr (Problem::enableInternalDirichletConstraints())
289 const auto internalDirichletConstraintsOwnElement = this->problem().hasInternalDirichletConstraint(this->element(), scv);
290 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
292 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
294 if (internalDirichletConstraintsOwnElement[eqIdx])
296 origResiduals[0][eqIdx] = origVolVars.priVars()[eqIdx] - dirichletValues[eqIdx];
297 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
300 A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx];
305 for (
const auto& dataJ : connectivityMap[globalI])
307 const auto& neighborElement = neighborElements[j-1];
308 const auto& neighborScv = fvGeometry.scv(dataJ.globalJ);
309 const auto internalDirichletConstraintsNeighbor = this->problem().hasInternalDirichletConstraint(neighborElement, neighborScv);
311 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
313 if (internalDirichletConstraintsNeighbor[eqIdx])
314 A[dataJ.globalJ][globalI][eqIdx][pvIdx] = 0.0;
316 A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j][eqIdx];
324 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
327 A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx];
331 for (
const auto& dataJ : connectivityMap[globalI])
332 A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j++][eqIdx];
343 if (enableGridFluxVarsCache)
344 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
347 return origResiduals[0];
357template<
class TypeTag,
class Assembler>
360 CCLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false>, false>
369 using Problem =
typename GridVariables::GridVolumeVariables::Problem;
384 if (this->assembler().isStationaryProblem())
385 DUNE_THROW(Dune::InvalidStateException,
"Using explicit jacobian assembler with stationary local residual");
388 auto residual = this->evalLocalResidual()[0];
389 const auto storageResidual = this->evalLocalStorageResidual()[0];
398 const auto& element = this->element();
399 const auto& fvGeometry = this->fvGeometry();
400 const auto& gridGeometry = this->assembler().gridGeometry();
401 auto&& curElemVolVars = this->curElemVolVars();
404 const auto globalI = gridGeometry.elementMapper().index(element);
405 const auto& scv = fvGeometry.scv(globalI);
406 auto& curVolVars = ParentType::getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
410 const auto& curSol = this->curSol();
411 const auto origPriVars = curSol[globalI];
412 const auto origVolVars = curVolVars;
417 NumEqVector partialDeriv;
420 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
425 auto evalStorage = [&](Scalar priVar)
429 elemSol[0][pvIdx] = priVar;
430 curVolVars.update(elemSol, this->problem(), element, scv);
431 return this->evalLocalStorageResidual()[0];
435 if (!this->elementIsGhost())
438 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(),
"Assembly.NumericDifferenceMethod");
440 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
446 else partialDeriv[pvIdx] = 1.0;
449 curVolVars = origVolVars;
452 elemSol[0][pvIdx] = origPriVars[pvIdx];
456 if constexpr (Problem::enableInternalDirichletConstraints())
459 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
460 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
462 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
464 if (internalDirichletConstraints[eqIdx])
466 residual[eqIdx] = origVolVars.priVars()[eqIdx] - dirichletValues[eqIdx];
467 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
470 A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx];
475 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
476 A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx];
490template<
class TypeTag,
class Assembler>
493 CCLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true>, true>
500 using Problem =
typename GridVariables::GridVolumeVariables::Problem;
516 if (this->elementIsGhost())
518 const auto globalI = this->assembler().gridGeometry().elementMapper().index(this->element());
519 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
520 A[globalI][globalI][pvIdx][pvIdx] = 1.0;
527 auto residual = this->evalLocalResidual()[0];
530 const auto& problem = this->problem();
531 const auto& element = this->element();
532 const auto& fvGeometry = this->fvGeometry();
533 const auto& curElemVolVars = this->curElemVolVars();
534 const auto& elemFluxVarsCache = this->elemFluxVarsCache();
537 const auto globalI = this->assembler().gridGeometry().elementMapper().index(element);
538 const auto& scv = fvGeometry.scv(globalI);
539 const auto& volVars = curElemVolVars[scv];
542 if (!this->assembler().isStationaryProblem())
543 this->localResidual().addStorageDerivatives(A[globalI][globalI], problem, element, fvGeometry, volVars, scv);
546 this->localResidual().addSourceDerivatives(A[globalI][globalI], problem, element, fvGeometry, volVars, scv);
549 for (
const auto& scvf : scvfs(fvGeometry))
552 if (!scvf.boundary())
553 this->localResidual().addFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
558 const auto& bcTypes = problem.boundaryTypes(element, scvf);
561 if (bcTypes.hasDirichlet() && !bcTypes.hasNeumann())
562 this->localResidual().addCCDirichletFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
565 else if (bcTypes.hasNeumann() && !bcTypes.hasDirichlet())
566 this->localResidual().addRobinFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
569 DUNE_THROW(Dune::NotImplemented,
"Mixed boundary conditions. Use pure boundary conditions by converting Dirichlet BCs to Robin BCs");
573 if constexpr (Problem::enableInternalDirichletConstraints())
576 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
577 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
579 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
581 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
583 if (internalDirichletConstraints[eqIdx])
585 residual[eqIdx] = volVars.priVars()[eqIdx] - dirichletValues[eqIdx];
586 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
589 for (
const auto& scvf : scvfs(fvGeometry))
590 if (!scvf.boundary())
591 A[globalI][fvGeometry.scv(scvf.outsideScvIdx()).dofIndex()][eqIdx][pvIdx] = 0.0;
607template<
class TypeTag,
class Assembler>
610 CCLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false>, false>
617 using Problem =
typename GridVariables::GridVolumeVariables::Problem;
633 if (this->elementIsGhost())
635 const auto globalI = this->assembler().gridGeometry().elementMapper().index(this->element());
636 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
637 A[globalI][globalI][pvIdx][pvIdx] = 1.0;
644 const auto residual = this->evalLocalResidual()[0];
647 const auto globalI = this->assembler().gridGeometry().elementMapper().index(this->element());
648 const auto& scv = this->fvGeometry().scv(globalI);
649 const auto& volVars = this->curElemVolVars()[scv];
652 this->localResidual().addStorageDerivatives(A[globalI][globalI], this->problem(), this->element(), this->fvGeometry(), volVars, scv);
654 if constexpr (Problem::enableInternalDirichletConstraints())
657 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
658 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
660 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
662 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
664 if (internalDirichletConstraints[eqIdx])
666 residual[eqIdx] = volVars.priVars()[eqIdx] - dirichletValues[eqIdx];
667 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
The local element solution class for cell-centered methods.
Cell-centered scheme local assembler using analytic (hand-coded) differentiation and explicit time di...
Definition: assembly/cclocalassembler.hh:611
NumEqVector assembleJacobianAndResidualImpl(JacobianMatrix &A, const GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: assembly/cclocalassembler.hh:630
Cell-centered scheme local assembler using analytic (hand-coded) differentiation and implicit time di...
Definition: assembly/cclocalassembler.hh:494
NumEqVector assembleJacobianAndResidualImpl(JacobianMatrix &A, const GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: assembly/cclocalassembler.hh:513
Cell-centered scheme local assembler using numeric differentiation and explicit time discretization.
Definition: assembly/cclocalassembler.hh:361
NumEqVector assembleJacobianAndResidualImpl(JacobianMatrix &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: assembly/cclocalassembler.hh:382
Cell-centered scheme local assembler using numeric differentiation and implicit time discretization.
Definition: assembly/cclocalassembler.hh:135
NumEqVector assembleJacobianAndResidualImpl(JacobianMatrix &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: assembly/cclocalassembler.hh:164
A base class for all local cell-centered assemblers.
Definition: assembly/cclocalassembler.hh:46
void assembleJacobian(JacobianMatrix &jac, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: assembly/cclocalassembler.hh:83
void assembleResidual(ResidualVector &res)
Assemble the residual only.
Definition: assembly/cclocalassembler.hh:93
void assembleJacobianAndResidual(JacobianMatrix &jac, ResidualVector &res, GridVariables &gridVariables, const PartialReassembler *partialReassembler)
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition: assembly/cclocalassembler.hh:63
An assembler for Jacobian and residual contribution per element (cell-centered methods)
Definition: assembly/cclocalassembler.hh:124
A base class for all local assemblers.
Definition: assembly/fvlocalassemblerbase.hh:36
ElementVolumeVariables & curElemVolVars()
The current element volume variables.
Definition: assembly/fvlocalassemblerbase.hh:253
Implementation & asImp_()
Definition: assembly/fvlocalassemblerbase.hh:297
const Assembler & assembler() const
The assembler.
Definition: assembly/fvlocalassemblerbase.hh:233
const Element & element() const
The current element.
Definition: assembly/fvlocalassemblerbase.hh:237
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
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
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:491
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
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 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.