3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
multidomain/boundary/freeflowporousmedium/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 * 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_BOUNDARY_FREEFLOW_POROUSMEDIUM_COUPLINGMANAGER_HH
26#define DUMUX_MULTIDOMAIN_BOUNDARY_FREEFLOW_POROUSMEDIUM_COUPLINGMANAGER_HH
27
28#include <utility>
29#include <memory>
30
31#include <dune/common/indices.hh>
37
38#include "couplingconditions.hh"
39
40namespace Dumux {
41
42namespace FreeFlowPorousMediumDetail {
43
44// global subdomain indices
45static constexpr auto freeFlowMomentumIndex = Dune::index_constant<0>();
46static constexpr auto freeFlowMassIndex = Dune::index_constant<1>();
47static constexpr auto porousMediumIndex = Dune::index_constant<2>();
48
49// coupling indices
50static constexpr auto freeFlowMassToFreeFlowMomentumIndex = Dune::index_constant<0>();
51static constexpr auto freeFlowMomentumToPorousMediumIndex = Dune::index_constant<1>();
52static constexpr auto freeFlowMassToPorousMediumIndex = Dune::index_constant<2>();
53static constexpr auto noCouplingIdx = Dune::index_constant<99>();
54
55constexpr auto makeCouplingManagerMap()
56{
57 auto map = std::array<std::array<std::size_t, 3>, 3>{};
58
59 // free flow (momentum-mass)
62
63 // free flow momentum - porous medium
66
67 // free flow mass - porous medium
70
71 return map;
72}
73
74template<std::size_t i>
75constexpr auto coupledDomains(Dune::index_constant<i> domainI)
76{
77 if constexpr (i == freeFlowMomentumIndex)
78 return std::make_tuple(freeFlowMassIndex, porousMediumIndex);
79 else if constexpr (i == freeFlowMassIndex)
80 return std::make_tuple(freeFlowMomentumIndex, porousMediumIndex);
81 else // i == porousMediumIndex
82 return std::make_tuple(freeFlowMomentumIndex, freeFlowMassIndex);
83}
84
85template<std::size_t i, std::size_t j>
86constexpr auto globalToLocalDomainIndices(Dune::index_constant<i>, Dune::index_constant<j>)
87{
88 static_assert(i <= 2 && j <= 2);
89 static_assert(i != j);
90
91 if constexpr (i < j)
92 return std::pair<Dune::index_constant<0>, Dune::index_constant<1>>{};
93 else
94 return std::pair<Dune::index_constant<1>, Dune::index_constant<0>>{};
95}
96
98{
99 static constexpr auto managerMap()
100 {
102 }
103
104 template<std::size_t i, std::size_t j>
105 static constexpr auto globalToLocal(Dune::index_constant<i> domainI, Dune::index_constant<j> domainJ)
106 {
108 }
109
110 template<std::size_t i>
111 static constexpr auto coupledDomains(Dune::index_constant<i> domainI)
112 {
114 }
115};
116
117template<class MDTraits>
119{
120 template<std::size_t id>
121 using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
122
125 >;
126
129 >;
130
133 >;
134
141};
142
143} // end namespace Detail
144
145
150template<class MDTraits>
153 MDTraits,
154 FreeFlowPorousMediumDetail::CouplingMaps,
155 typename FreeFlowPorousMediumDetail::CouplingManagers<MDTraits>::FreeFlowCouplingManager,
156 typename FreeFlowPorousMediumDetail::CouplingManagers<MDTraits>::FreeFlowMomentumPorousMediumCouplingManager,
157 typename FreeFlowPorousMediumDetail::CouplingManagers<MDTraits>::FreeFlowMassPorousMediumCouplingManager
158>
159{
161 MDTraits,
166 >;
167
168 using Scalar = typename MDTraits::Scalar;
169
170 // the sub domain type tags
171 template<std::size_t id>
172 using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
173
174 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
175 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
176 template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView;
177 template<std::size_t id> using SubControlVolumeFace = typename FVElementGeometry<id>::SubControlVolumeFace;
178 template<std::size_t id> using SubControlVolume = typename FVElementGeometry<id>::SubControlVolume;
179 template<std::size_t id> using ElementVolumeVariables = typename GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>::LocalView;
180 template<std::size_t id> using NumEqVector = typename Problem<id>::Traits::NumEqVector;
181
182 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
183 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
184 using SolutionVector = typename MDTraits::SolutionVector;
185
187
188public:
189
190 template<std::size_t i, std::size_t j>
191 using SubCouplingManager = typename ParentType::template SubCouplingManager<i, j>;
192
196
197public:
198 using ParentType::ParentType;
199
200 template<class GridVarsTuple>
201 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> freeFlowMomentumProblem,
202 std::shared_ptr<Problem<freeFlowMassIndex>> freeFlowMassProblem,
203 std::shared_ptr<Problem<porousMediumIndex>> porousMediumProblem,
204 GridVarsTuple&& gridVarsTuple,
205 const SolutionVector& curSol)
206 {
207 this->updateSolution(curSol); // generic coupling manager stores tuple of shared_ptr
208
209 // initialize the binary sub coupling managers
211 std::get<0>(ffSolVecTuple) = std::get<freeFlowMomentumIndex>(this->curSol());
212 std::get<1>(ffSolVecTuple) = std::get<freeFlowMassIndex>(this->curSol());
214 freeFlowMomentumProblem, freeFlowMassProblem,
215 std::make_tuple(std::get<freeFlowMomentumIndex>(gridVarsTuple), std::get<freeFlowMassIndex>(gridVarsTuple)),
216 ffSolVecTuple
217 );
218
220 std::get<0>(ffMassPmSolVecTuple) = std::get<freeFlowMassIndex>(this->curSol());
221 std::get<1>(ffMassPmSolVecTuple) = std::get<porousMediumIndex>(this->curSol());
223 freeFlowMassProblem, porousMediumProblem, ffMassPmSolVecTuple
224 );
225
227 std::get<0>(ffMomentumPmSolVecTuple) = std::get<freeFlowMomentumIndex>(this->curSol());
228 std::get<1>(ffMomentumPmSolVecTuple) = std::get<porousMediumIndex>(this->curSol());
230 freeFlowMomentumProblem, porousMediumProblem, ffMomentumPmSolVecTuple
231 );
232 }
233
237 template<std::size_t i, std::size_t j>
238 auto massCouplingCondition(Dune::index_constant<i> domainI, Dune::index_constant<j> domainJ,
239 const FVElementGeometry<i>& fvGeometry,
240 const typename FVElementGeometry<i>::SubControlVolumeFace& scvf,
241 const ElementVolumeVariables<i>& elemVolVars) const
242 {
243 static_assert(domainI != freeFlowMomentumIndex && domainJ != freeFlowMomentumIndex);
244
245 const auto& couplingContext = this->subApply(domainI, domainJ, [&](const auto& cm, auto&& ii, auto&& jj) -> const auto& {
246 return cm.couplingContext(ii, fvGeometry, scvf);
247 });
248
249 const auto& freeFlowElement = [&]
250 {
251 if constexpr (domainI == freeFlowMassIndex)
252 return fvGeometry.element();
253 else
254 return couplingContext.fvGeometry.element();
255 }();
256
257 const auto& freeFlowScvf = [&]
258 {
259 if constexpr (domainI == freeFlowMassIndex)
260 return scvf;
261 else
262 return couplingContext.fvGeometry.scvf(couplingContext.freeFlowMassScvfIdx);
263
264 }();
265
266 // todo revise velocity (see ff mom pm mgr)
267
268 couplingContext.velocity = this->subCouplingManager(freeFlowMomentumIndex, freeFlowMassIndex).faceVelocity(freeFlowElement, freeFlowScvf);
269 return CouplingConditions::massCouplingCondition(domainI, domainJ, fvGeometry, scvf, elemVolVars, couplingContext);
270 }
271
272
275
276 NumEqVector<freeFlowMomentumIndex> momentumCouplingCondition(Dune::index_constant<freeFlowMomentumIndex> domainI,
277 Dune::index_constant<porousMediumIndex> domainJ,
278 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
279 const typename FVElementGeometry<freeFlowMomentumIndex>::SubControlVolumeFace& scvf,
280 const ElementVolumeVariables<freeFlowMomentumIndex>& elemVolVars) const
281 {
282 if (scvf.isLateral())
283 return NumEqVector<freeFlowMomentumIndex>(0.0);
284
285 const auto& context = this->subCouplingManager(freeFlowMomentumIndex, porousMediumIndex).couplingContext(
286 domainI, fvGeometry, scvf
287 );
288
289 return CouplingConditions::momentumCouplingCondition(fvGeometry, scvf, elemVolVars, context);
290 }
291
295 auto darcyPermeability(const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
296 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const
297 {
298 if (scvf.isFrontal())
299 {
300 const auto& context = this->subCouplingManager(freeFlowMomentumIndex, porousMediumIndex).couplingContext(
301 Dune::index_constant<freeFlowMomentumIndex>(), fvGeometry, scvf
302 );
303
304 return CouplingConditions::darcyPermeability(fvGeometry, scvf, context);
305 }
306 else
307 {
308 const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf);
309 const auto& orthogonalScv = fvGeometry.scv(orthogonalScvf.insideScvIdx());
310 const auto& frontalScvfOnBoundary = fvGeometry.frontalScvfOnBoundary(orthogonalScv);
311 const auto& context = this->subCouplingManager(freeFlowMomentumIndex, porousMediumIndex).couplingContext(
312 Dune::index_constant<freeFlowMomentumIndex>(), fvGeometry, frontalScvfOnBoundary
313 );
314
315 return CouplingConditions::darcyPermeability(fvGeometry, frontalScvfOnBoundary, context);
316 }
317 }
318
321
325 Scalar pressure(const Element<freeFlowMomentumIndex>& element,
326 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
327 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const
328 {
330 element, fvGeometry, scvf
331 );
332 }
333
340 Scalar cellPressure(const Element<freeFlowMomentumIndex>& element,
341 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const
342 {
344 element, scvf
345 );
346 }
347
351 Scalar density(const Element<freeFlowMomentumIndex>& element,
352 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
353 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
354 const bool considerPreviousTimeStep = false) const
355 {
357 element, fvGeometry, scvf, considerPreviousTimeStep
358 );
359 }
360
361 auto insideAndOutsideDensity(const Element<freeFlowMomentumIndex>& element,
362 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
363 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
364 const bool considerPreviousTimeStep = false) const
365 {
366 return this->subCouplingManager(freeFlowMomentumIndex, freeFlowMassIndex).insideAndOutsideDensity(
367 element, fvGeometry, scvf, considerPreviousTimeStep
368 );
369 }
370
374 Scalar density(const Element<freeFlowMomentumIndex>& element,
375 const SubControlVolume<freeFlowMomentumIndex>& scv,
376 const bool considerPreviousTimeStep = false) const
377 {
379 element, scv, considerPreviousTimeStep
380 );
381 }
382
386 Scalar effectiveViscosity(const Element<freeFlowMomentumIndex>& element,
387 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
388 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const
389 {
390 return this->subCouplingManager(freeFlowMomentumIndex, freeFlowMassIndex).effectiveViscosity(
391 element, fvGeometry, scvf
392 );
393 }
394
398 auto faceVelocity(const Element<freeFlowMassIndex>& element,
399 const SubControlVolumeFace<freeFlowMassIndex>& scvf) const
400 {
402 element, scvf
403 );
404 }
405
406 template<std::size_t i>
407 const Problem<i>& problem(Dune::index_constant<i> domainI) const
408 {
409 return this->subApply(domainI, [&](const auto& cm, auto&& ii) -> const auto& {
410 return cm.problem(ii);
411 });
412 }
413
414 template<std::size_t i, std::size_t j>
415 bool isCoupled(Dune::index_constant<i> domainI,
416 Dune::index_constant<j> domainJ,
417 const SubControlVolumeFace<i>& scvf) const
418 {
419 return this->subApply(domainI, domainJ, [&](const auto& cm, auto&& ii, auto&& jj){
420 return cm.isCoupled(ii, scvf);
421 });
422 }
423
430 template<std::size_t i, std::size_t j>
431 bool isCoupled(Dune::index_constant<i> domainI,
432 Dune::index_constant<j> domainJ,
433 const SubControlVolume<i>& scv) const
434 {
435 return this->subApply(domainI, domainJ, [&](const auto& cm, auto&& ii, auto&& jj){
436 return cm.isCoupled(ii, scv);
437 });
438 }
439
443 bool isCoupledLateralScvf(Dune::index_constant<freeFlowMomentumIndex> domainI,
444 Dune::index_constant<porousMediumIndex> domainJ,
445 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const
446 {
447 return this->subApply(domainI, domainJ, [&](const auto& cm, auto&& ii, auto&& jj){
448 return cm.isCoupledLateralScvf(ii, scvf);
449 });
450 }
451
452
463 template<std::size_t j>
464 const auto& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
465 const Element<freeFlowMomentumIndex>& elementI,
466 const SubControlVolume<freeFlowMomentumIndex>& scvI,
467 Dune::index_constant<j> domainJ) const
468 {
469 static_assert(freeFlowMomentumIndex != j);
470 return this->subApply(domainI, domainJ, [&](const auto& cm, auto&& ii, auto&& jj) -> const auto& {
471 return cm.couplingStencil(ii, elementI, scvI, jj);
472 });
473 }
474};
475
476} // end namespace Dumux
477
478#endif
typename Detail::FreeFlowCouplingManagerSelector< Traits >::type FreeFlowCouplingManager
The interface of the coupling manager for free flow systems.
Definition: multidomain/freeflow/couplingmanager.hh:67
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
static constexpr auto freeFlowMomentumIndex
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:45
static constexpr auto freeFlowMassToFreeFlowMomentumIndex
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:50
constexpr auto globalToLocalDomainIndices(Dune::index_constant< i >, Dune::index_constant< j >)
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:86
static constexpr auto freeFlowMassToPorousMediumIndex
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:52
static constexpr auto porousMediumIndex
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:47
static constexpr auto noCouplingIdx
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:53
constexpr auto makeCouplingManagerMap()
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:55
static constexpr auto freeFlowMomentumToPorousMediumIndex
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:51
static constexpr auto freeFlowMassIndex
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:46
constexpr auto coupledDomains(Dune::index_constant< i > domainI)
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:75
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:55
Definition: common/properties.hh:100
The type for a global container for the volume variables.
Definition: common/properties.hh:107
Definition: freeflowporousmedium/couplingconditions.hh:84
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:98
static constexpr auto globalToLocal(Dune::index_constant< i > domainI, Dune::index_constant< j > domainJ)
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:105
static constexpr auto coupledDomains(Dune::index_constant< i > domainI)
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:111
static constexpr auto managerMap()
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:99
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:119
typename MDTraits::template SubDomain< id >::TypeTag SubDomainTypeTag
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:121
Dumux::FreeFlowCouplingManager< FreeFlowTraits > FreeFlowCouplingManager
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:136
Coupling manager for coupling freeflow and porous medium flow models.
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:159
static constexpr auto freeFlowMassIndex
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:194
auto faceVelocity(const Element< freeFlowMassIndex > &element, const SubControlVolumeFace< freeFlowMassIndex > &scvf) const
Returns the velocity at a given sub control volume face.
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:398
bool isCoupledLateralScvf(Dune::index_constant< freeFlowMomentumIndex > domainI, Dune::index_constant< porousMediumIndex > domainJ, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf) const
Returns whether a given scvf is coupled to the other domain.
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:443
Scalar density(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const bool considerPreviousTimeStep=false) const
Returns the density at a given sub control volume face.
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:351
typename ParentType::template SubCouplingManager< i, j > SubCouplingManager
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:191
const auto & couplingStencil(Dune::index_constant< freeFlowMomentumIndex > domainI, const Element< freeFlowMomentumIndex > &elementI, const SubControlVolume< freeFlowMomentumIndex > &scvI, 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/boundary/freeflowporousmedium/couplingmanager.hh:464
auto darcyPermeability(const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf) const
Returns the intrinsic permeability of the coupled Darcy element.
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:295
Scalar density(const Element< freeFlowMomentumIndex > &element, const SubControlVolume< freeFlowMomentumIndex > &scv, const bool considerPreviousTimeStep=false) const
Returns the density at a given sub control volume.
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:374
static constexpr auto porousMediumIndex
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:195
void init(std::shared_ptr< Problem< freeFlowMomentumIndex > > freeFlowMomentumProblem, std::shared_ptr< Problem< freeFlowMassIndex > > freeFlowMassProblem, std::shared_ptr< Problem< porousMediumIndex > > porousMediumProblem, GridVarsTuple &&gridVarsTuple, const SolutionVector &curSol)
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:201
static constexpr auto freeFlowMomentumIndex
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:193
NumEqVector< freeFlowMomentumIndex > momentumCouplingCondition(Dune::index_constant< freeFlowMomentumIndex > domainI, Dune::index_constant< porousMediumIndex > domainJ, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const typename FVElementGeometry< freeFlowMomentumIndex >::SubControlVolumeFace &scvf, const ElementVolumeVariables< freeFlowMomentumIndex > &elemVolVars) const
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:276
auto massCouplingCondition(Dune::index_constant< i > domainI, Dune::index_constant< j > domainJ, const FVElementGeometry< i > &fvGeometry, const typename FVElementGeometry< i >::SubControlVolumeFace &scvf, const ElementVolumeVariables< i > &elemVolVars) const
Returns the mass flux across the coupling boundary.
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:238
Scalar pressure(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf) const
Returns the pressure at a given sub control volume face.
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:325
bool isCoupled(Dune::index_constant< i > domainI, Dune::index_constant< j > domainJ, const SubControlVolumeFace< i > &scvf) const
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:415
Scalar cellPressure(const Element< freeFlowMomentumIndex > &element, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf) const
Returns the pressure at the center of a sub control volume corresponding to a given sub control volum...
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:340
bool isCoupled(Dune::index_constant< i > domainI, Dune::index_constant< j > domainJ, const SubControlVolume< i > &scv) const
If the boundary entity is on a coupling boundary.
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:431
Scalar effectiveViscosity(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf) const
Returns the pressure at a given sub control volume face.
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:386
const Problem< i > & problem(Dune::index_constant< i > domainI) const
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:407
auto insideAndOutsideDensity(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const bool considerPreviousTimeStep=false) const
Definition: multidomain/boundary/freeflowporousmedium/couplingmanager.hh:361
Coupling manager for Stokes and Darcy domains with equal dimension.
Definition: multidomain/boundary/freeflowporousmedium/ffmasspm/couplingmanager.hh:47
Coupling manager for Stokes and Darcy domains with equal dimension.
Definition: multidomain/boundary/freeflowporousmedium/ffmomentumpm/couplingmanager.hh:48
Coupling manager that combines an arbitrary number of binary coupling manager (coupling two domains e...
Definition: multibinarycouplingmanager.hh:62
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
Definition: multidomain/traits.hh:141
Declares all properties used in Dumux.
Freeflow coupling managers (Navier-Stokes mass-momentum coupling)