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>
66 static_assert(!Assembler::Problem::enableInternalDirichletConstraints(),
"Internal Dirichlet constraints are currently not implemented for cc-methods!");
70 using ParentType::ParentType;
76 template <
class PartialReassembler = DefaultPartialReassembler>
80 this->
asImp_().bindLocalViews();
81 const auto globalI = this->
assembler().gridGeometry().elementMapper().index(this->
element());
82 if (partialReassembler
85 res[globalI] = this->
asImp_().evalLocalResidual()[0];
89 res[globalI] = this->
asImp_().assembleJacobianAndResidualImpl(jac, gridVariables);
99 this->
asImp_().bindLocalViews();
100 this->
asImp_().assembleJacobianAndResidualImpl(jac, gridVariables);
108 this->
asImp_().bindLocalViews();
109 const auto globalI = this->
assembler().gridGeometry().elementMapper().index(this->
element());
110 res[globalI] = this->
asImp_().evalLocalResidual()[0];
122template<
class TypeTag,
class Assembler, DiffMethod diffMethod = DiffMethod::numeric,
bool implicit = true>
130template<
class TypeTag,
class Assembler>
133 CCLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true>, true >
141 using FVElementGeometry =
typename GridGeometry::LocalView;
149 static constexpr int maxElementStencilSize = GridGeometry::maxElementStencilSize;
154 using ParentType::ParentType;
173 const auto& gridGeometry = this->
assembler().gridGeometry();
178 const auto globalI = gridGeometry.elementMapper().index(
element);
179 const auto& connectivityMap = gridGeometry.connectivityMap();
180 const auto numNeighbors = connectivityMap[globalI].size();
183 Dune::ReservedVector<Element, maxElementStencilSize> neighborElements;
184 neighborElements.resize(numNeighbors);
188 Residuals origResiduals(numNeighbors + 1); origResiduals = 0.0;
194 auto evalNeighborFlux = [&] (
const auto& neighbor,
const auto& scvf)
196 if (neighbor.partitionType() == Dune::GhostEntity)
197 return NumEqVector(0.0);
209 for (
const auto& dataJ : connectivityMap[globalI])
211 neighborElements[j-1] = gridGeometry.element(dataJ.globalJ);
212 for (
const auto scvfIdx : dataJ.scvfsJ)
213 origResiduals[j] += evalNeighborFlux(neighborElements[j-1],
fvGeometry.scvf(scvfIdx));
225 const auto origPriVars =
curSol[globalI];
226 const auto origVolVars = curVolVars;
233 Residuals partialDerivs(numNeighbors + 1);
235 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
239 auto evalResiduals = [&](Scalar priVar)
241 Residuals partialDerivsTmp(numNeighbors + 1);
242 partialDerivsTmp = 0.0;
244 elemSol[0][pvIdx] = priVar;
247 if (enableGridFluxVarsCache)
254 for (std::size_t k = 0; k < numNeighbors; ++k)
255 for (
auto scvfIdx : connectivityMap[globalI][k].scvfsJ)
256 partialDerivsTmp[k+1] += evalNeighborFlux(neighborElements[k],
fvGeometry.scvf(scvfIdx));
258 return partialDerivsTmp;
265 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
272 partialDerivs[0] = 0.0;
273 partialDerivs[0][pvIdx] = 1.0;
278 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
281 A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx];
285 for (
const auto& dataJ : connectivityMap[globalI])
286 A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j++][eqIdx];
290 curVolVars = origVolVars;
293 elemSol[0][pvIdx] = origPriVars[pvIdx];
302 if (enableGridFluxVarsCache)
306 return origResiduals[0];
316template<
class TypeTag,
class Assembler>
319 CCLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false>, false>
332 using ParentType::ParentType;
342 if (this->
assembler().isStationaryProblem())
343 DUNE_THROW(Dune::InvalidStateException,
"Using explicit jacobian assembler with stationary local residual");
358 const auto& gridGeometry = this->
assembler().gridGeometry();
362 const auto globalI = gridGeometry.elementMapper().index(
element);
369 const auto origPriVars =
curSol[globalI];
370 const auto origVolVars = curVolVars;
375 NumEqVector partialDeriv;
378 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
383 auto evalStorage = [&](Scalar priVar)
387 elemSol[0][pvIdx] = priVar;
398 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
404 else partialDeriv[pvIdx] = 1.0;
408 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
409 A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx];
412 curVolVars = origVolVars;
415 elemSol[0][pvIdx] = origPriVars[pvIdx];
428template<
class TypeTag,
class Assembler>
431 CCLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true>, true>
442 using ParentType::ParentType;
455 const auto globalI = this->
assembler().gridGeometry().elementMapper().index(this->
element());
456 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
457 A[globalI][globalI][pvIdx][pvIdx] = 1.0;
460 return NumEqVector(0.0);
474 const auto globalI = this->
assembler().gridGeometry().elementMapper().index(
element);
479 if (!this->
assembler().isStationaryProblem())
489 if (!scvf.boundary())
498 if (bcTypes.hasDirichlet() && !bcTypes.hasNeumann())
502 else if (bcTypes.hasNeumann() && !bcTypes.hasDirichlet())
506 DUNE_THROW(Dune::NotImplemented,
"Mixed boundary conditions. Use pure boundary conditions by converting Dirichlet BCs to Robin BCs");
520template<
class TypeTag,
class Assembler>
523 CCLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false>, false>
534 using ParentType::ParentType;
547 const auto globalI = this->
assembler().gridGeometry().elementMapper().index(this->
element());
548 for (
int pvIdx = 0; pvIdx < numEq; ++pvIdx)
549 A[globalI][globalI][pvIdx][pvIdx] = 1.0;
552 return NumEqVector(0.0);
559 const auto globalI = this->
assembler().gridGeometry().elementMapper().index(this->
element());
560 const auto& scv = this->
fvGeometry().scv(globalI);
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.
A base class for all local assemblers.
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
@ analytic
Definition diffmethod.hh:38
@ numeric
Definition diffmethod.hh:38
@ green
does not need to be reassembled
Definition entitycolor.hh:52
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition parameters.hh:375
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:153
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:106
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:97
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:77
An assembler for Jacobian and residual contribution per element (cell-centered methods).
Definition cclocalassembler.hh:123
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:162
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:340
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:450
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:542
ElementVolumeVariables & curElemVolVars()
Definition fvlocalassemblerbase.hh:263
Implementation & asImp_()
Definition fvlocalassemblerbase.hh:307
ElementResidualVector evalLocalResidual() const
Definition fvlocalassemblerbase.hh:120
const Problem & problem() const
Definition fvlocalassemblerbase.hh:239
FVLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol)
Definition fvlocalassemblerbase.hh:73
FVElementGeometry & fvGeometry()
Definition fvlocalassemblerbase.hh:259
const Assembler & assembler() const
Definition fvlocalassemblerbase.hh:243
ElementFluxVariablesCache & elemFluxVarsCache()
Definition fvlocalassemblerbase.hh:271
bool elementIsGhost() const
Definition fvlocalassemblerbase.hh:251
LocalResidual & localResidual()
Definition fvlocalassemblerbase.hh:275
const Element & element() const
Definition fvlocalassemblerbase.hh:247
VolumeVariables & getVolVarAccess(GridVolumeVariables &gridVolVars, ElementVolumeVariables &elemVolVars, const SubControlVolume &scv)
Definition fvlocalassemblerbase.hh:314
ElementResidualVector evalLocalStorageResidual() const
Definition fvlocalassemblerbase.hh:176
const SolutionVector & curSol() const
Definition fvlocalassemblerbase.hh:255
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:424
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.