13#ifndef DUMUX_MULTIDOMAIN_COUPLING_MANAGER_HH
14#define DUMUX_MULTIDOMAIN_COUPLING_MANAGER_HH
19#include <dune/common/exceptions.hh>
20#include <dune/common/indices.hh>
21#include <dune/common/shared_ptr.hh>
22#include <dune/common/hybridutilities.hh>
23#include <dune/istl/multitypeblockvector.hh>
34template<
class... Args, std::size_t ...Is>
35auto toRef(
const std::tuple<Args...>& v, std::index_sequence<Is...> indices)
49 template<std::
size_t id>
using SubDomainTypeTag =
typename Traits::template SubDomain<id>::TypeTag;
52 template<std::
size_t id>
using Element =
typename GridView<id>::template Codim<0>::Entity;
54 template<std::
size_t id>
using ProblemPtr =
const Problem<id> *;
55 using ProblemPtrs =
typename Traits::template Tuple<ProblemPtr>;
57 template<std::
size_t id>
58 using SubSolutionVector
59 = std::decay_t<decltype(std::declval<typename Traits::SolutionVector>()[Dune::index_constant<id>()])>;
63 template<std::
size_t i, std::
size_t j>
84 using namespace Dune::Hybrid;
85 forEach(problems_, [](
auto&
problem){
89 forEach(curSols_, [](
auto& solutionVector){
90 solutionVector = std::make_shared<
typename std::decay_t<
decltype(solutionVector)>::element_type>();
113 template<std::
size_t i, std::
size_t j>
115 const Element<i>& elementI,
116 Dune::index_constant<j> domainJ)
const
118 static_assert(i != j,
"Domain i cannot be coupled to itself!");
120 "The coupling manager does not implement the couplingStencil() function");
131 template<std::
size_t id,
class JacobianPattern>
156 template<std::
size_t i,
class Assembler>
158 const Element<i>& elementI,
159 const Assembler& assembler)
182 template<std::
size_t i, std::
size_t j,
class LocalAssemblerI>
184 const LocalAssemblerI& localAssemblerI,
185 Dune::index_constant<j> domainJ,
186 std::size_t dofIdxGlobalJ,
187 const PrimaryVariables<j>& priVarsJ,
190 curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
207 template<std::
size_t i,
class LocalAssemblerI,
class UpdatableElementVolVars,
class UpdatableFluxVarCache>
209 const LocalAssemblerI& localAssemblerI,
210 UpdatableElementVolVars& elemVolVars,
211 UpdatableFluxVarCache& elemFluxVarsCache)
221 using namespace Dune::Hybrid;
222 forEach(integralRange(Dune::Hybrid::size(curSols_)), [&](
const auto id)
225 *std::get<id>(curSols_) =
curSol[id];
248 template<std::
size_t i, std::
size_t j,
class LocalAssemblerI>
250 const LocalAssemblerI& localAssemblerI,
251 Dune::index_constant<j> domainJ,
252 std::size_t dofIdxGlobalJ)
const
254 return localAssemblerI.evalLocalResidual();
263 template<std::
size_t i,
class LocalAssemblerI,
class JacobianMatrixDiagBlock,
class Gr
idVariables>
265 const LocalAssemblerI& localAssemblerI,
266 const typename LocalAssemblerI::LocalResidual::ElementResidualVector& origResiduals,
267 JacobianMatrixDiagBlock& A,
268 GridVariables& gridVariables)
274 template<std::
size_t i>
276 const std::string& paramGroup)
const
278 constexpr auto numEq = PrimaryVariables<i>::dimension;
286 template<
typename... SubProblems>
289 using namespace Dune::Hybrid;
290 forEach(integralRange(size(problems_)), [&](
const auto i)
299 template<
class SubProblem, std::
size_t i>
301 { std::get<i>(problems_) =
problem.get(); }
308 template<std::
size_t i>
309 const Problem<i>&
problem(Dune::index_constant<i> domainIdx)
const
311 const Problem<i>* p = std::get<i>(problems_);
312 assert(p &&
"The problem pointer is invalid. Use setSubProblems() before calling this function");
324 using namespace Dune::Hybrid;
325 forEach(integralRange(Dune::Hybrid::size(curSols_)), [&](
const auto id)
328 std::get<id>(curSols_) = Dune::stackobject_to_shared_ptr(*std::get<id>(
curSol));
337 template<std::
size_t i>
338 SubSolutionVector<i>&
curSol(Dune::index_constant<i> domainIdx)
339 {
return *std::get<i>(curSols_); }
346 template<std::
size_t i>
347 const SubSolutionVector<i>&
curSol(Dune::index_constant<i> domainIdx)
const
348 {
return *std::get<i>(curSols_); }
361 ProblemPtrs problems_;
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:48
void updateCoupledVariables(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, UpdatableElementVolVars &elemVolVars, UpdatableFluxVarCache &elemFluxVarsCache)
update variables of domain i that depend on variables in domain j after the coupling context has been...
Definition: multidomain/couplingmanager.hh:208
const SubSolutionVector< i > & curSol(Dune::index_constant< i > domainIdx) const
the solution vector of the subproblem
Definition: multidomain/couplingmanager.hh:347
void extendJacobianPattern(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: multidomain/couplingmanager.hh:132
void attachSolution(SolutionVectorStorage &curSol)
Attach a solution vector stored outside of this class.
Definition: multidomain/couplingmanager.hh:322
void setSubProblem(std::shared_ptr< SubProblem > problem, Dune::index_constant< i > domainIdx)
set a pointer to one of the sub problems
Definition: multidomain/couplingmanager.hh:300
void setSubProblems(const std::tuple< std::shared_ptr< SubProblems >... > &problems)
set the pointers to the sub problems
Definition: multidomain/couplingmanager.hh:287
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:309
void bindCouplingContext(Dune::index_constant< i > domainI, const Element< i > &elementI, const Assembler &assembler)
prepares all data and variables that are necessary to evaluate the residual of the element of domain ...
Definition: multidomain/couplingmanager.hh:157
std::vector< std::size_t > CouplingStencilType
default type used for coupling element stencils
Definition: multidomain/couplingmanager.hh:64
SubSolutionVector< i > & curSol(Dune::index_constant< i > domainIdx)
the solution vector of the subproblem
Definition: multidomain/couplingmanager.hh:338
void updateSolution(const SolutionVector &curSol)
Updates the entire solution vector, e.g. before assembly or after grid adaption Overload might want t...
Definition: multidomain/couplingmanager.hh:219
void evalAdditionalDomainDerivatives(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, const typename LocalAssemblerI::LocalResidual::ElementResidualVector &origResiduals, JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
evaluate additional derivatives of the element residual of a domain with respect to dofs in the same ...
Definition: multidomain/couplingmanager.hh:264
decltype(auto) numericEpsilon(Dune::index_constant< i >, const std::string ¶mGroup) const
return the numeric epsilon used for deflecting primary variables of coupled domain i
Definition: multidomain/couplingmanager.hh:275
CouplingManager()
Default constructor.
Definition: multidomain/couplingmanager.hh:82
const CouplingStencilType< i, j > & couplingStencil(Dune::index_constant< i > domainI, const Element< i > &elementI, Dune::index_constant< j > domainJ) const
returns an iterable container of all indices of degrees of freedom of domain j that couple with / inf...
Definition: multidomain/couplingmanager.hh:114
typename Traits::template TupleOfSharedPtr< SubSolutionVector > SolutionVectorStorage
the type in which the solution vector is stored in the manager
Definition: multidomain/couplingmanager.hh:71
typename Traits::SolutionVector SolutionVector
the type of the solution vector
Definition: multidomain/couplingmanager.hh:67
A helper class for local assemblers using numeric differentiation to determine the epsilon.
Definition: numericepsilon.hh:29
Definition: variablesbackend.hh:31
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:183
decltype(auto) evalCouplingResidual(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< j > domainJ, std::size_t dofIdxGlobalJ) const
evaluates the element residual of a coupled element of domain i which depends on the variables at the...
Definition: multidomain/couplingmanager.hh:249
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
auto toRef(const std::tuple< Args... > &v, std::index_sequence< Is... > indices)
Definition: multidomain/couplingmanager.hh:35
An adapter class for local assemblers using numeric differentiation.
Template which always yields a false value.
Definition: common/typetraits/typetraits.hh:24