25#ifndef DUMUX_MULTIDOMAIN_EMBEDDED_EXTENDEDSOURCESTENCIL_HH
26#define DUMUX_MULTIDOMAIN_EMBEDDED_EXTENDEDSOURCESTENCIL_HH
30#include <dune/common/indices.hh>
44template<
class CouplingManager>
47 using MDTraits =
typename CouplingManager::MultiDomainTraits;
48 using Scalar =
typename MDTraits::Scalar;
50 template<std::
size_t id>
using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
52 template<std::
size_t id>
using GridView =
typename GridGeometry<id>::GridView;
53 template<std::
size_t id>
using Element =
typename GridView<id>::template Codim<0>::Entity;
55 static constexpr auto bulkIdx =
typename MDTraits::template SubDomain<0>::Index();
56 static constexpr auto lowDimIdx =
typename MDTraits::template SubDomain<1>::Index();
58 template<std::
size_t id>
59 static constexpr bool isBox()
68 template<std::
size_t id,
class JacobianPattern>
72 for (
const auto& element : elements(couplingManager.gridView(domainI)))
74 const auto& dofs = extendedSourceStencil_(couplingManager, domainI, element);
75 if constexpr (isBox<domainI>())
77 for (
int i = 0; i < element.subEntities(GridView<domainI>::dimension); ++i)
78 for (
const auto globalJ : dofs)
79 pattern.add(couplingManager.
problem(domainI).gridGeometry().vertexMapper().subIndex(element, i, GridView<domainI>::dimension), globalJ);
83 const auto globalI = couplingManager.
problem(domainI).gridGeometry().elementMapper().index(element);
84 for (
const auto globalJ : dofs)
85 pattern.add(globalI, globalJ);
97 template<std::
size_t i,
class LocalAssemblerI,
class SolutionVector,
class JacobianMatrixDiagBlock,
class Gr
idVariables>
99 Dune::index_constant<i> domainI,
100 const LocalAssemblerI& localAssemblerI,
101 const SolutionVector& curSol,
102 JacobianMatrixDiagBlock& A,
103 GridVariables& gridVariables)
const
105 constexpr auto numEq = std::decay_t<
decltype(curSol[domainI][0])>::size();
106 const auto& elementI = localAssemblerI.element();
109 if (extendedSourceStencil_(couplingManager, domainI, elementI).empty())
113 const auto origResidual = localAssemblerI.evalLocalSourceResidual(elementI);
116 for (
const auto dofIndex : extendedSourceStencil_(couplingManager, domainI, elementI))
118 auto partialDerivs = origResidual;
119 const auto origPriVars = curSol[domainI][dofIndex];
122 for (
int pvIdx = 0; pvIdx < numEq; pvIdx++)
127 auto evalResiduals = [&](Scalar priVar)
130 auto priVars = origPriVars;
131 priVars[pvIdx] = priVar;
133 return localAssemblerI.evalLocalSourceResidual(elementI);
137 static const int numDiffMethod = getParam<int>(
"Assembly.NumericDifferenceMethod");
139 partialDerivs, origResidual, numDiffMethod);
142 for (
const auto& scvJ : scvs(localAssemblerI.fvGeometry()))
144 for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
150 A[scvJ.dofIndex()][dofIndex][eqIdx][pvIdx] += partialDerivs[scvJ.indexInElement()][eqIdx];
155 couplingManager.
updateCouplingContext(domainI, localAssemblerI, domainI, dofIndex, origPriVars, pvIdx);
161 void clear() { sourceStencils_.clear(); }
164 typename CouplingManager::template CouplingStencils<bulkIdx>&
stencil()
165 {
return sourceStencils_; }
169 template<std::
size_t id>
170 const auto& extendedSourceStencil_(
const CouplingManager& couplingManager, Dune::index_constant<id> domainId,
const Element<id>& element)
const
172 if constexpr (domainId == bulkIdx)
174 const auto bulkElementIdx = couplingManager.
problem(bulkIdx).gridGeometry().elementMapper().index(element);
175 if (sourceStencils_.count(bulkElementIdx))
176 return sourceStencils_.at(bulkElementIdx);
179 return couplingManager.emptyStencil(domainId);
183 typename CouplingManager::template CouplingStencils<bulkIdx> sourceStencils_;
A class for numeric differentiation.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
void updateCouplingContext(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< j > domainJ, std::size_t dofIdxGlobalJ, const PrimaryVariables< j > &priVarsJ, int pvIdxJ)
updates all data and variables that are necessary to evaluate the residual of the element of domain i...
Definition: multidomain/couplingmanager.hh:152
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: circlepoints.hh:36
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: common/properties.hh:101
Definition: multidomain/couplingmanager.hh:46
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:264
A class managing an extended source stencil.
Definition: extendedsourcestencil.hh:46
void extendJacobianPattern(const CouplingManager &couplingManager, Dune::index_constant< id > domainI, JacobianPattern &pattern) const
extend the jacobian pattern of the diagonal block of domain i by those entries that are not already i...
Definition: extendedsourcestencil.hh:69
void evalAdditionalDomainDerivatives(CouplingManager &couplingManager, Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, const SolutionVector &curSol, JacobianMatrixDiagBlock &A, GridVariables &gridVariables) const
evaluate additional derivatives of the element residual of a domain with respect to dofs in the same ...
Definition: extendedsourcestencil.hh:98
void clear()
clear the internal data
Definition: extendedsourcestencil.hh:161
CouplingManager::template CouplingStencils< bulkIdx > & stencil()
return a reference to the stencil
Definition: extendedsourcestencil.hh:164
Declares all properties used in Dumux.