3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
multidomain/facet/cellcentered/tpfa/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 *****************************************************************************/
24#ifndef DUMUX_CCTPFA_FACETCOUPLING_MANAGER_HH
25#define DUMUX_CCTPFA_FACETCOUPLING_MANAGER_HH
26
27#include <algorithm>
28#include <cassert>
29
36
37namespace Dumux {
38
50template<class MDTraits, class CouplingMapper, std::size_t bulkDomainId, std::size_t lowDimDomainId>
51class FacetCouplingManager<MDTraits, CouplingMapper, bulkDomainId, lowDimDomainId, DiscretizationMethod::cctpfa>
52: public virtual CouplingManager< MDTraits >
53{
55
56 // convenience aliases and instances of the two domain ids
57 using BulkIdType = typename MDTraits::template SubDomain<bulkDomainId>::Index;
58 using LowDimIdType = typename MDTraits::template SubDomain<lowDimDomainId>::Index;
59 static constexpr auto bulkId = BulkIdType();
60 static constexpr auto lowDimId = LowDimIdType();
61
62 // the sub-domain type tags
63 template<std::size_t id> using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
64
65 // further types specific to the sub-problems
66 template<std::size_t id> using PrimaryVariables = GetPropType<SubDomainTypeTag<id>, Properties::PrimaryVariables>;
67 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
68 template<std::size_t id> using NumEqVector = GetPropType<SubDomainTypeTag<id>, Properties::NumEqVector>;
69 template<std::size_t id> using LocalResidual = GetPropType<SubDomainTypeTag<id>, Properties::LocalResidual>;
70
71 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
72 template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView;
73 template<std::size_t id> using SubControlVolume = typename GridGeometry<id>::SubControlVolume;
74 template<std::size_t id> using SubControlVolumeFace = typename GridGeometry<id>::SubControlVolumeFace;
75 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
76 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
77 template<std::size_t id> using GridIndexType = typename IndexTraits< GridView<id> >::GridIndex;
78
79 template<std::size_t id> using GridVariables = GetPropType<SubDomainTypeTag<id>, Properties::GridVariables>;
80 template<std::size_t id> using GridVolumeVariables = typename GridVariables<id>::GridVolumeVariables;
81 template<std::size_t id> using ElementVolumeVariables = typename GridVolumeVariables<id>::LocalView;
82 template<std::size_t id> using VolumeVariables = typename ElementVolumeVariables<id>::VolumeVariables;
83 template<std::size_t id> using GridFluxVariablesCache = typename GridVariables<id>::GridFluxVariablesCache;
84 template<std::size_t id> using ElementFluxVariablesCache = typename GridFluxVariablesCache<id>::LocalView;
85
86 // this currently does not work for some grid-wide caches being active
87 static_assert(!getPropValue<SubDomainTypeTag<bulkId>, Properties::EnableGridFluxVariablesCache>(),
88 "Grid flux variables caching currently not supported in the bulk domain of cc-facet coupling models");
89 static_assert(!getPropValue<SubDomainTypeTag<lowDimId>, Properties::EnableGridVolumeVariablesCache>(),
90 "Grid volume variables caching currently not supported in the lower-dimensional domain of cc-facet coupling models");
91 static_assert(!getPropValue<SubDomainTypeTag<bulkId>, Properties::EnableGridVolumeVariablesCache>(),
92 "Grid volume variables caching currently not supported in the bulk domain of cc-facet coupling models");
93
94 // extract corresponding grid ids from the mapper
95 static constexpr int bulkDim = GridView<bulkDomainId>::dimension;
96 static constexpr int lowDimDim = GridView<lowDimDomainId>::dimension;
97 static constexpr auto bulkGridId = CouplingMapper::template gridId<bulkDim>();
98 static constexpr auto lowDimGridId = CouplingMapper::template gridId<lowDimDim>();
99
100 static constexpr bool lowDimUsesBox = GridGeometry<lowDimId>::discMethod == DiscretizationMethod::box;
101
108 struct BulkCouplingContext
109 {
110 bool isSet;
111 GridIndexType< bulkId > elementIdx;
112 std::vector< FVElementGeometry<lowDimId> > lowDimFvGeometries;
113 std::vector< VolumeVariables<lowDimId> > lowDimVolVars;
114
115 void reset()
116 {
117 lowDimFvGeometries.clear();
118 lowDimVolVars.clear();
119 isSet = false;
120 }
121 };
122
134 struct LowDimCouplingContext
135 {
136 bool isSet;
137 GridIndexType< lowDimId > elementIdx;
138 std::unique_ptr< FVElementGeometry<bulkId> > bulkFvGeometry;
139 std::unique_ptr< ElementVolumeVariables<bulkId> > bulkElemVolVars;
140 std::unique_ptr< ElementFluxVariablesCache<bulkId> > bulkElemFluxVarsCache;
141 std::unique_ptr< LocalResidual<bulkId> > bulkLocalResidual;
142
143 void reset()
144 {
145 bulkFvGeometry.reset(nullptr);
146 bulkElemVolVars.reset(nullptr);
147 bulkElemFluxVarsCache.reset(nullptr);
148 isSet = false;
149 }
150 };
151
152public:
153
155 template<std::size_t i, std::size_t j = (i == bulkId) ? lowDimId : bulkId>
156 using CouplingStencilType = typename CouplingMapper::template Stencil< CouplingMapper::template gridId<GridView<j>::dimension>() >;
157
159 using SolutionVector = typename MDTraits::SolutionVector;
160
169 void init(std::shared_ptr< Problem<bulkId> > bulkProblem,
170 std::shared_ptr< Problem<lowDimId> > lowDimProblem,
171 std::shared_ptr< CouplingMapper > couplingMapper,
172 const SolutionVector& curSol)
173 {
174 couplingMapperPtr_ = couplingMapper;
175
176 // set the sub problems
177 this->setSubProblem(bulkProblem, bulkId);
178 this->setSubProblem(lowDimProblem, lowDimId);
179
180 // copy the solution vector
181 ParentType::updateSolution(curSol);
182
183 // determine all bulk elements/scvfs that couple to low dim elements
184 bulkElemIsCoupled_.assign(bulkProblem->gridGeometry().gridView().size(0), false);
185 bulkScvfIsCoupled_.assign(bulkProblem->gridGeometry().numScvf(), false);
186
187 const auto& bulkMap = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
188 for (const auto& entry : bulkMap)
189 {
190 bulkElemIsCoupled_[entry.first] = true;
191 for (const auto& couplingEntry : entry.second.dofToCouplingScvfMap)
192 for (const auto& scvfIdx : couplingEntry.second)
193 bulkScvfIsCoupled_[scvfIdx] = true;
194 }
195 }
196
201 const Element<bulkId>& element,
202 LowDimIdType domainJ) const
203 {
204 const auto eIdx = this->problem(bulkId).gridGeometry().elementMapper().index(element);
205
206 if (bulkElemIsCoupled_[eIdx])
207 {
208 const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
209 auto it = map.find(eIdx);
210 assert(it != map.end());
211 return it->second.couplingStencil;
212 }
213
214 return getEmptyStencil(lowDimId);
215 }
216
221 const Element<lowDimId>& element,
222 BulkIdType domainJ) const
223 {
224 const auto eIdx = this->problem(lowDimId).gridGeometry().elementMapper().index(element);
225
226 const auto& map = couplingMapperPtr_->couplingMap(lowDimGridId, bulkGridId);
227 auto it = map.find(eIdx);
228 if (it != map.end()) return it->second.couplingStencil;
229 else return getEmptyStencil(bulkId);
230 }
231
235 bool isCoupled(const Element<bulkId>& element,
236 const SubControlVolumeFace<bulkId>& scvf) const
237 { return bulkScvfIsCoupled_[scvf.index()]; }
238
243 bool isOnInteriorBoundary(const Element<bulkId>& element,
244 const SubControlVolumeFace<bulkId>& scvf) const
245 { return isCoupled(element, scvf); }
246
250 const VolumeVariables<lowDimId>& getLowDimVolVars(const Element<bulkId>& element,
251 const SubControlVolumeFace<bulkId>& scvf) const
252 {
253 assert(bulkContext_.isSet);
254
255 const auto lowDimElemIdx = getLowDimElementIndex(element, scvf);
256 const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
257 const auto& s = map.find(bulkContext_.elementIdx)->second.couplingElementStencil;
258 const auto& idxInContext = std::distance( s.begin(), std::find(s.begin(), s.end(), lowDimElemIdx) );
259
260 assert(std::find(s.begin(), s.end(), lowDimElemIdx) != s.end());
261 return bulkContext_.lowDimVolVars[idxInContext];
262 }
263
267 const Element<lowDimId> getLowDimElement(const Element<bulkId>& element,
268 const SubControlVolumeFace<bulkId>& scvf) const
269 {
270 const auto lowDimElemIdx = getLowDimElementIndex(element, scvf);
271 return this->problem(lowDimId).gridGeometry().element(lowDimElemIdx);
272 }
273
277 const GridIndexType<lowDimId> getLowDimElementIndex(const Element<bulkId>& element,
278 const SubControlVolumeFace<bulkId>& scvf) const
279 {
280 assert(bulkScvfIsCoupled_[scvf.index()]);
281
282 const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
283 const auto& couplingData = map.at(scvf.insideScvIdx());
284
285 // search the low dim element idx this scvf is embedded in
286 auto it = std::find_if( couplingData.elementToScvfMap.begin(),
287 couplingData.elementToScvfMap.end(),
288 [&scvf] (auto& dataPair)
289 {
290 const auto& scvfs = dataPair.second;
291 return std::find(scvfs.begin(), scvfs.end(), scvf.index()) != scvfs.end();
292 } );
293
294 assert(it != couplingData.elementToScvfMap.end());
295 return it->first;
296 }
297
303 template< class BulkLocalAssembler >
304 typename LocalResidual<bulkId>::ElementResidualVector
306 const BulkLocalAssembler& bulkLocalAssembler,
307 LowDimIdType,
308 GridIndexType<lowDimId> dofIdxGlobalJ)
309 {
310 const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
311
312 assert(bulkContext_.isSet);
313 assert(bulkElemIsCoupled_[bulkContext_.elementIdx]);
314 assert(map.find(bulkContext_.elementIdx) != map.end());
315 assert(bulkContext_.elementIdx == this->problem(bulkId).gridGeometry().elementMapper().index(bulkLocalAssembler.element()));
316
317 typename LocalResidual<bulkId>::ElementResidualVector res(1);
318 res = 0.0;
319 res[0] = evalBulkFluxes(bulkLocalAssembler.element(),
320 bulkLocalAssembler.fvGeometry(),
321 bulkLocalAssembler.curElemVolVars(),
322 bulkLocalAssembler.elemFluxVarsCache(),
323 bulkLocalAssembler.localResidual(),
324 map.find(bulkContext_.elementIdx)->second.dofToCouplingScvfMap.at(dofIdxGlobalJ));
325 return res;
326 }
327
338 template< class LowDimLocalAssembler >
339 typename LocalResidual<lowDimId>::ElementResidualVector
341 const LowDimLocalAssembler& lowDimLocalAssembler,
342 BulkIdType,
343 GridIndexType<bulkId> dofIdxGlobalJ)
344 { return evalCouplingResidual(lowDimId, lowDimLocalAssembler, bulkId); }
345
351 template< class LowDimLocalAssembler >
352 typename LocalResidual<lowDimId>::ElementResidualVector
353 evalCouplingResidual(LowDimIdType, const LowDimLocalAssembler& lowDimLocalAssembler, BulkIdType)
354 {
355 // make sure this is called for the element for which the context was set
356 assert(lowDimContext_.isSet);
357 assert(this->problem(lowDimId).gridGeometry().elementMapper().index(lowDimLocalAssembler.element()) == lowDimContext_.elementIdx);
358
359 // evaluate sources for the first scv
360 // the sources are element-wise & scv-independent since we use tpfa in bulk domain
361 const auto source = evalSourcesFromBulk(lowDimLocalAssembler.element(),
362 lowDimLocalAssembler.fvGeometry(),
363 lowDimLocalAssembler.curElemVolVars(),
364 *scvs(lowDimLocalAssembler.fvGeometry()).begin());
365
366 // fill element residual vector with the sources
367 typename LocalResidual<lowDimId>::ElementResidualVector res(lowDimLocalAssembler.fvGeometry().numScv());
368 res = 0.0;
369 for (const auto& scv : scvs(lowDimLocalAssembler.fvGeometry()))
370 res[scv.localDofIndex()] -= source;
371
372 return res;
373 }
374
378 NumEqVector<lowDimId> evalSourcesFromBulk(const Element<lowDimId>& element,
379 const FVElementGeometry<lowDimId>& fvGeometry,
380 const ElementVolumeVariables<lowDimId>& elemVolVars,
381 const SubControlVolume<lowDimId>& scv)
382 {
383 // make sure the this is called for the element of the context
384 assert(this->problem(lowDimId).gridGeometry().elementMapper().index(element) == lowDimContext_.elementIdx);
385
386 NumEqVector<lowDimId> sources(0.0);
387
388 const auto& map = couplingMapperPtr_->couplingMap(lowDimGridId, bulkGridId);
389 auto it = map.find(lowDimContext_.elementIdx);
390 if (it == map.end())
391 return sources;
392
393 assert(lowDimContext_.isSet);
394 for (const auto& embedment : it->second.embedments)
395 sources += evalBulkFluxes(this->problem(bulkId).gridGeometry().element(embedment.first),
396 *lowDimContext_.bulkFvGeometry,
397 *lowDimContext_.bulkElemVolVars,
398 *lowDimContext_.bulkElemFluxVarsCache,
399 *lowDimContext_.bulkLocalResidual,
400 embedment.second);
401
402 // if lowdim domain uses box, we distribute the sources equally among the scvs
403 if (lowDimUsesBox)
404 sources /= fvGeometry.numScv();
405
406 return sources;
407 }
408
414 template< class Assembler >
415 void bindCouplingContext(BulkIdType, const Element<bulkId>& element, const Assembler& assembler)
416 {
417 // clear context
418 bulkContext_.reset();
419
420 // set index in context in any case
421 const auto bulkElemIdx = this->problem(bulkId).gridGeometry().elementMapper().index(element);
422 bulkContext_.elementIdx = bulkElemIdx;
423
424 // if element is coupled, actually set the context
425 if (bulkElemIsCoupled_[bulkElemIdx])
426 {
427 const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
428
429 auto it = map.find(bulkElemIdx); assert(it != map.end());
430 const auto& elementStencil = it->second.couplingElementStencil;
431 bulkContext_.lowDimFvGeometries.reserve(elementStencil.size());
432 bulkContext_.lowDimVolVars.reserve(elementStencil.size());
433
434 for (const auto lowDimElemIdx : elementStencil)
435 {
436 const auto& ldSol = Assembler::isImplicit() ? this->curSol()[lowDimId] : assembler.prevSol()[lowDimId];
437 const auto& ldProblem = this->problem(lowDimId);
438 const auto& ldGridGeometry = this->problem(lowDimId).gridGeometry();
439
440 const auto elemJ = ldGridGeometry.element(lowDimElemIdx);
441 auto fvGeom = localView(ldGridGeometry);
442 fvGeom.bindElement(elemJ);
443
444 VolumeVariables<lowDimId> volVars;
445
446 // if low dim domain uses the box scheme, we have to create interpolated vol vars
447 if (lowDimUsesBox)
448 {
449 const auto elemGeom = elemJ.geometry();
450 FacetCoupling::makeInterpolatedVolVars(volVars, ldProblem, ldSol, fvGeom, elemJ, elemGeom, elemGeom.center());
451 }
452 // if low dim domain uses a cc scheme we can directly update the vol vars
453 else
454 volVars.update( elementSolution(elemJ, ldSol, ldGridGeometry),
455 ldProblem,
456 elemJ,
457 fvGeom.scv(lowDimElemIdx) );
458
459 bulkContext_.isSet = true;
460 bulkContext_.lowDimFvGeometries.emplace_back( std::move(fvGeom) );
461 bulkContext_.lowDimVolVars.emplace_back( std::move(volVars) );
462 }
463 }
464 }
465
475 template< class Assembler >
476 void bindCouplingContext(LowDimIdType, const Element<lowDimId>& element, const Assembler& assembler)
477 {
478 // reset contexts
479 bulkContext_.reset();
480 lowDimContext_.reset();
481
482 // set index in context in any case
483 const auto lowDimElemIdx = this->problem(lowDimId).gridGeometry().elementMapper().index(element);
484 lowDimContext_.elementIdx = lowDimElemIdx;
485
486 const auto& map = couplingMapperPtr_->couplingMap(lowDimGridId, bulkGridId);
487 auto it = map.find(lowDimElemIdx);
488
489 // if element is coupled, actually set the context
490 if (it != map.end())
491 {
492 // first bind the low dim context for the first neighboring bulk element
493 const auto& bulkGridGeom = this->problem(bulkId).gridGeometry();
494 const auto bulkElem = bulkGridGeom.element(it->second.embedments[0].first);
495 bindCouplingContext(bulkId, bulkElem, assembler);
496
497 // then simply bind the local views of that first neighbor
498 auto bulkFvGeom = localView(bulkGridGeom);
499 auto bulkElemVolVars = Assembler::isImplicit() ? localView(assembler.gridVariables(bulkId).curGridVolVars())
500 : localView(assembler.gridVariables(bulkId).prevGridVolVars());
501 auto bulkElemFluxVarsCache = localView(assembler.gridVariables(bulkId).gridFluxVarsCache());
502
503 // evaluate variables on old/new time level depending on time disc scheme
504 const auto& bulkSol = Assembler::isImplicit() ? this->curSol()[bulkId] : assembler.prevSol()[bulkId];
505 bulkFvGeom.bind(bulkElem);
506 bulkElemVolVars.bind(bulkElem, bulkFvGeom, bulkSol);
507 bulkElemFluxVarsCache.bind(bulkElem, bulkFvGeom, bulkElemVolVars);
508
509 lowDimContext_.isSet = true;
510 lowDimContext_.bulkFvGeometry = std::make_unique< FVElementGeometry<bulkId> >( std::move(bulkFvGeom) );
511 lowDimContext_.bulkElemVolVars = std::make_unique< ElementVolumeVariables<bulkId> >( std::move(bulkElemVolVars) );
512 lowDimContext_.bulkElemFluxVarsCache = std::make_unique< ElementFluxVariablesCache<bulkId> >( std::move(bulkElemFluxVarsCache) );
513 lowDimContext_.bulkLocalResidual = std::make_unique< LocalResidual<bulkId> >(assembler.localResidual(bulkId));
514 }
515 }
516
521 template< class BulkLocalAssembler >
522 void updateCouplingContext(BulkIdType domainI,
523 const BulkLocalAssembler& bulkLocalAssembler,
524 LowDimIdType domainJ,
525 GridIndexType<lowDimId> dofIdxGlobalJ,
526 const PrimaryVariables<lowDimId>& priVarsJ,
527 unsigned int pvIdxJ)
528 {
529 // communicate deflected solution
530 ParentType::updateCouplingContext(domainI, bulkLocalAssembler, domainJ, dofIdxGlobalJ, priVarsJ, pvIdxJ);
531
532 // Since coupling only occurs via the fluxes, the context does not
533 // have to be updated in explicit time discretization schemes, where
534 // they are strictly evaluated on the old time level
535 if (!BulkLocalAssembler::isImplicit())
536 return;
537
538 // skip the rest if context is empty
539 if (bulkContext_.isSet)
540 {
541 const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
542 const auto& couplingElemStencil = map.find(bulkContext_.elementIdx)->second.couplingElementStencil;
543 const auto& ldSol = this->curSol()[lowDimId];
544 const auto& ldProblem = this->problem(lowDimId);
545 const auto& ldGridGeometry = this->problem(lowDimId).gridGeometry();
546
547 // find the low-dim elements in coupling stencil, where this dof is contained in
548 const auto couplingElements = [&] ()
549 {
550 if (lowDimUsesBox)
551 {
552 std::vector< Element<lowDimId> > lowDimElems;
553 std::for_each( couplingElemStencil.begin(), couplingElemStencil.end(),
554 [&] (auto lowDimElemIdx)
555 {
556 auto element = ldGridGeometry.element(lowDimElemIdx);
557 for (int i = 0; i < element.geometry().corners(); ++i)
558 {
559 const auto dofIdx = ldGridGeometry.vertexMapper().subIndex(element, i, lowDimDim);
560 if (dofIdxGlobalJ == dofIdx) { lowDimElems.emplace_back( std::move(element) ); break; }
561 }
562 } );
563 return lowDimElems;
564 }
565 // dof index = element index for cc schemes
566 else
567 return std::vector<Element<lowDimId>>( {ldGridGeometry.element(dofIdxGlobalJ)} );
568 } ();
569
570 // update all necessary vol vars in context
571 for (const auto& element : couplingElements)
572 {
573 // find index in coupling context
574 const auto eIdxGlobal = ldGridGeometry.elementMapper().index(element);
575 auto it = std::find(couplingElemStencil.begin(), couplingElemStencil.end(), eIdxGlobal);
576 const auto idxInContext = std::distance(couplingElemStencil.begin(), it);
577 assert(it != couplingElemStencil.end());
578
579 auto& volVars = bulkContext_.lowDimVolVars[idxInContext];
580 const auto& fvGeom = bulkContext_.lowDimFvGeometries[idxInContext];
581 // if low dim domain uses the box scheme, we have to create interpolated vol vars
582 if (lowDimUsesBox)
583 {
584 const auto elemGeom = element.geometry();
585 FacetCoupling::makeInterpolatedVolVars(volVars, ldProblem, ldSol, fvGeom, element, elemGeom, elemGeom.center());
586 }
587 // if low dim domain uses a cc scheme we can directly update the vol vars
588 else
589 volVars.update( elementSolution(element, ldSol, ldGridGeometry),
590 ldProblem,
591 element,
592 fvGeom.scv(eIdxGlobal) );
593 }
594 }
595 }
596
601 template< class BulkLocalAssembler >
602 void updateCouplingContext(BulkIdType domainI,
603 const BulkLocalAssembler& bulkLocalAssembler,
604 BulkIdType domainJ,
605 GridIndexType<bulkId> dofIdxGlobalJ,
606 const PrimaryVariables<bulkId>& priVarsJ,
607 unsigned int pvIdxJ)
608 {
609 // communicate deflected solution
610 ParentType::updateCouplingContext(domainI, bulkLocalAssembler, domainJ, dofIdxGlobalJ, priVarsJ, pvIdxJ);
611 }
612
618 template< class LowDimLocalAssembler >
619 void updateCouplingContext(LowDimIdType domainI,
620 const LowDimLocalAssembler& lowDimLocalAssembler,
621 BulkIdType domainJ,
622 GridIndexType<bulkId> dofIdxGlobalJ,
623 const PrimaryVariables<bulkId>& priVarsJ,
624 unsigned int pvIdxJ)
625 {
626 // communicate deflected solution
627 ParentType::updateCouplingContext(domainI, lowDimLocalAssembler, domainJ, dofIdxGlobalJ, priVarsJ, pvIdxJ);
628
629 // Since coupling only occurs via the fluxes, the context does not
630 // have to be updated in explicit time discretization schemes, where
631 // they are strictly evaluated on the old time level
632 if (!LowDimLocalAssembler::isImplicit())
633 return;
634
635 // skip the rest if context is empty
636 if (lowDimContext_.isSet)
637 {
638 assert(lowDimContext_.elementIdx == this->problem(lowDimId).gridGeometry().elementMapper().index(lowDimLocalAssembler.element()));
639
640 // since we use cc scheme in bulk domain: dof index = element index
641 const auto& bulkGridGeom = this->problem(bulkId).gridGeometry();
642 const auto elementJ = bulkGridGeom.element(dofIdxGlobalJ);
643
644 // update corresponding vol vars in context
645 const auto& scv = lowDimContext_.bulkFvGeometry->scv(dofIdxGlobalJ);
646 const auto elemSol = elementSolution(elementJ, this->curSol()[bulkId], bulkGridGeom);
647 (*lowDimContext_.bulkElemVolVars)[dofIdxGlobalJ].update(elemSol, this->problem(bulkId), elementJ, scv);
648
649 // update the element flux variables cache (tij might be solution-dependent)
650 if (dofIdxGlobalJ == bulkContext_.elementIdx)
651 lowDimContext_.bulkElemFluxVarsCache->update( elementJ, *lowDimContext_.bulkFvGeometry, *lowDimContext_.bulkElemVolVars);
652 else
653 lowDimContext_.bulkElemFluxVarsCache->update( this->problem(bulkId).gridGeometry().element(bulkContext_.elementIdx),
654 *lowDimContext_.bulkFvGeometry,
655 *lowDimContext_.bulkElemVolVars );
656 }
657 }
658
665 template< class LowDimLocalAssembler >
666 void updateCouplingContext(LowDimIdType domainI,
667 const LowDimLocalAssembler& lowDimLocalAssembler,
668 LowDimIdType domainJ,
669 GridIndexType<lowDimId> dofIdxGlobalJ,
670 const PrimaryVariables<lowDimId>& priVarsJ,
671 unsigned int pvIdxJ)
672 {
673 // communicate deflected solution
674 ParentType::updateCouplingContext(domainI, lowDimLocalAssembler, domainJ, dofIdxGlobalJ, priVarsJ, pvIdxJ);
675
676 // Since coupling only occurs via the fluxes, the context does not
677 // have to be updated in explicit time discretization schemes, where
678 // they are strictly evaluated on the old time level
679 if (!LowDimLocalAssembler::isImplicit())
680 return;
681
682 // skip the rest if context is empty
683 if (lowDimContext_.isSet)
684 {
685 const auto& ldSol = this->curSol()[lowDimId];
686 const auto& ldProblem = this->problem(lowDimId);
687 const auto& ldGridGeometry = this->problem(lowDimId).gridGeometry();
688
689 assert(bulkContext_.isSet);
690 assert(lowDimContext_.elementIdx == ldGridGeometry.elementMapper().index(lowDimLocalAssembler.element()));
691
692 // update the corresponding vol vars in the bulk context
693 const auto& bulkMap = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
694 const auto& couplingElementStencil = bulkMap.find(bulkContext_.elementIdx)->second.couplingElementStencil;
695 auto it = std::find(couplingElementStencil.begin(), couplingElementStencil.end(), lowDimContext_.elementIdx);
696 assert(it != couplingElementStencil.end());
697 const auto idxInContext = std::distance(couplingElementStencil.begin(), it);
698
699 auto& volVars = bulkContext_.lowDimVolVars[idxInContext];
700 const auto& fvGeom = bulkContext_.lowDimFvGeometries[idxInContext];
701 const auto& element = lowDimLocalAssembler.element();
702 // if low dim domain uses the box scheme, we have to create interpolated vol vars
703 if (lowDimUsesBox)
704 {
705 const auto elemGeom = element.geometry();
706 FacetCoupling::makeInterpolatedVolVars(volVars, ldProblem, ldSol, fvGeom, element, elemGeom, elemGeom.center());
707 }
708 // if low dim domain uses a cc scheme we can directly update the vol vars
709 else
710 volVars.update( elementSolution(element, ldSol, ldGridGeometry),
711 ldProblem,
712 element,
713 fvGeom.scv(lowDimContext_.elementIdx) );
714
715 // update the element flux variables cache (tij depend on low dim values in context)
716 const auto contextElem = this->problem(bulkId).gridGeometry().element(bulkContext_.elementIdx);
717 lowDimContext_.bulkElemFluxVarsCache->update(contextElem, *lowDimContext_.bulkFvGeometry, *lowDimContext_.bulkElemVolVars);
718 }
719 }
720
722 using ParentType::updateCoupledVariables;
723
728 template< class BulkLocalAssembler, class UpdatableFluxVarCache >
729 void updateCoupledVariables(BulkIdType domainI,
730 const BulkLocalAssembler& bulkLocalAssembler,
731 ElementVolumeVariables<bulkId>& elemVolVars,
732 UpdatableFluxVarCache& fluxVarsCache)
733 {
734 // update transmissibilities after low dim context has changed (implicit only)
735 if (BulkLocalAssembler::isImplicit())
736 fluxVarsCache.update(bulkLocalAssembler.element(),
737 bulkLocalAssembler.fvGeometry(),
738 bulkLocalAssembler.curElemVolVars());
739 }
740
745 template< class BulkLocalAssembler, class UpdatableFluxVarCache >
746 void updateCoupledVariables(BulkIdType domainI,
747 const BulkLocalAssembler& bulkLocalAssembler,
748 GridVolumeVariables<bulkId>& gridVolVars,
749 UpdatableFluxVarCache& fluxVarsCache)
750 {
751 // update transmissibilities after low dim context has changed (implicit only)
752 if (BulkLocalAssembler::isImplicit())
753 {
754 auto elemVolVars = localView(gridVolVars);
755 elemVolVars.bind(bulkLocalAssembler.element(), bulkLocalAssembler.fvGeometry(), this->curSol()[bulkId]);
756 fluxVarsCache.update(bulkLocalAssembler.element(), bulkLocalAssembler.fvGeometry(), elemVolVars);
757 }
758 }
759
761 template<std::size_t id, std::enable_if_t<(id == bulkId || id == lowDimId), int> = 0>
762 const typename CouplingMapper::template Stencil<id>&
763 getEmptyStencil(Dune::index_constant<id>) const
764 { return std::get<(id == bulkId ? 0 : 1)>(emptyStencilTuple_); }
765
766protected:
768 const BulkCouplingContext& bulkCouplingContext() const { return bulkContext_; }
769 const LowDimCouplingContext& lowDimCouplingContext() const { return lowDimContext_; }
770
772 BulkCouplingContext& bulkCouplingContext() { return bulkContext_; }
773 LowDimCouplingContext& lowDimCouplingContext() { return lowDimContext_; }
774
776 template<class BulkScvfIndices>
777 NumEqVector<bulkId> evalBulkFluxes(const Element<bulkId>& elementI,
778 const FVElementGeometry<bulkId>& fvGeometry,
779 const ElementVolumeVariables<bulkId>& elemVolVars,
780 const ElementFluxVariablesCache<bulkId>& elemFluxVarsCache,
781 const LocalResidual<bulkId>& localResidual,
782 const BulkScvfIndices& scvfIndices) const
783 {
784
785 NumEqVector<bulkId> coupledFluxes(0.0);
786 for (const auto& scvfIdx : scvfIndices)
787 coupledFluxes += localResidual.evalFlux(this->problem(bulkId),
788 elementI,
789 fvGeometry,
790 elemVolVars,
791 elemFluxVarsCache,
792 fvGeometry.scvf(scvfIdx));
793 return coupledFluxes;
794 }
795
796private:
797 std::shared_ptr<CouplingMapper> couplingMapperPtr_;
798
801 std::vector<bool> bulkElemIsCoupled_;
802 std::vector<bool> bulkScvfIsCoupled_;
803
805 using BulkStencil = typename CouplingMapper::template Stencil<bulkId>;
806 using LowDimStencil = typename CouplingMapper::template Stencil<lowDimId>;
807 std::tuple<BulkStencil, LowDimStencil> emptyStencilTuple_;
808
810 BulkCouplingContext bulkContext_;
811 LowDimCouplingContext lowDimContext_;
812};
813
814} // end namespace Dumux
815
816#endif
Defines the index types used for grid and local indices.
Element solution classes and factory functions.
The available discretization methods in Dumux.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethod::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:115
void makeInterpolatedVolVars(VolumeVariables &volVars, const Problem &problem, const SolutionVector &sol, const FVGeometry &fvGeometry, const typename FVGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const typename FVGeometry::GridGeometry::GridView::template Codim< 0 >::Entity::Geometry &elemGeom, const typename FVGeometry::GridGeometry::GridView::template Codim< 0 >::Entity::Geometry::GlobalCoordinate &pos)
Free function that allows the creation of a volume variables object interpolated to a given position ...
Definition: multidomain/facet/couplingmanager.hh:52
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
Struture to define the index types used for grid and local indices.
Definition: indextraits.hh:38
A vector of primary variables.
Definition: common/properties.hh:59
A vector of size number equations that can be used for Neumann fluxes, sources, residuals,...
Definition: common/properties.hh:61
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:69
Definition: common/properties.hh:91
Definition: common/properties.hh:150
If disabled, the volume variables are not stored (reduces memory, but is slower)
Definition: common/properties.hh:178
specifies if data on flux vars should be saved (faster, but more memory consuming)
Definition: common/properties.hh:188
The grid variables object managing variable data on the grid (volvars/fluxvars cache)
Definition: common/properties.hh:190
Definition: multidomain/couplingmanager.hh:46
void updateCouplingContext(BulkIdType domainI, const BulkLocalAssembler &bulkLocalAssembler, LowDimIdType domainJ, GridIndexType< lowDimId > dofIdxGlobalJ, const PrimaryVariables< lowDimId > &priVarsJ, unsigned int pvIdxJ)
After deflecting the solution of the lower-dimensional domain, we have to update the element volume v...
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:522
typename CouplingMapper::template Stencil< CouplingMapper::template gridId< GridView< j >::dimension >() > CouplingStencilType
types used for coupling stencils
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:156
typename MDTraits::SolutionVector SolutionVector
the type of the solution vector
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:159
void updateCoupledVariables(BulkIdType domainI, const BulkLocalAssembler &bulkLocalAssembler, ElementVolumeVariables< bulkId > &elemVolVars, UpdatableFluxVarCache &fluxVarsCache)
Update the transmissibilities in the bulk domain after the coupling context changed.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:729
const LowDimCouplingContext & lowDimCouplingContext() const
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:769
const CouplingStencilType< lowDimId > & couplingStencil(LowDimIdType domainI, const Element< lowDimId > &element, BulkIdType domainJ) const
The coupling stencil of the lower-dimensional domain with the bulk domain.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:220
LocalResidual< lowDimId >::ElementResidualVector evalCouplingResidual(LowDimIdType, const LowDimLocalAssembler &lowDimLocalAssembler, BulkIdType)
Evaluates the coupling element residual of a lower-dimensional domain element with respect to a dof i...
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:353
void bindCouplingContext(BulkIdType, const Element< bulkId > &element, const Assembler &assembler)
For the assembly of the element residual of a bulk domain element we need to prepare all variables of...
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:415
void updateCouplingContext(BulkIdType domainI, const BulkLocalAssembler &bulkLocalAssembler, BulkIdType domainJ, GridIndexType< bulkId > dofIdxGlobalJ, const PrimaryVariables< bulkId > &priVarsJ, unsigned int pvIdxJ)
Update the coupling context for a derivative bulk -> bulk. Here, we simply have to update the solutio...
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:602
BulkCouplingContext & bulkCouplingContext()
Return references to the bulk coupling contexts.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:772
bool isOnInteriorBoundary(const Element< bulkId > &element, const SubControlVolumeFace< bulkId > &scvf) const
returns true if a bulk scvf coincides with a facet element.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:243
const VolumeVariables< lowDimId > & getLowDimVolVars(const Element< bulkId > &element, const SubControlVolumeFace< bulkId > &scvf) const
returns the vol vars of a lower-dimensional element coinciding with a bulk scvf.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:250
LowDimCouplingContext & lowDimCouplingContext()
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:773
NumEqVector< bulkId > evalBulkFluxes(const Element< bulkId > &elementI, const FVElementGeometry< bulkId > &fvGeometry, const ElementVolumeVariables< bulkId > &elemVolVars, const ElementFluxVariablesCache< bulkId > &elemFluxVarsCache, const LocalResidual< bulkId > &localResidual, const BulkScvfIndices &scvfIndices) const
evaluates the bulk-facet exchange fluxes for a given facet element
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:777
const CouplingStencilType< bulkId > & couplingStencil(BulkIdType domainI, const Element< bulkId > &element, LowDimIdType domainJ) const
The coupling stencil of a given bulk domain element.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:200
const BulkCouplingContext & bulkCouplingContext() const
Return const references to the bulk coupling contexts.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:768
void updateCoupledVariables(BulkIdType domainI, const BulkLocalAssembler &bulkLocalAssembler, GridVolumeVariables< bulkId > &gridVolVars, UpdatableFluxVarCache &fluxVarsCache)
Update the transmissibilities in the bulk domain after the coupling context changed.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:746
const Element< lowDimId > getLowDimElement(const Element< bulkId > &element, const SubControlVolumeFace< bulkId > &scvf) const
returns the lower-dimensional element coinciding with a bulk scvf.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:267
const CouplingMapper::template Stencil< id > & getEmptyStencil(Dune::index_constant< id >) const
Empty stencil to be returned for elements that aren't coupled.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:763
void updateCouplingContext(LowDimIdType domainI, const LowDimLocalAssembler &lowDimLocalAssembler, BulkIdType domainJ, GridIndexType< bulkId > dofIdxGlobalJ, const PrimaryVariables< bulkId > &priVarsJ, unsigned int pvIdxJ)
After deflecting the solution of the bulk domain, we have to update the element volume variables and ...
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:619
const GridIndexType< lowDimId > getLowDimElementIndex(const Element< bulkId > &element, const SubControlVolumeFace< bulkId > &scvf) const
returns the index of the lower-dimensional element coinciding with a bulk scvf.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:277
void init(std::shared_ptr< Problem< bulkId > > bulkProblem, std::shared_ptr< Problem< lowDimId > > lowDimProblem, std::shared_ptr< CouplingMapper > couplingMapper, const SolutionVector &curSol)
Initialize the coupling manager.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:169
bool isCoupled(const Element< bulkId > &element, const SubControlVolumeFace< bulkId > &scvf) const
returns true if a bulk scvf flux depends on data in the facet domain.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:235
LocalResidual< lowDimId >::ElementResidualVector evalCouplingResidual(LowDimIdType, const LowDimLocalAssembler &lowDimLocalAssembler, BulkIdType, GridIndexType< bulkId > dofIdxGlobalJ)
Evaluates the coupling element residual of a lower-dimensional domain element with respect to a dof i...
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:340
LocalResidual< bulkId >::ElementResidualVector evalCouplingResidual(BulkIdType, const BulkLocalAssembler &bulkLocalAssembler, LowDimIdType, GridIndexType< lowDimId > dofIdxGlobalJ)
Evaluates the coupling element residual of a bulk domain element with respect to a dof in the lower-d...
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:305
void updateCouplingContext(LowDimIdType domainI, const LowDimLocalAssembler &lowDimLocalAssembler, LowDimIdType domainJ, GridIndexType< lowDimId > dofIdxGlobalJ, const PrimaryVariables< lowDimId > &priVarsJ, unsigned int pvIdxJ)
After deflecting the solution of the lower-dimensional domain has been deflected during the assembly ...
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:666
NumEqVector< lowDimId > evalSourcesFromBulk(const Element< lowDimId > &element, const FVElementGeometry< lowDimId > &fvGeometry, const ElementVolumeVariables< lowDimId > &elemVolVars, const SubControlVolume< lowDimId > &scv)
Computes the sources in a lower-dimensional element stemming from the bulk domain.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:378
void bindCouplingContext(LowDimIdType, const Element< lowDimId > &element, const Assembler &assembler)
For the assembly of the element residual of a bulk domain element we need to prepare the local views ...
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:476
Implementation for the coupling manager between two domains of dimension d and (d-1) for models consi...
Definition: multidomain/facet/couplingmanager.hh:95
Declares all properties used in Dumux.
The interface of the coupling manager for multi domain problems.