version 3.11-dev
multidomain/couplingmanager.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_MULTIDOMAIN_COUPLING_MANAGER_HH
14#define DUMUX_MULTIDOMAIN_COUPLING_MANAGER_HH
15
16#include <memory>
17#include <tuple>
18#include <vector>
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>
24
28
29namespace Dumux {
30
35template<class Traits>
37{
38 template<std::size_t id> using SubDomainTypeTag = typename Traits::template SubDomain<id>::TypeTag;
39 template<std::size_t id> using PrimaryVariables = GetPropType<SubDomainTypeTag<id>, Properties::PrimaryVariables>;
40 template<std::size_t id> using GridView = typename GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>::GridView;
41 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
42 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
43 template<std::size_t id> using ProblemPtr = const Problem<id> *;
44 using ProblemPtrs = typename Traits::template Tuple<ProblemPtr>;
45
46 template<std::size_t id>
47 using SubSolutionVector
48 = std::decay_t<decltype(std::declval<typename Traits::SolutionVector>()[Dune::index_constant<id>()])>;
49
50public:
52 template<std::size_t i, std::size_t j>
53 using CouplingStencilType = std::vector<std::size_t>;
54
56 using SolutionVector = typename Traits::SolutionVector;
57
58protected:
60 using SolutionVectorStorage = typename Traits::template TupleOfSharedPtr<SubSolutionVector>;
61
62public:
72 {
73 using namespace Dune::Hybrid;
74 forEach(problems_, [](auto& problem){
75 problem = nullptr;
76 });
77
78 forEach(curSols_, [](auto& solutionVector){
79 solutionVector = std::make_shared<typename std::decay_t<decltype(solutionVector)>::element_type>();
80 });
81 }
82
86 // \{
87
102 template<std::size_t i, std::size_t j>
103 const CouplingStencilType<i, j>& couplingStencil(Dune::index_constant<i> domainI,
104 const Element<i>& elementI,
105 Dune::index_constant<j> domainJ) const
106 {
107 static_assert(i != j, "Domain i cannot be coupled to itself!");
108 static_assert(AlwaysFalse<Dune::index_constant<i>>::value,
109 "The coupling manager does not implement the couplingStencil() function");
110 }
111
120 template<std::size_t id, class JacobianPattern>
121 void extendJacobianPattern(Dune::index_constant<id> domainI, JacobianPattern& pattern) const
122 {}
123
124 // \}
125
129 // \{
130
145 template<std::size_t i, class Assembler>
146 void bindCouplingContext(Dune::index_constant<i> domainI,
147 const Element<i>& elementI,
148 const Assembler& assembler)
149 {}
150
151
171 template<std::size_t i, std::size_t j, class LocalAssemblerI>
172 void updateCouplingContext(Dune::index_constant<i> domainI,
173 const LocalAssemblerI& localAssemblerI,
174 Dune::index_constant<j> domainJ,
175 std::size_t dofIdxGlobalJ,
176 const PrimaryVariables<j>& priVarsJ,
177 int pvIdxJ)
178 {
179 curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
180 }
181
196 template<std::size_t i, class LocalAssemblerI, class UpdatableElementVars, class UpdatableFluxVarCache>
197 void updateCoupledVariables(Dune::index_constant<i> domainI,
198 const LocalAssemblerI& localAssemblerI,
199 UpdatableElementVars& elemVars,
200 UpdatableFluxVarCache& elemFluxVarsCache)
201 {}
202
218 template<std::size_t i, class LocalAssemblerI, class UpdatableElementVars>
219 void updateCoupledVariables(Dune::index_constant<i> domainI,
220 const LocalAssemblerI& localAssemblerI,
221 UpdatableElementVars& elemVars)
222 {}
223
230 {
231 using namespace Dune::Hybrid;
232 forEach(integralRange(Dune::Hybrid::size(curSols_)), [&](const auto id)
233 {
234 // copy external solution into object stored in this class
235 *std::get<id>(curSols_) = curSol[id];
236 });
237 }
238
239 // \}
240
258 template<std::size_t i, std::size_t j, class LocalAssemblerI>
259 decltype(auto) evalCouplingResidual(Dune::index_constant<i> domainI,
260 const LocalAssemblerI& localAssemblerI,
261 Dune::index_constant<j> domainJ,
262 std::size_t dofIdxGlobalJ) const
263 {
264 return localAssemblerI.evalLocalResidual();
265 }
266
273 template<std::size_t i, class LocalAssemblerI, class JacobianMatrixDiagBlock, class GridVariables>
274 void evalAdditionalDomainDerivatives(Dune::index_constant<i> domainI,
275 const LocalAssemblerI& localAssemblerI,
276 const typename LocalAssemblerI::LocalResidual::ElementResidualVector& origResiduals,
277 JacobianMatrixDiagBlock& A,
278 GridVariables& gridVariables)
279 {}
280
284 template<std::size_t i>
285 decltype(auto) numericEpsilon(Dune::index_constant<i>,
286 const std::string& paramGroup) const
287 {
288 constexpr auto numEq = PrimaryVariables<i>::dimension;
290 }
291
296 template<typename... SubProblems>
297 void setSubProblems(const std::tuple<std::shared_ptr<SubProblems>...>& problems)
298 {
299 using namespace Dune::Hybrid;
300 forEach(integralRange(size(problems_)), [&](const auto i)
301 { setSubProblem(std::get<i>(problems), i); });
302 }
303
309 template<class SubProblem, std::size_t i>
310 void setSubProblem(std::shared_ptr<SubProblem> problem, Dune::index_constant<i> domainIdx)
311 { std::get<i>(problems_) = problem.get(); }
312
318 template<std::size_t i>
319 const Problem<i>& problem(Dune::index_constant<i> domainIdx) const
320 {
321 const Problem<i>* p = std::get<i>(problems_);
322 assert(p && "The problem pointer is invalid. Use setSubProblems() before calling this function");
323 return *p;
324 }
325
326protected:
333 {
334 using namespace Dune::Hybrid;
335 forEach(integralRange(Dune::Hybrid::size(curSols_)), [&](const auto id)
336 {
337 // do not take ownership of the external pointer's object
338 std::get<id>(curSols_) = Dune::stackobject_to_shared_ptr(*std::get<id>(curSol));
339 });
340 }
341
347 template<std::size_t i>
348 SubSolutionVector<i>& curSol(Dune::index_constant<i> domainIdx)
349 { return *std::get<i>(curSols_); }
350
356 template<std::size_t i>
357 const SubSolutionVector<i>& curSol(Dune::index_constant<i> domainIdx) const
358 { return *std::get<i>(curSols_); }
359
360private:
365 SolutionVectorStorage curSols_;
366
371 ProblemPtrs problems_;
372};
373
374} // end namespace Dumux
375
376#endif
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:37
const SubSolutionVector< i > & curSol(Dune::index_constant< i > domainIdx) const
the solution vector of the subproblem
Definition: multidomain/couplingmanager.hh:357
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:332
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:310
void setSubProblems(const std::tuple< std::shared_ptr< SubProblems >... > &problems)
set the pointers to the sub problems
Definition: multidomain/couplingmanager.hh:297
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:319
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:348
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:229
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:274
decltype(auto) numericEpsilon(Dune::index_constant< i >, const std::string &paramGroup) const
return the numeric epsilon used for deflecting primary variables of coupled domain i
Definition: multidomain/couplingmanager.hh:285
void updateCoupledVariables(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, UpdatableElementVars &elemVars, 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
void updateCoupledVariables(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, UpdatableElementVars &elemVars)
update variables of domain i that depend on variables in domain j after the coupling context has been...
Definition: multidomain/couplingmanager.hh:219
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:32
Defines all properties used in Dumux.
Type traits.
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:259
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Definition: adapt.hh:17
An adapter class for local assemblers using numeric differentiation.
Template which always yields a false value.
Definition: common/typetraits/typetraits.hh:24