GCC Code Coverage Report


Directory: ../../../builds/dumux-repositories/
File: dumux/dumux/multidomain/boundary/freeflowporenetwork/couplingmanager.hh
Date: 2025-04-12 19:19:20
Exec Total Coverage
Lines: 91 96 94.8%
Functions: 28 28 100.0%
Branches: 65 138 47.1%

Line Branch Exec Source
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 //
7 /*!
8 * \file
9 * \ingroup FreeFlowPoreNetworkCoupling
10 * \copydoc Dumux::FreeFlowPoreNetworkCouplingManager
11 */
12
13 #ifndef DUMUX_MULTIDOMAIN_BOUNDARY_FREEFLOW_PORENETWORK_COUPLINGMANAGER_HH
14 #define DUMUX_MULTIDOMAIN_BOUNDARY_FREEFLOW_PORENETWORK_COUPLINGMANAGER_HH
15
16 #include <utility>
17 #include <memory>
18
19 #include <dune/common/indices.hh>
20 #include <dumux/common/properties.hh>
21 #include <dumux/multidomain/boundary/freeflowporenetwork/ffmassporenetwork/couplingmanager.hh>
22 #include <dumux/multidomain/boundary/freeflowporenetwork/ffmomentumporenetwork/couplingmanager.hh>
23 #include <dumux/multidomain/freeflow/couplingmanager.hh>
24 #include <dumux/multidomain/multibinarycouplingmanager.hh>
25
26 #include "couplingconditions.hh"
27 #include "couplingmapper.hh"
28
29 namespace Dumux {
30
31 namespace FreeFlowPoreNetworkDetail {
32
33 // global subdomain indices
34 static constexpr auto freeFlowMomentumIndex = Dune::index_constant<0>();
35 static constexpr auto freeFlowMassIndex = Dune::index_constant<1>();
36 static constexpr auto poreNetworkIndex = Dune::index_constant<2>();
37
38 // coupling indices
39 static constexpr auto freeFlowMassToFreeFlowMomentumIndex = Dune::index_constant<0>();
40 static constexpr auto freeFlowMomentumToPoreNetworkIndex = Dune::index_constant<1>();
41 static constexpr auto freeFlowMassToPoreNetworkIndex = Dune::index_constant<2>();
42 static constexpr auto noCouplingIdx = Dune::index_constant<99>();
43
44 constexpr auto makeCouplingManagerMap()
45 {
46 auto map = std::array<std::array<std::size_t, 3>, 3>{};
47
48 // free flow (momentum-mass)
49 map[freeFlowMomentumIndex][freeFlowMassIndex] = freeFlowMassToFreeFlowMomentumIndex;
50 map[freeFlowMassIndex][freeFlowMomentumIndex] = freeFlowMassToFreeFlowMomentumIndex;
51
52 // free flow momentum - porous medium
53 map[freeFlowMomentumIndex][poreNetworkIndex] = freeFlowMomentumToPoreNetworkIndex;
54 map[poreNetworkIndex][freeFlowMomentumIndex] = freeFlowMomentumToPoreNetworkIndex;
55
56 // free flow mass - porous medium
57 map[freeFlowMassIndex][poreNetworkIndex] = freeFlowMassToPoreNetworkIndex;
58 map[poreNetworkIndex][freeFlowMassIndex] = freeFlowMassToPoreNetworkIndex;
59
60 return map;
61 }
62
63 template<std::size_t i>
64 8218617 constexpr auto coupledDomains(Dune::index_constant<i> domainI)
65 {
66 if constexpr (i == freeFlowMomentumIndex)
67 7215534 return std::make_tuple(freeFlowMassIndex, poreNetworkIndex);
68 else if constexpr (i == freeFlowMassIndex)
69 998456 return std::make_tuple(freeFlowMomentumIndex, poreNetworkIndex);
70 else // i == poreNetworkIndex
71 4627 return std::make_tuple(freeFlowMomentumIndex, freeFlowMassIndex);
72 }
73
74 template<std::size_t i, std::size_t j>
75 constexpr auto globalToLocalDomainIndices(Dune::index_constant<i>, Dune::index_constant<j>)
76 {
77 static_assert(i <= 2 && j <= 2);
78 static_assert(i != j);
79
80 if constexpr (i < j)
81 return std::pair<Dune::index_constant<0>, Dune::index_constant<1>>{};
82 else
83 return std::pair<Dune::index_constant<1>, Dune::index_constant<0>>{};
84 }
85
86 struct CouplingMaps
87 {
88 static constexpr auto managerMap()
89 {
90 return FreeFlowPoreNetworkDetail::makeCouplingManagerMap();
91 }
92
93 template<std::size_t i, std::size_t j>
94 static constexpr auto globalToLocal(Dune::index_constant<i> domainI, Dune::index_constant<j> domainJ)
95 {
96 return FreeFlowPoreNetworkDetail::globalToLocalDomainIndices(domainI, domainJ);
97 }
98
99 template<std::size_t i>
100 8218617 static constexpr auto coupledDomains(Dune::index_constant<i> domainI)
101 {
102 return FreeFlowPoreNetworkDetail::coupledDomains(domainI);
103 }
104 };
105
106 template<class MDTraits>
107 struct CouplingManagers
108 {
109 template<std::size_t id>
110 using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
111
112 using FreeFlowTraits = MultiDomainTraits<
113 SubDomainTypeTag<freeFlowMomentumIndex>, SubDomainTypeTag<freeFlowMassIndex>
114 >;
115
116 using FreeFlowMomentumPoreNetworkTraits = MultiDomainTraits<
117 SubDomainTypeTag<freeFlowMomentumIndex>, SubDomainTypeTag<poreNetworkIndex>
118 >;
119
120 using FreeFlowMassPoreNetworkTraits = MultiDomainTraits<
121 SubDomainTypeTag<freeFlowMassIndex>, SubDomainTypeTag<poreNetworkIndex>
122 >;
123
124 using FreeFlowCouplingManager
125 = Dumux::FreeFlowCouplingManager<FreeFlowTraits>;
126 using FreeFlowMomentumPoreNetworkCouplingManager
127 = Dumux::FreeFlowMomentumPoreNetworkCouplingManager<FreeFlowMomentumPoreNetworkTraits>;
128 using FreeFlowMassPoreNetworkCouplingManager
129 = Dumux::FreeFlowMassPoreNetworkCouplingManager<FreeFlowMassPoreNetworkTraits>;
130 };
131
132 } // end namespace Detail
133
134 /*!
135 * \ingroup FreeFlowPoreNetworkCoupling
136 * \brief Coupling manager for coupling freeflow and pore-network models
137 */
138 template<class MDTraits>
139 12 class FreeFlowPoreNetworkCouplingManager
140 : public MultiBinaryCouplingManager<
141 MDTraits,
142 FreeFlowPoreNetworkDetail::CouplingMaps,
143 typename FreeFlowPoreNetworkDetail::CouplingManagers<MDTraits>::FreeFlowCouplingManager,
144 typename FreeFlowPoreNetworkDetail::CouplingManagers<MDTraits>::FreeFlowMomentumPoreNetworkCouplingManager,
145 typename FreeFlowPoreNetworkDetail::CouplingManagers<MDTraits>::FreeFlowMassPoreNetworkCouplingManager
146 >
147 {
148 using ParentType = MultiBinaryCouplingManager<
149 MDTraits,
150 FreeFlowPoreNetworkDetail::CouplingMaps,
151 typename FreeFlowPoreNetworkDetail::CouplingManagers<MDTraits>::FreeFlowCouplingManager,
152 typename FreeFlowPoreNetworkDetail::CouplingManagers<MDTraits>::FreeFlowMomentumPoreNetworkCouplingManager,
153 typename FreeFlowPoreNetworkDetail::CouplingManagers<MDTraits>::FreeFlowMassPoreNetworkCouplingManager
154 >;
155
156 using ThisType = FreeFlowPoreNetworkCouplingManager<MDTraits>;
157
158 using Scalar = typename MDTraits::Scalar;
159
160 // the sub domain type tags
161 template<std::size_t id>
162 using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
163
164 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
165 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
166 template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView;
167 template<std::size_t id> using SubControlVolumeFace = typename FVElementGeometry<id>::SubControlVolumeFace;
168 template<std::size_t id> using SubControlVolume = typename FVElementGeometry<id>::SubControlVolume;
169 template<std::size_t id> using ElementVolumeVariables = typename GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>::LocalView;
170 template<std::size_t id> using NumEqVector = typename Problem<id>::Traits::NumEqVector;
171
172 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
173 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
174 using SolutionVector = typename MDTraits::SolutionVector;
175
176 template<std::size_t id>
177 using SubSolutionVector
178 = std::decay_t<decltype(std::declval<SolutionVector>()[Dune::index_constant<id>()])>;
179
180 using CouplingConditions = FreeFlowPoreNetworkCouplingConditions<MDTraits, FreeFlowPoreNetworkCouplingManager<MDTraits>>;
181 using CouplingMapper = StaggeredFreeFlowPoreNetworkCouplingMapper;
182
183 public:
184
185 template<std::size_t i, std::size_t j>
186 using SubCouplingManager = typename ParentType::template SubCouplingManager<i, j>;
187
188 static constexpr auto freeFlowMomentumIndex = FreeFlowPoreNetworkDetail::freeFlowMomentumIndex;
189 static constexpr auto freeFlowMassIndex = FreeFlowPoreNetworkDetail::freeFlowMassIndex;
190 static constexpr auto poreNetworkIndex = FreeFlowPoreNetworkDetail::poreNetworkIndex;
191
192 public:
193 using ParentType::ParentType;
194
195 template<class GridVarsTuple>
196 2 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> freeFlowMomentumProblem,
197 std::shared_ptr<Problem<freeFlowMassIndex>> freeFlowMassProblem,
198 std::shared_ptr<Problem<poreNetworkIndex>> poreNetworkProblem,
199 GridVarsTuple&& gridVarsTuple,
200 const SolutionVector& curSol)
201 {
202 // initialize sub coupling manager that are not stationary or transient problem specific
203
3/10
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
6 this->init_(freeFlowMomentumProblem, freeFlowMassProblem, poreNetworkProblem, std::forward<GridVarsTuple>(gridVarsTuple), curSol);
204
205 // initialize stationary-specific sub coupling manager for free-flow
206 using FFSol = typename SubCouplingManager<freeFlowMomentumIndex, freeFlowMassIndex>::SolutionVectorStorage;
207
2/6
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✗ Branch 8 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
10 this->subCouplingManager(freeFlowMomentumIndex, freeFlowMassIndex).init(
208 freeFlowMomentumProblem, freeFlowMassProblem,
209 2 std::make_tuple(std::get<freeFlowMomentumIndex>(gridVarsTuple), std::get<freeFlowMassIndex>(gridVarsTuple)),
210 2 FFSol{ std::get<freeFlowMomentumIndex>(this->curSol()), std::get<freeFlowMassIndex>(this->curSol()) }
211 );
212 2 }
213
214 template<class GridVarsTuple>
215 4 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> freeFlowMomentumProblem,
216 std::shared_ptr<Problem<freeFlowMassIndex>> freeFlowMassProblem,
217 std::shared_ptr<Problem<poreNetworkIndex>> poreNetworkProblem,
218 GridVarsTuple&& gridVarsTuple,
219 const SolutionVector& curSol,
220 const SolutionVector& prevSol)
221 {
222 // initialize sub coupling manager that are not stationary or transient problem specific
223
3/10
✓ Branch 3 taken 4 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 4 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 4 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
12 this->init_(freeFlowMomentumProblem, freeFlowMassProblem, poreNetworkProblem, std::forward<GridVarsTuple>(gridVarsTuple), curSol);
224
225 using FFSol = typename SubCouplingManager<freeFlowMomentumIndex, freeFlowMassIndex>::SolutionVectorStorage;
226 using FFPrevSol = std::tuple<const SubSolutionVector<freeFlowMomentumIndex>*, const SubSolutionVector<freeFlowMassIndex>*>;
227
2/6
✓ Branch 5 taken 4 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 4 times.
✗ Branch 8 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
20 this->subCouplingManager(freeFlowMomentumIndex, freeFlowMassIndex).init(
228 freeFlowMomentumProblem, freeFlowMassProblem,
229 4 std::make_tuple(std::get<freeFlowMomentumIndex>(gridVarsTuple), std::get<freeFlowMassIndex>(gridVarsTuple)),
230 4 FFSol{ std::get<freeFlowMomentumIndex>(this->curSol()), std::get<freeFlowMassIndex>(this->curSol()) },
231 4 FFPrevSol{ &prevSol[freeFlowMomentumIndex], &prevSol[freeFlowMassIndex] }
232 );
233 4 }
234
235 /*!
236 * \brief Returns the mass flux across the coupling boundary.
237 */
238 3791 auto massCouplingCondition(Dune::index_constant<poreNetworkIndex> domainI,
239 Dune::index_constant<freeFlowMassIndex> domainJ,
240 const FVElementGeometry<poreNetworkIndex>& fvGeometry,
241 const typename FVElementGeometry<poreNetworkIndex>::SubControlVolume& scv,
242 const ElementVolumeVariables<poreNetworkIndex>& elemVolVars) const
243 {
244 3791 const auto& couplingContext = this->subApply(domainI, domainJ, [&](const auto& cm, auto&& ii, auto&& jj) -> const auto& {
245 3791 return cm.couplingContext(ii, fvGeometry, scv);
246 });
247
248
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 3791 times.
3791 const auto& freeFlowMassGridGeometry = this->subApply(domainI, domainJ, [&](const auto& cm, auto&& ii, auto&& jj) -> const auto& {
249
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 3791 times.
3791 return cm.problem(jj).gridGeometry();
250 });
251
252
2/2
✓ Branch 0 taken 12883 times.
✓ Branch 1 taken 3791 times.
16674 for (auto& c : couplingContext)
253 {
254 12883 const auto& freeFlowElement = freeFlowMassGridGeometry.element(c.scv.elementIndex());
255 12883 c.velocity = this->subCouplingManager(freeFlowMomentumIndex, freeFlowMassIndex).faceVelocity(freeFlowElement, c.scvf);
256 }
257
258 3791 return CouplingConditions::massCouplingCondition(domainI, domainJ, fvGeometry, scv, elemVolVars, couplingContext);
259 }
260
261 /*!
262 * \brief Returns the mass flux across the coupling boundary.
263 */
264 7981 auto massCouplingCondition(Dune::index_constant<freeFlowMassIndex> domainI,
265 Dune::index_constant<poreNetworkIndex> domainJ,
266 const FVElementGeometry<freeFlowMassIndex>& fvGeometry,
267 const typename FVElementGeometry<freeFlowMassIndex>::SubControlVolumeFace& scvf,
268 const ElementVolumeVariables<freeFlowMassIndex>& elemVolVars) const
269 {
270 7981 const auto& couplingContext = this->subApply(domainI, domainJ, [&](const auto& cm, auto&& ii, auto&& jj) -> const auto& {
271
1/2
✓ Branch 0 taken 7981 times.
✗ Branch 1 not taken.
7981 return cm.couplingContext(ii, fvGeometry, scvf);
272 });
273
274 7981 couplingContext.velocity = this->subCouplingManager(freeFlowMomentumIndex, freeFlowMassIndex).faceVelocity(fvGeometry.element(), scvf);
275 7981 return CouplingConditions::massCouplingCondition(domainI, domainJ, fvGeometry, scvf, elemVolVars, couplingContext);
276 }
277
278 /*!
279 * \brief Returns the energy flux across the coupling boundary.
280 */
281 1835 auto energyCouplingCondition(Dune::index_constant<poreNetworkIndex> domainI,
282 Dune::index_constant<freeFlowMassIndex> domainJ,
283 const FVElementGeometry<poreNetworkIndex>& fvGeometry,
284 const typename FVElementGeometry<poreNetworkIndex>::SubControlVolume& scv,
285 const ElementVolumeVariables<poreNetworkIndex>& elemVolVars) const
286 {
287 1835 const auto& couplingContext = this->subApply(domainI, domainJ, [&](const auto& cm, auto&& ii, auto&& jj) -> const auto& {
288 1835 return cm.couplingContext(ii, fvGeometry, scv);
289 });
290
291
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 1835 times.
1835 const auto& freeFlowMassGridGeometry = this->subApply(domainI, domainJ, [&](const auto& cm, auto&& ii, auto&& jj) -> const auto& {
292
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 1835 times.
1835 return cm.problem(jj).gridGeometry();
293 });
294
295
2/2
✓ Branch 0 taken 6825 times.
✓ Branch 1 taken 1835 times.
8660 for (auto& c : couplingContext)
296 {
297 6825 const auto& freeFlowElement = freeFlowMassGridGeometry.element(c.scv.elementIndex());
298 6825 c.velocity = this->subCouplingManager(freeFlowMomentumIndex, freeFlowMassIndex).faceVelocity(freeFlowElement, c.scvf);
299 }
300
301 1835 return CouplingConditions::energyCouplingCondition(domainI, domainJ, fvGeometry, scv, elemVolVars, couplingContext);
302 }
303
304 /*!
305 * \brief Returns the energy flux across the coupling boundary.
306 */
307 3796 auto energyCouplingCondition(Dune::index_constant<freeFlowMassIndex> domainI,
308 Dune::index_constant<poreNetworkIndex> domainJ,
309 const FVElementGeometry<freeFlowMassIndex>& fvGeometry,
310 const typename FVElementGeometry<freeFlowMassIndex>::SubControlVolumeFace& scvf,
311 const ElementVolumeVariables<freeFlowMassIndex>& elemVolVars) const
312 {
313 3796 const auto& couplingContext = this->subApply(domainI, domainJ, [&](const auto& cm, auto&& ii, auto&& jj) -> const auto& {
314
1/2
✓ Branch 0 taken 3796 times.
✗ Branch 1 not taken.
3796 return cm.couplingContext(ii, fvGeometry, scvf);
315 });
316
317 3796 couplingContext.velocity = this->subCouplingManager(freeFlowMomentumIndex, freeFlowMassIndex).faceVelocity(fvGeometry.element(), scvf);
318 3796 return CouplingConditions::energyCouplingCondition(domainI, domainJ, fvGeometry, scvf, elemVolVars, couplingContext);
319 }
320
321 //////////////////////// Conditions for FreeFlowMomentum - PoreNetwork coupling //////////
322 ///////////////////////////////////////////////////////////////////////////////////////////
323
324
2/2
✓ Branch 0 taken 93628 times.
✓ Branch 1 taken 23798 times.
117426 NumEqVector<freeFlowMomentumIndex> momentumCouplingCondition(Dune::index_constant<freeFlowMomentumIndex> domainI,
325 Dune::index_constant<poreNetworkIndex> domainJ,
326 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
327 const typename FVElementGeometry<freeFlowMomentumIndex>::SubControlVolumeFace& scvf,
328 const ElementVolumeVariables<freeFlowMomentumIndex>& elemVolVars) const
329 {
330
2/2
✓ Branch 0 taken 93628 times.
✓ Branch 1 taken 23798 times.
117426 if (scvf.isLateral())
331 93628 return NumEqVector<freeFlowMomentumIndex>(0.0);
332
333
1/2
✓ Branch 0 taken 23798 times.
✗ Branch 1 not taken.
23798 const auto& context = this->subCouplingManager(freeFlowMomentumIndex, poreNetworkIndex).couplingContext(
334 fvGeometry, scvf
335 );
336
337 23798 return CouplingConditions::momentumCouplingCondition(fvGeometry, scvf, elemVolVars, context);
338 }
339
340 93628 Scalar coupledPoreInscribedRadius(const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
341 const typename FVElementGeometry<freeFlowMomentumIndex>::SubControlVolumeFace& scvf) const
342 {
343
1/2
✓ Branch 0 taken 93628 times.
✗ Branch 1 not taken.
93628 const auto& context = this->subCouplingManager(freeFlowMomentumIndex, poreNetworkIndex).couplingContext(
344 fvGeometry, scvf
345 );
346
347 280884 const auto& pnmScv = [&]
348 {
349
3/4
✓ Branch 0 taken 93628 times.
✓ Branch 1 taken 93628 times.
✓ Branch 2 taken 187256 times.
✗ Branch 3 not taken.
187256 for (const auto& scv : scvs(context.fvGeometry))
350
2/2
✓ Branch 0 taken 93628 times.
✓ Branch 1 taken 93628 times.
187256 if (scv.dofIndex() == context.poreNetworkDofIdx)
351 93628 return scv;
352
353 DUNE_THROW(Dune::InvalidStateException, "No scv found");
354 93628 }();
355
356 93628 return context.elemVolVars[pnmScv].poreInscribedRadius();
357 }
358
359 93628 auto interfaceThroatVelocity(const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
360 const typename FVElementGeometry<freeFlowMomentumIndex>::SubControlVolumeFace& scvf) const
361 {
362
1/2
✓ Branch 0 taken 93628 times.
✗ Branch 1 not taken.
93628 const auto& context = this->subCouplingManager(freeFlowMomentumIndex, poreNetworkIndex).couplingContext(
363 fvGeometry, scvf
364 );
365
366 93628 return CouplingConditions::interfaceThroatVelocity(fvGeometry, scvf, context);
367 }
368
369 // //////////////////////// Conditions for FreeFlowMomentum - FreeFlowMass coupling //////////
370 // ///////////////////////////////////////////////////////////////////////////////////////////
371
372 /*!
373 * \brief Returns the pressure at a given sub control volume face.
374 */
375 5077382 Scalar pressure(const Element<freeFlowMomentumIndex>& element,
376 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
377 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const
378 {
379 5077382 return this->subCouplingManager(freeFlowMomentumIndex, freeFlowMassIndex).pressure(
380 element, fvGeometry, scvf
381 );
382 }
383
384 /*!
385 * \brief Returns the pressure at the center of a sub control volume corresponding to a given sub control volume face.
386 * This is used for setting a Dirichlet pressure for the mass model when a fixed pressure for the momentum balance is set at another
387 * boundary. Since the the pressure at the given scvf is solution-dependent and thus unknown a priori, we just use the value
388 * of the interior cell here.
389 */
390 Scalar cellPressure(const Element<freeFlowMomentumIndex>& element,
391 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const
392 {
393 return this->subCouplingManager(freeFlowMomentumIndex, freeFlowMassIndex).cellPressure(
394 element, scvf
395 );
396 }
397
398 /*!
399 * \brief Returns the density at a given sub control volume face.
400 */
401 Scalar density(const Element<freeFlowMomentumIndex>& element,
402 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
403 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
404 const bool considerPreviousTimeStep = false) const
405 {
406 return this->subCouplingManager(freeFlowMomentumIndex, freeFlowMassIndex).density(
407 element, fvGeometry, scvf, considerPreviousTimeStep
408 );
409 }
410
411 auto insideAndOutsideDensity(const Element<freeFlowMomentumIndex>& element,
412 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
413 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
414 const bool considerPreviousTimeStep = false) const
415 {
416 return this->subCouplingManager(freeFlowMomentumIndex, freeFlowMassIndex).insideAndOutsideDensity(
417 element, fvGeometry, scvf, considerPreviousTimeStep
418 );
419 }
420
421 /*!
422 * \brief Returns the density at a given sub control volume.
423 */
424 15448448 Scalar density(const Element<freeFlowMomentumIndex>& element,
425 const SubControlVolume<freeFlowMomentumIndex>& scv,
426 const bool considerPreviousTimeStep = false) const
427 {
428 15448448 return this->subCouplingManager(freeFlowMomentumIndex, freeFlowMassIndex).density(
429 element, scv, considerPreviousTimeStep
430 );
431 }
432
433 /*!
434 * \brief Returns the pressure at a given sub control volume face.
435 */
436 15232146 Scalar effectiveViscosity(const Element<freeFlowMomentumIndex>& element,
437 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
438 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const
439 {
440 15232146 return this->subCouplingManager(freeFlowMomentumIndex, freeFlowMassIndex).effectiveViscosity(
441 element, fvGeometry, scvf
442 );
443 }
444
445 /*!
446 * \brief Returns the velocity at a given sub control volume face.
447 */
448 9645569 auto faceVelocity(const Element<freeFlowMassIndex>& element,
449 const SubControlVolumeFace<freeFlowMassIndex>& scvf) const
450 {
451 9645569 return this->subCouplingManager(freeFlowMomentumIndex, freeFlowMassIndex).faceVelocity(
452 element, scvf
453 );
454 }
455
456 template<std::size_t i>
457 const Problem<i>& problem(Dune::index_constant<i> domainI) const
458 {
459 return this->subApply(domainI, [&](const auto& cm, auto&& ii) -> const auto& {
460 return cm.problem(ii);
461 });
462 }
463
464 template<std::size_t i, std::size_t j>
465 1262208 bool isCoupled(Dune::index_constant<i> domainI,
466 Dune::index_constant<j> domainJ,
467 const SubControlVolumeFace<i>& scvf) const
468 {
469
10/10
✓ Branch 1 taken 189488 times.
✓ Branch 2 taken 93628 times.
✓ Branch 4 taken 117426 times.
✓ Branch 5 taken 446224 times.
✓ Branch 8 taken 89795 times.
✓ Branch 9 taken 39016 times.
✓ Branch 11 taken 118930 times.
✓ Branch 12 taken 77650 times.
✓ Branch 14 taken 2892 times.
✓ Branch 15 taken 69288 times.
1605470 return this->subApply(domainI, domainJ, [&](const auto& cm, auto&& ii, auto&& jj){
470
4/4
✓ Branch 4 taken 89795 times.
✓ Branch 5 taken 39016 times.
✓ Branch 7 taken 118930 times.
✓ Branch 8 taken 77650 times.
1605470 return cm.isCoupled(ii, scvf);
471 });
472 }
473
474 /*!
475 * \brief If the boundary entity is on a coupling boundary
476 * \param domainI the domain index of domain i for which to compute the flux
477 * \param domainJ the domain index of domain j for which to compute the flux
478 * \param scv the sub control volume
479 */
480 template<std::size_t i, std::size_t j>
481 8092 bool isCoupled(Dune::index_constant<i> domainI,
482 Dune::index_constant<j> domainJ,
483 const SubControlVolume<i>& scv) const
484 {
485
4/4
✓ Branch 0 taken 3791 times.
✓ Branch 1 taken 3927 times.
✓ Branch 2 taken 201 times.
✓ Branch 3 taken 173 times.
8092 return this->subApply(domainI, domainJ, [&](const auto& cm, auto&& ii, auto&& jj){
486
4/4
✓ Branch 0 taken 3791 times.
✓ Branch 1 taken 3927 times.
✓ Branch 2 taken 201 times.
✓ Branch 3 taken 173 times.
8092 return cm.isCoupled(ii, scv);
487 });
488 }
489
490 using ParentType::couplingStencil;
491 /*!
492 * \brief returns an iterable container of all indices of degrees of freedom of domain j
493 * that couple with / influence the residual of the given sub-control volume of domain i
494 *
495 * \param domainI the domain index of domain i
496 * \param elementI the coupled element of domain í
497 * \param scvI the sub-control volume of domain i
498 * \param domainJ the domain index of domain j
499 */
500 template<std::size_t j>
501 606736 const auto& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
502 const Element<freeFlowMomentumIndex>& elementI,
503 const SubControlVolume<freeFlowMomentumIndex>& scvI,
504 Dune::index_constant<j> domainJ) const
505 {
506 static_assert(freeFlowMomentumIndex != j);
507
1/2
✓ Branch 2 taken 9152 times.
✗ Branch 3 not taken.
606736 return this->subApply(domainI, domainJ, [&](const auto& cm, auto&& ii, auto&& jj) -> const auto& {
508
1/2
✓ Branch 2 taken 9152 times.
✗ Branch 3 not taken.
606736 return cm.couplingStencil(ii, elementI, scvI, jj);
509 });
510 }
511
512 private:
513 /*
514 * \brief Initializes sub-coupling managers for stationary and transient problems
515 */
516 template<class GridVarsTuple>
517 6 void init_(std::shared_ptr<Problem<freeFlowMomentumIndex>> freeFlowMomentumProblem,
518 std::shared_ptr<Problem<freeFlowMassIndex>> freeFlowMassProblem,
519 std::shared_ptr<Problem<poreNetworkIndex>> poreNetworkProblem,
520 GridVarsTuple&& gridVarsTuple,
521 const SolutionVector& curSol)
522 {
523 6 this->updateSolution(curSol); // generic coupling manager stores tuple of shared_ptr
524
525 auto couplingMapper = std::make_shared<CouplingMapper>();
526 6 couplingMapper->update(freeFlowMomentumProblem->gridGeometry(),
527
1/2
✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
6 freeFlowMassProblem->gridGeometry(),
528
1/2
✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
6 poreNetworkProblem->gridGeometry()
529 );
530
531 // initialize the binary sub coupling managers for stationary and transient problems
532 using FFMassPNMSol = typename SubCouplingManager<freeFlowMassIndex, poreNetworkIndex>::SolutionVectorStorage;
533
3/10
✓ Branch 5 taken 6 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 6 times.
✗ Branch 8 not taken.
✓ Branch 9 taken 6 times.
✗ Branch 10 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
30 this->subCouplingManager(freeFlowMassIndex, poreNetworkIndex).init(
534 freeFlowMassProblem, poreNetworkProblem, couplingMapper,
535 6 FFMassPNMSol{ std::get<freeFlowMassIndex>(this->curSol()), std::get<poreNetworkIndex>(this->curSol()) }
536 );
537
538 using FFMomPNMSol = typename SubCouplingManager<freeFlowMomentumIndex, poreNetworkIndex>::SolutionVectorStorage;
539
4/14
✓ Branch 6 taken 6 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 6 times.
✗ Branch 9 not taken.
✓ Branch 11 taken 6 times.
✗ Branch 12 not taken.
✓ Branch 14 taken 6 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
36 this->subCouplingManager(freeFlowMomentumIndex, poreNetworkIndex).init(
540 freeFlowMomentumProblem, poreNetworkProblem,
541 6 std::make_tuple(std::get<freeFlowMomentumIndex>(gridVarsTuple), std::get<poreNetworkIndex>(gridVarsTuple)),
542 6 couplingMapper, FFMomPNMSol{ std::get<freeFlowMomentumIndex>(this->curSol()), std::get<poreNetworkIndex>(this->curSol()) }
543 );
544 6 }
545
546 };
547
548 } // end namespace Dumux
549
550 #endif
551