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>
38 template<std::
size_t id>
using SubDomainTypeTag =
typename Traits::template SubDomain<id>::TypeTag;
41 template<std::
size_t id>
using Element =
typename GridView<id>::template Codim<0>::Entity;
43 template<std::
size_t id>
using ProblemPtr =
const Problem<id> *;
44 using ProblemPtrs =
typename Traits::template Tuple<ProblemPtr>;
46 template<std::
size_t id>
47 using SubSolutionVector
48 = std::decay_t<decltype(std::declval<typename Traits::SolutionVector>()[Dune::index_constant<id>()])>;
52 template<std::
size_t i, std::
size_t j>
73 using namespace Dune::Hybrid;
74 forEach(problems_, [](
auto&
problem){
78 forEach(curSols_, [](
auto& solutionVector){
79 solutionVector = std::make_shared<
typename std::decay_t<
decltype(solutionVector)>::element_type>();
102 template<std::
size_t i, std::
size_t j>
104 const Element<i>& elementI,
105 Dune::index_constant<j> domainJ)
const
107 static_assert(i != j,
"Domain i cannot be coupled to itself!");
109 "The coupling manager does not implement the couplingStencil() function");
120 template<std::
size_t id,
class JacobianPattern>
145 template<std::
size_t i,
class Assembler>
147 const Element<i>& elementI,
148 const Assembler& assembler)
171 template<std::
size_t i, std::
size_t j,
class LocalAssemblerI>
173 const LocalAssemblerI& localAssemblerI,
174 Dune::index_constant<j> domainJ,
175 std::size_t dofIdxGlobalJ,
176 const PrimaryVariables<j>& priVarsJ,
179 curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
196 template<std::
size_t i,
class LocalAssemblerI,
class UpdatableElementVolVars,
class UpdatableFluxVarCache>
198 const LocalAssemblerI& localAssemblerI,
199 UpdatableElementVolVars& elemVolVars,
200 UpdatableFluxVarCache& elemFluxVarsCache)
210 using namespace Dune::Hybrid;
211 forEach(integralRange(Dune::Hybrid::size(curSols_)), [&](
const auto id)
214 *std::get<id>(curSols_) =
curSol[id];
237 template<std::
size_t i, std::
size_t j,
class LocalAssemblerI>
239 const LocalAssemblerI& localAssemblerI,
240 Dune::index_constant<j> domainJ,
241 std::size_t dofIdxGlobalJ)
const
243 return localAssemblerI.evalLocalResidual();
252 template<std::
size_t i,
class LocalAssemblerI,
class JacobianMatrixDiagBlock,
class Gr
idVariables>
254 const LocalAssemblerI& localAssemblerI,
255 const typename LocalAssemblerI::LocalResidual::ElementResidualVector& origResiduals,
256 JacobianMatrixDiagBlock& A,
257 GridVariables& gridVariables)
263 template<std::
size_t i>
265 const std::string& paramGroup)
const
267 constexpr auto numEq = PrimaryVariables<i>::dimension;
275 template<
typename... SubProblems>
278 using namespace Dune::Hybrid;
279 forEach(integralRange(size(problems_)), [&](
const auto i)
288 template<
class SubProblem, std::
size_t i>
290 { std::get<i>(problems_) =
problem.get(); }
297 template<std::
size_t i>
298 const Problem<i>&
problem(Dune::index_constant<i> domainIdx)
const
300 const Problem<i>* p = std::get<i>(problems_);
301 assert(p &&
"The problem pointer is invalid. Use setSubProblems() before calling this function");
313 using namespace Dune::Hybrid;
314 forEach(integralRange(Dune::Hybrid::size(curSols_)), [&](
const auto id)
317 std::get<id>(curSols_) = Dune::stackobject_to_shared_ptr(*std::get<id>(
curSol));
326 template<std::
size_t i>
327 SubSolutionVector<i>&
curSol(Dune::index_constant<i> domainIdx)
328 {
return *std::get<i>(curSols_); }
335 template<std::
size_t i>
336 const SubSolutionVector<i>&
curSol(Dune::index_constant<i> domainIdx)
const
337 {
return *std::get<i>(curSols_); }
350 ProblemPtrs problems_;
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:37
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:197
const SubSolutionVector< i > & curSol(Dune::index_constant< i > domainIdx) const
the solution vector of the subproblem
Definition: multidomain/couplingmanager.hh:336
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:121
void attachSolution(const SolutionVectorStorage &curSol)
Attach a solution vector stored outside of this class.
Definition: multidomain/couplingmanager.hh:311
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:289
void setSubProblems(const std::tuple< std::shared_ptr< SubProblems >... > &problems)
set the pointers to the sub problems
Definition: multidomain/couplingmanager.hh:276
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:298
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:146
std::vector< std::size_t > CouplingStencilType
default type used for coupling element stencils
Definition: multidomain/couplingmanager.hh:53
SubSolutionVector< i > & curSol(Dune::index_constant< i > domainIdx)
the solution vector of the subproblem
Definition: multidomain/couplingmanager.hh:327
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:208
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:253
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:264
CouplingManager()
Default constructor.
Definition: multidomain/couplingmanager.hh:71
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:103
typename Traits::template TupleOfSharedPtr< SubSolutionVector > SolutionVectorStorage
the type in which the solution vector is stored in the manager
Definition: multidomain/couplingmanager.hh:60
typename Traits::SolutionVector SolutionVector
the type of the solution vector
Definition: multidomain/couplingmanager.hh:56
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:172
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:238
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
An adapter class for local assemblers using numeric differentiation.
Template which always yields a false value.
Definition: common/typetraits/typetraits.hh:24