25#ifndef DUMUX_BOX_LOCAL_ASSEMBLER_HH
26#define DUMUX_BOX_LOCAL_ASSEMBLER_HH
28#include <dune/istl/matrixindexset.hh>
29#include <dune/istl/bvector.hh>
52template<
class TypeTag,
class Assembler,
class Implementation,
bool implicit>
65 using ParentType::ParentType;
78 template <
class PartialReassembler = DefaultPartialReassembler>
82 this->
asImp_().bindLocalViews();
83 const auto eIdxGlobal = this->
assembler().gridGeometry().elementMapper().index(this->
element());
84 if (partialReassembler
87 const auto residual = this->
asImp_().evalLocalResidual();
88 for (
const auto& scv : scvs(this->
fvGeometry()))
89 res[scv.dofIndex()] += residual[scv.localDofIndex()];
93 const auto residual = this->
asImp_().assembleJacobianAndResidualImpl(jac, gridVariables, partialReassembler);
94 for (
const auto& scv : scvs(this->
fvGeometry()))
95 res[scv.dofIndex()] += residual[scv.localDofIndex()];
99 using GridGeometry =
typename GridVariables::GridGeometry;
100 using GridView =
typename GridGeometry::GridView;
101 static constexpr auto dim = GridView::dimension;
103 int numVerticesLocal = this->
element().subEntities(dim);
105 for (
int i = 0; i < numVerticesLocal; ++i) {
106 auto vertex = this->
element().template subEntity<dim>(i);
108 if (vertex.partitionType() == Dune::InteriorEntity ||
109 vertex.partitionType() == Dune::BorderEntity)
116 int vIdx = this->
assembler().gridGeometry().vertexMapper().index(vertex);
118 typedef typename JacobianMatrix::block_type BlockType;
119 BlockType &J = jac[vIdx][vIdx];
120 for (
int j = 0; j < BlockType::rows; ++j)
128 auto applyDirichlet = [&] (
const auto& scvI,
129 const auto& dirichletValues,
133 res[scvI.dofIndex()][eqIdx] = this->
curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
135 auto& row = jac[scvI.dofIndex()];
136 for (
auto col = row.begin(); col != row.end(); ++col)
137 row[col.index()][eqIdx] = 0.0;
139 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
142 if (this->
assembler().gridGeometry().dofOnPeriodicBoundary(scvI.dofIndex()))
144 const auto periodicDof = this->
assembler().gridGeometry().periodicallyMappedDof(scvI.dofIndex());
145 res[periodicDof][eqIdx] = this->
curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
146 const auto end = jac[periodicDof].end();
147 for (
auto it = jac[periodicDof].begin(); it != end; ++it)
148 (*it) = periodicDof != it.index() ? 0.0 : 1.0;
152 this->
asImp_().enforceDirichletConstraints(applyDirichlet);
161 this->
asImp_().bindLocalViews();
162 this->
asImp_().assembleJacobianAndResidualImpl(jac, gridVariables);
164 auto applyDirichlet = [&] (
const auto& scvI,
165 const auto& dirichletValues,
169 auto& row = jac[scvI.dofIndex()];
170 for (
auto col = row.begin(); col != row.end(); ++col)
171 row[col.index()][eqIdx] = 0.0;
173 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
176 this->
asImp_().enforceDirichletConstraints(applyDirichlet);
184 this->
asImp_().bindLocalViews();
187 for (
const auto& scv : scvs(this->
fvGeometry()))
188 res[scv.dofIndex()] += residual[scv.localDofIndex()];
190 auto applyDirichlet = [&] (
const auto& scvI,
191 const auto& dirichletValues,
195 res[scvI.dofIndex()][eqIdx] = this->
curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
198 this->
asImp_().enforceDirichletConstraints(applyDirichlet);
202 template<
typename ApplyFunction>
206 this->
asImp_().evalDirichletBoundaries(applyDirichlet);
208 this->
asImp_().enforceInternalDirichletConstraints(applyDirichlet);
214 template<
typename ApplyDirichletFunctionType >
221 for (
const auto& scvI : scvs(this->
fvGeometry()))
223 const auto bcTypes = this->
elemBcTypes()[scvI.localDofIndex()];
224 if (bcTypes.hasDirichlet())
226 const auto dirichletValues = this->
problem().dirichlet(this->
element(), scvI);
229 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
231 if (bcTypes.isDirichlet(eqIdx))
233 const auto pvIdx = bcTypes.eqToDirichletIndex(eqIdx);
234 assert(0 <= pvIdx && pvIdx < numEq);
235 applyDirichlet(scvI, dirichletValues, eqIdx, pvIdx);
253template<
class TypeTag,
class Assembler, DiffMethod diffMethod = DiffMethod::numeric,
bool implicit = true>
261template<
class TypeTag,
class Assembler>
264 BoxLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true>, true>
273 using ElementResidualVector =
typename LocalResidual::ElementResidualVector;
282 using ParentType::ParentType;
290 template <
class PartialReassembler = DefaultPartialReassembler>
315 ElementResidualVector partialDerivs(
element.subEntities(dim));
321 const auto dofIdx = scv.dofIndex();
323 const VolumeVariables origVolVars(curVolVars);
326 for (
int pvIdx = 0; pvIdx < numEq; pvIdx++)
330 auto evalResiduals = [&](Scalar priVar)
333 elemSol[scv.localDofIndex()][pvIdx] = priVar;
342 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
348 if (!partialReassembler
351 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
357 A[scvJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[scvJ.localDofIndex()][eqIdx];
363 curVolVars = origVolVars;
366 elemSol[scv.localDofIndex()][pvIdx] =
curSol[scv.dofIndex()][pvIdx];
370 return origResiduals;
380template<
class TypeTag,
class Assembler>
383 BoxLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false>, false>
392 using ElementResidualVector =
typename LocalResidual::ElementResidualVector;
399 using ParentType::ParentType;
407 template <
class PartialReassembler = DefaultPartialReassembler>
411 if (partialReassembler)
412 DUNE_THROW(Dune::NotImplemented,
"partial reassembly for explicit time discretization");
436 ElementResidualVector partialDerivs(
element.subEntities(dim));
442 const auto dofIdx = scv.dofIndex();
444 const VolumeVariables origVolVars(curVolVars);
447 for (
int pvIdx = 0; pvIdx < numEq; pvIdx++)
451 auto evalStorage = [&](Scalar priVar)
454 elemSol[scv.localDofIndex()][pvIdx] = priVar;
463 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
466 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
472 A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
476 curVolVars = origVolVars;
479 elemSol[scv.localDofIndex()][pvIdx] =
curSol[scv.dofIndex()][pvIdx];
483 return origResiduals;
492template<
class TypeTag,
class Assembler>
495 BoxLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true>, true>
502 using ElementResidualVector =
typename LocalResidual::ElementResidualVector;
506 using ParentType::ParentType;
514 template <
class PartialReassembler = DefaultPartialReassembler>
518 if (partialReassembler)
519 DUNE_THROW(Dune::NotImplemented,
"partial reassembly for analytic differentiation");
543 const auto dofIdx = scv.dofIndex();
549 if (!this->
assembler().isStationaryProblem())
550 this->
localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
559 this->
localResidual().addSourceDerivatives(A[dofIdx][dofIdx],
570 if (!scvf.boundary())
586 const auto& insideScv =
fvGeometry.scv(scvf.insideScvIdx());
587 if (this->
elemBcTypes()[insideScv.localDofIndex()].hasNeumann())
590 this->
localResidual().addRobinFluxDerivatives(A[insideScv.dofIndex()],
601 return origResiduals;
611template<
class TypeTag,
class Assembler>
614 BoxLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false>, false>
621 using ElementResidualVector =
typename LocalResidual::ElementResidualVector;
625 using ParentType::ParentType;
633 template <
class PartialReassembler = DefaultPartialReassembler>
637 if (partialReassembler)
638 DUNE_THROW(Dune::NotImplemented,
"partial reassembly for explicit time discretization");
661 const auto dofIdx = scv.dofIndex();
666 this->
localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
674 return origResiduals;
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.
A base class for all local assemblers.
An adapter class for local assemblers using numeric differentiation.
Detects which entries in the Jacobian have to be recomputed.
A class for numeric differentiation.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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:438
make the local view function available whenever we use the grid geometry
Definition adapt.hh:29
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 box assemblers.
Definition boxlocalassembler.hh:54
void assembleJacobian(JacobianMatrix &jac, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition boxlocalassembler.hh:159
void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
Evaluates Dirichlet boundaries.
Definition boxlocalassembler.hh:215
void enforceDirichletConstraints(const ApplyFunction &applyDirichlet)
Enforce Dirichlet constraints.
Definition boxlocalassembler.hh:203
void bindLocalViews()
Definition boxlocalassembler.hh:67
void assembleResidual(SolutionVector &res)
Assemble the residual only.
Definition boxlocalassembler.hh:182
void assembleJacobianAndResidual(JacobianMatrix &jac, SolutionVector &res, GridVariables &gridVariables, const PartialReassembler *partialReassembler=nullptr)
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition boxlocalassembler.hh:79
An assembler for Jacobian and residual contribution per element (box methods).
Definition boxlocalassembler.hh:254
ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix &A, GridVariables &gridVariables, const PartialReassembler *partialReassembler=nullptr)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition boxlocalassembler.hh:291
ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix &A, GridVariables &gridVariables, const PartialReassembler *partialReassembler=nullptr)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition boxlocalassembler.hh:408
ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix &A, GridVariables &gridVariables, const PartialReassembler *partialReassembler=nullptr)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition boxlocalassembler.hh:515
ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix &A, GridVariables &gridVariables, const PartialReassembler *partialReassembler=nullptr)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition boxlocalassembler.hh:634
void bindLocalViews()
Definition fvlocalassemblerbase.hh:185
ElementVolumeVariables & curElemVolVars()
Definition fvlocalassemblerbase.hh:263
ElementBoundaryTypes & elemBcTypes()
Definition fvlocalassemblerbase.hh:279
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
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
Declares all properties used in Dumux.
The local element solution class for the box method.