3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
multibinarycouplingmanager.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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_MULTIDOMAIN_MULTIBINARY_COUPLINGMANAGER_HH
26#define DUMUX_MULTIDOMAIN_MULTIBINARY_COUPLINGMANAGER_HH
27
28#include <utility>
29#include <memory>
30#include <dune/common/hybridutilities.hh>
33
34namespace Dumux {
35
36namespace Detail {
37
38template <std::size_t, typename Tuple>
39struct HasIndex;
40
41template <std::size_t i, typename... Indices>
42struct HasIndex<i, std::tuple<Indices...>>
43: std::disjunction<std::is_same<Dune::index_constant<i>, Indices>...>
44{};
45
46} // end namespace Detail
47
60template<class MDTraits, class CouplingMap, class ...CouplingMgrs>
62{
63 // the sub domain type tags
64 template<std::size_t id>
65 using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
66
67 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
68 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
69 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
70 template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView;
71 template<std::size_t id> using SubControlVolumeFace = typename FVElementGeometry<id>::SubControlVolumeFace;
72 template<std::size_t id> using SubControlVolume = typename FVElementGeometry<id>::SubControlVolume;
73 using CouplingStencil = std::vector<std::size_t>;
74
75 template<std::size_t id>
76 using SubCouplingManagerT = typename std::tuple_element_t<id, std::tuple<CouplingMgrs...>>;
77
78 using CMIndices = std::make_index_sequence<sizeof...(CouplingMgrs)>;
80
81 template<std::size_t id>
82 using SubSolutionVector = std::decay_t<decltype(std::declval<typename MDTraits::SolutionVector>()[Dune::index_constant<id>()])>;
83 using SolutionVectors = typename MDTraits::template TupleOfSharedPtr<SubSolutionVector>;
84
85 static constexpr auto couplingManagerMap_ = CouplingMap::managerMap();
86
88 template<std::size_t i, std::size_t j>
89 static constexpr auto globalToLocal_(Dune::index_constant<i> domainI, Dune::index_constant<j> domainJ)
90 {
91 static_assert(i <= MDTraits::numSubDomains && j <= MDTraits::numSubDomains);
92 return CouplingMap::globalToLocal(domainI, domainJ);
93 }
94
96 template<class Map, std::size_t i, std::size_t j>
97 static constexpr bool isCoupled_(Dune::index_constant<i> domainI, Dune::index_constant<j>)
98 { return Detail::HasIndex<j, std::decay_t<decltype(Map::coupledDomains(domainI))>>::value; }
99
101 template<std::size_t i, std::size_t j>
102 static constexpr auto subCouplingManagerIndex_(Dune::index_constant<i> domainI, Dune::index_constant<j> domainJ)
103 {
104 static_assert(
105 isCoupled_<CouplingMap>(domainI, domainJ),
106 "Sub-coupling manager only exists for coupled domains."
107 );
108 return couplingManagerMap_[i][j];
109 }
110
111public:
112 template<std::size_t i, std::size_t j>
113 using SubCouplingManager = SubCouplingManagerT<couplingManagerMap_[i][j]>;
114
116 {
117 using namespace Dune::Hybrid;
118 forEach(couplingManagers_, [&](auto&& couplingManager)
119 {
120 couplingManager = std::make_shared<typename std::decay_t<decltype(couplingManager)>::element_type>();
121 });
122
123 forEach(solutionVectors_, [&](auto&& solutionVector)
124 {
125 solutionVector = std::make_shared<typename std::decay_t<decltype(solutionVector)>::element_type>();
126 });
127 }
128
130 template<std::size_t i, std::size_t j>
131 auto& subCouplingManager(Dune::index_constant<i> domainI,
132 Dune::index_constant<j> domainJ)
133 {
134 constexpr auto idx = subCouplingManagerIndex_(domainI, domainJ);
135 return *std::get<idx>(couplingManagers_);
136 }
137
139 template<std::size_t i, std::size_t j>
140 const auto& subCouplingManager(Dune::index_constant<i> domainI,
141 Dune::index_constant<j> domainJ) const
142 {
143 constexpr auto idx = subCouplingManagerIndex_(domainI, domainJ);
144 return *std::get<idx>(couplingManagers_);
145 }
146
148 template<std::size_t i, std::size_t j, class Apply>
149 decltype(auto) subApply(Dune::index_constant<i> domainI,
150 Dune::index_constant<j> domainJ,
151 Apply&& apply)
152 {
153 constexpr auto localIndices = globalToLocal_(domainI, domainJ);
154 return apply(subCouplingManager(domainI, domainJ), localIndices.first, localIndices.second);
155 }
156
158 template<std::size_t i, std::size_t j, class Apply>
159 decltype(auto) subApply(Dune::index_constant<i> domainI,
160 Dune::index_constant<j> domainJ,
161 const Apply& apply) const
162 {
163 constexpr auto localIndices = globalToLocal_(domainI, domainJ);
164 return apply(subCouplingManager(domainI, domainJ), localIndices.first, localIndices.second);
165 }
166
168 template<std::size_t i, class Apply>
169 decltype(auto) subApply(Dune::index_constant<i> domainI, Apply&& apply)
170 {
171 constexpr auto dm = CouplingMap::coupledDomains(domainI);
172 static_assert(std::tuple_size_v<std::decay_t<decltype(dm)>> != 0, "subApply on uncoupled domain!");
173 constexpr auto domainJ = std::get<0>(dm);
174 constexpr auto localIndices = globalToLocal_(domainI, domainJ);
175 return apply(subCouplingManager(domainI, domainJ), localIndices.first);
176 }
177
179 template<std::size_t i, class Apply>
180 decltype(auto) subApply(Dune::index_constant<i> domainI, const Apply& apply) const
181 {
182 constexpr auto dm = CouplingMap::coupledDomains(domainI);
183 static_assert(std::tuple_size_v<std::decay_t<decltype(dm)>> != 0, "subApply on uncoupled domain!");
184 constexpr auto domainJ = std::get<0>(dm);
185 constexpr auto localIndices = globalToLocal_(domainI, domainJ);
186 return apply(subCouplingManager(domainI, domainJ), localIndices.first);
187 }
188
190 void updateSolution(const typename MDTraits::SolutionVector& curSol)
191 {
192 using namespace Dune::Hybrid;
193 forEach(integralRange(Dune::Hybrid::size(solutionVectors_)), [&](const auto id)
194 {
195 *std::get<id>(solutionVectors_) = curSol[id];
196 });
197 }
198
207 template<std::size_t id, class JacobianPattern>
208 void extendJacobianPattern(Dune::index_constant<id> domainI, JacobianPattern& pattern) const
209 {}
210
214 template<std::size_t i, class Entity, std::size_t j>
215 const auto& couplingStencil(Dune::index_constant<i> domainI,
216 const Entity& entity,
217 Dune::index_constant<j> domainJ) const
218 {
219 // if the domains are coupled according to the map, forward to sub-coupling manager
220 if constexpr (isCoupled_<CouplingMap>(domainI, domainJ))
221 return subApply(domainI, domainJ, [&](const auto& cm, auto&& ii, auto&& jj) -> const auto& {
222 return cm.couplingStencil(ii, entity, jj);
223 });
224 else
225 return emptyStencil_;
226 }
227
228 // ! evaluate coupling residual for the derivative low dim DOF with respect to bulk DOF
229 // ! we only need to evaluate the part of the residual that will be influence by the bulk DOF
230 template<std::size_t i, class LocalAssemblerI, std::size_t j>
231 decltype(auto) evalCouplingResidual(Dune::index_constant<i> domainI,
232 const SubControlVolumeFace<i>& scvfI,
233 const LocalAssemblerI& localAssemblerI,
234 Dune::index_constant<j> domainJ,
235 std::size_t dofIdxGlobalJ) const
236 {
237 if constexpr (i == j || isCoupled_<CouplingMap>(domainI, domainJ))
238 return subApply(domainI, domainJ, [&](const auto& cm, auto&& ii, auto&& jj) -> decltype(auto) {
239 return cm.evalCouplingResidual(ii, scvfI, localAssemblerI, jj, dofIdxGlobalJ);
240 });
241 else
242 {
243 DUNE_THROW(Dune::InvalidStateException,
244 "Calling evalCouplingResidual for uncoupled domains " << i << " and " << j
245 );
246 return localAssemblerI.evalLocalResidual();
247 }
248 }
249
252 template<std::size_t i, class LocalAssemblerI, std::size_t j>
253 decltype(auto) evalCouplingResidual(Dune::index_constant<i> domainI,
254 const LocalAssemblerI& localAssemblerI,
255 Dune::index_constant<j> domainJ,
256 std::size_t dofIdxGlobalJ) const
257 {
258 if constexpr (i == j || isCoupled_<CouplingMap>(domainI, domainJ))
259 return subApply(domainI, domainJ, [&](const auto& cm, auto&& ii, auto&& jj) -> decltype(auto) {
260 return cm.evalCouplingResidual(ii, localAssemblerI, jj, dofIdxGlobalJ);
261 });
262 else
263 {
264 DUNE_THROW(Dune::InvalidStateException,
265 "Calling evalCouplingResidual for uncoupled domains " << i << " and " << j
266 );
267 return localAssemblerI.evalLocalResidual();
268 }
269 }
270
273 template<std::size_t i, class LocalAssemblerI, std::size_t j>
274 decltype(auto) evalCouplingResidual(Dune::index_constant<i> domainI,
275 const LocalAssemblerI& localAssemblerI,
276 const SubControlVolume<i>& scvI,
277 Dune::index_constant<j> domainJ,
278 std::size_t dofIdxGlobalJ) const
279 {
280 if constexpr (i == j || isCoupled_<CouplingMap>(domainI, domainJ))
281 return subApply(domainI, domainJ, [&](const auto& cm, auto&& ii, auto&& jj) -> decltype(auto) {
282 return cm.evalCouplingResidual(ii, localAssemblerI, scvI, jj, dofIdxGlobalJ);
283 });
284 else
285 {
286 DUNE_THROW(Dune::InvalidStateException,
287 "Calling evalCouplingResidual for uncoupled domains " << i << " and " << j
288 );
289 return localAssemblerI.evalLocalResidual();
290 }
291 }
292
296 template<std::size_t i, class LocalAssemblerI, std::size_t j, class PrimaryVariables>
297 void updateCouplingContext(Dune::index_constant<i> domainI,
298 const LocalAssemblerI& localAssemblerI,
299 Dune::index_constant<j> domainJ,
300 const std::size_t dofIdxGlobalJ,
301 const PrimaryVariables& priVars,
302 int pvIdxJ)
303 {
304 // only one other manager needs to take action if i != j
305 if constexpr (i != j)
306 {
307 if constexpr (isCoupled_<CouplingMap>(domainI, domainJ))
308 subApply(domainI, domainJ, [&](auto& cm, auto&& ii, auto&& jj){
309 cm.updateCouplingContext(ii, localAssemblerI, jj, dofIdxGlobalJ, priVars, pvIdxJ);
310 });
311 }
312 else
313 {
314 // for i == j, we need to call all relevant managers where domainI is involved and make sure that the
315 // context is updated with respect to its own domain (I-I coupling context)
316 using namespace Dune::Hybrid;
317 forEach(CouplingMap::coupledDomains(domainI), [&](const auto domainJ){
318 static_assert(domainI != domainJ);
319 subApply(domainI, domainJ, [&](auto& cm, auto&& ii, auto&& jj){
320 cm.updateCouplingContext(ii, localAssemblerI, ii, dofIdxGlobalJ, priVars, pvIdxJ);
321 });
322 });
323 }
324 }
325
327 template<std::size_t i, class Assembler = int>
328 void bindCouplingContext(Dune::index_constant<i> domainI, const Element<i>& element, const Assembler& assembler = 0)
329 {
330 // for the coupling blocks
331 using namespace Dune::Hybrid;
332 forEach(CouplingMap::coupledDomains(domainI), [&](const auto domainJ){
333 subApply(domainI, domainJ, [&](auto& cm, auto&& ii, auto&& jj) -> void {
334 cm.bindCouplingContext(ii, element, assembler);
335 });
336 });
337 }
338
343 template<std::size_t i>
344 decltype(auto) numericEpsilon(Dune::index_constant<i> domainI,
345 const std::string& paramGroup) const
346 {
347 return subApply(domainI, [&](const auto& cm, auto&& ii) -> decltype(auto) {
348 return cm.numericEpsilon(ii, paramGroup);
349 });
350 }
351
358 template<std::size_t i, class LocalAssemblerI, class JacobianMatrixDiagBlock, class GridVariables>
359 void evalAdditionalDomainDerivatives(Dune::index_constant<i> domainI,
360 const LocalAssemblerI& localAssemblerI,
361 const typename LocalAssemblerI::LocalResidual::ElementResidualVector& origResiduals,
362 JacobianMatrixDiagBlock& A,
363 GridVariables& gridVariables)
364 {}
365
366 template<std::size_t i, class LocalAssemblerI, class UpdatableElementVolVars, class UpdatableFluxVarCache>
367 void updateCoupledVariables(Dune::index_constant<i> domainI,
368 const LocalAssemblerI& localAssemblerI,
369 UpdatableElementVolVars& elemVolVars,
370 UpdatableFluxVarCache& elemFluxVarsCache)
371 {
372 // for the coupling blocks
373 using namespace Dune::Hybrid;
374 forEach(CouplingMap::coupledDomains(domainI), [&](const auto domainJ){
375 subApply(domainI, domainJ, [&](auto& cm, auto&& ii, auto&& jj){
376 cm.updateCoupledVariables(ii, localAssemblerI, elemVolVars, elemFluxVarsCache);
377 });
378 });
379 }
380
381protected:
382 SolutionVectors& curSol()
383 { return solutionVectors_; }
384
385 const SolutionVectors& curSol() const
386 { return solutionVectors_; }
387
388private:
389 CouplingManagers couplingManagers_;
390 SolutionVectors solutionVectors_;
391
392 CouplingStencil emptyStencil_;
393};
394
395} // end namespace Dumux
396
397#endif
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
constexpr auto coupledDomains(Dune::index_constant< i > domainI)
Definition: multidomain/boundary/freeflowporenetwork/couplingmanager.hh:76
Definition: common/properties.hh:100
Definition: multibinarycouplingmanager.hh:39
Coupling manager that combines an arbitrary number of binary coupling manager (coupling two domains e...
Definition: multibinarycouplingmanager.hh:62
decltype(auto) evalCouplingResidual(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, const SubControlVolume< i > &scvI, Dune::index_constant< j > domainJ, std::size_t dofIdxGlobalJ) const
Definition: multibinarycouplingmanager.hh:274
decltype(auto) evalCouplingResidual(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< j > domainJ, std::size_t dofIdxGlobalJ) const
Definition: multibinarycouplingmanager.hh:253
SubCouplingManagerT< couplingManagerMap_[i][j]> SubCouplingManager
Definition: multibinarycouplingmanager.hh:113
const auto & couplingStencil(Dune::index_constant< i > domainI, const Entity &entity, Dune::index_constant< j > domainJ) const
Return the coupling element stencil for a given bulk domain element.
Definition: multibinarycouplingmanager.hh:215
MultiBinaryCouplingManager()
Definition: multibinarycouplingmanager.hh:115
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: multibinarycouplingmanager.hh:359
decltype(auto) numericEpsilon(Dune::index_constant< i > domainI, const std::string &paramGroup) const
return the numeric epsilon used for deflecting primary variables of coupled domain i.
Definition: multibinarycouplingmanager.hh:344
decltype(auto) subApply(Dune::index_constant< i > domainI, Dune::index_constant< j > domainJ, Apply &&apply)
apply a function to the domainI-domainJ sub coupling manager using its local indices
Definition: multibinarycouplingmanager.hh:149
auto & subCouplingManager(Dune::index_constant< i > domainI, Dune::index_constant< j > domainJ)
return the binary sub-coupling manager
Definition: multibinarycouplingmanager.hh:131
void updateCoupledVariables(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, UpdatableElementVolVars &elemVolVars, UpdatableFluxVarCache &elemFluxVarsCache)
Definition: multibinarycouplingmanager.hh:367
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: multibinarycouplingmanager.hh:208
decltype(auto) subApply(Dune::index_constant< i > domainI, Apply &&apply)
apply a function to a sub coupling manager containing this domain
Definition: multibinarycouplingmanager.hh:169
void bindCouplingContext(Dune::index_constant< i > domainI, const Element< i > &element, const Assembler &assembler=0)
Bind the coupling context for a low dim element TODO remove Assembler.
Definition: multibinarycouplingmanager.hh:328
void updateCouplingContext(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< j > domainJ, const std::size_t dofIdxGlobalJ, const PrimaryVariables &priVars, int pvIdxJ)
Update the coupling context for the bulk face residual w.r.t to the lowDim dofs.
Definition: multibinarycouplingmanager.hh:297
decltype(auto) evalCouplingResidual(Dune::index_constant< i > domainI, const SubControlVolumeFace< i > &scvfI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< j > domainJ, std::size_t dofIdxGlobalJ) const
Definition: multibinarycouplingmanager.hh:231
const SolutionVectors & curSol() const
Definition: multibinarycouplingmanager.hh:385
decltype(auto) subApply(Dune::index_constant< i > domainI, const Apply &apply) const
apply a function to a sub coupling manager containing this domain
Definition: multibinarycouplingmanager.hh:180
SolutionVectors & curSol()
Definition: multibinarycouplingmanager.hh:382
void updateSolution(const typename MDTraits::SolutionVector &curSol)
Update the solution vector before assembly.
Definition: multibinarycouplingmanager.hh:190
const auto & subCouplingManager(Dune::index_constant< i > domainI, Dune::index_constant< j > domainJ) const
return the binary sub-coupling manager
Definition: multibinarycouplingmanager.hh:140
decltype(auto) subApply(Dune::index_constant< i > domainI, Dune::index_constant< j > domainJ, const Apply &apply) const
apply a function to the domainI-domainJ sub coupling manager using its local indices
Definition: multibinarycouplingmanager.hh:159
typename makeFromIndexedType< std::tuple, PtrType, Indices >::type type
Definition: multidomain/traits.hh:89
Declares all properties used in Dumux.
Traits for multidomain problems.