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;
156 using ParentType::ParentType;
175 const auto& gridGeometry = this->
assembler().gridGeometry();
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;
196 auto evalNeighborFlux = [&] (
const auto& neighbor,
const auto& scvf)
198 if (neighbor.partitionType() == Dune::GhostEntity)
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));
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;
249 if (enableGridFluxVarsCache)
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;
267 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
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)
347 return origResiduals[0];
357template<
class TypeTag,
class Assembler>
360 CCLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false>, false>
369 using Problem =
typename GridVariables::GridVolumeVariables::Problem;
374 using ParentType::ParentType;
384 if (this->
assembler().isStationaryProblem())
385 DUNE_THROW(Dune::InvalidStateException,
"Using explicit jacobian assembler with stationary local residual");
400 const auto& gridGeometry = this->
assembler().gridGeometry();
404 const auto globalI = gridGeometry.elementMapper().index(
element);
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;
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;
505 using ParentType::ParentType;
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;
537 const auto globalI = this->
assembler().gridGeometry().elementMapper().index(
element);
542 if (!this->
assembler().isStationaryProblem())
552 if (!scvf.boundary())
561 if (bcTypes.hasDirichlet() && !bcTypes.hasNeumann())
565 else if (bcTypes.hasNeumann() && !bcTypes.hasDirichlet())
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;
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;
622 using ParentType::ParentType;
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;
647 const auto globalI = this->
assembler().gridGeometry().elementMapper().index(this->
element());
648 const auto& scv = this->
fvGeometry().scv(globalI);
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;
A base class for all local assemblers.
The local element solution class for cell-centered methods.
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
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
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
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
ElementVolumeVariables & curElemVolVars()
Definition assembly/fvlocalassemblerbase.hh:253
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
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
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: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
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
@ analytic
Definition diffmethod.hh:26
@ numeric
Definition diffmethod.hh:26
@ green
does not need to be reassembled
Definition entitycolor.hh:40
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
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 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.