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