13#ifndef DUMUX_MULTIDOMAIN_EMBEDDED_EXTENDEDSOURCESTENCIL_HH 
   14#define DUMUX_MULTIDOMAIN_EMBEDDED_EXTENDEDSOURCESTENCIL_HH 
   18#include <dune/common/indices.hh> 
   34template<
class CouplingManager>
 
   37    using MDTraits = 
typename CouplingManager::MultiDomainTraits;
 
   38    using Scalar = 
typename MDTraits::Scalar;
 
   40    template<std::
size_t id> 
using SubDomainTypeTag = 
typename MDTraits::template SubDomain<id>::TypeTag;
 
   42    template<std::
size_t id> 
using GridView = 
typename GridGeometry<id>::GridView;
 
   43    template<std::
size_t id> 
using Element = 
typename GridView<id>::template Codim<0>::Entity;
 
   45    static constexpr auto bulkIdx = 
typename MDTraits::template SubDomain<0>::Index();
 
   46    static constexpr auto lowDimIdx = 
typename MDTraits::template SubDomain<1>::Index();
 
   48    template<std::
size_t id>
 
   49    static constexpr bool isBox()
 
   58    template<std::
size_t id, 
class JacobianPattern>
 
   62        for (
const auto& element : elements(couplingManager.gridView(domainI)))
 
   64            const auto& dofs = extendedSourceStencil_(couplingManager, domainI, element);
 
   65            if constexpr (isBox<domainI>())
 
   67                for (
int i = 0; i < element.subEntities(GridView<domainI>::dimension); ++i)
 
   68                    for (
const auto globalJ : dofs)
 
   69                        pattern.add(couplingManager.
problem(domainI).gridGeometry().vertexMapper().subIndex(element, i, GridView<domainI>::dimension), globalJ);
 
   73                const auto globalI = couplingManager.
problem(domainI).gridGeometry().elementMapper().index(element);
 
   74                for (
const auto globalJ : dofs)
 
   75                    pattern.add(globalI, globalJ);
 
 
   87    template<std::
size_t i, 
class LocalAssemblerI, 
class JacobianMatrixDiagBlock, 
class Gr
idVariables>
 
   89                                         Dune::index_constant<i> domainI,
 
   90                                         const LocalAssemblerI& localAssemblerI,
 
   91                                         JacobianMatrixDiagBlock& A,
 
   92                                         GridVariables& gridVariables)
 const 
   94        const auto& curSolI = couplingManager.
curSol(domainI);
 
   95        constexpr auto numEq = std::decay_t<
decltype(curSolI[0])>::size();
 
   96        const auto& elementI = localAssemblerI.element();
 
   99        if (extendedSourceStencil_(couplingManager, domainI, elementI).empty())
 
  103        const auto origResidual = localAssemblerI.evalLocalSourceResidual(elementI);
 
  106        for (
const auto dofIndex : extendedSourceStencil_(couplingManager, domainI, elementI))
 
  108            auto partialDerivs = origResidual;
 
  109            const auto origPriVars = curSolI[dofIndex];
 
  110            auto priVars = origPriVars;
 
  113            for (
int pvIdx = 0; pvIdx < numEq; pvIdx++)
 
  118                const auto evalResiduals = [&](
const Scalar priVar)
 
  121                    priVars[pvIdx] = priVar;
 
  123                    return localAssemblerI.evalLocalSourceResidual(elementI);
 
  128                static const int numDiffMethod = 
getParamFromGroup<int>(localAssemblerI.problem().paramGroup(), 
"Assembly.NumericDifferenceMethod");
 
  130                    evalResiduals, priVars[pvIdx], partialDerivs, origResidual, eps_(priVars[pvIdx], pvIdx), numDiffMethod
 
  134                for (
const auto& scvJ : scvs(localAssemblerI.fvGeometry()))
 
  136                    for (
int eqIdx = 0; eqIdx < numEq; eqIdx++)
 
  142                        A[scvJ.dofIndex()][dofIndex][eqIdx][pvIdx] += partialDerivs[scvJ.indexInElement()][eqIdx];
 
  147                priVars[pvIdx] = origPriVars[pvIdx];
 
  150                couplingManager.
updateCouplingContext(domainI, localAssemblerI, domainI, dofIndex, origPriVars, pvIdx);
 
 
  156    void clear() { sourceStencils_.clear(); }
 
  159    typename CouplingManager::template CouplingStencils<bulkIdx>& 
stencil()
 
  160    { 
return sourceStencils_; }
 
 
  163    const typename CouplingManager::template CouplingStencils<bulkIdx>& 
stencil()
 const 
  164    { 
return sourceStencils_; }
 
 
  168    template<std::
size_t id>
 
  169    const auto& extendedSourceStencil_(
const CouplingManager& couplingManager, Dune::index_constant<id> domainId, 
const Element<id>& element)
 const 
  171        if constexpr (domainId == bulkIdx)
 
  173            const auto bulkElementIdx = couplingManager.
problem(bulkIdx).gridGeometry().elementMapper().index(element);
 
  174            if (
auto stencil = sourceStencils_.find(bulkElementIdx); 
stencil != sourceStencils_.end())
 
  178        return couplingManager.emptyStencil(domainId);
 
  182    typename CouplingManager::template CouplingStencils<bulkIdx> sourceStencils_;
 
 
The interface of the coupling manager for multi domain problems.
Definition multidomain/couplingmanager.hh:37
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition multidomain/couplingmanager.hh:297
SubSolutionVector< i > & curSol(Dune::index_constant< i > domainIdx)
the solution vector of the subproblem
Definition multidomain/couplingmanager.hh:326
A class managing an extended source stencil.
Definition embedded/extendedsourcestencil.hh:36
void evalAdditionalDomainDerivatives(CouplingManager &couplingManager, Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, JacobianMatrixDiagBlock &A, GridVariables &gridVariables) const
evaluate additional derivatives of the element residual of a domain with respect to dofs in the same ...
Definition embedded/extendedsourcestencil.hh:88
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 embedded/extendedsourcestencil.hh:59
const CouplingManager::template CouplingStencils< bulkIdx > & stencil() const
return a const reference to the stencil
Definition embedded/extendedsourcestencil.hh:163
void clear()
clear the internal data
Definition embedded/extendedsourcestencil.hh:156
CouplingManager::template CouplingStencils< bulkIdx > & stencil()
return a reference to the stencil
Definition embedded/extendedsourcestencil.hh:159
static void partialDerivative(const Function &function, Scalar x0, FunctionEvalType &derivative, const FunctionEvalType &fx0, const int numericDifferenceMethod=1)
Computes the derivative of a function with respect to a function parameter.
Definition numericdifferentiation.hh:50
A helper class for local assemblers using numeric differentiation to determine the epsilon.
Definition numericepsilon.hh:29
Defines all properties used 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:171
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition parameters.hh:149
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition propertysystem.hh:296
The available discretization methods in Dumux.
constexpr Box box
Definition method.hh:147
Definition circlepoints.hh:24
A class for numeric differentiation.
An adapter class for local assemblers using numeric differentiation.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.