3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
multidomain/facet/box/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_BOX_FACETCOUPLING_MANAGER_HH
25#define DUMUX_BOX_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, DiscretizationMethods::Box>
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 = Dumux::NumEqVector<PrimaryVariables<id>>;
69 template<std::size_t id> using ElementBoundaryTypes = GetPropType<SubDomainTypeTag<id>, Properties::ElementBoundaryTypes>;
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 GridView<id>::IndexSet::IndexType;
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 // extract corresponding grid ids from the mapper
88 static constexpr int bulkDim = GridView<bulkDomainId>::dimension;
89 static constexpr int lowDimDim = GridView<lowDimDomainId>::dimension;
90 static constexpr auto bulkGridId = CouplingMapper::template gridId<bulkDim>();
91 static constexpr auto lowDimGridId = CouplingMapper::template gridId<lowDimDim>();
92
93 static constexpr bool lowDimUsesBox = GridGeometry<lowDimId>::discMethod == DiscretizationMethods::box;
94
101 struct BulkCouplingContext
102 {
103 bool isSet;
104 GridIndexType< bulkId > elementIdx;
105 std::vector< FVElementGeometry<lowDimId> > lowDimFvGeometries;
106 std::vector< std::vector<VolumeVariables<lowDimId>> > lowDimElemVolVars;
107
108 void reset()
109 {
110 lowDimFvGeometries.clear();
111 lowDimElemVolVars.clear();
112 isSet = false;
113 }
114 };
115
126 struct LowDimCouplingContext
127 {
128 private:
129 // For dim == dimworld in bulk domain, we know there is always going
130 // to be two bulk neighbors for a lower-dimensional element. Choose
131 // the container type depending on this
132 static constexpr int bulkDimWorld = GridView<bulkId>::dimensionworld;
133 static constexpr bool bulkIsSurfaceGrid = bulkDim != bulkDimWorld;
134
135 template< class T >
136 using ContainerType = typename std::conditional< bulkIsSurfaceGrid,
137 std::vector< T >,
138 std::array< T, 2 > >::type;
139
140 public:
141 bool isSet;
142 GridIndexType< lowDimId > elementIdx;
143 ContainerType< std::unique_ptr<FVElementGeometry<bulkId>> > bulkFvGeometries;
144 ContainerType< std::unique_ptr<ElementVolumeVariables<bulkId>> > bulkElemVolVars;
145 ContainerType< std::unique_ptr<ElementFluxVariablesCache<bulkId>> > bulkElemFluxVarsCache;
146 std::unique_ptr< LocalResidual<bulkId> > bulkLocalResidual;
147 ContainerType< ElementBoundaryTypes<bulkId> > bulkElemBcTypes;
148
149 void reset()
150 {
151 isSet = false;
152 auto resetPtr = [&] (auto&& uniquePtr) { uniquePtr.reset(); };
153 std::for_each(bulkFvGeometries.begin(), bulkFvGeometries.end(), resetPtr);
154 std::for_each(bulkElemVolVars.begin(), bulkElemVolVars.end(), resetPtr);
155 std::for_each(bulkElemFluxVarsCache.begin(), bulkElemFluxVarsCache.end(), resetPtr);
156 bulkLocalResidual.reset();
157 }
158
159 // resize containers for surface grids
160 template<bool s = bulkIsSurfaceGrid, std::enable_if_t<s, int> = 0>
161 void resize(std::size_t numEmbedments)
162 {
163 bulkFvGeometries.resize(numEmbedments);
164 bulkElemVolVars.resize(numEmbedments);
165 bulkElemFluxVarsCache.resize(numEmbedments);
166 bulkElemBcTypes.resize(numEmbedments);
167 }
168
169 // no memory reservation necessary for non-surface grids
170 template<bool s = bulkIsSurfaceGrid, std::enable_if_t<!s, int> = 0>
171 void resize(std::size_t numEmbedments)
172 {}
173 };
174
175public:
176
178 template<std::size_t i, std::size_t j = (i == bulkId) ? lowDimId : bulkId>
179 using CouplingStencilType = typename CouplingMapper::template Stencil< CouplingMapper::template gridId<GridView<j>::dimension>() >;
180
182 using SolutionVector = typename MDTraits::SolutionVector;
183
192 void init(std::shared_ptr< Problem<bulkId> > bulkProblem,
193 std::shared_ptr< Problem<lowDimId> > lowDimProblem,
194 std::shared_ptr< CouplingMapper > couplingMapper,
195 const SolutionVector& curSol)
196 {
197 couplingMapperPtr_ = couplingMapper;
198
199 // set the sub problems
200 this->setSubProblem(bulkProblem, bulkId);
201 this->setSubProblem(lowDimProblem, lowDimId);
202
203 // copy the solution vector
204 ParentType::updateSolution(curSol);
205
206 // determine all bulk elements that couple to low dim elements
207 const auto& bulkMap = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
208 bulkElemIsCoupled_.assign(bulkProblem->gridGeometry().gridView().size(0), false);
209 std::for_each( bulkMap.begin(),
210 bulkMap.end(),
211 [&] (const auto& entry) { bulkElemIsCoupled_[entry.first] = true; });
212 }
213
218 const Element<bulkId>& element,
219 LowDimIdType domainJ) const
220 {
221 const auto eIdx = this->problem(domainI).gridGeometry().elementMapper().index(element);
222
223 if (bulkElemIsCoupled_[eIdx])
224 {
225 const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
226 auto it = map.find(eIdx);
227 assert(it != map.end());
228 return it->second.couplingStencil;
229 }
230
231 return getEmptyStencil(lowDimId);
232 }
233
238 const Element<lowDimId>& element,
239 BulkIdType domainJ) const
240 {
241 const auto eIdx = this->problem(lowDimId).gridGeometry().elementMapper().index(element);
242
243 const auto& map = couplingMapperPtr_->couplingMap(lowDimGridId, bulkGridId);
244 auto it = map.find(eIdx);
245 if (it != map.end()) return it->second.couplingStencil;
246 else return getEmptyStencil(bulkId);
247 }
248
253 bool isCoupled(const Element<bulkId>& element,
254 const SubControlVolumeFace<bulkId>& scvf) const
255 { return scvf.interiorBoundary(); }
256
261 bool isOnInteriorBoundary(const Element<bulkId>& element,
262 const SubControlVolumeFace<bulkId>& scvf) const
263 { return isCoupled(element, scvf); }
264
268 const VolumeVariables<lowDimId>& getLowDimVolVars(const Element<bulkId>& element,
269 const SubControlVolumeFace<bulkId>& scvf) const
270 {
271 const auto eIdx = this->problem(bulkId).gridGeometry().elementMapper().index(element);
272 assert(bulkContext_.isSet);
273 assert(bulkElemIsCoupled_[eIdx]);
274
275 const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
276 const auto& couplingData = map.find(eIdx)->second;
277
278 // search the low dim element idx this scvf is embedded in
279 // and find out the local index of the scv it couples to
280 unsigned int coupledScvIdx;
281 auto it = std::find_if( couplingData.elementToScvfMap.begin(),
282 couplingData.elementToScvfMap.end(),
283 [&] (auto& dataPair)
284 {
285 const auto& scvfList = dataPair.second;
286 auto it = std::find(scvfList.begin(), scvfList.end(), scvf.index());
287 coupledScvIdx = std::distance(scvfList.begin(), it);
288 return it != scvfList.end();
289 } );
290
291 assert(it != couplingData.elementToScvfMap.end());
292 const auto lowDimElemIdx = it->first;
293 const auto& s = map.find(bulkContext_.elementIdx)->second.couplingElementStencil;
294 const auto& idxInContext = std::distance( s.begin(), std::find(s.begin(), s.end(), lowDimElemIdx) );
295 assert(std::find(s.begin(), s.end(), lowDimElemIdx) != s.end());
296
297 if (lowDimUsesBox)
298 return bulkContext_.lowDimElemVolVars[idxInContext][coupledScvIdx];
299 else
300 return bulkContext_.lowDimElemVolVars[idxInContext][0];
301 }
302
306 const Element<lowDimId> getLowDimElement(const Element<bulkId>& element,
307 const SubControlVolumeFace<bulkId>& scvf) const
308 {
309 const auto lowDimElemIdx = getLowDimElementIndex(element, scvf);
310 return this->problem(lowDimId).gridGeometry().element(lowDimElemIdx);
311 }
312
316 const GridIndexType<lowDimId> getLowDimElementIndex(const Element<bulkId>& element,
317 const SubControlVolumeFace<bulkId>& scvf) const
318 {
319 const auto eIdx = this->problem(bulkId).gridGeometry().elementMapper().index(element);
320
321 assert(bulkElemIsCoupled_[eIdx]);
322 const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
323 const auto& couplingData = map.at(eIdx);
324
325 // search the low dim element idx this scvf is embedded in
326 auto it = std::find_if( couplingData.elementToScvfMap.begin(),
327 couplingData.elementToScvfMap.end(),
328 [&] (auto& dataPair)
329 {
330 const auto& scvfList = dataPair.second;
331 auto it = std::find(scvfList.begin(), scvfList.end(), scvf.index());
332 return it != scvfList.end();
333 } );
334
335 assert(it != couplingData.elementToScvfMap.end());
336 return it->first;
337 }
338
345 template< class BulkLocalAssembler >
346 typename LocalResidual<bulkId>::ElementResidualVector
348 const BulkLocalAssembler& bulkLocalAssembler,
349 LowDimIdType,
350 GridIndexType<lowDimId> dofIdxGlobalJ)
351 {
352 const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
353
354 assert(bulkContext_.isSet);
355 assert(bulkElemIsCoupled_[bulkContext_.elementIdx]);
356 assert(map.find(bulkContext_.elementIdx) != map.end());
357 assert(bulkContext_.elementIdx == this->problem(bulkId).gridGeometry().elementMapper().index(bulkLocalAssembler.element()));
358
359 const auto& fvGeometry = bulkLocalAssembler.fvGeometry();
360 typename LocalResidual<bulkId>::ElementResidualVector res(fvGeometry.numScv());
361 res = 0.0;
362
363 // compute fluxes across the coupling scvfs
364 const auto& couplingScvfs = map.find(bulkContext_.elementIdx)->second.dofToCouplingScvfMap.at(dofIdxGlobalJ);
365 for (auto scvfIdx : couplingScvfs)
366 {
367 const auto& scvf = fvGeometry.scvf(scvfIdx);
368 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
369 res[insideScv.localDofIndex()] += evalBulkFluxes_( bulkLocalAssembler.element(),
370 bulkLocalAssembler.fvGeometry(),
371 bulkLocalAssembler.curElemVolVars(),
372 bulkLocalAssembler.elemBcTypes(),
373 bulkLocalAssembler.elemFluxVarsCache(),
374 bulkLocalAssembler.localResidual(),
375 std::array<GridIndexType<bulkId>, 1>({scvfIdx}) );
376 }
377
378 return res;
379 }
380
387 template< class LowDimLocalAssembler >
388 typename LocalResidual<lowDimId>::ElementResidualVector
390 const LowDimLocalAssembler& lowDimLocalAssembler,
391 BulkIdType,
392 GridIndexType<bulkId> dofIdxGlobalJ)
393 {
394 // make sure this is called for the element for which the context was set
395 assert(lowDimContext_.isSet);
396 assert(this->problem(lowDimId).gridGeometry().elementMapper().index(lowDimLocalAssembler.element()) == lowDimContext_.elementIdx);
397
398 // evaluate sources for all scvs in lower-dimensional element
399 typename LocalResidual<lowDimId>::ElementResidualVector res(lowDimLocalAssembler.fvGeometry().numScv());
400 res = 0.0;
401 for (const auto& scv : scvs(lowDimLocalAssembler.fvGeometry()))
402 res[scv.localDofIndex()] -= evalSourcesFromBulk(lowDimLocalAssembler.element(),
403 lowDimLocalAssembler.fvGeometry(),
404 lowDimLocalAssembler.curElemVolVars(),
405 scv);
406 return res;
407 }
408
412 NumEqVector<lowDimId> evalSourcesFromBulk(const Element<lowDimId>& element,
413 const FVElementGeometry<lowDimId>& fvGeometry,
414 const ElementVolumeVariables<lowDimId>& elemVolVars,
415 const SubControlVolume<lowDimId>& scv)
416 {
417 // make sure the this is called for the element of the context
418 assert(this->problem(lowDimId).gridGeometry().elementMapper().index(element) == lowDimContext_.elementIdx);
419
420 NumEqVector<lowDimId> sources(0.0);
421 const auto& map = couplingMapperPtr_->couplingMap(lowDimGridId, bulkGridId);
422 auto it = map.find(lowDimContext_.elementIdx);
423 if (it == map.end())
424 return sources;
425
426 assert(lowDimContext_.isSet);
427 for (unsigned int i = 0; i < it->second.embedments.size(); ++i)
428 {
429 const auto& embedment = it->second.embedments[i];
430
431 // list of scvfs in the bulk domain whose fluxes enter this scv
432 // if low dim domain uses tpfa, this is all scvfs lying on this element
433 // if it uses box, it is the one scvf coinciding with the given scv
434 const auto& coincidingScvfs = embedment.second;
435 const auto& scvfList = lowDimUsesBox ? std::vector<GridIndexType<lowDimId>>{ coincidingScvfs[scv.localDofIndex()] }
436 : coincidingScvfs;
437
438 sources += evalBulkFluxes_(this->problem(bulkId).gridGeometry().element(embedment.first),
439 *(lowDimContext_.bulkFvGeometries[i]),
440 *(lowDimContext_.bulkElemVolVars[i]),
441 lowDimContext_.bulkElemBcTypes[i],
442 *(lowDimContext_.bulkElemFluxVarsCache[i]),
443 *lowDimContext_.bulkLocalResidual,
444 scvfList);
445 }
446
447 return sources;
448 }
449
455 template< class Assembler >
456 void bindCouplingContext(BulkIdType, const Element<bulkId>& element, const Assembler& assembler)
457 {
458 // clear context
459 bulkContext_.reset();
460
461 // set index in context in any case
462 const auto bulkElemIdx = this->problem(bulkId).gridGeometry().elementMapper().index(element);
463 bulkContext_.elementIdx = bulkElemIdx;
464
465 // if element is coupled, actually set the context
466 if (bulkElemIsCoupled_[bulkElemIdx])
467 {
468 const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
469
470 auto it = map.find(bulkElemIdx); assert(it != map.end());
471 const auto& elementStencil = it->second.couplingElementStencil;
472 bulkContext_.lowDimFvGeometries.reserve(elementStencil.size());
473 bulkContext_.lowDimElemVolVars.reserve(elementStencil.size());
474
475 // local view on the bulk fv geometry
476 const auto bulkFvGeometry = localView(this->problem(bulkId).gridGeometry()).bindElement(element);
477
478 const auto& ldGridGeometry = this->problem(lowDimId).gridGeometry();
479 auto fvGeom = localView(ldGridGeometry);
480 for (const auto lowDimElemIdx : elementStencil)
481 {
482 const auto& ldSol = Assembler::isImplicit() ? this->curSol(lowDimId) : assembler.prevSol()[lowDimId];
483 const auto elemJ = ldGridGeometry.element(lowDimElemIdx);
484 fvGeom.bindElement(elemJ);
485
486 std::vector<VolumeVariables<lowDimId>> elemVolVars(fvGeom.numScv());
487 const auto& coupledScvfIndices = it->second.elementToScvfMap.at(lowDimElemIdx);
488 makeCoupledLowDimElemVolVars_(element, bulkFvGeometry, elemJ, fvGeom, ldSol, coupledScvfIndices, elemVolVars);
489
490 bulkContext_.isSet = true;
491 bulkContext_.lowDimFvGeometries.emplace_back( std::move(fvGeom) );
492 bulkContext_.lowDimElemVolVars.emplace_back( std::move(elemVolVars) );
493 }
494 }
495 }
496
506 template< class Assembler >
507 void bindCouplingContext(LowDimIdType, const Element<lowDimId>& element, const Assembler& assembler)
508 {
509 // reset contexts
510 bulkContext_.reset();
511 lowDimContext_.reset();
512
513 // set index in context in any case
514 const auto lowDimElemIdx = this->problem(lowDimId).gridGeometry().elementMapper().index(element);
515 lowDimContext_.elementIdx = lowDimElemIdx;
516
517 const auto& map = couplingMapperPtr_->couplingMap(lowDimGridId, bulkGridId);
518 auto it = map.find(lowDimElemIdx);
519
520 // if element is coupled, actually set the context
521 if (it != map.end())
522 {
523 // first bind the low dim context for the first neighboring bulk element
524 // this will set up the lower-dimensional data on this current lowdim element
525 // which will be enough to compute the fluxes in all neighboring elements
526 const auto& bulkGridGeom = this->problem(bulkId).gridGeometry();
527 const auto bulkElem = bulkGridGeom.element(it->second.embedments[0].first);
528 bindCouplingContext(bulkId, bulkElem, assembler);
529
530 // reserve space in the context and bind local views of neighbors
531 const auto& embedments = it->second.embedments;
532 const auto numEmbedments = embedments.size();
533 lowDimContext_.resize(numEmbedments);
534
535 auto bulkFvGeom = localView(bulkGridGeom);
536 auto bulkElemVolVars = localView(assembler.gridVariables(bulkId).curGridVolVars());
537 auto bulkElemFluxVarsCache = localView(assembler.gridVariables(bulkId).gridFluxVarsCache());
538
539 for (unsigned int i = 0; i < numEmbedments; ++i)
540 {
541 const auto& bulkSol = Assembler::isImplicit() ? this->curSol(bulkId) : assembler.prevSol()[bulkId];
542 const auto curBulkElem = bulkGridGeom.element(embedments[i].first);
543
544 bulkFvGeom.bind(curBulkElem);
545 bulkElemVolVars.bind(curBulkElem, bulkFvGeom, bulkSol);
546 bulkElemFluxVarsCache.bind(curBulkElem, bulkFvGeom, bulkElemVolVars);
547
548 lowDimContext_.isSet = true;
549 lowDimContext_.bulkElemBcTypes[i].update(this->problem(bulkId), curBulkElem, bulkFvGeom);
550 lowDimContext_.bulkFvGeometries[i] = std::make_unique< FVElementGeometry<bulkId> >( std::move(bulkFvGeom) );
551 lowDimContext_.bulkElemVolVars[i] = std::make_unique< ElementVolumeVariables<bulkId> >( std::move(bulkElemVolVars) );
552 lowDimContext_.bulkElemFluxVarsCache[i] = std::make_unique< ElementFluxVariablesCache<bulkId> >( std::move(bulkElemFluxVarsCache) );
553 }
554
555 // finally, set the local residual
556 lowDimContext_.bulkLocalResidual = std::make_unique< LocalResidual<bulkId> >(assembler.localResidual(bulkId));
557 }
558 }
559
564 template< class BulkLocalAssembler >
565 void updateCouplingContext(BulkIdType domainI,
566 const BulkLocalAssembler& bulkLocalAssembler,
567 LowDimIdType domainJ,
568 GridIndexType<lowDimId> dofIdxGlobalJ,
569 const PrimaryVariables<lowDimId>& priVarsJ,
570 unsigned int pvIdxJ)
571 {
572 // communicate deflected solution
573 ParentType::updateCouplingContext(domainI, bulkLocalAssembler, domainJ, dofIdxGlobalJ, priVarsJ, pvIdxJ);
574
575 // Since coupling only occurs via the fluxes, the context does not
576 // have to be updated in explicit time discretization schemes, where
577 // they are strictly evaluated on the old time level
578 if (!BulkLocalAssembler::isImplicit())
579 return;
580
581 // skip the rest if context is empty
582 if (bulkContext_.isSet)
583 {
584 const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
585 const auto& couplingEntry = map.at(bulkContext_.elementIdx);
586 const auto& couplingElemStencil = couplingEntry.couplingElementStencil;
587 const auto& ldGridGeometry = this->problem(lowDimId).gridGeometry();
588
589 // find the low-dim elements in coupling stencil, where this dof is contained in
590 const auto couplingElements = [&] ()
591 {
592 if (lowDimUsesBox)
593 {
594 std::vector< Element<lowDimId> > lowDimElems;
595 std::for_each( couplingElemStencil.begin(), couplingElemStencil.end(),
596 [&] (auto lowDimElemIdx)
597 {
598 auto element = ldGridGeometry.element(lowDimElemIdx);
599 for (unsigned int i = 0; i < element.geometry().corners(); ++i)
600 {
601 const auto dofIdx = ldGridGeometry.vertexMapper().subIndex(element, i, lowDimDim);
602 if (dofIdxGlobalJ == dofIdx) { lowDimElems.emplace_back( std::move(element) ); break; }
603 }
604 } );
605 return lowDimElems;
606 }
607 // dof index = element index for cc schemes
608 else
609 return std::vector<Element<lowDimId>>( {ldGridGeometry.element(dofIdxGlobalJ)} );
610 } ();
611
612 // update all necessary vol vars in context
613 for (const auto& element : couplingElements)
614 {
615 // find index in coupling context
616 const auto eIdxGlobal = ldGridGeometry.elementMapper().index(element);
617 auto it = std::find(couplingElemStencil.begin(), couplingElemStencil.end(), eIdxGlobal);
618 const auto idxInContext = std::distance(couplingElemStencil.begin(), it);
619 assert(it != couplingElemStencil.end());
620
621 auto& elemVolVars = bulkContext_.lowDimElemVolVars[idxInContext];
622 const auto& fvGeom = bulkContext_.lowDimFvGeometries[idxInContext];
623 const auto& coupledScvfIndices = couplingEntry.elementToScvfMap.at(eIdxGlobal);
624 makeCoupledLowDimElemVolVars_(bulkLocalAssembler.element(), bulkLocalAssembler.fvGeometry(),
625 element, fvGeom, this->curSol(lowDimId), coupledScvfIndices, elemVolVars);
626 }
627 }
628 }
629
634 template< class BulkLocalAssembler >
635 void updateCouplingContext(BulkIdType domainI,
636 const BulkLocalAssembler& bulkLocalAssembler,
637 BulkIdType domainJ,
638 GridIndexType<bulkId> dofIdxGlobalJ,
639 const PrimaryVariables<bulkId>& priVarsJ,
640 unsigned int pvIdxJ)
641 {
642 // communicate deflected solution
643 ParentType::updateCouplingContext(domainI, bulkLocalAssembler, domainJ, dofIdxGlobalJ, priVarsJ, pvIdxJ);
644 }
645
651 template< class LowDimLocalAssembler >
652 void updateCouplingContext(LowDimIdType domainI,
653 const LowDimLocalAssembler& lowDimLocalAssembler,
654 BulkIdType domainJ,
655 GridIndexType<bulkId> dofIdxGlobalJ,
656 const PrimaryVariables<bulkId>& priVarsJ,
657 unsigned int pvIdxJ)
658 {
659 // communicate deflected solution
660 ParentType::updateCouplingContext(domainI, lowDimLocalAssembler, domainJ, dofIdxGlobalJ, priVarsJ, pvIdxJ);
661
662 // Since coupling only occurs via the fluxes, the context does not
663 // have to be updated in explicit time discretization schemes, where
664 // they are strictly evaluated on the old time level
665 if (!LowDimLocalAssembler::isImplicit())
666 return;
667
668 // skip the rest if context is empty
669 if (lowDimContext_.isSet)
670 {
671 assert(lowDimContext_.elementIdx == this->problem(lowDimId).gridGeometry().elementMapper().index(lowDimLocalAssembler.element()));
672
673 const auto& map = couplingMapperPtr_->couplingMap(lowDimGridId, bulkGridId);
674 auto it = map.find(lowDimContext_.elementIdx);
675
676 assert(it != map.end());
677 const auto& embedments = it->second.embedments;
678 const auto& bulkGridGeometry = this->problem(bulkId).gridGeometry();
679
680 // TODO more efficient (only one dof update per embedment)
681 // update the elem volvars in context of those elements that contain globalJ
682 unsigned int embedmentIdx = 0;
683 for (const auto& embedment : embedments)
684 {
685 const auto elementJ = bulkGridGeometry.element(embedment.first);
686
687 for (unsigned int i = 0; i < elementJ.subEntities(bulkDim); ++i)
688 {
689 const auto dofIdx = bulkGridGeometry.vertexMapper().subIndex(elementJ, i, bulkDim);
690 if (dofIdx == dofIdxGlobalJ)
691 {
692 // element contains the deflected dof
693 const auto& fvGeom = *lowDimContext_.bulkFvGeometries[embedmentIdx];
694 (*lowDimContext_.bulkElemVolVars[embedmentIdx]).bindElement(elementJ, fvGeom, this->curSol(bulkId));
695 }
696 }
697
698 // keep track of embedments
699 embedmentIdx++;
700 }
701 }
702 }
703
710 template< class LowDimLocalAssembler >
711 void updateCouplingContext(LowDimIdType domainI,
712 const LowDimLocalAssembler& lowDimLocalAssembler,
713 LowDimIdType domainJ,
714 GridIndexType<lowDimId> dofIdxGlobalJ,
715 const PrimaryVariables<lowDimId>& priVarsJ,
716 unsigned int pvIdxJ)
717 {
718 // communicate deflected solution
719 ParentType::updateCouplingContext(domainI, lowDimLocalAssembler, domainJ, dofIdxGlobalJ, priVarsJ, pvIdxJ);
720
721 // Since coupling only occurs via the fluxes, the context does not
722 // have to be updated in explicit time discretization schemes, where
723 // they are strictly evaluated on the old time level
724 if (!LowDimLocalAssembler::isImplicit())
725 return;
726
727 // skip the rest if context is empty
728 if (lowDimContext_.isSet)
729 {
730 assert(bulkContext_.isSet);
731 assert(lowDimContext_.elementIdx == this->problem(lowDimId).gridGeometry().elementMapper().index(lowDimLocalAssembler.element()));
732
733 // update the corresponding vol vars in the bulk context
734 const auto& bulkMap = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
735 const auto& bulkCouplingEntry = bulkMap.at(bulkContext_.elementIdx);
736 const auto& couplingElementStencil = bulkCouplingEntry.couplingElementStencil;
737 auto it = std::find(couplingElementStencil.begin(), couplingElementStencil.end(), lowDimContext_.elementIdx);
738 assert(it != couplingElementStencil.end());
739 const auto idxInContext = std::distance(couplingElementStencil.begin(), it);
740
741 // get neighboring bulk element from the bulk context (is the same elemet as first entry in low dim context)
742 const auto& bulkElement = this->problem(bulkId).gridGeometry().element(bulkContext_.elementIdx);
743 const auto& bulkFvGeometry = *lowDimContext_.bulkFvGeometries[0];
744
745 auto& elemVolVars = bulkContext_.lowDimElemVolVars[idxInContext];
746 const auto& fvGeom = bulkContext_.lowDimFvGeometries[idxInContext];
747 const auto& element = lowDimLocalAssembler.element();
748 const auto& scvfIndices = bulkCouplingEntry.elementToScvfMap.at(lowDimContext_.elementIdx);
749
750 makeCoupledLowDimElemVolVars_(bulkElement, bulkFvGeometry, element, fvGeom,
751 this->curSol(lowDimId), scvfIndices, elemVolVars);
752 }
753 }
754
756 template<std::size_t id, std::enable_if_t<(id == bulkId || id == lowDimId), int> = 0>
757 const typename CouplingMapper::template Stencil<id>&
758 getEmptyStencil(Dune::index_constant<id>) const
759 { return std::get<(id == bulkId ? 0 : 1)>(emptyStencilTuple_); }
760
761private:
762
767 template<class LowDimSolution, class ScvfIdxStorage>
768 void makeCoupledLowDimElemVolVars_(const Element<bulkId>& bulkElement,
769 const FVElementGeometry<bulkId>& bulkFvGeometry,
770 const Element<lowDimId>& lowDimElement,
771 const FVElementGeometry<lowDimId>& lowDimFvGeometry,
772 const LowDimSolution& lowDimSol,
773 const ScvfIdxStorage& coupledBulkScvfIndices,
774 std::vector<VolumeVariables<lowDimId>>& elemVolVars)
775 {
776 const auto lowDimElemSol = elementSolution(lowDimElement, lowDimSol, lowDimFvGeometry.gridGeometry());
777
778 // for tpfa, make only one volume variables object for the element
779 // for the update, use the first (and only - for tpfa) scv of the element
780 if (!lowDimUsesBox)
781 elemVolVars[0].update(lowDimElemSol, this->problem(lowDimId), lowDimElement, *scvs(lowDimFvGeometry).begin());
782 // for box, make vol vars either:
783 // - in each scv center, which geometrically coincides with the scvf integration point (Neumann coupling)
784 // - at the vertex position (Dirichlet coupling)
785 else
786 {
787 // recall that the scvfs in the coupling map were ordered such
788 // that the i-th scvf in the map couples to the i-th lowdim scv
789 for (unsigned int i = 0; i < coupledBulkScvfIndices.size(); ++i)
790 {
791 const auto& scvf = bulkFvGeometry.scvf(coupledBulkScvfIndices[i]);
792 const auto bcTypes = this->problem(bulkId).interiorBoundaryTypes(bulkElement, scvf);
793 if (bcTypes.hasOnlyNeumann())
794 FacetCoupling::makeInterpolatedVolVars(elemVolVars[i], this->problem(lowDimId), lowDimSol,
795 lowDimFvGeometry, lowDimElement, lowDimElement.geometry(),
796 scvf.ipGlobal());
797 else
798 elemVolVars[i].update( lowDimElemSol, this->problem(lowDimId), lowDimElement, lowDimFvGeometry.scv(i) );
799 }
800 }
801 }
802
804 template<class BulkScvfIndices>
805 NumEqVector<bulkId> evalBulkFluxes_(const Element<bulkId>& elementI,
806 const FVElementGeometry<bulkId>& fvGeometry,
807 const ElementVolumeVariables<bulkId>& elemVolVars,
808 const ElementBoundaryTypes<bulkId>& elemBcTypes,
809 const ElementFluxVariablesCache<bulkId>& elemFluxVarsCache,
810 const LocalResidual<bulkId>& localResidual,
811 const BulkScvfIndices& scvfIndices) const
812 {
813
814 NumEqVector<bulkId> coupledFluxes(0.0);
815 for (const auto& scvfIdx : scvfIndices)
816 coupledFluxes += localResidual.evalFlux(this->problem(bulkId),
817 elementI,
818 fvGeometry,
819 elemVolVars,
820 elemBcTypes,
821 elemFluxVarsCache,
822 fvGeometry.scvf(scvfIdx));
823 return coupledFluxes;
824 }
825
826 std::shared_ptr<CouplingMapper> couplingMapperPtr_;
827
830 std::vector<bool> bulkElemIsCoupled_;
831
833 using BulkStencil = typename CouplingMapper::template Stencil<bulkId>;
834 using LowDimStencil = typename CouplingMapper::template Stencil<lowDimId>;
835 std::tuple<BulkStencil, LowDimStencil> emptyStencilTuple_;
836
838 BulkCouplingContext bulkContext_;
839 LowDimCouplingContext lowDimContext_;
840};
841
842} // end namespace Dumux
843
844#endif
A helper to deduce a vector with the same size as numbers of equations.
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:292
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
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
constexpr Box box
Definition: method.hh:139
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
Stores the boundary types on an element.
Definition: common/properties.hh:99
Definition: common/properties.hh:102
The grid variables object managing variable data on the grid (volvars/fluxvars cache)
Definition: common/properties.hh:123
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:60
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/box/couplingmanager.hh:456
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/box/couplingmanager.hh:565
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/box/couplingmanager.hh:507
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/box/couplingmanager.hh:261
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/box/couplingmanager.hh:758
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/box/couplingmanager.hh:192
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/box/couplingmanager.hh:237
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/box/couplingmanager.hh:316
typename CouplingMapper::template Stencil< CouplingMapper::template gridId< GridView< j >::dimension >() > CouplingStencilType
types used for coupling stencils
Definition: multidomain/facet/box/couplingmanager.hh:179
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/box/couplingmanager.hh:347
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/box/couplingmanager.hh:268
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/box/couplingmanager.hh:635
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 of t...
Definition: multidomain/facet/box/couplingmanager.hh:652
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/box/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/box/couplingmanager.hh:711
typename MDTraits::SolutionVector SolutionVector
the type of the solution vector
Definition: multidomain/facet/box/couplingmanager.hh:182
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/box/couplingmanager.hh:389
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/box/couplingmanager.hh:253
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/box/couplingmanager.hh:412
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/box/couplingmanager.hh:217
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.