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>
56template<std::
size_t id,
class TypeTag,
class Assembler,
class Implementation,
bool implicit>
65 using SolutionVector =
typename Assembler::SolutionVector;
70 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
71 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
72 using ElementFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache::LocalView;
73 using Scalar =
typename GridVariables::Scalar;
75 using GridGeometry =
typename GridVariables::GridGeometry;
76 using FVElementGeometry =
typename GridGeometry::LocalView;
77 using SubControlVolume =
typename GridGeometry::SubControlVolume;
78 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
80 using Element =
typename GridView::template Codim<0>::Entity;
82 using CouplingManager =
typename Assembler::CouplingManager;
88 static constexpr auto domainId =
typename Dune::index_constant<id>();
90 using ParentType::ParentType;
95 const SolutionVector&
curSol,
113 template<
class JacobianMatrixRow,
class Gr
idVariablesTuple>
116 this->
asImp_().bindLocalViews();
121 const auto residual = this->
asImp_().assembleJacobianAndResidualImpl(jacRow[
domainId], *std::get<domainId>(gridVariables));
124 for (
const auto& scv : scvs(this->
fvGeometry()))
125 res[scv.dofIndex()] += residual[scv.localDofIndex()];
128 using namespace Dune::Hybrid;
129 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](
auto&& i)
136 auto applyDirichlet = [&] (
const auto& scvI,
137 const auto& dirichletValues,
141 res[scvI.dofIndex()][eqIdx] = this->
curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
147 for (
const auto& scvJ : scvs(this->
fvGeometry()))
148 jacRow[
domainId][scvI.dofIndex()][scvJ.dofIndex()][eqIdx] = 0.0;
151 jacRow[
domainId][scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
155 this->
asImp_().enforceDirichletConstraints(applyDirichlet);
162 template<std::size_t otherId,
class JacRow,
class GridVariables,
163 typename std::enable_if_t<(otherId == id),
int> = 0>
165 const ElementResidualVector& res, GridVariables& gridVariables)
171 template<std::size_t otherId,
class JacRow,
class GridVariables,
172 typename std::enable_if_t<(otherId != id),
int> = 0>
174 const ElementResidualVector& res, GridVariables& gridVariables)
176 this->
asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
184 this->
asImp_().bindLocalViews();
188 for (
const auto& scv : scvs(this->
fvGeometry()))
189 res[scv.dofIndex()] += residual[scv.localDofIndex()];
191 auto applyDirichlet = [&] (
const auto& scvI,
192 const auto& dirichletValues,
196 res[scvI.dofIndex()][eqIdx] = this->
curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
199 this->
asImp_().enforceDirichletConstraints(applyDirichlet);
208 ElementResidualVector residual(this->
fvGeometry().numScv());
214 const auto& curVolVars = elemVolVars[scv];
216 source *= -scv.volume()*curVolVars.extrusionFactor();
217 residual[scv.localDofIndex()] = std::move(source);
230 template<
typename ApplyFunction>
234 this->
asImp_().evalDirichletBoundaries(applyDirichlet);
236 this->
asImp_().enforceInternalDirichletConstraints(applyDirichlet);
243 template<
typename ApplyDirichletFunctionType >
250 for (
const auto& scvI : scvs(this->
fvGeometry()))
252 const auto bcTypes = this->
elemBcTypes()[scvI.localDofIndex()];
253 if (bcTypes.hasDirichlet())
255 const auto dirichletValues = this->
problem().dirichlet(this->
element(), scvI);
258 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
260 if (bcTypes.isDirichlet(eqIdx))
262 const auto pvIdx = bcTypes.eqToDirichletIndex(eqIdx);
263 assert(0 <= pvIdx && pvIdx < numEq);
264 applyDirichlet(scvI, dirichletValues, eqIdx, pvIdx);
292 if (!this->
assembler().isStationaryProblem())
312 {
return couplingManager_; }
329template<std::
size_t id,
class TypeTag,
class Assembler, DiffMethod DM = DiffMethod::numeric,
bool implicit = true>
338template<std::
size_t id,
class TypeTag,
class Assembler>
347 using ElementResidualVector =
typename ParentType::LocalResidual::ElementResidualVector;
351 using FVElementGeometry =
typename GridGeometry::LocalView;
352 using GridView =
typename GridGeometry::GridView;
353 using Element =
typename GridView::template Codim<0>::Entity;
358 static constexpr bool enableGridFluxVarsCache = getPropValue<TypeTag, Properties::EnableGridFluxVariablesCache>();
359 static constexpr bool enableGridVolVarsCache = getPropValue<TypeTag, Properties::EnableGridVolumeVariablesCache>();
360 static constexpr auto domainI = Dune::index_constant<id>();
371 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
381 const auto origResiduals = this->evalLocalResidual();
386 const auto& element = this->element();
387 const auto& fvGeometry = this->fvGeometry();
388 const auto& curSol = this->curSol()[domainI];
389 auto&& curElemVolVars = this->curElemVolVars();
392 auto elemSol =
elementSolution(element, curSol, fvGeometry.gridGeometry());
394 auto partialDerivs = origResiduals;
397 for (
auto&& scv : scvs(fvGeometry))
400 const auto dofIdx = scv.dofIndex();
401 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
402 const VolumeVariables origVolVars(curVolVars);
405 for (
int pvIdx = 0; pvIdx < numEq; pvIdx++)
409 auto evalResiduals = [&](Scalar priVar)
412 const auto localDofIndex = scv.localDofIndex();
413 elemSol[localDofIndex][pvIdx] = priVar;
414 curVolVars.update(elemSol, this->problem(), element, scv);
415 this->couplingManager().updateCouplingContext(domainI, *
this, domainI, scv.dofIndex(), elemSol[localDofIndex], pvIdx);
416 return this->evalLocalResidual();
420 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(),
"Assembly.NumericDifferenceMethod");
423 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
426 for (
auto&& scvJ : scvs(fvGeometry))
428 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
434 A[scvJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[scvJ.localDofIndex()][eqIdx];
439 curVolVars = origVolVars;
442 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
443 this->couplingManager().updateCouplingContext(domainI, *
this, domainI, scv.dofIndex(), elemSol[scv.localDofIndex()], pvIdx);
448 this->couplingManager().evalAdditionalDomainDerivatives(domainI, *
this, origResiduals, A, gridVariables);
450 return origResiduals;
459 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
461 const ElementResidualVector& res, GridVariables& gridVariables)
468 const auto& element = this->element();
469 const auto& fvGeometry = this->fvGeometry();
470 auto&& curElemVolVars = this->curElemVolVars();
471 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
474 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
477 auto updateSelf = [&] ()
481 if (enableGridFluxVarsCache)
483 if (enableGridVolVarsCache)
484 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
486 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, gridVariables.gridFluxVarsCache());
490 if (enableGridVolVarsCache)
491 this->couplingManager().updateCoupledVariables(domainI, *
this, gridVariables.curGridVolVars(), elemFluxVarsCache);
493 this->couplingManager().updateCoupledVariables(domainI, *
this, curElemVolVars, elemFluxVarsCache);
497 const auto& curSolJ = this->curSol()[domainJ];
498 for (
const auto globalJ : stencil)
501 const auto origPriVarsJ = curSolJ[globalJ];
502 auto priVarsJ = origPriVarsJ;
505 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *
this, domainJ, globalJ);
507 for (
int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
509 auto evalCouplingResidual = [&](Scalar priVar)
511 priVarsJ[pvIdx] = priVar;
512 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
514 return this->couplingManager().evalCouplingResidual(domainI, *
this, domainJ, globalJ);
518 ElementResidualVector partialDerivs(element.subEntities(dim));
520 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
521 static const int numDiffMethod = getParamFromGroup<int>(paramGroup,
"Assembly.NumericDifferenceMethod");
522 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
525 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
528 for (
auto&& scv : scvs(fvGeometry))
530 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
540 const auto bcTypes = this->elemBcTypes()[scv.localDofIndex()];
541 if (bcTypes.isCouplingDirichlet(eqIdx))
542 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = partialDerivs[scv.localDofIndex()][eqIdx];
543 else if (bcTypes.isDirichlet(eqIdx))
544 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
546 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
551 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
554 this->couplingManager().updateCouplingContext(domainI, *
this, domainJ, globalJ, priVarsJ, pvIdx);
574template<std::
size_t id,
class TypeTag,
class Assembler>
577 SubDomainBoxLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, false>, false >
583 using ElementResidualVector =
typename ParentType::LocalResidual::ElementResidualVector;
586 using FVElementGeometry =
typename GridGeometry::LocalView;
587 using GridView =
typename GridGeometry::GridView;
588 using Element =
typename GridView::template Codim<0>::Entity;
593 static constexpr bool enableGridFluxVarsCache = getPropValue<TypeTag, Properties::EnableGridFluxVariablesCache>();
594 static constexpr bool enableGridVolVarsCache = getPropValue<TypeTag, Properties::EnableGridVolumeVariablesCache>();
595 static constexpr auto domainI = Dune::index_constant<id>();
606 template<
class JacobianMatrixDiagBlock,
class Gr
idVariables>
609 if (this->assembler().isStationaryProblem())
610 DUNE_THROW(Dune::InvalidStateException,
"Using explicit jacobian assembler with stationary local residual");
613 const auto& element = this->element();
614 const auto& fvGeometry = this->fvGeometry();
615 const auto& curSol = this->curSol()[domainI];
616 auto&& curElemVolVars = this->curElemVolVars();
619 const auto origResiduals = this->evalLocalResidual();
620 const auto origStorageResiduals = this->evalLocalStorageResidual();
629 auto elemSol =
elementSolution(element, curSol, fvGeometry.gridGeometry());
632 ElementResidualVector partialDerivs(element.subEntities(dim));
635 for (
auto&& scv : scvs(fvGeometry))
638 const auto dofIdx = scv.dofIndex();
639 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
640 const VolumeVariables origVolVars(curVolVars);
643 for (
int pvIdx = 0; pvIdx < numEq; pvIdx++)
647 auto evalStorage = [&](Scalar priVar)
650 elemSol[scv.localDofIndex()][pvIdx] = priVar;
651 curVolVars.update(elemSol, this->problem(), element, scv);
652 return this->evalLocalStorageResidual();
657 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(),
"Assembly.NumericDifferenceMethod");
659 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
662 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
668 A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
672 curVolVars = origVolVars;
675 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
681 return origResiduals;
690 template<std::
size_t otherId,
class JacobianBlock,
class Gr
idVariables>
692 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.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A arithmetic block vector type based on DUNE's reserved vector.
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:115
DiffMethod
Differentiation methods in order to compute the derivatives of the residual i.e. the entries in the j...
Definition: diffmethod.hh:37
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
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
Definition: common/properties/model.hh:34
A base class for all local assemblers.
Definition: fvlocalassemblerbase.hh:48
ElementVolumeVariables & curElemVolVars()
The current element volume variables.
Definition: fvlocalassemblerbase.hh:263
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: fvlocalassemblerbase.hh:68
ElementBoundaryTypes & elemBcTypes()
The element's boundary types.
Definition: fvlocalassemblerbase.hh:279
Implementation & asImp_()
Definition: fvlocalassemblerbase.hh:307
ElementVolumeVariables & prevElemVolVars()
The element volume variables of the provious time step.
Definition: fvlocalassemblerbase.hh:267
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:259
const Assembler & assembler() const
The assembler.
Definition: fvlocalassemblerbase.hh:243
ElementFluxVariablesCache & elemFluxVarsCache()
The element flux variables cache.
Definition: fvlocalassemblerbase.hh:271
LocalResidual & localResidual()
The local residual for the current element.
Definition: fvlocalassemblerbase.hh:275
const Element & element() const
The current element.
Definition: fvlocalassemblerbase.hh:247
const SolutionVector & curSol() const
The current solution.
Definition: fvlocalassemblerbase.hh:255
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:58
void enforceDirichletConstraints(const ApplyFunction &applyDirichlet)
Enforce Dirichlet constraints.
Definition: subdomainboxlocalassembler.hh:231
void bindLocalViews()
Prepares all local views necessary for local assembly.
Definition: subdomainboxlocalassembler.hh:275
static constexpr auto domainId
export the domain id of this sub-domain
Definition: subdomainboxlocalassembler.hh:88
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:205
CouplingManager & couplingManager()
return reference to the coupling manager
Definition: subdomainboxlocalassembler.hh:311
ElementResidualVector evalLocalSourceResidual(const Element &neighbor) const
Evaluates the local source term depending on time discretization scheme.
Definition: subdomainboxlocalassembler.hh:226
void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
Incorporate Dirichlet boundary conditions.
Definition: subdomainboxlocalassembler.hh:244
void assembleResidual(SubSolutionVector &res)
Assemble the residual vector entries only.
Definition: subdomainboxlocalassembler.hh:182
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:164
const Problem & problem() const
return reference to the underlying problem
Definition: subdomainboxlocalassembler.hh:307
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:114
SubDomainBoxLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, CouplingManager &couplingManager)
Definition: subdomainboxlocalassembler.hh:93
The box scheme multidomain local assembler.
Definition: subdomainboxlocalassembler.hh:330
Box scheme multi domain local assembler using numeric differentiation and implicit time discretizatio...
Definition: subdomainboxlocalassembler.hh:342
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:372
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:460
Box scheme multi domain local assembler using numeric differentiation and explicit time discretizatio...
Definition: subdomainboxlocalassembler.hh:578
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:607
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:691
Declares all properties used in Dumux.