3.6-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
37
38namespace Dumux {
39
51template<class MDTraits, class CouplingMapper, std::size_t bulkDomainId, std::size_t lowDimDomainId>
52class FacetCouplingManager<MDTraits, CouplingMapper, bulkDomainId, lowDimDomainId, DiscretizationMethods::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(!GetPropType<SubDomainTypeTag<bulkId>, Properties::GridVariables>::GridFluxVariablesCache::cachingEnabled,
89 "Grid flux variables caching currently not supported in the bulk domain of cc-facet coupling models");
90 static_assert(!GetPropType<SubDomainTypeTag<lowDimId>, Properties::GridVariables>::GridVolumeVariables::cachingEnabled,
91 "Grid volume variables caching currently not supported in the lower-dimensional domain of cc-facet coupling models");
92 static_assert(!GetPropType<SubDomainTypeTag<bulkId>, Properties::GridVariables>::GridVolumeVariables::cachingEnabled,
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 == DiscretizationMethods::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 const auto& ldGridGeometry = this->problem(lowDimId).gridGeometry();
436 auto fvGeom = localView(ldGridGeometry);
437 for (const auto lowDimElemIdx : elementStencil)
438 {
439 const auto& ldSol = Assembler::isImplicit() ? this->curSol(lowDimId) : assembler.prevSol()[lowDimId];
440 const auto& ldProblem = this->problem(lowDimId);
441
442 const auto elemJ = ldGridGeometry.element(lowDimElemIdx);
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 // evaluate variables on old/new time level depending on time disc scheme
499 const auto& bulkSol = Assembler::isImplicit() ? this->curSol(bulkId) : assembler.prevSol()[bulkId];
500
501 // then simply bind the local views of that first neighbor
502 auto bulkFvGeom = localView(bulkGridGeom).bind(bulkElem);
503 auto bulkElemVolVars = Assembler::isImplicit() ? localView(assembler.gridVariables(bulkId).curGridVolVars()).bind(bulkElem, bulkFvGeom, bulkSol)
504 : localView(assembler.gridVariables(bulkId).prevGridVolVars()).bind(bulkElem, bulkFvGeom, bulkSol);
505 auto bulkElemFluxVarsCache = localView(assembler.gridVariables(bulkId).gridFluxVarsCache()).bind(bulkElem, bulkFvGeom, bulkElemVolVars);
506
507 lowDimContext_.isSet = true;
508 lowDimContext_.bulkFvGeometry = std::make_unique< FVElementGeometry<bulkId> >( std::move(bulkFvGeom) );
509 lowDimContext_.bulkElemVolVars = std::make_unique< ElementVolumeVariables<bulkId> >( std::move(bulkElemVolVars) );
510 lowDimContext_.bulkElemFluxVarsCache = std::make_unique< ElementFluxVariablesCache<bulkId> >( std::move(bulkElemFluxVarsCache) );
511 lowDimContext_.bulkLocalResidual = std::make_unique< LocalResidual<bulkId> >(assembler.localResidual(bulkId));
512 }
513 }
514
519 template< class BulkLocalAssembler >
520 void updateCouplingContext(BulkIdType domainI,
521 const BulkLocalAssembler& bulkLocalAssembler,
522 LowDimIdType domainJ,
523 GridIndexType<lowDimId> dofIdxGlobalJ,
524 const PrimaryVariables<lowDimId>& priVarsJ,
525 unsigned int pvIdxJ)
526 {
527 // communicate deflected solution
528 ParentType::updateCouplingContext(domainI, bulkLocalAssembler, domainJ, dofIdxGlobalJ, priVarsJ, pvIdxJ);
529
530 // Since coupling only occurs via the fluxes, the context does not
531 // have to be updated in explicit time discretization schemes, where
532 // they are strictly evaluated on the old time level
533 if (!BulkLocalAssembler::isImplicit())
534 return;
535
536 // skip the rest if context is empty
537 if (bulkContext_.isSet)
538 {
539 const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
540 const auto& couplingElemStencil = map.find(bulkContext_.elementIdx)->second.couplingElementStencil;
541 const auto& ldSol = this->curSol(lowDimId);
542 const auto& ldProblem = this->problem(lowDimId);
543 const auto& ldGridGeometry = this->problem(lowDimId).gridGeometry();
544
545 // find the low-dim elements in coupling stencil, where this dof is contained in
546 const auto couplingElements = [&] ()
547 {
548 if (lowDimUsesBox)
549 {
550 std::vector< Element<lowDimId> > lowDimElems;
551 std::for_each( couplingElemStencil.begin(), couplingElemStencil.end(),
552 [&] (auto lowDimElemIdx)
553 {
554 auto element = ldGridGeometry.element(lowDimElemIdx);
555 for (int i = 0; i < element.geometry().corners(); ++i)
556 {
557 const auto dofIdx = ldGridGeometry.vertexMapper().subIndex(element, i, lowDimDim);
558 if (dofIdxGlobalJ == dofIdx) { lowDimElems.emplace_back( std::move(element) ); break; }
559 }
560 } );
561 return lowDimElems;
562 }
563 // dof index = element index for cc schemes
564 else
565 return std::vector<Element<lowDimId>>( {ldGridGeometry.element(dofIdxGlobalJ)} );
566 } ();
567
568 // update all necessary vol vars in context
569 for (const auto& element : couplingElements)
570 {
571 // find index in coupling context
572 const auto eIdxGlobal = ldGridGeometry.elementMapper().index(element);
573 auto it = std::find(couplingElemStencil.begin(), couplingElemStencil.end(), eIdxGlobal);
574 const auto idxInContext = std::distance(couplingElemStencil.begin(), it);
575 assert(it != couplingElemStencil.end());
576
577 auto& volVars = bulkContext_.lowDimVolVars[idxInContext];
578 const auto& fvGeom = bulkContext_.lowDimFvGeometries[idxInContext];
579 // if low dim domain uses the box scheme, we have to create interpolated vol vars
580 if (lowDimUsesBox)
581 {
582 const auto elemGeom = element.geometry();
583 FacetCoupling::makeInterpolatedVolVars(volVars, ldProblem, ldSol, fvGeom, element, elemGeom, elemGeom.center());
584 }
585 // if low dim domain uses a cc scheme we can directly update the vol vars
586 else
587 volVars.update( elementSolution(element, ldSol, ldGridGeometry),
588 ldProblem,
589 element,
590 fvGeom.scv(eIdxGlobal) );
591 }
592 }
593 }
594
599 template< class BulkLocalAssembler >
600 void updateCouplingContext(BulkIdType domainI,
601 const BulkLocalAssembler& bulkLocalAssembler,
602 BulkIdType domainJ,
603 GridIndexType<bulkId> dofIdxGlobalJ,
604 const PrimaryVariables<bulkId>& priVarsJ,
605 unsigned int pvIdxJ)
606 {
607 // communicate deflected solution
608 ParentType::updateCouplingContext(domainI, bulkLocalAssembler, domainJ, dofIdxGlobalJ, priVarsJ, pvIdxJ);
609 }
610
616 template< class LowDimLocalAssembler >
617 void updateCouplingContext(LowDimIdType domainI,
618 const LowDimLocalAssembler& lowDimLocalAssembler,
619 BulkIdType domainJ,
620 GridIndexType<bulkId> dofIdxGlobalJ,
621 const PrimaryVariables<bulkId>& priVarsJ,
622 unsigned int pvIdxJ)
623 {
624 // communicate deflected solution
625 ParentType::updateCouplingContext(domainI, lowDimLocalAssembler, domainJ, dofIdxGlobalJ, priVarsJ, pvIdxJ);
626
627 // Since coupling only occurs via the fluxes, the context does not
628 // have to be updated in explicit time discretization schemes, where
629 // they are strictly evaluated on the old time level
630 if (!LowDimLocalAssembler::isImplicit())
631 return;
632
633 // skip the rest if context is empty
634 if (lowDimContext_.isSet)
635 {
636 assert(lowDimContext_.elementIdx == this->problem(lowDimId).gridGeometry().elementMapper().index(lowDimLocalAssembler.element()));
637
638 // since we use cc scheme in bulk domain: dof index = element index
639 const auto& bulkGridGeom = this->problem(bulkId).gridGeometry();
640 const auto elementJ = bulkGridGeom.element(dofIdxGlobalJ);
641
642 // update corresponding vol vars in context
643 const auto& scv = lowDimContext_.bulkFvGeometry->scv(dofIdxGlobalJ);
644 const auto elemSol = elementSolution(elementJ, this->curSol(bulkId), bulkGridGeom);
645 (*lowDimContext_.bulkElemVolVars)[dofIdxGlobalJ].update(elemSol, this->problem(bulkId), elementJ, scv);
646
647 // update the element flux variables cache (tij might be solution-dependent)
648 if (dofIdxGlobalJ == bulkContext_.elementIdx)
649 lowDimContext_.bulkElemFluxVarsCache->update( elementJ, *lowDimContext_.bulkFvGeometry, *lowDimContext_.bulkElemVolVars);
650 else
651 lowDimContext_.bulkElemFluxVarsCache->update( this->problem(bulkId).gridGeometry().element(bulkContext_.elementIdx),
652 *lowDimContext_.bulkFvGeometry,
653 *lowDimContext_.bulkElemVolVars );
654 }
655 }
656
663 template< class LowDimLocalAssembler >
664 void updateCouplingContext(LowDimIdType domainI,
665 const LowDimLocalAssembler& lowDimLocalAssembler,
666 LowDimIdType domainJ,
667 GridIndexType<lowDimId> dofIdxGlobalJ,
668 const PrimaryVariables<lowDimId>& priVarsJ,
669 unsigned int pvIdxJ)
670 {
671 // communicate deflected solution
672 ParentType::updateCouplingContext(domainI, lowDimLocalAssembler, domainJ, dofIdxGlobalJ, priVarsJ, pvIdxJ);
673
674 // Since coupling only occurs via the fluxes, the context does not
675 // have to be updated in explicit time discretization schemes, where
676 // they are strictly evaluated on the old time level
677 if (!LowDimLocalAssembler::isImplicit())
678 return;
679
680 // skip the rest if context is empty
681 if (lowDimContext_.isSet)
682 {
683 const auto& ldSol = this->curSol(lowDimId);
684 const auto& ldProblem = this->problem(lowDimId);
685 const auto& ldGridGeometry = this->problem(lowDimId).gridGeometry();
686
687 assert(bulkContext_.isSet);
688 assert(lowDimContext_.elementIdx == ldGridGeometry.elementMapper().index(lowDimLocalAssembler.element()));
689
690 // update the corresponding vol vars in the bulk context
691 const auto& bulkMap = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
692 const auto& couplingElementStencil = bulkMap.find(bulkContext_.elementIdx)->second.couplingElementStencil;
693 auto it = std::find(couplingElementStencil.begin(), couplingElementStencil.end(), lowDimContext_.elementIdx);
694 assert(it != couplingElementStencil.end());
695 const auto idxInContext = std::distance(couplingElementStencil.begin(), it);
696
697 auto& volVars = bulkContext_.lowDimVolVars[idxInContext];
698 const auto& fvGeom = bulkContext_.lowDimFvGeometries[idxInContext];
699 const auto& element = lowDimLocalAssembler.element();
700 // if low dim domain uses the box scheme, we have to create interpolated vol vars
701 if (lowDimUsesBox)
702 {
703 const auto elemGeom = element.geometry();
704 FacetCoupling::makeInterpolatedVolVars(volVars, ldProblem, ldSol, fvGeom, element, elemGeom, elemGeom.center());
705 }
706 // if low dim domain uses a cc scheme we can directly update the vol vars
707 else
708 volVars.update( elementSolution(element, ldSol, ldGridGeometry),
709 ldProblem,
710 element,
711 fvGeom.scv(lowDimContext_.elementIdx) );
712
713 // update the element flux variables cache (tij depend on low dim values in context)
714 const auto contextElem = this->problem(bulkId).gridGeometry().element(bulkContext_.elementIdx);
715 lowDimContext_.bulkElemFluxVarsCache->update(contextElem, *lowDimContext_.bulkFvGeometry, *lowDimContext_.bulkElemVolVars);
716 }
717 }
718
720 using ParentType::updateCoupledVariables;
721
726 template< class BulkLocalAssembler, class UpdatableFluxVarCache >
727 void updateCoupledVariables(BulkIdType domainI,
728 const BulkLocalAssembler& bulkLocalAssembler,
729 ElementVolumeVariables<bulkId>& elemVolVars,
730 UpdatableFluxVarCache& fluxVarsCache)
731 {
732 // update transmissibilities after low dim context has changed (implicit only)
733 if (BulkLocalAssembler::isImplicit())
734 fluxVarsCache.update(bulkLocalAssembler.element(),
735 bulkLocalAssembler.fvGeometry(),
736 bulkLocalAssembler.curElemVolVars());
737 }
738
743 template< class BulkLocalAssembler, class UpdatableFluxVarCache >
744 void updateCoupledVariables(BulkIdType domainI,
745 const BulkLocalAssembler& bulkLocalAssembler,
746 GridVolumeVariables<bulkId>& gridVolVars,
747 UpdatableFluxVarCache& fluxVarsCache)
748 {
749 // update transmissibilities after low dim context has changed (implicit only)
750 if (BulkLocalAssembler::isImplicit())
751 {
752 const auto elemVolVars = localView(gridVolVars).bind(bulkLocalAssembler.element(), bulkLocalAssembler.fvGeometry(), this->curSol(bulkId));
753 fluxVarsCache.update(bulkLocalAssembler.element(), bulkLocalAssembler.fvGeometry(), elemVolVars);
754 }
755 }
756
758 template<std::size_t id, std::enable_if_t<(id == bulkId || id == lowDimId), int> = 0>
759 const typename CouplingMapper::template Stencil<id>&
760 getEmptyStencil(Dune::index_constant<id>) const
761 { return std::get<(id == bulkId ? 0 : 1)>(emptyStencilTuple_); }
762
763protected:
765 const BulkCouplingContext& bulkCouplingContext() const { return bulkContext_; }
766 const LowDimCouplingContext& lowDimCouplingContext() const { return lowDimContext_; }
767
769 BulkCouplingContext& bulkCouplingContext() { return bulkContext_; }
770 LowDimCouplingContext& lowDimCouplingContext() { return lowDimContext_; }
771
773 template<class BulkScvfIndices>
774 NumEqVector<bulkId> evalBulkFluxes(const Element<bulkId>& elementI,
775 const FVElementGeometry<bulkId>& fvGeometry,
776 const ElementVolumeVariables<bulkId>& elemVolVars,
777 const ElementFluxVariablesCache<bulkId>& elemFluxVarsCache,
778 const LocalResidual<bulkId>& localResidual,
779 const BulkScvfIndices& scvfIndices) const
780 {
781
782 NumEqVector<bulkId> coupledFluxes(0.0);
783 for (const auto& scvfIdx : scvfIndices)
784 coupledFluxes += localResidual.evalFlux(this->problem(bulkId),
785 elementI,
786 fvGeometry,
787 elemVolVars,
788 elemFluxVarsCache,
789 fvGeometry.scvf(scvfIdx));
790 return coupledFluxes;
791 }
792
793private:
794 std::shared_ptr<CouplingMapper> couplingMapperPtr_;
795
798 std::vector<bool> bulkElemIsCoupled_;
799 std::vector<bool> bulkScvfIsCoupled_;
800
802 using BulkStencil = typename CouplingMapper::template Stencil<bulkId>;
803 using LowDimStencil = typename CouplingMapper::template Stencil<lowDimId>;
804 std::tuple<BulkStencil, LowDimStencil> emptyStencilTuple_;
805
807 BulkCouplingContext bulkContext_;
808 LowDimCouplingContext lowDimContext_;
809};
810
811} // end namespace Dumux
812
813#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.
Element solution classes and factory functions.
The available discretization methods in Dumux.
static ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:294
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
constexpr Box box
Definition: method.hh:136
Structure 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:55
Definition: common/properties.hh:72
Definition: common/properties.hh:100
The grid variables object managing variable data on the grid (volvars/fluxvars cache)
Definition: common/properties.hh:121
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:60
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
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
typename CouplingMapper::template Stencil< CouplingMapper::template gridId< GridView< j >::dimension >() > CouplingStencilType
types used for coupling stencils
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:157
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
BulkCouplingContext & bulkCouplingContext()
Return references to the bulk coupling contexts.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:769
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
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:774
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:727
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
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:770
typename MDTraits::SolutionVector SolutionVector
the type of the solution vector
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:160
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
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
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 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:744
const LowDimCouplingContext & lowDimCouplingContext() const
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:766
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:520
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
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:664
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:617
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
const BulkCouplingContext & bulkCouplingContext() const
Return const references to the bulk coupling contexts.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:765
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:760
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
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
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:600
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.