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-FileCopyrightInfo: 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 FreeFlowPorousMediumCoupling | ||
10 | * \copydoc Dumux::FreeFlowPorousMediumCouplingManagerStaggeredCCTpfa | ||
11 | */ | ||
12 | |||
13 | #ifndef DUMUX_MULTIDOMAIN_BOUNDARY_FREEFLOW_POROUSMEDIUM_COUPLINGMANAGER_STAGGERED_CCTPFA_HH | ||
14 | #define DUMUX_MULTIDOMAIN_BOUNDARY_FREEFLOW_POROUSMEDIUM_COUPLINGMANAGER_STAGGERED_CCTPFA_HH | ||
15 | |||
16 | #include "couplingmanager_base.hh" | ||
17 | #include "couplingconditions_staggered_cctpfa.hh" | ||
18 | |||
19 | namespace Dumux { | ||
20 | |||
21 | /*! | ||
22 | * \ingroup FreeFlowPorousMediumCoupling | ||
23 | * \brief Coupling manager for coupling freeflow and porous medium flow models | ||
24 | * specialization for staggered-cctpfa coupling. | ||
25 | */ | ||
26 | template<class MDTraits> | ||
27 |
1/4✗ Branch 4 not taken.
✗ Branch 5 not taken.
✓ Branch 8 taken 14 times.
✗ Branch 9 not taken.
|
28 | class FreeFlowPorousMediumCouplingManagerStaggeredCCTpfa |
28 | : public FreeFlowPorousMediumCouplingManagerBase<MDTraits> | ||
29 | { | ||
30 | using ParentType = FreeFlowPorousMediumCouplingManagerBase<MDTraits>; | ||
31 | |||
32 | using Scalar = typename MDTraits::Scalar; | ||
33 | |||
34 | // the sub domain type tags | ||
35 | template<std::size_t id> | ||
36 | using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag; | ||
37 | |||
38 | template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>; | ||
39 | template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>; | ||
40 | template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView; | ||
41 | template<std::size_t id> using SubControlVolumeFace = typename FVElementGeometry<id>::SubControlVolumeFace; | ||
42 | template<std::size_t id> using SubControlVolume = typename FVElementGeometry<id>::SubControlVolume; | ||
43 | template<std::size_t id> using ElementVolumeVariables = typename GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>::LocalView; | ||
44 | template<std::size_t id> using NumEqVector = typename Problem<id>::Traits::NumEqVector; | ||
45 | |||
46 | template<std::size_t id> using GridView = typename GridGeometry<id>::GridView; | ||
47 | template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity; | ||
48 | using SolutionVector = typename MDTraits::SolutionVector; | ||
49 | |||
50 | using CouplingConditions = FFPMCouplingConditionsStaggeredCCTpfa<MDTraits, FreeFlowPorousMediumCouplingManagerStaggeredCCTpfa<MDTraits>>; | ||
51 | |||
52 | public: | ||
53 | static constexpr auto freeFlowMomentumIndex = ParentType::freeFlowMomentumIndex; | ||
54 | static constexpr auto freeFlowMassIndex = ParentType::freeFlowMassIndex; | ||
55 | static constexpr auto porousMediumIndex = ParentType::porousMediumIndex; | ||
56 | |||
57 | public: | ||
58 | /*! | ||
59 | * \brief Returns the mass flux across the coupling boundary. | ||
60 | */ | ||
61 | template<std::size_t i, std::size_t j> | ||
62 | 118320 | auto massCouplingCondition(Dune::index_constant<i> domainI, Dune::index_constant<j> domainJ, | |
63 | const FVElementGeometry<i>& fvGeometry, | ||
64 | const typename FVElementGeometry<i>::SubControlVolumeFace& scvf, | ||
65 | const ElementVolumeVariables<i>& elemVolVars) const | ||
66 | { | ||
67 | static_assert(domainI != freeFlowMomentumIndex && domainJ != freeFlowMomentumIndex); | ||
68 | |||
69 | 118320 | const auto& couplingContext = this->subApply(domainI, domainJ, [&](const auto& cm, auto&& ii, auto&& jj) -> const auto& { | |
70 |
0/4✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
|
118320 | return cm.couplingContext(ii, fvGeometry, scvf); |
71 | }); | ||
72 | |||
73 | 236640 | const auto& freeFlowElement = [&] | |
74 | { | ||
75 | if constexpr (domainI == freeFlowMassIndex) | ||
76 | 151680 | return fvGeometry.element(); | |
77 | else | ||
78 | 84960 | return couplingContext.fvGeometry.element(); | |
79 | }(); | ||
80 | |||
81 | 118320 | const auto& freeFlowScvf = [&] | |
82 | { | ||
83 | if constexpr (domainI == freeFlowMassIndex) | ||
84 | 75840 | return scvf; | |
85 | else | ||
86 | 84960 | return couplingContext.fvGeometry.scvf(couplingContext.freeFlowMassScvfIdx); | |
87 | |||
88 | 118320 | }(); | |
89 | |||
90 | // todo revise velocity (see ff mom pm mgr) | ||
91 | |||
92 | 236640 | couplingContext.velocity = this->subCouplingManager(freeFlowMomentumIndex, freeFlowMassIndex).faceVelocity(freeFlowElement, freeFlowScvf); | |
93 | 118320 | return CouplingConditions::massCouplingCondition(domainI, domainJ, fvGeometry, scvf, elemVolVars, couplingContext); | |
94 | } | ||
95 | |||
96 | |||
97 | //////////////////////// Conditions for FreeFlowMomentum - PorousMedium coupling ////////// | ||
98 | /////////////////////////////////////////////////////////////////////////////////////////// | ||
99 | |||
100 | 178310 | NumEqVector<freeFlowMomentumIndex> momentumCouplingCondition(Dune::index_constant<freeFlowMomentumIndex> domainI, | |
101 | Dune::index_constant<porousMediumIndex> domainJ, | ||
102 | const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry, | ||
103 | const typename FVElementGeometry<freeFlowMomentumIndex>::SubControlVolumeFace& scvf, | ||
104 | const ElementVolumeVariables<freeFlowMomentumIndex>& elemVolVars) const | ||
105 | { | ||
106 |
2/2✓ Branch 0 taken 143624 times.
✓ Branch 1 taken 34686 times.
|
178310 | if (scvf.isLateral()) |
107 | 143624 | return NumEqVector<freeFlowMomentumIndex>(0.0); | |
108 | |||
109 | 69372 | const auto& context = this->subCouplingManager(freeFlowMomentumIndex, porousMediumIndex).couplingContext( | |
110 | domainI, fvGeometry, scvf | ||
111 | ); | ||
112 | |||
113 |
1/2✓ Branch 0 taken 34686 times.
✗ Branch 1 not taken.
|
34686 | return CouplingConditions::momentumCouplingCondition(fvGeometry, scvf, elemVolVars, context); |
114 | } | ||
115 | |||
116 | /*! | ||
117 | * \brief Returns the intrinsic permeability of the coupled Darcy element. | ||
118 | */ | ||
119 | 143704 | auto darcyPermeability(const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry, | |
120 | const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const | ||
121 | { | ||
122 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 143704 times.
|
143704 | if (scvf.isFrontal()) |
123 | { | ||
124 | ✗ | const auto& context = this->subCouplingManager(freeFlowMomentumIndex, porousMediumIndex).couplingContext( | |
125 | Dune::index_constant<freeFlowMomentumIndex>(), fvGeometry, scvf | ||
126 | ); | ||
127 | |||
128 | ✗ | return CouplingConditions::darcyPermeability(fvGeometry, scvf, context); | |
129 | } | ||
130 | else | ||
131 | { | ||
132 | 143704 | const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf); | |
133 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 143704 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 143704 times.
|
287408 | const auto& orthogonalScv = fvGeometry.scv(orthogonalScvf.insideScvIdx()); |
134 | 143704 | const auto& frontalScvfOnBoundary = fvGeometry.frontalScvfOnBoundary(orthogonalScv); | |
135 | 287408 | const auto& context = this->subCouplingManager(freeFlowMomentumIndex, porousMediumIndex).couplingContext( | |
136 | Dune::index_constant<freeFlowMomentumIndex>(), fvGeometry, frontalScvfOnBoundary | ||
137 | ); | ||
138 | |||
139 | 147312 | return CouplingConditions::darcyPermeability(fvGeometry, frontalScvfOnBoundary, context); | |
140 | } | ||
141 | } | ||
142 | |||
143 | //////////////////////// Conditions for FreeFlowMomentum - FreeFlowMass coupling ////////// | ||
144 | /////////////////////////////////////////////////////////////////////////////////////////// | ||
145 | |||
146 | /*! | ||
147 | * \brief Returns the pressure at a given sub control volume face. | ||
148 | */ | ||
149 | ✗ | Scalar pressure(const Element<freeFlowMomentumIndex>& element, | |
150 | const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry, | ||
151 | const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const | ||
152 | { | ||
153 | 16255742 | return this->subCouplingManager(freeFlowMomentumIndex, freeFlowMassIndex).pressure( | |
154 | element, fvGeometry, scvf | ||
155 | 16255742 | ); | |
156 | } | ||
157 | |||
158 | /*! | ||
159 | * \brief Returns the pressure at the center of a sub control volume corresponding to a given sub control volume face. | ||
160 | * This is used for setting a Dirichlet pressure for the mass model when a fixed pressure for the momentum balance is set at another | ||
161 | * boundary. Since the the pressure at the given scvf is solution-dependent and thus unknown a priori, we just use the value | ||
162 | * of the interior cell here. | ||
163 | */ | ||
164 | Scalar cellPressure(const Element<freeFlowMomentumIndex>& element, | ||
165 | const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const | ||
166 | { | ||
167 | return this->subCouplingManager(freeFlowMomentumIndex, freeFlowMassIndex).cellPressure( | ||
168 | element, scvf | ||
169 | ); | ||
170 | } | ||
171 | |||
172 | /*! | ||
173 | * \brief Returns the density at a given sub control volume face. | ||
174 | */ | ||
175 | Scalar density(const Element<freeFlowMomentumIndex>& element, | ||
176 | const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry, | ||
177 | const SubControlVolumeFace<freeFlowMomentumIndex>& scvf, | ||
178 | const bool considerPreviousTimeStep = false) const | ||
179 | { | ||
180 | 7355630 | return this->subCouplingManager(freeFlowMomentumIndex, freeFlowMassIndex).density( | |
181 | element, fvGeometry, scvf, considerPreviousTimeStep | ||
182 | 7355630 | ); | |
183 | } | ||
184 | |||
185 | auto insideAndOutsideDensity(const Element<freeFlowMomentumIndex>& element, | ||
186 | const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry, | ||
187 | const SubControlVolumeFace<freeFlowMomentumIndex>& scvf, | ||
188 | const bool considerPreviousTimeStep = false) const | ||
189 | { | ||
190 | 14602600 | return this->subCouplingManager(freeFlowMomentumIndex, freeFlowMassIndex).insideAndOutsideDensity( | |
191 | element, fvGeometry, scvf, considerPreviousTimeStep | ||
192 | 14602600 | ); | |
193 | } | ||
194 | |||
195 | /*! | ||
196 | * \brief Returns the density at a given sub control volume. | ||
197 | */ | ||
198 | Scalar density(const Element<freeFlowMomentumIndex>& element, | ||
199 | const SubControlVolume<freeFlowMomentumIndex>& scv, | ||
200 | const bool considerPreviousTimeStep = false) const | ||
201 | { | ||
202 |
0/2✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
125190 | return this->subCouplingManager(freeFlowMomentumIndex, freeFlowMassIndex).density( |
203 | element, scv, considerPreviousTimeStep | ||
204 |
0/2✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
125190 | ); |
205 | } | ||
206 | |||
207 | /*! | ||
208 | * \brief Returns the pressure at a given sub control volume face. | ||
209 | */ | ||
210 | Scalar effectiveViscosity(const Element<freeFlowMomentumIndex>& element, | ||
211 | const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry, | ||
212 | const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const | ||
213 | { | ||
214 | 48628984 | return this->subCouplingManager(freeFlowMomentumIndex, freeFlowMassIndex).effectiveViscosity( | |
215 | element, fvGeometry, scvf | ||
216 | 48628984 | ); | |
217 | } | ||
218 | |||
219 | /*! | ||
220 | * \brief Returns the velocity at a given sub control volume face. | ||
221 | */ | ||
222 | auto faceVelocity(const Element<freeFlowMassIndex>& element, | ||
223 | const SubControlVolumeFace<freeFlowMassIndex>& scvf) const | ||
224 | { | ||
225 | 18858080 | return this->subCouplingManager(freeFlowMomentumIndex, freeFlowMassIndex).faceVelocity( | |
226 | element, scvf | ||
227 | 18858080 | ); | |
228 | } | ||
229 | |||
230 | /*! | ||
231 | * \brief Returns whether a given scvf is coupled to the other domain | ||
232 | */ | ||
233 | bool isCoupledLateralScvf(Dune::index_constant<freeFlowMomentumIndex> domainI, | ||
234 | Dune::index_constant<porousMediumIndex> domainJ, | ||
235 | const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const | ||
236 | { | ||
237 | return this->subApply(domainI, domainJ, [&](const auto& cm, auto&& ii, auto&& jj){ | ||
238 | return cm.isCoupledLateralScvf(ii, scvf); | ||
239 | }); | ||
240 | } | ||
241 | }; | ||
242 | |||
243 | } // end namespace Dumux | ||
244 | |||
245 | #endif | ||
246 |