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 FacetCoupling | ||
10 | * \copydoc Dumux::FacetCouplingManager | ||
11 | */ | ||
12 | #ifndef DUMUX_CCTPFA_FACETCOUPLING_MANAGER_HH | ||
13 | #define DUMUX_CCTPFA_FACETCOUPLING_MANAGER_HH | ||
14 | |||
15 | #include <algorithm> | ||
16 | #include <cassert> | ||
17 | |||
18 | #include <dumux/common/properties.hh> | ||
19 | #include <dumux/common/indextraits.hh> | ||
20 | #include <dumux/common/numeqvector.hh> | ||
21 | #include <dumux/discretization/method.hh> | ||
22 | #include <dumux/discretization/elementsolution.hh> | ||
23 | #include <dumux/multidomain/couplingmanager.hh> | ||
24 | #include <dumux/multidomain/facet/couplingmanager.hh> | ||
25 | |||
26 | namespace Dumux { | ||
27 | |||
28 | /*! | ||
29 | * \ingroup FacetCoupling | ||
30 | * \brief Manages the coupling between bulk elements and lower dimensional elements | ||
31 | * where the coupling occurs across the facets of the bulk grid. This implementation | ||
32 | * is to be used in conjunction with models using the cell-centered tpfa scheme. | ||
33 | * | ||
34 | * \tparam MDTraits The multidomain traits containing the types on all sub-domains | ||
35 | * \tparam CouplingMapper Class containing maps on the coupling between dofs of different grids | ||
36 | * \tparam bulkDomainId The domain id of the bulk problem | ||
37 | * \tparam lowDimDomainId The domain id of the lower-dimensional problem | ||
38 | */ | ||
39 | template<class MDTraits, class CouplingMapper, std::size_t bulkDomainId, std::size_t lowDimDomainId> | ||
40 | class FacetCouplingManager<MDTraits, CouplingMapper, bulkDomainId, lowDimDomainId, DiscretizationMethods::CCTpfa> | ||
41 | : public virtual CouplingManager< MDTraits > | ||
42 | { | ||
43 | using ParentType = CouplingManager< MDTraits >; | ||
44 | |||
45 | // convenience aliases and instances of the two domain ids | ||
46 | using BulkIdType = typename MDTraits::template SubDomain<bulkDomainId>::Index; | ||
47 | using LowDimIdType = typename MDTraits::template SubDomain<lowDimDomainId>::Index; | ||
48 | static constexpr auto bulkId = BulkIdType(); | ||
49 | static constexpr auto lowDimId = LowDimIdType(); | ||
50 | |||
51 | // the sub-domain type tags | ||
52 | template<std::size_t id> using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag; | ||
53 | |||
54 | // further types specific to the sub-problems | ||
55 | template<std::size_t id> using PrimaryVariables = GetPropType<SubDomainTypeTag<id>, Properties::PrimaryVariables>; | ||
56 | template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>; | ||
57 | template<std::size_t id> using NumEqVector = Dumux::NumEqVector<PrimaryVariables<id>>; | ||
58 | template<std::size_t id> using LocalResidual = typename MDTraits::template SubDomain<id>::LocalResidual; | ||
59 | |||
60 | template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>; | ||
61 | template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView; | ||
62 | template<std::size_t id> using SubControlVolume = typename GridGeometry<id>::SubControlVolume; | ||
63 | template<std::size_t id> using SubControlVolumeFace = typename GridGeometry<id>::SubControlVolumeFace; | ||
64 | template<std::size_t id> using GridView = typename GridGeometry<id>::GridView; | ||
65 | template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity; | ||
66 | template<std::size_t id> using GridIndexType = typename IndexTraits< GridView<id> >::GridIndex; | ||
67 | |||
68 | template<std::size_t id> using GridVariables = GetPropType<SubDomainTypeTag<id>, Properties::GridVariables>; | ||
69 | template<std::size_t id> using GridVolumeVariables = typename GridVariables<id>::GridVolumeVariables; | ||
70 | template<std::size_t id> using ElementVolumeVariables = typename GridVolumeVariables<id>::LocalView; | ||
71 | template<std::size_t id> using VolumeVariables = typename ElementVolumeVariables<id>::VolumeVariables; | ||
72 | template<std::size_t id> using GridFluxVariablesCache = typename GridVariables<id>::GridFluxVariablesCache; | ||
73 | template<std::size_t id> using ElementFluxVariablesCache = typename GridFluxVariablesCache<id>::LocalView; | ||
74 | |||
75 | // this currently does not work for some grid-wide caches being active | ||
76 | static_assert(!GetPropType<SubDomainTypeTag<bulkId>, Properties::GridVariables>::GridFluxVariablesCache::cachingEnabled, | ||
77 | "Grid flux variables caching currently not supported in the bulk domain of cc-facet coupling models"); | ||
78 | static_assert(!GetPropType<SubDomainTypeTag<lowDimId>, Properties::GridVariables>::GridVolumeVariables::cachingEnabled, | ||
79 | "Grid volume variables caching currently not supported in the lower-dimensional domain of cc-facet coupling models"); | ||
80 | static_assert(!GetPropType<SubDomainTypeTag<bulkId>, Properties::GridVariables>::GridVolumeVariables::cachingEnabled, | ||
81 | "Grid volume variables caching currently not supported in the bulk domain of cc-facet coupling models"); | ||
82 | |||
83 | // extract corresponding grid ids from the mapper | ||
84 | static constexpr int bulkDim = GridView<bulkDomainId>::dimension; | ||
85 | static constexpr int lowDimDim = GridView<lowDimDomainId>::dimension; | ||
86 | static constexpr auto bulkGridId = CouplingMapper::template gridId<bulkDim>(); | ||
87 | static constexpr auto lowDimGridId = CouplingMapper::template gridId<lowDimDim>(); | ||
88 | |||
89 | static constexpr bool lowDimUsesBox = GridGeometry<lowDimId>::discMethod == DiscretizationMethods::box; | ||
90 | |||
91 | /*! | ||
92 | * \brief The coupling context of the bulk domain. Contains all data of the lower- | ||
93 | * dimensional domain which is required for the computation of a bulk element | ||
94 | * residual. This boils down to the geometries and volume variables of all | ||
95 | * lower-dimensional elements connected to a given bulk element. | ||
96 | */ | ||
97 | struct BulkCouplingContext | ||
98 | { | ||
99 | bool isSet; | ||
100 | GridIndexType< bulkId > elementIdx; | ||
101 | std::vector< FVElementGeometry<lowDimId> > lowDimFvGeometries; | ||
102 | std::vector< VolumeVariables<lowDimId> > lowDimVolVars; | ||
103 | |||
104 | 4783626 | void reset() | |
105 | { | ||
106 |
2/2✓ Branch 0 taken 432551 times.
✓ Branch 1 taken 2199454 times.
|
4783626 | lowDimFvGeometries.clear(); |
107 |
2/2✓ Branch 0 taken 432551 times.
✓ Branch 1 taken 2199454 times.
|
4783626 | lowDimVolVars.clear(); |
108 | 4783626 | isSet = false; | |
109 | 4783626 | } | |
110 | }; | ||
111 | |||
112 | /*! | ||
113 | * \brief The coupling context of the lower-dimensional (codim 1) domain. Contains | ||
114 | * all data of the bulk domain which is required for computation of element | ||
115 | * residuals in the lower-dimensional domain. This is essentially everything | ||
116 | * that is necessary to compute the fluxes in bulk domain entering a given | ||
117 | * lower-dimensional element. Thus, we store and bind the local views of one | ||
118 | * of the neighboring elements, which will be enough to compute the fluxes | ||
119 | * stemming from all neighboring elements. | ||
120 | * | ||
121 | * \note We need unique ptrs here because the local views have no default constructor. | ||
122 | */ | ||
123 | struct LowDimCouplingContext | ||
124 | { | ||
125 | bool isSet; | ||
126 | GridIndexType< lowDimId > elementIdx; | ||
127 | std::unique_ptr< FVElementGeometry<bulkId> > bulkFvGeometry; | ||
128 | std::unique_ptr< ElementVolumeVariables<bulkId> > bulkElemVolVars; | ||
129 | std::unique_ptr< ElementFluxVariablesCache<bulkId> > bulkElemFluxVarsCache; | ||
130 | std::unique_ptr< LocalResidual<bulkId> > bulkLocalResidual; | ||
131 | |||
132 | 209602 | void reset() | |
133 | { | ||
134 |
2/2✓ Branch 0 taken 108488 times.
✓ Branch 1 taken 38 times.
|
209602 | bulkFvGeometry.reset(nullptr); |
135 |
2/2✓ Branch 0 taken 108488 times.
✓ Branch 1 taken 38 times.
|
209602 | bulkElemVolVars.reset(nullptr); |
136 |
2/2✓ Branch 0 taken 108488 times.
✓ Branch 1 taken 38 times.
|
209602 | bulkElemFluxVarsCache.reset(nullptr); |
137 | 209602 | isSet = false; | |
138 | 209602 | } | |
139 | }; | ||
140 | |||
141 | public: | ||
142 | |||
143 | //! types used for coupling stencils | ||
144 | template<std::size_t i, std::size_t j = (i == bulkId) ? lowDimId : bulkId> | ||
145 | using CouplingStencilType = typename CouplingMapper::template Stencil< CouplingMapper::template gridId<GridView<j>::dimension>() >; | ||
146 | |||
147 | //! the type of the solution vector | ||
148 | using SolutionVector = typename MDTraits::SolutionVector; | ||
149 | |||
150 | /*! | ||
151 | * \brief Initialize the coupling manager. | ||
152 | * | ||
153 | * \param bulkProblem The problem to be solved on the bulk domain | ||
154 | * \param lowDimProblem The problem to be solved on the lower-dimensional domain | ||
155 | * \param couplingMapper The mapper object containing the connectivity between the domains | ||
156 | * \param curSol The current solution | ||
157 | */ | ||
158 | 44 | void init(std::shared_ptr< Problem<bulkId> > bulkProblem, | |
159 | std::shared_ptr< Problem<lowDimId> > lowDimProblem, | ||
160 | std::shared_ptr< CouplingMapper > couplingMapper, | ||
161 | const SolutionVector& curSol) | ||
162 | { | ||
163 | 44 | couplingMapperPtr_ = couplingMapper; | |
164 | |||
165 | // set the sub problems | ||
166 |
2/4✓ Branch 1 taken 38 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 38 times.
✗ Branch 4 not taken.
|
44 | this->setSubProblem(bulkProblem, bulkId); |
167 |
2/4✓ Branch 1 taken 38 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 38 times.
✗ Branch 4 not taken.
|
44 | this->setSubProblem(lowDimProblem, lowDimId); |
168 | |||
169 | // copy the solution vector | ||
170 | 44 | ParentType::updateSolution(curSol); | |
171 | |||
172 | // determine all bulk elements/scvfs that couple to low dim elements | ||
173 | 176 | bulkElemIsCoupled_.assign(bulkProblem->gridGeometry().gridView().size(0), false); | |
174 | 147 | bulkScvfIsCoupled_.assign(bulkProblem->gridGeometry().numScvf(), false); | |
175 | |||
176 | 88 | const auto& bulkMap = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId); | |
177 |
2/2✓ Branch 0 taken 3615 times.
✓ Branch 1 taken 38 times.
|
5145 | for (const auto& entry : bulkMap) |
178 | { | ||
179 | 10026 | bulkElemIsCoupled_[entry.first] = true; | |
180 |
2/2✓ Branch 0 taken 5934 times.
✓ Branch 1 taken 3615 times.
|
23633 | for (const auto& couplingEntry : entry.second.dofToCouplingScvfMap) |
181 |
4/4✓ Branch 0 taken 11516 times.
✓ Branch 1 taken 5934 times.
✓ Branch 2 taken 11516 times.
✓ Branch 3 taken 5934 times.
|
59870 | for (const auto& scvfIdx : couplingEntry.second) |
182 | 51132 | bulkScvfIsCoupled_[scvfIdx] = true; | |
183 | } | ||
184 | 44 | } | |
185 | |||
186 | /*! | ||
187 | * \brief The coupling stencil of a given bulk domain element. | ||
188 | */ | ||
189 | 586164 | const CouplingStencilType<bulkId>& couplingStencil(BulkIdType domainI, | |
190 | const Element<bulkId>& element, | ||
191 | LowDimIdType domainJ) const | ||
192 | { | ||
193 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 578104 times.
|
586164 | const auto eIdx = this->problem(bulkId).gridGeometry().elementMapper().index(element); |
194 | |||
195 |
4/4✓ Branch 0 taken 20956 times.
✓ Branch 1 taken 557148 times.
✓ Branch 2 taken 20956 times.
✓ Branch 3 taken 557148 times.
|
1172328 | if (bulkElemIsCoupled_[eIdx]) |
196 | { | ||
197 | 43320 | const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId); | |
198 | 21660 | auto it = map.find(eIdx); | |
199 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 20956 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 20956 times.
|
43320 | assert(it != map.end()); |
200 | 43320 | return it->second.couplingStencil; | |
201 | } | ||
202 | |||
203 | 1129008 | return getEmptyStencil(lowDimId); | |
204 | } | ||
205 | |||
206 | /*! | ||
207 | * \brief The coupling stencil of the lower-dimensional domain with the bulk domain. | ||
208 | */ | ||
209 | 10096 | const CouplingStencilType<lowDimId>& couplingStencil(LowDimIdType domainI, | |
210 | const Element<lowDimId>& element, | ||
211 | BulkIdType domainJ) const | ||
212 | { | ||
213 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 9728 times.
|
10096 | const auto eIdx = this->problem(lowDimId).gridGeometry().elementMapper().index(element); |
214 | |||
215 | 20192 | const auto& map = couplingMapperPtr_->couplingMap(lowDimGridId, bulkGridId); | |
216 | 10096 | auto it = map.find(eIdx); | |
217 |
2/4✓ Branch 0 taken 9728 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 9728 times.
✗ Branch 3 not taken.
|
20192 | if (it != map.end()) return it->second.couplingStencil; |
218 | ✗ | else return getEmptyStencil(bulkId); | |
219 | } | ||
220 | |||
221 | /*! | ||
222 | * \brief returns true if a bulk scvf flux depends on data in the facet domain. | ||
223 | */ | ||
224 | ✗ | bool isCoupled(const Element<bulkId>& element, | |
225 | const SubControlVolumeFace<bulkId>& scvf) const | ||
226 |
12/12✓ Branch 0 taken 8577748 times.
✓ Branch 1 taken 237934 times.
✓ Branch 2 taken 8577748 times.
✓ Branch 3 taken 237934 times.
✓ Branch 4 taken 5307252 times.
✓ Branch 5 taken 179676 times.
✓ Branch 6 taken 5307252 times.
✓ Branch 7 taken 179676 times.
✓ Branch 8 taken 3301344 times.
✓ Branch 9 taken 111888 times.
✓ Branch 10 taken 3301344 times.
✓ Branch 11 taken 111888 times.
|
35431684 | { return bulkScvfIsCoupled_[scvf.index()]; } |
227 | |||
228 | /*! | ||
229 | * \brief returns true if a bulk scvf coincides with a facet element. | ||
230 | * \note for tpfa, this is always true for coupled scvfs | ||
231 | */ | ||
232 | ✗ | bool isOnInteriorBoundary(const Element<bulkId>& element, | |
233 | const SubControlVolumeFace<bulkId>& scvf) const | ||
234 |
24/24✓ Branch 0 taken 2738696 times.
✓ Branch 1 taken 7492972 times.
✓ Branch 2 taken 2738696 times.
✓ Branch 3 taken 7492972 times.
✓ Branch 4 taken 3480372 times.
✓ Branch 5 taken 2567178 times.
✓ Branch 6 taken 3480372 times.
✓ Branch 7 taken 2567178 times.
✓ Branch 8 taken 1652768 times.
✓ Branch 9 taken 9900860 times.
✓ Branch 10 taken 1652768 times.
✓ Branch 11 taken 9900860 times.
✓ Branch 12 taken 2550600 times.
✓ Branch 13 taken 1419012 times.
✓ Branch 14 taken 2550600 times.
✓ Branch 15 taken 1419012 times.
✓ Branch 16 taken 2423816 times.
✓ Branch 17 taken 80891 times.
✓ Branch 18 taken 2423816 times.
✓ Branch 19 taken 80891 times.
✓ Branch 20 taken 150725 times.
✓ Branch 21 taken 2538440 times.
✓ Branch 22 taken 150725 times.
✓ Branch 23 taken 2538440 times.
|
73992660 | { return isCoupled(element, scvf); } |
235 | |||
236 | /*! | ||
237 | * \brief returns the vol vars of a lower-dimensional element coinciding with a bulk scvf. | ||
238 | */ | ||
239 | 11123284 | const VolumeVariables<lowDimId>& getLowDimVolVars(const Element<bulkId>& element, | |
240 | const SubControlVolumeFace<bulkId>& scvf) const | ||
241 | { | ||
242 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10207266 times.
|
11123284 | assert(bulkContext_.isSet); |
243 | |||
244 | 11123284 | const auto lowDimElemIdx = getLowDimElementIndex(element, scvf); | |
245 | 22246568 | const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId); | |
246 | 22246568 | const auto& s = map.find(bulkContext_.elementIdx)->second.couplingElementStencil; | |
247 | 33369852 | const auto& idxInContext = std::distance( s.begin(), std::find(s.begin(), s.end(), lowDimElemIdx) ); | |
248 | |||
249 |
2/4✗ Branch 4 not taken.
✓ Branch 5 taken 10207266 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 10207266 times.
|
44493136 | assert(std::find(s.begin(), s.end(), lowDimElemIdx) != s.end()); |
250 | 22246568 | return bulkContext_.lowDimVolVars[idxInContext]; | |
251 | } | ||
252 | |||
253 | /*! | ||
254 | * \brief returns the lower-dimensional element coinciding with a bulk scvf. | ||
255 | */ | ||
256 | const Element<lowDimId> getLowDimElement(const Element<bulkId>& element, | ||
257 | const SubControlVolumeFace<bulkId>& scvf) const | ||
258 | { | ||
259 | const auto lowDimElemIdx = getLowDimElementIndex(element, scvf); | ||
260 | return this->problem(lowDimId).gridGeometry().element(lowDimElemIdx); | ||
261 | } | ||
262 | |||
263 | /*! | ||
264 | * \brief returns the index of the lower-dimensional element coinciding with a bulk scvf. | ||
265 | */ | ||
266 | 11123284 | const GridIndexType<lowDimId> getLowDimElementIndex(const Element<bulkId>& element, | |
267 | const SubControlVolumeFace<bulkId>& scvf) const | ||
268 | { | ||
269 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 10207266 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10207266 times.
|
22246568 | assert(bulkScvfIsCoupled_[scvf.index()]); |
270 | |||
271 | 22246568 | const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId); | |
272 | 12801166 | const auto& couplingData = map.at(scvf.insideScvIdx()); | |
273 | |||
274 | // search the low dim element idx this scvf is embedded in | ||
275 | 34303260 | auto it = std::find_if( couplingData.elementToScvfMap.begin(), | |
276 | couplingData.elementToScvfMap.end(), | ||
277 | 19515904 | [&scvf] (auto& dataPair) | |
278 | { | ||
279 | 10224656 | const auto& scvfs = dataPair.second; | |
280 | 40898624 | return std::find(scvfs.begin(), scvfs.end(), scvf.index()) != scvfs.end(); | |
281 | } ); | ||
282 | |||
283 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 10207266 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10207266 times.
|
22246568 | assert(it != couplingData.elementToScvfMap.end()); |
284 | 22246568 | return it->first; | |
285 | } | ||
286 | |||
287 | /*! | ||
288 | * \brief Evaluates the coupling element residual of a bulk domain element with respect | ||
289 | * to a dof in the lower-dimensional domain (dofIdxGlobalJ). This is essentially | ||
290 | * the fluxes across the bulk element facets that depend on dofIdxGlobalJ. | ||
291 | */ | ||
292 | template< class BulkLocalAssembler > | ||
293 | typename LocalResidual<bulkId>::ElementResidualVector | ||
294 | 98024 | evalCouplingResidual(BulkIdType, | |
295 | const BulkLocalAssembler& bulkLocalAssembler, | ||
296 | LowDimIdType, | ||
297 | GridIndexType<lowDimId> dofIdxGlobalJ) | ||
298 | { | ||
299 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 96872 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 96872 times.
|
196048 | const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId); |
300 | |||
301 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 96872 times.
|
98024 | assert(bulkContext_.isSet); |
302 |
3/6✗ Branch 0 not taken.
✓ Branch 1 taken 96872 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 96872 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 96872 times.
|
294072 | assert(bulkElemIsCoupled_[bulkContext_.elementIdx]); |
303 |
1/2✗ Branch 2 not taken.
✓ Branch 3 taken 96872 times.
|
294072 | assert(map.find(bulkContext_.elementIdx) != map.end()); |
304 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 96872 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 96872 times.
|
98024 | assert(bulkContext_.elementIdx == this->problem(bulkId).gridGeometry().elementMapper().index(bulkLocalAssembler.element())); |
305 | |||
306 | 98024 | typename LocalResidual<bulkId>::ElementResidualVector res(1); | |
307 | 98024 | res = 0.0; | |
308 | 98024 | res[0] = evalBulkFluxes(bulkLocalAssembler.element(), | |
309 | bulkLocalAssembler.fvGeometry(), | ||
310 | bulkLocalAssembler.curElemVolVars(), | ||
311 | bulkLocalAssembler.elemFluxVarsCache(), | ||
312 | bulkLocalAssembler.localResidual(), | ||
313 | 196048 | map.find(bulkContext_.elementIdx)->second.dofToCouplingScvfMap.at(dofIdxGlobalJ)); | |
314 | 98024 | return res; | |
315 | } | ||
316 | |||
317 | /*! | ||
318 | * \brief Evaluates the coupling element residual of a lower-dimensional domain element | ||
319 | * with respect to a dof in the bulk domain (dofIdxGlobalJ). This is essentially | ||
320 | * the fluxes across the facets of the neighboring bulk elements that coincide with | ||
321 | * the given element. | ||
322 | * | ||
323 | * \note The coupling residual in this case is always the entire transfer flux from bulk | ||
324 | * to the lowDim domain. It is therefore independent of the given dof index in the | ||
325 | * bulk domain, which is why we directly forward to the index-independent function. | ||
326 | */ | ||
327 | template< class LowDimLocalAssembler > | ||
328 | typename LocalResidual<lowDimId>::ElementResidualVector | ||
329 | ✗ | evalCouplingResidual(LowDimIdType, | |
330 | const LowDimLocalAssembler& lowDimLocalAssembler, | ||
331 | BulkIdType, | ||
332 | GridIndexType<bulkId> dofIdxGlobalJ) | ||
333 |
3/5✓ Branch 1 taken 48 times.
✓ Branch 2 taken 9034 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 528 times.
✗ Branch 6 not taken.
|
27368 | { return evalCouplingResidual(lowDimId, lowDimLocalAssembler, bulkId); } |
334 | |||
335 | /*! | ||
336 | * \brief Evaluates the coupling element residual of a lower-dimensional domain element | ||
337 | * with respect to a dof in the bulk domain (dofIdxGlobalJ). This is essentially | ||
338 | * the fluxes across the facets of the neighboring bulk elements. | ||
339 | */ | ||
340 | template< class LowDimLocalAssembler > | ||
341 | typename LocalResidual<lowDimId>::ElementResidualVector | ||
342 | 28520 | evalCouplingResidual(LowDimIdType, const LowDimLocalAssembler& lowDimLocalAssembler, BulkIdType) | |
343 | { | ||
344 | // make sure this is called for the element for which the context was set | ||
345 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 27368 times.
|
28520 | assert(lowDimContext_.isSet); |
346 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 27368 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 27368 times.
|
28520 | assert(this->problem(lowDimId).gridGeometry().elementMapper().index(lowDimLocalAssembler.element()) == lowDimContext_.elementIdx); |
347 | |||
348 | // evaluate sources for the first scv | ||
349 | // the sources are element-wise & scv-independent since we use tpfa in bulk domain | ||
350 | 57040 | const auto source = evalSourcesFromBulk(lowDimLocalAssembler.element(), | |
351 | lowDimLocalAssembler.fvGeometry(), | ||
352 | lowDimLocalAssembler.curElemVolVars(), | ||
353 | 85560 | *scvs(lowDimLocalAssembler.fvGeometry()).begin()); | |
354 | |||
355 | // fill element residual vector with the sources | ||
356 | 57040 | typename LocalResidual<lowDimId>::ElementResidualVector res(lowDimLocalAssembler.fvGeometry().numScv()); | |
357 | 28520 | res = 0.0; | |
358 |
2/2✓ Branch 0 taken 27368 times.
✓ Branch 1 taken 27368 times.
|
57040 | for (const auto& scv : scvs(lowDimLocalAssembler.fvGeometry())) |
359 | 48768 | res[scv.localDofIndex()] -= source; | |
360 | |||
361 | 28520 | return res; | |
362 | } | ||
363 | |||
364 | /*! | ||
365 | * \brief Computes the sources in a lower-dimensional element stemming from the bulk domain. | ||
366 | */ | ||
367 | 144592 | NumEqVector<lowDimId> evalSourcesFromBulk(const Element<lowDimId>& element, | |
368 | const FVElementGeometry<lowDimId>& fvGeometry, | ||
369 | const ElementVolumeVariables<lowDimId>& elemVolVars, | ||
370 | const SubControlVolume<lowDimId>& scv) const | ||
371 | { | ||
372 | // make sure the this is called for the element of the context | ||
373 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 91088 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 91088 times.
|
144592 | assert(this->problem(lowDimId).gridGeometry().elementMapper().index(element) == lowDimContext_.elementIdx); |
374 | |||
375 | 144592 | NumEqVector<lowDimId> sources(0.0); | |
376 | |||
377 | 289184 | const auto& map = couplingMapperPtr_->couplingMap(lowDimGridId, bulkGridId); | |
378 | 144592 | auto it = map.find(lowDimContext_.elementIdx); | |
379 |
4/4✓ Branch 0 taken 18144 times.
✓ Branch 1 taken 72944 times.
✓ Branch 2 taken 18144 times.
✓ Branch 3 taken 72944 times.
|
289184 | if (it == map.end()) |
380 | ✗ | return sources; | |
381 | |||
382 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 91088 times.
|
144592 | assert(lowDimContext_.isSet); |
383 |
4/4✓ Branch 0 taken 182176 times.
✓ Branch 1 taken 91088 times.
✓ Branch 2 taken 182176 times.
✓ Branch 3 taken 91088 times.
|
1157216 | for (const auto& embedment : it->second.embedments) |
384 |
1/2✓ Branch 3 taken 181696 times.
✗ Branch 4 not taken.
|
634996 | sources += evalBulkFluxes(this->problem(bulkId).gridGeometry().element(embedment.first), |
385 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 182176 times.
|
289424 | *lowDimContext_.bulkFvGeometry, |
386 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 182176 times.
|
289424 | *lowDimContext_.bulkElemVolVars, |
387 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 182176 times.
|
289424 | *lowDimContext_.bulkElemFluxVarsCache, |
388 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 182176 times.
|
289424 | *lowDimContext_.bulkLocalResidual, |
389 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 182176 times.
|
289424 | embedment.second); |
390 | |||
391 | // if lowdim domain uses box, we distribute the sources equally among the scvs | ||
392 | if (lowDimUsesBox) | ||
393 | sources /= fvGeometry.numScv(); | ||
394 | |||
395 | 126448 | return sources; | |
396 | } | ||
397 | |||
398 | /*! | ||
399 | * \brief For the assembly of the element residual of a bulk domain element | ||
400 | * we need to prepare all variables of lower-dimensional domain elements | ||
401 | * that are coupled to the given bulk element | ||
402 | */ | ||
403 | template< class Assembler > | ||
404 | 5046958 | void bindCouplingContext(BulkIdType, const Element<bulkId>& element, const Assembler& assembler) | |
405 | { | ||
406 | // clear context | ||
407 | 5046958 | bulkContext_.reset(); | |
408 | |||
409 | // set index in context in any case | ||
410 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2523479 times.
|
5046958 | const auto bulkElemIdx = this->problem(bulkId).gridGeometry().elementMapper().index(element); |
411 | 5046958 | bulkContext_.elementIdx = bulkElemIdx; | |
412 | |||
413 | // if element is coupled, actually set the context | ||
414 |
6/6✓ Branch 0 taken 432589 times.
✓ Branch 1 taken 2090890 times.
✓ Branch 2 taken 432589 times.
✓ Branch 3 taken 2090890 times.
✓ Branch 4 taken 432589 times.
✓ Branch 5 taken 2090890 times.
|
15140874 | if (bulkElemIsCoupled_[bulkElemIdx]) |
415 | { | ||
416 | 1730356 | const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId); | |
417 | |||
418 |
2/4✗ Branch 1 not taken.
✓ Branch 2 taken 432589 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 432589 times.
|
2595534 | auto it = map.find(bulkElemIdx); assert(it != map.end()); |
419 | 865178 | const auto& elementStencil = it->second.couplingElementStencil; | |
420 | 1730356 | bulkContext_.lowDimFvGeometries.reserve(elementStencil.size()); | |
421 | 1730356 | bulkContext_.lowDimVolVars.reserve(elementStencil.size()); | |
422 | |||
423 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 432589 times.
|
865178 | const auto& ldGridGeometry = this->problem(lowDimId).gridGeometry(); |
424 | 1730356 | auto fvGeom = localView(ldGridGeometry); | |
425 |
4/4✓ Branch 0 taken 877504 times.
✓ Branch 1 taken 432589 times.
✓ Branch 2 taken 877504 times.
✓ Branch 3 taken 432589 times.
|
4350542 | for (const auto lowDimElemIdx : elementStencil) |
426 | { | ||
427 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 877504 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 829000 times.
|
3413008 | const auto& ldSol = assembler.isImplicit() ? this->curSol(lowDimId) : assembler.prevSol()[lowDimId]; |
428 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 877504 times.
|
1755008 | const auto& ldProblem = this->problem(lowDimId); |
429 | |||
430 |
1/2✓ Branch 1 taken 877504 times.
✗ Branch 2 not taken.
|
1755008 | const auto elemJ = ldGridGeometry.element(lowDimElemIdx); |
431 |
1/2✓ Branch 1 taken 877504 times.
✗ Branch 2 not taken.
|
1755008 | fvGeom.bindElement(elemJ); |
432 | |||
433 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 877504 times.
|
1755008 | VolumeVariables<lowDimId> volVars; |
434 | |||
435 | // if low dim domain uses the box scheme, we have to create interpolated vol vars | ||
436 | if (lowDimUsesBox) | ||
437 | { | ||
438 | const auto elemGeom = elemJ.geometry(); | ||
439 | FacetCoupling::makeInterpolatedVolVars(volVars, ldProblem, ldSol, fvGeom, elemJ, elemGeom, elemGeom.center()); | ||
440 | } | ||
441 | // if low dim domain uses a cc scheme we can directly update the vol vars | ||
442 | else | ||
443 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 877504 times.
✓ Branch 4 taken 877504 times.
✗ Branch 5 not taken.
|
1755008 | volVars.update( elementSolution(elemJ, ldSol, ldGridGeometry), |
444 | ldProblem, | ||
445 | elemJ, | ||
446 | fvGeom.scv(lowDimElemIdx) ); | ||
447 | |||
448 | 1755008 | bulkContext_.isSet = true; | |
449 |
1/2✓ Branch 1 taken 877504 times.
✗ Branch 2 not taken.
|
1755008 | bulkContext_.lowDimFvGeometries.emplace_back( std::move(fvGeom) ); |
450 |
1/2✓ Branch 1 taken 877504 times.
✗ Branch 2 not taken.
|
1755008 | bulkContext_.lowDimVolVars.emplace_back( std::move(volVars) ); |
451 | } | ||
452 | } | ||
453 | 5046958 | } | |
454 | |||
455 | /*! | ||
456 | * \brief For the assembly of the element residual of a bulk domain element | ||
457 | * we need to prepare the local views of one of the neighboring bulk | ||
458 | * domain elements. These are used later to compute the fluxes across | ||
459 | * the faces over which the coupling occurs | ||
460 | * \note Since the low-dim coupling residua are fluxes stemming from | ||
461 | * the bulk domain, we have to prepare the bulk coupling context | ||
462 | * for the neighboring element (where fluxes are calculated) as well. | ||
463 | */ | ||
464 | template< class Assembler > | ||
465 | 209602 | void bindCouplingContext(LowDimIdType, const Element<lowDimId>& element, const Assembler& assembler) | |
466 | { | ||
467 | // reset contexts | ||
468 | 209602 | bulkContext_.reset(); | |
469 | 209602 | lowDimContext_.reset(); | |
470 | |||
471 | // set index in context in any case | ||
472 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 108526 times.
|
209602 | const auto lowDimElemIdx = this->problem(lowDimId).gridGeometry().elementMapper().index(element); |
473 | 209602 | lowDimContext_.elementIdx = lowDimElemIdx; | |
474 | |||
475 | 419204 | const auto& map = couplingMapperPtr_->couplingMap(lowDimGridId, bulkGridId); | |
476 | 209602 | auto it = map.find(lowDimElemIdx); | |
477 | |||
478 | // if element is coupled, actually set the context | ||
479 |
2/4✓ Branch 0 taken 108526 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 108526 times.
✗ Branch 3 not taken.
|
419204 | if (it != map.end()) |
480 | { | ||
481 | // first bind the low dim context for the first neighboring bulk element | ||
482 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 108526 times.
|
209602 | const auto& bulkGridGeom = this->problem(bulkId).gridGeometry(); |
483 | 628782 | const auto bulkElem = bulkGridGeom.element(it->second.embedments[0].first); | |
484 |
1/2✓ Branch 1 taken 108514 times.
✗ Branch 2 not taken.
|
209602 | bindCouplingContext(bulkId, bulkElem, assembler); |
485 | |||
486 | // evaluate variables on old/new time level depending on time disc scheme | ||
487 |
2/4✓ Branch 1 taken 108526 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 100000 times.
✗ Branch 5 not taken.
|
409602 | const auto& bulkSol = assembler.isImplicit() ? this->curSol(bulkId) : assembler.prevSol()[bulkId]; |
488 | |||
489 | // then simply bind the local views of that first neighbor | ||
490 |
2/4✓ Branch 1 taken 108526 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 108526 times.
✗ Branch 5 not taken.
|
838408 | auto bulkFvGeom = localView(bulkGridGeom).bind(bulkElem); |
491 |
4/8✓ Branch 1 taken 108526 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 8526 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 8526 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 8526 times.
✗ Branch 11 not taken.
|
448010 | auto bulkElemVolVars = assembler.isImplicit() ? localView(assembler.gridVariables(bulkId).curGridVolVars()).bind(bulkElem, bulkFvGeom, bulkSol) |
492 |
3/6✓ Branch 1 taken 100000 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 100000 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 100000 times.
✗ Branch 8 not taken.
|
600000 | : localView(assembler.gridVariables(bulkId).prevGridVolVars()).bind(bulkElem, bulkFvGeom, bulkSol); |
493 |
4/8✓ Branch 1 taken 108526 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 108526 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 108526 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 108526 times.
✗ Branch 11 not taken.
|
942411 | auto bulkElemFluxVarsCache = localView(assembler.gridVariables(bulkId).gridFluxVarsCache()).bind(bulkElem, bulkFvGeom, bulkElemVolVars); |
494 | |||
495 | 209602 | lowDimContext_.isSet = true; | |
496 |
3/6✓ Branch 1 taken 108526 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 108526 times.
✓ Branch 6 taken 108526 times.
✗ Branch 7 not taken.
|
419204 | lowDimContext_.bulkFvGeometry = std::make_unique< FVElementGeometry<bulkId> >( std::move(bulkFvGeom) ); |
497 |
3/6✓ Branch 1 taken 108526 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 108526 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 108526 times.
|
209602 | lowDimContext_.bulkElemVolVars = std::make_unique< ElementVolumeVariables<bulkId> >( std::move(bulkElemVolVars) ); |
498 |
3/7✓ Branch 1 taken 108526 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 108526 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 108526 times.
✗ Branch 7 not taken.
|
313605 | lowDimContext_.bulkElemFluxVarsCache = std::make_unique< ElementFluxVariablesCache<bulkId> >( std::move(bulkElemFluxVarsCache) ); |
499 |
4/6✓ Branch 1 taken 108526 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 108526 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 108488 times.
✓ Branch 7 taken 38 times.
|
838320 | lowDimContext_.bulkLocalResidual = std::make_unique< LocalResidual<bulkId> >(assembler.localResidual(bulkId)); |
500 | } | ||
501 | 209602 | } | |
502 | |||
503 | /*! | ||
504 | * \brief After deflecting the solution of the lower-dimensional domain, | ||
505 | * we have to update the element volume variables object if the context. | ||
506 | */ | ||
507 | template< class BulkLocalAssembler > | ||
508 | 131288 | void updateCouplingContext(BulkIdType domainI, | |
509 | const BulkLocalAssembler& bulkLocalAssembler, | ||
510 | LowDimIdType domainJ, | ||
511 | GridIndexType<lowDimId> dofIdxGlobalJ, | ||
512 | const PrimaryVariables<lowDimId>& priVarsJ, | ||
513 | unsigned int pvIdxJ) | ||
514 | { | ||
515 | // communicate deflected solution | ||
516 |
1/2✓ Branch 0 taken 130136 times.
✗ Branch 1 not taken.
|
131288 | ParentType::updateCouplingContext(domainI, bulkLocalAssembler, domainJ, dofIdxGlobalJ, priVarsJ, pvIdxJ); |
517 | |||
518 | // Since coupling only occurs via the fluxes, the context does not | ||
519 | // have to be updated in explicit time discretization schemes, where | ||
520 | // they are strictly evaluated on the old time level | ||
521 | if (!bulkLocalAssembler.isImplicit()) | ||
522 | return; | ||
523 | |||
524 | // skip the rest if context is empty | ||
525 |
1/2✓ Branch 0 taken 130136 times.
✗ Branch 1 not taken.
|
131288 | if (bulkContext_.isSet) |
526 | { | ||
527 | 262576 | const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId); | |
528 |
1/2✗ Branch 1 not taken.
✓ Branch 2 taken 130136 times.
|
262576 | const auto& couplingElemStencil = map.find(bulkContext_.elementIdx)->second.couplingElementStencil; |
529 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 130136 times.
|
131288 | const auto& ldSol = this->curSol(lowDimId); |
530 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 130136 times.
|
131288 | const auto& ldProblem = this->problem(lowDimId); |
531 | 262576 | const auto& ldGridGeometry = this->problem(lowDimId).gridGeometry(); | |
532 | |||
533 | // find the low-dim elements in coupling stencil, where this dof is contained in | ||
534 |
1/2✓ Branch 1 taken 130136 times.
✗ Branch 2 not taken.
|
652984 | const auto couplingElements = [&] () |
535 | { | ||
536 | if (lowDimUsesBox) | ||
537 | { | ||
538 | std::vector< Element<lowDimId> > lowDimElems; | ||
539 | std::for_each( couplingElemStencil.begin(), couplingElemStencil.end(), | ||
540 | [&] (auto lowDimElemIdx) | ||
541 | { | ||
542 | auto element = ldGridGeometry.element(lowDimElemIdx); | ||
543 | for (int i = 0; i < element.geometry().corners(); ++i) | ||
544 | { | ||
545 | const auto dofIdx = ldGridGeometry.vertexMapper().subIndex(element, i, lowDimDim); | ||
546 | if (dofIdxGlobalJ == dofIdx) { lowDimElems.emplace_back( std::move(element) ); break; } | ||
547 | } | ||
548 | } ); | ||
549 | return lowDimElems; | ||
550 | } | ||
551 | // dof index = element index for cc schemes | ||
552 | else | ||
553 |
0/6✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
|
777360 | return std::vector<Element<lowDimId>>( {ldGridGeometry.element(dofIdxGlobalJ)} ); |
554 | } (); | ||
555 | |||
556 | // update all necessary vol vars in context | ||
557 |
4/4✓ Branch 0 taken 130136 times.
✓ Branch 1 taken 130136 times.
✓ Branch 2 taken 130136 times.
✓ Branch 3 taken 130136 times.
|
656440 | for (const auto& element : couplingElements) |
558 | { | ||
559 | // find index in coupling context | ||
560 | 262576 | const auto eIdxGlobal = ldGridGeometry.elementMapper().index(element); | |
561 | 393864 | auto it = std::find(couplingElemStencil.begin(), couplingElemStencil.end(), eIdxGlobal); | |
562 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 130136 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 130136 times.
|
262576 | const auto idxInContext = std::distance(couplingElemStencil.begin(), it); |
563 |
3/6✗ Branch 0 not taken.
✓ Branch 1 taken 130136 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 130136 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 130136 times.
|
393864 | assert(it != couplingElemStencil.end()); |
564 | |||
565 |
2/2✓ Branch 0 taken 108192 times.
✓ Branch 1 taken 21944 times.
|
131288 | auto& volVars = bulkContext_.lowDimVolVars[idxInContext]; |
566 |
2/2✓ Branch 0 taken 108192 times.
✓ Branch 1 taken 21944 times.
|
131288 | const auto& fvGeom = bulkContext_.lowDimFvGeometries[idxInContext]; |
567 | // if low dim domain uses the box scheme, we have to create interpolated vol vars | ||
568 | if (lowDimUsesBox) | ||
569 | { | ||
570 | const auto elemGeom = element.geometry(); | ||
571 | FacetCoupling::makeInterpolatedVolVars(volVars, ldProblem, ldSol, fvGeom, element, elemGeom, elemGeom.center()); | ||
572 | } | ||
573 | // if low dim domain uses a cc scheme we can directly update the vol vars | ||
574 | else | ||
575 |
3/6✓ Branch 0 taken 108192 times.
✓ Branch 1 taken 21944 times.
✓ Branch 4 taken 130136 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
239480 | volVars.update( elementSolution(element, ldSol, ldGridGeometry), |
576 | ldProblem, | ||
577 | element, | ||
578 | fvGeom.scv(eIdxGlobal) ); | ||
579 | } | ||
580 | } | ||
581 | } | ||
582 | |||
583 | /*! | ||
584 | * \brief Update the coupling context for a derivative bulk -> bulk. | ||
585 | * Here, we simply have to update the solution. | ||
586 | */ | ||
587 | template< class BulkLocalAssembler > | ||
588 | ✗ | void updateCouplingContext(BulkIdType domainI, | |
589 | const BulkLocalAssembler& bulkLocalAssembler, | ||
590 | BulkIdType domainJ, | ||
591 | GridIndexType<bulkId> dofIdxGlobalJ, | ||
592 | const PrimaryVariables<bulkId>& priVarsJ, | ||
593 | unsigned int pvIdxJ) | ||
594 | { | ||
595 | // communicate deflected solution | ||
596 | 3649412 | ParentType::updateCouplingContext(domainI, bulkLocalAssembler, domainJ, dofIdxGlobalJ, priVarsJ, pvIdxJ); | |
597 | ✗ | } | |
598 | |||
599 | /*! | ||
600 | * \brief After deflecting the solution of the bulk domain, we have to update | ||
601 | * the element volume variables and transmissibilities of the neighboring | ||
602 | * bulk element stored in the context. | ||
603 | */ | ||
604 | template< class LowDimLocalAssembler > | ||
605 | 131288 | void updateCouplingContext(LowDimIdType domainI, | |
606 | const LowDimLocalAssembler& lowDimLocalAssembler, | ||
607 | BulkIdType domainJ, | ||
608 | GridIndexType<bulkId> dofIdxGlobalJ, | ||
609 | const PrimaryVariables<bulkId>& priVarsJ, | ||
610 | unsigned int pvIdxJ) | ||
611 | { | ||
612 | // communicate deflected solution | ||
613 |
1/2✓ Branch 0 taken 130136 times.
✗ Branch 1 not taken.
|
131288 | ParentType::updateCouplingContext(domainI, lowDimLocalAssembler, domainJ, dofIdxGlobalJ, priVarsJ, pvIdxJ); |
614 | |||
615 | // Since coupling only occurs via the fluxes, the context does not | ||
616 | // have to be updated in explicit time discretization schemes, where | ||
617 | // they are strictly evaluated on the old time level | ||
618 | if (!lowDimLocalAssembler.isImplicit()) | ||
619 | return; | ||
620 | |||
621 | // skip the rest if context is empty | ||
622 |
1/2✓ Branch 0 taken 130136 times.
✗ Branch 1 not taken.
|
131288 | if (lowDimContext_.isSet) |
623 | { | ||
624 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 130136 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 130136 times.
|
131288 | assert(lowDimContext_.elementIdx == this->problem(lowDimId).gridGeometry().elementMapper().index(lowDimLocalAssembler.element())); |
625 | |||
626 | // since we use cc scheme in bulk domain: dof index = element index | ||
627 |
3/4✗ Branch 0 not taken.
✓ Branch 1 taken 130136 times.
✓ Branch 2 taken 72 times.
✓ Branch 3 taken 24 times.
|
131288 | const auto& bulkGridGeom = this->problem(bulkId).gridGeometry(); |
628 |
2/2✓ Branch 0 taken 72 times.
✓ Branch 1 taken 24 times.
|
262384 | const auto elementJ = bulkGridGeom.element(dofIdxGlobalJ); |
629 | |||
630 | // update corresponding vol vars in context | ||
631 |
4/4✓ Branch 0 taken 42532 times.
✓ Branch 1 taken 87604 times.
✓ Branch 2 taken 42532 times.
✓ Branch 3 taken 87604 times.
|
262576 | const auto& scv = lowDimContext_.bulkFvGeometry->scv(dofIdxGlobalJ); |
632 |
2/4✓ Branch 1 taken 130040 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 130040 times.
✗ Branch 5 not taken.
|
262576 | const auto elemSol = elementSolution(elementJ, this->curSol(bulkId), bulkGridGeom); |
633 |
2/4✗ Branch 2 not taken.
✓ Branch 3 taken 130136 times.
✓ Branch 5 taken 130040 times.
✗ Branch 6 not taken.
|
262576 | (*lowDimContext_.bulkElemVolVars)[dofIdxGlobalJ].update(elemSol, this->problem(bulkId), elementJ, scv); |
634 | |||
635 | // update the element flux variables cache (tij might be solution-dependent) | ||
636 |
2/2✓ Branch 0 taken 33172 times.
✓ Branch 1 taken 96964 times.
|
131288 | if (dofIdxGlobalJ == bulkContext_.elementIdx) |
637 |
4/8✓ Branch 1 taken 33148 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 33148 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 33148 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 33148 times.
✗ Branch 11 not taken.
|
134896 | lowDimContext_.bulkElemFluxVarsCache->update( elementJ, *lowDimContext_.bulkFvGeometry, *lowDimContext_.bulkElemVolVars); |
638 | else | ||
639 |
4/8✗ Branch 0 not taken.
✓ Branch 1 taken 96964 times.
✓ Branch 3 taken 96892 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 96892 times.
✗ Branch 7 not taken.
✓ Branch 9 taken 96892 times.
✗ Branch 10 not taken.
|
97564 | lowDimContext_.bulkElemFluxVarsCache->update( this->problem(bulkId).gridGeometry().element(bulkContext_.elementIdx), |
640 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 96964 times.
|
97564 | *lowDimContext_.bulkFvGeometry, |
641 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 96964 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 96964 times.
|
195128 | *lowDimContext_.bulkElemVolVars ); |
642 | } | ||
643 | } | ||
644 | |||
645 | /*! | ||
646 | * \brief After deflecting the solution of the lower-dimensional domain has been deflected | ||
647 | * during the assembly of the element residual of a lower-dimensional element, we | ||
648 | * have to communicate this to the volume variables stored in the context as well | ||
649 | * as the transmissibilities. | ||
650 | */ | ||
651 | template< class LowDimLocalAssembler > | ||
652 | 33724 | void updateCouplingContext(LowDimIdType domainI, | |
653 | const LowDimLocalAssembler& lowDimLocalAssembler, | ||
654 | LowDimIdType domainJ, | ||
655 | GridIndexType<lowDimId> dofIdxGlobalJ, | ||
656 | const PrimaryVariables<lowDimId>& priVarsJ, | ||
657 | unsigned int pvIdxJ) | ||
658 | { | ||
659 | // communicate deflected solution | ||
660 |
1/2✓ Branch 0 taken 33172 times.
✗ Branch 1 not taken.
|
33724 | ParentType::updateCouplingContext(domainI, lowDimLocalAssembler, domainJ, dofIdxGlobalJ, priVarsJ, pvIdxJ); |
661 | |||
662 | // Since coupling only occurs via the fluxes, the context does not | ||
663 | // have to be updated in explicit time discretization schemes, where | ||
664 | // they are strictly evaluated on the old time level | ||
665 | if (!lowDimLocalAssembler.isImplicit()) | ||
666 | return; | ||
667 | |||
668 | // skip the rest if context is empty | ||
669 |
1/2✓ Branch 0 taken 33172 times.
✗ Branch 1 not taken.
|
33724 | if (lowDimContext_.isSet) |
670 | { | ||
671 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 33172 times.
|
33724 | const auto& ldSol = this->curSol(lowDimId); |
672 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 33172 times.
|
33724 | const auto& ldProblem = this->problem(lowDimId); |
673 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 33172 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 33172 times.
|
67448 | const auto& ldGridGeometry = this->problem(lowDimId).gridGeometry(); |
674 | |||
675 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 33172 times.
|
33724 | assert(bulkContext_.isSet); |
676 |
1/2✗ Branch 2 not taken.
✓ Branch 3 taken 33172 times.
|
67448 | assert(lowDimContext_.elementIdx == ldGridGeometry.elementMapper().index(lowDimLocalAssembler.element())); |
677 | |||
678 | // update the corresponding vol vars in the bulk context | ||
679 | 67448 | const auto& bulkMap = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId); | |
680 | 67448 | const auto& couplingElementStencil = bulkMap.find(bulkContext_.elementIdx)->second.couplingElementStencil; | |
681 | 101172 | auto it = std::find(couplingElementStencil.begin(), couplingElementStencil.end(), lowDimContext_.elementIdx); | |
682 |
3/6✗ Branch 0 not taken.
✓ Branch 1 taken 33172 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 33172 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 33172 times.
|
101172 | assert(it != couplingElementStencil.end()); |
683 |
4/4✓ Branch 0 taken 26880 times.
✓ Branch 1 taken 6292 times.
✓ Branch 2 taken 26880 times.
✓ Branch 3 taken 6292 times.
|
67448 | const auto idxInContext = std::distance(couplingElementStencil.begin(), it); |
684 | |||
685 |
2/2✓ Branch 0 taken 26880 times.
✓ Branch 1 taken 6292 times.
|
33724 | auto& volVars = bulkContext_.lowDimVolVars[idxInContext]; |
686 |
2/2✓ Branch 0 taken 26880 times.
✓ Branch 1 taken 6292 times.
|
33724 | const auto& fvGeom = bulkContext_.lowDimFvGeometries[idxInContext]; |
687 | 33724 | const auto& element = lowDimLocalAssembler.element(); | |
688 | // if low dim domain uses the box scheme, we have to create interpolated vol vars | ||
689 | if (lowDimUsesBox) | ||
690 | { | ||
691 | const auto elemGeom = element.geometry(); | ||
692 | FacetCoupling::makeInterpolatedVolVars(volVars, ldProblem, ldSol, fvGeom, element, elemGeom, elemGeom.center()); | ||
693 | } | ||
694 | // if low dim domain uses a cc scheme we can directly update the vol vars | ||
695 | else | ||
696 |
2/2✓ Branch 0 taken 26880 times.
✓ Branch 1 taken 6292 times.
|
60604 | volVars.update( elementSolution(element, ldSol, ldGridGeometry), |
697 | ldProblem, | ||
698 | element, | ||
699 | fvGeom.scv(lowDimContext_.elementIdx) ); | ||
700 | |||
701 | // update the element flux variables cache (tij depend on low dim values in context) | ||
702 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 33172 times.
|
67400 | const auto contextElem = this->problem(bulkId).gridGeometry().element(bulkContext_.elementIdx); |
703 |
4/8✓ Branch 1 taken 33148 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 33148 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 33148 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 33148 times.
✗ Branch 11 not taken.
|
134896 | lowDimContext_.bulkElemFluxVarsCache->update(contextElem, *lowDimContext_.bulkFvGeometry, *lowDimContext_.bulkElemVolVars); |
704 | } | ||
705 | } | ||
706 | |||
707 | //! pull up empty update function to be used for low-dim domain | ||
708 | using ParentType::updateCoupledVariables; | ||
709 | |||
710 | /*! | ||
711 | * \brief Update the transmissibilities in the bulk domain after the coupling context changed | ||
712 | * \note Specialization of the function for deactivated grid-wide volume variables caching | ||
713 | */ | ||
714 | template< class BulkLocalAssembler, class UpdatableFluxVarCache > | ||
715 | ✗ | void updateCoupledVariables(BulkIdType domainI, | |
716 | const BulkLocalAssembler& bulkLocalAssembler, | ||
717 | ElementVolumeVariables<bulkId>& elemVolVars, | ||
718 | UpdatableFluxVarCache& fluxVarsCache) | ||
719 | { | ||
720 | // update transmissibilities after low dim context has changed (implicit only) | ||
721 | if (bulkLocalAssembler.isImplicit()) | ||
722 |
10/26✓ Branch 1 taken 48 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 31276 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 31276 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 31228 times.
✗ Branch 11 not taken.
✓ Branch 16 taken 528 times.
✗ Branch 17 not taken.
✓ Branch 19 taken 528 times.
✗ Branch 20 not taken.
✓ Branch 22 taken 528 times.
✗ Branch 23 not taken.
✗ Branch 31 not taken.
✗ Branch 32 not taken.
✗ Branch 34 not taken.
✗ Branch 35 not taken.
✗ Branch 37 not taken.
✗ Branch 38 not taken.
✓ Branch 46 taken 528 times.
✗ Branch 47 not taken.
✓ Branch 49 taken 528 times.
✗ Branch 50 not taken.
✓ Branch 52 taken 528 times.
✗ Branch 53 not taken.
|
293784 | fluxVarsCache.update(bulkLocalAssembler.element(), |
723 | bulkLocalAssembler.fvGeometry(), | ||
724 | bulkLocalAssembler.curElemVolVars()); | ||
725 | ✗ | } | |
726 | |||
727 | /*! | ||
728 | * \brief Update the transmissibilities in the bulk domain after the coupling context changed | ||
729 | * \note Specialization of the function for enabled grid-wide volume variables caching | ||
730 | */ | ||
731 | template< class BulkLocalAssembler, class UpdatableFluxVarCache > | ||
732 | void updateCoupledVariables(BulkIdType domainI, | ||
733 | const BulkLocalAssembler& bulkLocalAssembler, | ||
734 | GridVolumeVariables<bulkId>& gridVolVars, | ||
735 | UpdatableFluxVarCache& fluxVarsCache) | ||
736 | { | ||
737 | // update transmissibilities after low dim context has changed (implicit only) | ||
738 | if (bulkLocalAssembler.isImplicit()) | ||
739 | { | ||
740 | const auto elemVolVars = localView(gridVolVars).bind(bulkLocalAssembler.element(), bulkLocalAssembler.fvGeometry(), this->curSol(bulkId)); | ||
741 | fluxVarsCache.update(bulkLocalAssembler.element(), bulkLocalAssembler.fvGeometry(), elemVolVars); | ||
742 | } | ||
743 | } | ||
744 | |||
745 | //! Empty stencil to be returned for elements that aren't coupled | ||
746 | template<std::size_t id, std::enable_if_t<(id == bulkId || id == lowDimId), int> = 0> | ||
747 | const typename CouplingMapper::template Stencil<id>& | ||
748 | ✗ | getEmptyStencil(Dune::index_constant<id>) const | |
749 |
2/4✓ Branch 7 taken 1927 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 1927 times.
✗ Branch 11 not taken.
|
1129744 | { return std::get<(id == bulkId ? 0 : 1)>(emptyStencilTuple_); } |
750 | |||
751 | protected: | ||
752 | //! Return const references to the bulk coupling contexts | ||
753 | const BulkCouplingContext& bulkCouplingContext() const { return bulkContext_; } | ||
754 | const LowDimCouplingContext& lowDimCouplingContext() const { return lowDimContext_; } | ||
755 | |||
756 | //! Return references to the bulk coupling contexts | ||
757 | BulkCouplingContext& bulkCouplingContext() { return bulkContext_; } | ||
758 | LowDimCouplingContext& lowDimCouplingContext() { return lowDimContext_; } | ||
759 | |||
760 | //! evaluates the bulk-facet exchange fluxes for a given facet element | ||
761 | template<class BulkScvfIndices> | ||
762 | 575424 | NumEqVector<bulkId> evalBulkFluxes(const Element<bulkId>& elementI, | |
763 | const FVElementGeometry<bulkId>& fvGeometry, | ||
764 | const ElementVolumeVariables<bulkId>& elemVolVars, | ||
765 | const ElementFluxVariablesCache<bulkId>& elemFluxVarsCache, | ||
766 | const LocalResidual<bulkId>& localResidual, | ||
767 | const BulkScvfIndices& scvfIndices) const | ||
768 | { | ||
769 | |||
770 | 575424 | NumEqVector<bulkId> coupledFluxes(0.0); | |
771 |
4/4✓ Branch 0 taken 630368 times.
✓ Branch 1 taken 344092 times.
✓ Branch 2 taken 630368 times.
✓ Branch 3 taken 344092 times.
|
3697936 | for (const auto& scvfIdx : scvfIndices) |
772 |
1/2✗ Branch 1 not taken.
✓ Branch 2 taken 630368 times.
|
1535000 | coupledFluxes += localResidual.evalFlux(this->problem(bulkId), |
773 | elementI, | ||
774 | fvGeometry, | ||
775 | elemVolVars, | ||
776 | elemFluxVarsCache, | ||
777 | fvGeometry.scvf(scvfIdx)); | ||
778 | 575424 | return coupledFluxes; | |
779 | } | ||
780 | |||
781 | private: | ||
782 | std::shared_ptr<CouplingMapper> couplingMapperPtr_; | ||
783 | |||
784 | //! store bools for all bulk elements/scvfs that indicate if they | ||
785 | //! are coupled, so that we don't have to search in the map every time | ||
786 | std::vector<bool> bulkElemIsCoupled_; | ||
787 | std::vector<bool> bulkScvfIsCoupled_; | ||
788 | |||
789 | //! empty stencil to return for non-coupled elements | ||
790 | using BulkStencil = typename CouplingMapper::template Stencil<bulkId>; | ||
791 | using LowDimStencil = typename CouplingMapper::template Stencil<lowDimId>; | ||
792 | std::tuple<BulkStencil, LowDimStencil> emptyStencilTuple_; | ||
793 | |||
794 | //! The coupling contexts of the two domains | ||
795 | BulkCouplingContext bulkContext_; | ||
796 | LowDimCouplingContext lowDimContext_; | ||
797 | }; | ||
798 | |||
799 | } // end namespace Dumux | ||
800 | |||
801 | #endif | ||
802 |