26#ifndef DUMUX_MULTIDOMAIN_BOX_LOCAL_ASSEMBLER_HH
27#define DUMUX_MULTIDOMAIN_BOX_LOCAL_ASSEMBLER_HH
29#include <dune/common/reservedvector.hh>
30#include <dune/common/indices.hh>
31#include <dune/common/hybridutilities.hh>
32#include <dune/grid/common/gridenums.hh>
33#include <dune/istl/matrixindexset.hh>
58template<std::
size_t id,
class TypeTag,
class Assembler,
class Implementation,
bool implicit>
67 using SolutionVector =
typename Assembler::SolutionVector;
72 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
73 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
74 using ElementFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache::LocalView;
75 using Scalar =
typename GridVariables::Scalar;
77 using GridGeometry =
typename GridVariables::GridGeometry;
78 using FVElementGeometry =
typename GridGeometry::LocalView;
79 using SubControlVolume =
typename GridGeometry::SubControlVolume;
80 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
82 using GridView =
typename GridGeometry::GridView;
83 using Element =
typename GridView::template Codim<0>::Entity;
85 using CouplingManager =
typename Assembler::CouplingManager;
91 static constexpr auto domainId =
typename Dune::index_constant<id>();
93 using ParentType::ParentType;
98 const SolutionVector&
curSol,
116 template<
class JacobianMatrixRow,
class Gr
idVariablesTuple>
119 this->
asImp_().bindLocalViews();
126 const auto residual = this->
asImp_().assembleJacobianAndResidualImpl(jacRow[
domainId], *std::get<domainId>(gridVariables));
129 for (
const auto& scv : scvs(this->
fvGeometry()))
130 res[scv.dofIndex()] += residual[scv.localDofIndex()];
133 using namespace Dune::Hybrid;
134 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](
auto&& i)
142 using GridGeometry =
typename GridVariables::GridGeometry;
143 using GridView =
typename GridGeometry::GridView;
144 static constexpr auto dim = GridView::dimension;
146 int numVerticesLocal = this->
element().subEntities(dim);
148 for (
int i = 0; i < numVerticesLocal; ++i)
150 const auto vertex = this->
element().template subEntity<dim>(i);
152 if (vertex.partitionType() == Dune::InteriorEntity ||
153 vertex.partitionType() == Dune::BorderEntity)
160 const auto vIdx = this->
assembler().gridGeometry(
domainId).vertexMapper().index(vertex);
162 typedef typename JacobianMatrix::block_type
BlockType;
164 for (
int j = 0; j < BlockType::rows; ++j)
173 auto applyDirichlet = [&] (
const auto& scvI,
174 const auto& dirichletValues,
178 res[scvI.dofIndex()][eqIdx] = this->
curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
184 for (
const auto& scvJ : scvs(this->
fvGeometry()))
185 jacRow[
domainId][scvI.dofIndex()][scvJ.dofIndex()][eqIdx] = 0.0;
188 jacRow[
domainId][scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
192 this->
asImp_().enforceDirichletConstraints(applyDirichlet);
199 template<std::size_t otherId,
class JacRow,
class GridVariables,
200 typename std::enable_if_t<(otherId == id),
int> = 0>
202 const ElementResidualVector& res, GridVariables& gridVariables)
208 template<std::size_t otherId,
class JacRow,
class GridVariables,
209 typename std::enable_if_t<(otherId != id),
int> = 0>
211 const ElementResidualVector& res, GridVariables& gridVariables)
213 this->
asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
221 this->
asImp_().bindLocalViews();
225 for (
const auto& scv : scvs(this->
fvGeometry()))
226 res[scv.dofIndex()] += residual[scv.localDofIndex()];
228 auto applyDirichlet = [&] (
const auto& scvI,
229 const auto& dirichletValues,
233 res[scvI.dofIndex()][eqIdx] = this->
curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
236 this->
asImp_().enforceDirichletConstraints(applyDirichlet);
245 ElementResidualVector residual(this->
fvGeometry().numScv());
251 const auto& curVolVars = elemVolVars[scv];
254 residual[scv.localDofIndex()] = std::move(source);
267 template<
typename ApplyFunction>
271 this->
asImp_().evalDirichletBoundaries(applyDirichlet);
273 this->
asImp_().enforceInternalDirichletConstraints(applyDirichlet);
280 template<
typename ApplyDirichletFunctionType >
287 for (
const auto& scvI : scvs(this->
fvGeometry()))
289 const auto bcTypes = this->
elemBcTypes()[scvI.localDofIndex()];
290 if (bcTypes.hasDirichlet())
292 const auto dirichletValues = this->
problem().dirichlet(this->
element(), scvI);
295 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
297 if (bcTypes.isDirichlet(eqIdx))
299 const auto pvIdx = bcTypes.eqToDirichletIndex(eqIdx);
300 assert(0 <= pvIdx && pvIdx < numEq);
301 applyDirichlet(scvI, dirichletValues, eqIdx, pvIdx);
329 if (!this->
assembler().isStationaryProblem())
349 {
return couplingManager_; }
366template<std::
size_t id,
class TypeTag,
class Assembler, DiffMethod DM = DiffMethod::numeric,
bool implicit = true>
375template<std::
size_t id,
class TypeTag,
class Assembler>
384 using ElementResidualVector =
typename ParentType::LocalResidual::ElementResidualVector;
388 using FVElementGeometry =
typename GridGeometry::LocalView;
389 using GridView =
typename GridGeometry::GridView;
390 using Element =
typename GridView::template Codim<0>::Entity;
394 enum { dim = GridView::dimension };
396 static constexpr bool enableGridFluxVarsCache = getPropValue<TypeTag, Properties::EnableGridFluxVariablesCache>();
397 static constexpr bool enableGridVolVarsCache = getPropValue<TypeTag, Properties::EnableGridVolumeVariablesCache>();
398 static constexpr auto domainI = Dune::index_constant<id>();
409 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
419 const auto origResiduals = this->evalLocalResidual();
424 const auto& element = this->element();
425 const auto& fvGeometry = this->fvGeometry();
426 const auto& curSol = this->curSol()[domainI];
427 auto&& curElemVolVars = this->curElemVolVars();
430 auto elemSol =
elementSolution(element, curSol, fvGeometry.gridGeometry());
432 auto partialDerivs = origResiduals;
435 for (
auto&& scv : scvs(fvGeometry))
438 const auto dofIdx = scv.dofIndex();
439 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
440 const VolumeVariables origVolVars(curVolVars);
443 for (
int pvIdx = 0; pvIdx < numEq; pvIdx++)
447 auto evalResiduals = [&](Scalar priVar)
450 const auto localDofIndex = scv.localDofIndex();
451 elemSol[localDofIndex][pvIdx] = priVar;
452 curVolVars.update(elemSol, this->problem(), element, scv);
453 this->couplingManager().updateCouplingContext(domainI, *
this, domainI, scv.dofIndex(), elemSol[localDofIndex], pvIdx);
454 return this->evalLocalResidual();
458 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(),
"Assembly.NumericDifferenceMethod");
461 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
464 for (
auto&& scvJ : scvs(fvGeometry))
466 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
472 A[scvJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[scvJ.localDofIndex()][eqIdx];
477 curVolVars = origVolVars;
480 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
481 this->couplingManager().updateCouplingContext(domainI, *
this, domainI, scv.dofIndex(), elemSol[scv.localDofIndex()], pvIdx);
486 this->couplingManager().evalAdditionalDomainDerivatives(domainI, *
this, origResiduals, A, gridVariables);
488 return origResiduals;
497 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
499 const ElementResidualVector& res, GridVariables& gridVariables)
506 const auto& element = this->element();
507 const auto& fvGeometry = this->fvGeometry();
508 auto&& curElemVolVars = this->curElemVolVars();
509 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
512 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
515 auto updateCoupledVariables = [&] ()
519 if (enableGridFluxVarsCache)
521 if (enableGridVolVarsCache)
522 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
524 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, gridVariables.gridFluxVarsCache());
528 if (enableGridVolVarsCache)
529 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), elemFluxVarsCache);
531 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, elemFluxVarsCache);
535 const auto& curSolJ = this->curSol()[domainJ];
536 for (
const auto globalJ : stencil)
539 const auto origPriVarsJ = curSolJ[globalJ];
540 auto priVarsJ = origPriVarsJ;
543 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *
this, domainJ, globalJ);
545 for (
int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
547 auto evalCouplingResidual = [&](Scalar priVar)
549 priVarsJ[pvIdx] = priVar;
550 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
551 updateCoupledVariables();
552 return this->couplingManager().evalCouplingResidual(domainI, *
this, domainJ, globalJ);
556 ElementResidualVector partialDerivs(element.subEntities(dim));
558 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
559 static const int numDiffMethod = getParamFromGroup<int>(paramGroup,
"Assembly.NumericDifferenceMethod");
560 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
563 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
566 for (
auto&& scv : scvs(fvGeometry))
568 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
578 const auto bcTypes = this->elemBcTypes()[scv.localDofIndex()];
579 if (bcTypes.isCouplingDirichlet(eqIdx))
580 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = partialDerivs[scv.localDofIndex()][eqIdx];
581 else if (bcTypes.isDirichlet(eqIdx))
582 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
584 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
587 if constexpr (Problem::enableInternalDirichletConstraints())
589 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
590 if (internalDirichletConstraints[eqIdx])
591 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
598 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
601 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
608 updateCoupledVariables();
619template<std::
size_t id,
class TypeTag,
class Assembler>
622 SubDomainBoxLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, false>, false >
628 using ElementResidualVector =
typename ParentType::LocalResidual::ElementResidualVector;
631 using FVElementGeometry =
typename GridGeometry::LocalView;
632 using GridView =
typename GridGeometry::GridView;
633 using Element =
typename GridView::template Codim<0>::Entity;
636 enum { dim = GridView::dimension };
638 static constexpr bool enableGridFluxVarsCache = getPropValue<TypeTag, Properties::EnableGridFluxVariablesCache>();
639 static constexpr bool enableGridVolVarsCache = getPropValue<TypeTag, Properties::EnableGridVolumeVariablesCache>();
640 static constexpr auto domainI = Dune::index_constant<id>();
651 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
654 if (this->assembler().isStationaryProblem())
655 DUNE_THROW(Dune::InvalidStateException,
"Using explicit jacobian assembler with stationary local residual");
658 const auto& element = this->element();
659 const auto& fvGeometry = this->fvGeometry();
660 const auto& curSol = this->curSol()[domainI];
661 auto&& curElemVolVars = this->curElemVolVars();
664 const auto origResiduals = this->evalLocalResidual();
665 const auto origStorageResiduals = this->evalLocalStorageResidual();
674 auto elemSol =
elementSolution(element, curSol, fvGeometry.gridGeometry());
677 ElementResidualVector partialDerivs(element.subEntities(dim));
680 for (
auto&& scv : scvs(fvGeometry))
683 const auto dofIdx = scv.dofIndex();
684 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
685 const VolumeVariables origVolVars(curVolVars);
688 for (
int pvIdx = 0; pvIdx < numEq; pvIdx++)
692 auto evalStorage = [&](Scalar priVar)
695 elemSol[scv.localDofIndex()][pvIdx] = priVar;
696 curVolVars.update(elemSol, this->problem(), element, scv);
697 return this->evalLocalStorageResidual();
702 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(),
"Assembly.NumericDifferenceMethod");
704 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
707 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
713 A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
717 curVolVars = origVolVars;
720 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
726 return origResiduals;
735 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
737 const ElementResidualVector& res, GridVariables& gridVariables)
An enum class to define various differentiation methods available in order to compute the derivatives...
An adapter class for local assemblers using numeric differentiation.
A class for numeric differentiation.
A arithmetic block vector type based on DUNE's reserved vector.
A helper to deduce a vector with the same size as numbers of equations.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Helper classes to compute the integration elements.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
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:118
DiffMethod
Differentiation methods in order to compute the derivatives of the residual i.e. the entries in the j...
Definition: diffmethod.hh:37
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:46
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
typename BlockTypeHelper< SolutionVector, Dune::IsNumber< SolutionVector >::value >::type BlockType
Definition: nonlinear/newtonsolver.hh:198
Definition: common/pdesolver.hh:36
Scalar volume(Shape shape, Scalar inscribedRadius)
Returns the volume of a given geometry based on the inscribed radius.
Definition: poreproperties.hh:73
A base class for all local assemblers.
Definition: fvlocalassemblerbase.hh:48
ElementVolumeVariables & curElemVolVars()
The current element volume variables.
Definition: fvlocalassemblerbase.hh:265
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: fvlocalassemblerbase.hh:68
ElementBoundaryTypes & elemBcTypes()
The element's boundary types.
Definition: fvlocalassemblerbase.hh:281
Implementation & asImp_()
Definition: fvlocalassemblerbase.hh:309
ElementVolumeVariables & prevElemVolVars()
The element volume variables of the provious time step.
Definition: fvlocalassemblerbase.hh:269
ElementResidualVector evalLocalResidual() const
Convenience function to evaluate the complete local residual for the current element....
Definition: fvlocalassemblerbase.hh:120
FVElementGeometry & fvGeometry()
The global finite volume geometry.
Definition: fvlocalassemblerbase.hh:261
const Assembler & assembler() const
The assembler.
Definition: fvlocalassemblerbase.hh:245
ElementFluxVariablesCache & elemFluxVarsCache()
The element flux variables cache.
Definition: fvlocalassemblerbase.hh:273
bool elementIsGhost() const
Returns if element is a ghost entity.
Definition: fvlocalassemblerbase.hh:253
LocalResidual & localResidual()
The local residual for the current element.
Definition: fvlocalassemblerbase.hh:277
const Element & element() const
The current element.
Definition: fvlocalassemblerbase.hh:249
const SolutionVector & curSol() const
The current solution.
Definition: fvlocalassemblerbase.hh:257
A helper class for local assemblers using numeric differentiation to determine the epsilon.
Definition: numericepsilon.hh:41
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
Definition: multidomain/couplingmanager.hh:46
A base class for all box local assemblers.
Definition: subdomainboxlocalassembler.hh:60
void enforceDirichletConstraints(const ApplyFunction &applyDirichlet)
Enforce Dirichlet constraints.
Definition: subdomainboxlocalassembler.hh:268
void bindLocalViews()
Prepares all local views necessary for local assembly.
Definition: subdomainboxlocalassembler.hh:312
static constexpr auto domainId
export the domain id of this sub-domain
Definition: subdomainboxlocalassembler.hh:91
ElementResidualVector evalLocalSourceResidual(const Element &element, const ElementVolumeVariables &elemVolVars) const
Evaluates the local source term for an element and given element volume variables.
Definition: subdomainboxlocalassembler.hh:242
CouplingManager & couplingManager()
return reference to the coupling manager
Definition: subdomainboxlocalassembler.hh:348
ElementResidualVector evalLocalSourceResidual(const Element &neighbor) const
Evaluates the local source term depending on time discretization scheme.
Definition: subdomainboxlocalassembler.hh:263
void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
Incorporate Dirichlet boundary conditions.
Definition: subdomainboxlocalassembler.hh:281
void assembleResidual(SubSolutionVector &res)
Assemble the residual vector entries only.
Definition: subdomainboxlocalassembler.hh:219
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacRow &jacRow, const ElementResidualVector &res, GridVariables &gridVariables)
Assemble the entries in a coupling block of the jacobian. There is no coupling block between a domain...
Definition: subdomainboxlocalassembler.hh:201
const Problem & problem() const
return reference to the underlying problem
Definition: subdomainboxlocalassembler.hh:344
void assembleJacobianAndResidual(JacobianMatrixRow &jacRow, SubSolutionVector &res, GridVariablesTuple &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition: subdomainboxlocalassembler.hh:117
SubDomainBoxLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, CouplingManager &couplingManager)
Definition: subdomainboxlocalassembler.hh:96
The box scheme multidomain local assembler.
Definition: subdomainboxlocalassembler.hh:367
Box scheme multi domain local assembler using numeric differentiation and implicit time discretizatio...
Definition: subdomainboxlocalassembler.hh:379
ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomainboxlocalassembler.hh:410
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const ElementResidualVector &res, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomainboxlocalassembler.hh:498
Box scheme multi domain local assembler using numeric differentiation and explicit time discretizatio...
Definition: subdomainboxlocalassembler.hh:623
ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomainboxlocalassembler.hh:652
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const ElementResidualVector &res, GridVariables &gridVariables)
Computes the coupling derivatives with respect to the given element and adds them to the global matri...
Definition: subdomainboxlocalassembler.hh:736
Declares all properties used in Dumux.