25#ifndef DUMUX_CC_LOCAL_ASSEMBLER_HH
26#define DUMUX_CC_LOCAL_ASSEMBLER_HH
28#include <dune/common/reservedvector.hh>
29#include <dune/grid/common/gridenums.hh>
30#include <dune/istl/matrixindexset.hh>
55template<
class TypeTag,
class Assembler,
class Implementation,
bool implicit>
68 using ParentType::ParentType;
74 template <
class PartialReassembler = DefaultPartialReassembler>
78 this->
asImp_().bindLocalViews();
79 const auto globalI = this->
assembler().gridGeometry().elementMapper().index(this->
element());
80 if (partialReassembler
83 res[globalI] = this->
asImp_().evalLocalResidual()[0];
87 res[globalI] = this->
asImp_().assembleJacobianAndResidualImpl(jac, gridVariables);
97 this->
asImp_().bindLocalViews();
98 this->
asImp_().assembleJacobianAndResidualImpl(jac, gridVariables);
106 this->
asImp_().bindLocalViews();
107 const auto globalI = this->
assembler().gridGeometry().elementMapper().index(this->
element());
108 res[globalI] = this->
asImp_().evalLocalResidual()[0];
120template<
class TypeTag,
class Assembler, DiffMethod diffMethod = DiffMethod::numeric,
bool implicit = true>
128template<
class TypeTag,
class Assembler>
131 CCLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true>, true >
139 using FVElementGeometry =
typename GridGeometry::LocalView;
142 using Problem =
typename GridVariables::GridVolumeVariables::Problem;
148 static constexpr int maxElementStencilSize = GridGeometry::maxElementStencilSize;
149 static constexpr bool enableGridFluxVarsCache = getPropValue<TypeTag, Properties::EnableGridFluxVariablesCache>();
170 const auto& element = this->element();
171 const auto& fvGeometry = this->fvGeometry();
172 const auto& gridGeometry = this->assembler().gridGeometry();
173 auto&& curElemVolVars = this->curElemVolVars();
174 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
177 const auto globalI = gridGeometry.elementMapper().index(element);
178 const auto& connectivityMap = gridGeometry.connectivityMap();
179 const auto numNeighbors = connectivityMap[globalI].size();
182 Dune::ReservedVector<Element, maxElementStencilSize> neighborElements;
183 neighborElements.resize(numNeighbors);
187 Residuals origResiduals(numNeighbors + 1); origResiduals = 0.0;
188 origResiduals[0] = this->evalLocalResidual()[0];
193 auto evalNeighborFlux = [&] (
const auto& neighbor,
const auto& scvf)
195 if (neighbor.partitionType() == Dune::GhostEntity)
196 return NumEqVector(0.0);
198 return this->localResidual().evalFlux(this->problem(),
201 this->curElemVolVars(),
202 this->elemFluxVarsCache(), scvf);
208 for (
const auto& dataJ : connectivityMap[globalI])
210 neighborElements[j-1] = gridGeometry.element(dataJ.globalJ);
211 for (
const auto scvfIdx : dataJ.scvfsJ)
212 origResiduals[j] += evalNeighborFlux(neighborElements[j-1], fvGeometry.scvf(scvfIdx));
218 const auto& scv = fvGeometry.scv(globalI);
219 auto& curVolVars = ParentType::getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
223 const auto& curSol = this->curSol();
224 const auto origPriVars = curSol[globalI];
225 const auto origVolVars = curVolVars;
232 Residuals partialDerivs(numNeighbors + 1);
234 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
238 auto evalResiduals = [&](Scalar priVar)
240 Residuals partialDerivsTmp(numNeighbors + 1);
241 partialDerivsTmp = 0.0;
243 elemSol[0][pvIdx] = priVar;
244 curVolVars.update(elemSol, this->problem(), element, scv);
245 elemFluxVarsCache.update(element, fvGeometry, curElemVolVars);
246 if (enableGridFluxVarsCache)
247 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
250 partialDerivsTmp[0] = this->evalLocalResidual()[0];
253 for (std::size_t k = 0; k < numNeighbors; ++k)
254 for (
auto scvfIdx : connectivityMap[globalI][k].scvfsJ)
255 partialDerivsTmp[k+1] += evalNeighborFlux(neighborElements[k], fvGeometry.scvf(scvfIdx));
257 return partialDerivsTmp;
262 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(),
"Assembly.NumericDifferenceMethod");
264 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
269 if (this->elementIsGhost())
271 partialDerivs[0] = 0.0;
272 partialDerivs[0][pvIdx] = 1.0;
277 if constexpr (Problem::enableInternalDirichletConstraints())
280 const auto internalDirichletConstraintsOwnElement = this->problem().hasInternalDirichletConstraint(this->element(), scv);
281 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
283 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
285 if (internalDirichletConstraintsOwnElement[eqIdx])
287 origResiduals[0][eqIdx] = origVolVars.priVars()[eqIdx] - dirichletValues[eqIdx];
288 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
291 A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx];
296 for (
const auto& dataJ : connectivityMap[globalI])
298 const auto& neighborElement = neighborElements[j-1];
299 const auto& neighborScv = fvGeometry.scv(dataJ.globalJ);
300 const auto internalDirichletConstraintsNeighbor = this->problem().hasInternalDirichletConstraint(neighborElement, neighborScv);
302 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
304 if (internalDirichletConstraintsNeighbor[eqIdx])
305 A[dataJ.globalJ][globalI][eqIdx][pvIdx] = 0.0;
307 A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j][eqIdx];
315 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
318 A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx];
322 for (
const auto& dataJ : connectivityMap[globalI])
323 A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j++][eqIdx];
328 curVolVars = origVolVars;
331 elemSol[0][pvIdx] = origPriVars[pvIdx];
340 if (enableGridFluxVarsCache)
341 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
344 return origResiduals[0];
354template<
class TypeTag,
class Assembler>
357 CCLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false>, false>
366 using Problem =
typename GridVariables::GridVolumeVariables::Problem;
381 if (this->assembler().isStationaryProblem())
382 DUNE_THROW(Dune::InvalidStateException,
"Using explicit jacobian assembler with stationary local residual");
385 auto residual = this->evalLocalResidual()[0];
386 const auto storageResidual = this->evalLocalStorageResidual()[0];
395 const auto& element = this->element();
396 const auto& fvGeometry = this->fvGeometry();
397 const auto& gridGeometry = this->assembler().gridGeometry();
398 auto&& curElemVolVars = this->curElemVolVars();
401 const auto globalI = gridGeometry.elementMapper().index(element);
402 const auto& scv = fvGeometry.scv(globalI);
403 auto& curVolVars = ParentType::getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
407 const auto& curSol = this->curSol();
408 const auto origPriVars = curSol[globalI];
409 const auto origVolVars = curVolVars;
414 NumEqVector partialDeriv;
417 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
422 auto evalStorage = [&](Scalar priVar)
426 elemSol[0][pvIdx] = priVar;
427 curVolVars.update(elemSol, this->problem(), element, scv);
428 return this->evalLocalStorageResidual()[0];
432 if (!this->elementIsGhost())
435 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(),
"Assembly.NumericDifferenceMethod");
437 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
443 else partialDeriv[pvIdx] = 1.0;
447 if constexpr (Problem::enableInternalDirichletConstraints())
450 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
451 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
453 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
455 if (internalDirichletConstraints[eqIdx])
457 residual[eqIdx] = origVolVars.priVars()[eqIdx] - dirichletValues[eqIdx];
458 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
461 A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx];
466 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
467 A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx];
471 curVolVars = origVolVars;
474 elemSol[0][pvIdx] = origPriVars[pvIdx];
487template<
class TypeTag,
class Assembler>
490 CCLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true>, true>
497 using Problem =
typename GridVariables::GridVolumeVariables::Problem;
513 if (this->elementIsGhost())
515 const auto globalI = this->assembler().gridGeometry().elementMapper().index(this->element());
516 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
517 A[globalI][globalI][pvIdx][pvIdx] = 1.0;
520 return NumEqVector(0.0);
524 auto residual = this->evalLocalResidual()[0];
527 const auto& problem = this->problem();
528 const auto& element = this->element();
529 const auto& fvGeometry = this->fvGeometry();
530 const auto& curElemVolVars = this->curElemVolVars();
531 const auto& elemFluxVarsCache = this->elemFluxVarsCache();
534 const auto globalI = this->assembler().gridGeometry().elementMapper().index(element);
535 const auto& scv = fvGeometry.scv(globalI);
536 const auto& volVars = curElemVolVars[scv];
539 if (!this->assembler().isStationaryProblem())
540 this->localResidual().addStorageDerivatives(A[globalI][globalI], problem, element, fvGeometry, volVars, scv);
543 this->localResidual().addSourceDerivatives(A[globalI][globalI], problem, element, fvGeometry, volVars, scv);
546 for (
const auto& scvf : scvfs(fvGeometry))
549 if (!scvf.boundary())
550 this->localResidual().addFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
555 const auto& bcTypes = problem.boundaryTypes(element, scvf);
558 if (bcTypes.hasDirichlet() && !bcTypes.hasNeumann())
559 this->localResidual().addCCDirichletFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
562 else if (bcTypes.hasNeumann() && !bcTypes.hasDirichlet())
563 this->localResidual().addRobinFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
566 DUNE_THROW(Dune::NotImplemented,
"Mixed boundary conditions. Use pure boundary conditions by converting Dirichlet BCs to Robin BCs");
570 if constexpr (Problem::enableInternalDirichletConstraints())
573 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
574 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
576 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
578 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
580 if (internalDirichletConstraints[eqIdx])
582 residual[eqIdx] = volVars.priVars()[eqIdx] - dirichletValues[eqIdx];
583 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
586 for (
const auto& scvf : scvfs(fvGeometry))
587 if (!scvf.boundary())
588 A[globalI][fvGeometry.scv(scvf.outsideScvIdx()).dofIndex()][eqIdx][pvIdx] = 0.0;
604template<
class TypeTag,
class Assembler>
607 CCLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false>, false>
614 using Problem =
typename GridVariables::GridVolumeVariables::Problem;
630 if (this->elementIsGhost())
632 const auto globalI = this->assembler().gridGeometry().elementMapper().index(this->element());
633 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
634 A[globalI][globalI][pvIdx][pvIdx] = 1.0;
637 return NumEqVector(0.0);
641 const auto residual = this->evalLocalResidual()[0];
644 const auto globalI = this->assembler().gridGeometry().elementMapper().index(this->element());
645 const auto& scv = this->fvGeometry().scv(globalI);
646 const auto& volVars = this->curElemVolVars()[scv];
649 this->localResidual().addStorageDerivatives(A[globalI][globalI], this->problem(), this->element(), this->fvGeometry(), volVars, scv);
651 if constexpr (Problem::enableInternalDirichletConstraints())
654 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
655 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
657 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
659 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
661 if (internalDirichletConstraints[eqIdx])
663 residual[eqIdx] = volVars.priVars()[eqIdx] - dirichletValues[eqIdx];
664 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
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.
An adapter class for local assemblers using numeric differentiation.
Detects which entries in the Jacobian have to be recomputed.
A class for numeric differentiation.
A arithmetic block vector type based on DUNE's reserved vector.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The flux stencil specialized for different discretization schemes.
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethod::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:115
DiffMethod
Differentiation methods in order to compute the derivatives of the residual i.e. the entries in the j...
Definition: diffmethod.hh:37
@ green
does not need to be reassembled
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
A base class for all local cell-centered assemblers.
Definition: cclocalassembler.hh:57
void assembleResidual(SolutionVector &res)
Assemble the residual only.
Definition: cclocalassembler.hh:104
void assembleJacobian(JacobianMatrix &jac, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: cclocalassembler.hh:95
void assembleJacobianAndResidual(JacobianMatrix &jac, SolutionVector &res, GridVariables &gridVariables, const PartialReassembler *partialReassembler)
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition: cclocalassembler.hh:75
An assembler for Jacobian and residual contribution per element (cell-centered methods)
Definition: cclocalassembler.hh:121
Cell-centered scheme local assembler using numeric differentiation and implicit time discretization.
Definition: cclocalassembler.hh:132
NumEqVector assembleJacobianAndResidualImpl(JacobianMatrix &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: cclocalassembler.hh:161
Cell-centered scheme local assembler using numeric differentiation and explicit time discretization.
Definition: cclocalassembler.hh:358
NumEqVector assembleJacobianAndResidualImpl(JacobianMatrix &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: cclocalassembler.hh:379
Cell-centered scheme local assembler using analytic (hand-coded) differentiation and implicit time di...
Definition: cclocalassembler.hh:491
NumEqVector assembleJacobianAndResidualImpl(JacobianMatrix &A, const GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: cclocalassembler.hh:510
Cell-centered scheme local assembler using analytic (hand-coded) differentiation and explicit time di...
Definition: cclocalassembler.hh:608
NumEqVector assembleJacobianAndResidualImpl(JacobianMatrix &A, const GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: cclocalassembler.hh:627
A base class for all local assemblers.
Definition: fvlocalassemblerbase.hh:48
Implementation & asImp_()
Definition: fvlocalassemblerbase.hh:309
const Assembler & assembler() const
The assembler.
Definition: fvlocalassemblerbase.hh:245
const Element & element() const
The current element.
Definition: fvlocalassemblerbase.hh:249
A helper class for local assemblers using numeric differentiation to determine the epsilon.
Definition: numericepsilon.hh:41
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:425
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:496
static void partialDerivative(const Function &function, Scalar x0, FunctionEvalType &derivative, const FunctionEvalType &fx0, const int numericDifferenceMethod=1)
Computes the derivative of a function with repect to a function parameter.
Definition: numericdifferentiation.hh:61
A arithmetic block vector type based on DUNE's reserved vector.
Definition: reservedblockvector.hh:38
The flux stencil specialized for different discretization schemes.
Definition: fluxstencil.hh:45
Declares all properties used in Dumux.
The local element solution class for cell-centered methods.