24#ifndef DUMUX_CCTPFA_FACETCOUPLING_MANAGER_HH
25#define DUMUX_CCTPFA_FACETCOUPLING_MANAGER_HH
50template<
class MDTraits,
class CouplingMapper, std::
size_t bulkDomainId, std::
size_t lowDimDomainId>
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();
63 template<std::
size_t id>
using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
72 template<std::
size_t id>
using FVElementGeometry =
typename GridGeometry<id>::LocalView;
73 template<std::
size_t id>
using SubControlVolume =
typename GridGeometry<id>::SubControlVolume;
74 template<std::
size_t id>
using SubControlVolumeFace =
typename GridGeometry<id>::SubControlVolumeFace;
75 template<std::
size_t id>
using GridView =
typename GridGeometry<id>::GridView;
76 template<std::
size_t id>
using Element =
typename GridView<id>::template Codim<0>::Entity;
80 template<std::
size_t id>
using GridVolumeVariables =
typename GridVariables<id>::GridVolumeVariables;
81 template<std::
size_t id>
using ElementVolumeVariables =
typename GridVolumeVariables<id>::LocalView;
82 template<std::
size_t id>
using VolumeVariables =
typename ElementVolumeVariables<id>::VolumeVariables;
83 template<std::
size_t id>
using GridFluxVariablesCache =
typename GridVariables<id>::GridFluxVariablesCache;
84 template<std::
size_t id>
using ElementFluxVariablesCache =
typename GridFluxVariablesCache<id>::LocalView;
88 "Grid flux variables caching currently not supported in the bulk domain of cc-facet coupling models");
90 "Grid volume variables caching currently not supported in the lower-dimensional domain of cc-facet coupling models");
92 "Grid volume variables caching currently not supported in the bulk domain of cc-facet coupling models");
95 static constexpr int bulkDim = GridView<bulkDomainId>::dimension;
96 static constexpr int lowDimDim = GridView<lowDimDomainId>::dimension;
97 static constexpr auto bulkGridId = CouplingMapper::template gridId<bulkDim>();
98 static constexpr auto lowDimGridId = CouplingMapper::template gridId<lowDimDim>();
108 struct BulkCouplingContext
111 GridIndexType< bulkId > elementIdx;
112 std::vector< FVElementGeometry<lowDimId> > lowDimFvGeometries;
113 std::vector< VolumeVariables<lowDimId> > lowDimVolVars;
117 lowDimFvGeometries.clear();
118 lowDimVolVars.clear();
134 struct LowDimCouplingContext
137 GridIndexType< lowDimId > elementIdx;
138 std::unique_ptr< FVElementGeometry<bulkId> > bulkFvGeometry;
139 std::unique_ptr< ElementVolumeVariables<bulkId> > bulkElemVolVars;
140 std::unique_ptr< ElementFluxVariablesCache<bulkId> > bulkElemFluxVarsCache;
141 std::unique_ptr< LocalResidual<bulkId> > bulkLocalResidual;
145 bulkFvGeometry.reset(
nullptr);
146 bulkElemVolVars.reset(
nullptr);
147 bulkElemFluxVarsCache.reset(
nullptr);
155 template<std::
size_t i, std::
size_t j = (i == bulkId) ? lowDimId : bulkId>
156 using CouplingStencilType =
typename CouplingMapper::template Stencil< CouplingMapper::template gridId<GridView<j>::dimension>() >;
169 void init(std::shared_ptr< Problem<bulkId> > bulkProblem,
170 std::shared_ptr< Problem<lowDimId> > lowDimProblem,
171 std::shared_ptr< CouplingMapper > couplingMapper,
174 couplingMapperPtr_ = couplingMapper;
177 this->setSubProblem(bulkProblem, bulkId);
178 this->setSubProblem(lowDimProblem, lowDimId);
181 ParentType::updateSolution(curSol);
184 bulkElemIsCoupled_.assign(bulkProblem->gridGeometry().gridView().size(0),
false);
185 bulkScvfIsCoupled_.assign(bulkProblem->gridGeometry().numScvf(),
false);
187 const auto& bulkMap = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
188 for (
const auto& entry : bulkMap)
190 bulkElemIsCoupled_[entry.first] =
true;
191 for (
const auto& couplingEntry : entry.second.dofToCouplingScvfMap)
192 for (
const auto& scvfIdx : couplingEntry.second)
193 bulkScvfIsCoupled_[scvfIdx] =
true;
201 const Element<bulkId>& element,
202 LowDimIdType domainJ)
const
204 const auto eIdx = this->problem(bulkId).gridGeometry().elementMapper().index(element);
206 if (bulkElemIsCoupled_[eIdx])
208 const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
209 auto it = map.find(eIdx);
210 assert(it != map.end());
211 return it->second.couplingStencil;
214 return getEmptyStencil(lowDimId);
221 const Element<lowDimId>& element,
222 BulkIdType domainJ)
const
224 const auto eIdx = this->problem(lowDimId).gridGeometry().elementMapper().index(element);
226 const auto& map = couplingMapperPtr_->couplingMap(lowDimGridId, bulkGridId);
227 auto it = map.find(eIdx);
228 if (it != map.end())
return it->second.couplingStencil;
229 else return getEmptyStencil(bulkId);
236 const SubControlVolumeFace<bulkId>& scvf)
const
237 {
return bulkScvfIsCoupled_[scvf.index()]; }
244 const SubControlVolumeFace<bulkId>& scvf)
const
245 {
return isCoupled(element, scvf); }
251 const SubControlVolumeFace<bulkId>& scvf)
const
253 assert(bulkContext_.isSet);
255 const auto lowDimElemIdx = getLowDimElementIndex(element, scvf);
256 const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
257 const auto& s = map.find(bulkContext_.elementIdx)->second.couplingElementStencil;
258 const auto& idxInContext = std::distance( s.begin(), std::find(s.begin(), s.end(), lowDimElemIdx) );
260 assert(std::find(s.begin(), s.end(), lowDimElemIdx) != s.end());
261 return bulkContext_.lowDimVolVars[idxInContext];
268 const SubControlVolumeFace<bulkId>& scvf)
const
270 const auto lowDimElemIdx = getLowDimElementIndex(element, scvf);
271 return this->problem(lowDimId).gridGeometry().element(lowDimElemIdx);
278 const SubControlVolumeFace<bulkId>& scvf)
const
280 assert(bulkScvfIsCoupled_[scvf.index()]);
282 const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
283 const auto& couplingData = map.at(scvf.insideScvIdx());
286 auto it = std::find_if( couplingData.elementToScvfMap.begin(),
287 couplingData.elementToScvfMap.end(),
288 [&scvf] (
auto& dataPair)
290 const auto& scvfs = dataPair.second;
291 return std::find(scvfs.begin(), scvfs.end(), scvf.index()) != scvfs.end();
294 assert(it != couplingData.elementToScvfMap.end());
303 template<
class BulkLocalAssembler >
304 typename LocalResidual<bulkId>::ElementResidualVector
306 const BulkLocalAssembler& bulkLocalAssembler,
308 GridIndexType<lowDimId> dofIdxGlobalJ)
310 const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
312 assert(bulkContext_.isSet);
313 assert(bulkElemIsCoupled_[bulkContext_.elementIdx]);
314 assert(map.find(bulkContext_.elementIdx) != map.end());
315 assert(bulkContext_.elementIdx == this->problem(bulkId).gridGeometry().elementMapper().index(bulkLocalAssembler.element()));
317 typename LocalResidual<bulkId>::ElementResidualVector res(1);
319 res[0] = evalBulkFluxes(bulkLocalAssembler.element(),
320 bulkLocalAssembler.fvGeometry(),
321 bulkLocalAssembler.curElemVolVars(),
322 bulkLocalAssembler.elemFluxVarsCache(),
323 bulkLocalAssembler.localResidual(),
324 map.find(bulkContext_.elementIdx)->second.dofToCouplingScvfMap.at(dofIdxGlobalJ));
338 template<
class LowDimLocalAssembler >
339 typename LocalResidual<lowDimId>::ElementResidualVector
341 const LowDimLocalAssembler& lowDimLocalAssembler,
343 GridIndexType<bulkId> dofIdxGlobalJ)
344 {
return evalCouplingResidual(lowDimId, lowDimLocalAssembler, bulkId); }
351 template<
class LowDimLocalAssembler >
352 typename LocalResidual<lowDimId>::ElementResidualVector
356 assert(lowDimContext_.isSet);
357 assert(this->problem(lowDimId).gridGeometry().elementMapper().index(lowDimLocalAssembler.element()) == lowDimContext_.elementIdx);
361 const auto source = evalSourcesFromBulk(lowDimLocalAssembler.element(),
362 lowDimLocalAssembler.fvGeometry(),
363 lowDimLocalAssembler.curElemVolVars(),
364 *scvs(lowDimLocalAssembler.fvGeometry()).begin());
367 typename LocalResidual<lowDimId>::ElementResidualVector res(lowDimLocalAssembler.fvGeometry().numScv());
369 for (
const auto& scv : scvs(lowDimLocalAssembler.fvGeometry()))
370 res[scv.localDofIndex()] -= source;
379 const FVElementGeometry<lowDimId>& fvGeometry,
380 const ElementVolumeVariables<lowDimId>& elemVolVars,
381 const SubControlVolume<lowDimId>& scv)
384 assert(this->problem(lowDimId).gridGeometry().elementMapper().index(element) == lowDimContext_.elementIdx);
386 NumEqVector<lowDimId> sources(0.0);
388 const auto& map = couplingMapperPtr_->couplingMap(lowDimGridId, bulkGridId);
389 auto it = map.find(lowDimContext_.elementIdx);
393 assert(lowDimContext_.isSet);
394 for (
const auto& embedment : it->second.embedments)
395 sources += evalBulkFluxes(this->problem(bulkId).gridGeometry().element(embedment.first),
396 *lowDimContext_.bulkFvGeometry,
397 *lowDimContext_.bulkElemVolVars,
398 *lowDimContext_.bulkElemFluxVarsCache,
399 *lowDimContext_.bulkLocalResidual,
404 sources /= fvGeometry.numScv();
414 template<
class Assembler >
418 bulkContext_.reset();
421 const auto bulkElemIdx = this->problem(bulkId).gridGeometry().elementMapper().index(element);
422 bulkContext_.elementIdx = bulkElemIdx;
425 if (bulkElemIsCoupled_[bulkElemIdx])
427 const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
429 auto it = map.find(bulkElemIdx); assert(it != map.end());
430 const auto& elementStencil = it->second.couplingElementStencil;
431 bulkContext_.lowDimFvGeometries.reserve(elementStencil.size());
432 bulkContext_.lowDimVolVars.reserve(elementStencil.size());
434 for (
const auto lowDimElemIdx : elementStencil)
436 const auto& ldSol = Assembler::isImplicit() ? this->curSol()[lowDimId] : assembler.prevSol()[lowDimId];
437 const auto& ldProblem = this->problem(lowDimId);
438 const auto& ldGridGeometry = this->problem(lowDimId).gridGeometry();
440 const auto elemJ = ldGridGeometry.element(lowDimElemIdx);
442 fvGeom.bindElement(elemJ);
444 VolumeVariables<lowDimId> volVars;
449 const auto elemGeom = elemJ.geometry();
457 fvGeom.scv(lowDimElemIdx) );
459 bulkContext_.isSet =
true;
460 bulkContext_.lowDimFvGeometries.emplace_back( std::move(fvGeom) );
461 bulkContext_.lowDimVolVars.emplace_back( std::move(volVars) );
475 template<
class Assembler >
479 bulkContext_.reset();
480 lowDimContext_.reset();
483 const auto lowDimElemIdx = this->problem(lowDimId).gridGeometry().elementMapper().index(element);
484 lowDimContext_.elementIdx = lowDimElemIdx;
486 const auto& map = couplingMapperPtr_->couplingMap(lowDimGridId, bulkGridId);
487 auto it = map.find(lowDimElemIdx);
493 const auto& bulkGridGeom = this->problem(bulkId).gridGeometry();
494 const auto bulkElem = bulkGridGeom.element(it->second.embedments[0].first);
495 bindCouplingContext(bulkId, bulkElem, assembler);
498 auto bulkFvGeom =
localView(bulkGridGeom);
499 auto bulkElemVolVars = Assembler::isImplicit() ?
localView(assembler.gridVariables(bulkId).curGridVolVars())
500 :
localView(assembler.gridVariables(bulkId).prevGridVolVars());
501 auto bulkElemFluxVarsCache =
localView(assembler.gridVariables(bulkId).gridFluxVarsCache());
504 const auto& bulkSol = Assembler::isImplicit() ? this->curSol()[bulkId] : assembler.prevSol()[bulkId];
505 bulkFvGeom.bind(bulkElem);
506 bulkElemVolVars.bind(bulkElem, bulkFvGeom, bulkSol);
507 bulkElemFluxVarsCache.bind(bulkElem, bulkFvGeom, bulkElemVolVars);
509 lowDimContext_.isSet =
true;
510 lowDimContext_.bulkFvGeometry = std::make_unique< FVElementGeometry<bulkId> >( std::move(bulkFvGeom) );
511 lowDimContext_.bulkElemVolVars = std::make_unique< ElementVolumeVariables<bulkId> >( std::move(bulkElemVolVars) );
512 lowDimContext_.bulkElemFluxVarsCache = std::make_unique< ElementFluxVariablesCache<bulkId> >( std::move(bulkElemFluxVarsCache) );
513 lowDimContext_.bulkLocalResidual = std::make_unique< LocalResidual<bulkId> >(assembler.localResidual(bulkId));
521 template<
class BulkLocalAssembler >
523 const BulkLocalAssembler& bulkLocalAssembler,
524 LowDimIdType domainJ,
525 GridIndexType<lowDimId> dofIdxGlobalJ,
526 const PrimaryVariables<lowDimId>& priVarsJ,
530 ParentType::updateCouplingContext(domainI, bulkLocalAssembler, domainJ, dofIdxGlobalJ, priVarsJ, pvIdxJ);
535 if (!BulkLocalAssembler::isImplicit())
539 if (bulkContext_.isSet)
541 const auto& map = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
542 const auto& couplingElemStencil = map.find(bulkContext_.elementIdx)->second.couplingElementStencil;
543 const auto& ldSol = this->curSol()[lowDimId];
544 const auto& ldProblem = this->problem(lowDimId);
545 const auto& ldGridGeometry = this->problem(lowDimId).gridGeometry();
548 const auto couplingElements = [&] ()
552 std::vector< Element<lowDimId> > lowDimElems;
553 std::for_each( couplingElemStencil.begin(), couplingElemStencil.end(),
554 [&] (
auto lowDimElemIdx)
556 auto element = ldGridGeometry.element(lowDimElemIdx);
557 for (int i = 0; i < element.geometry().corners(); ++i)
559 const auto dofIdx = ldGridGeometry.vertexMapper().subIndex(element, i, lowDimDim);
560 if (dofIdxGlobalJ == dofIdx) { lowDimElems.emplace_back( std::move(element) ); break; }
567 return std::vector<Element<lowDimId>>( {ldGridGeometry.element(dofIdxGlobalJ)} );
571 for (
const auto& element : couplingElements)
574 const auto eIdxGlobal = ldGridGeometry.elementMapper().index(element);
575 auto it = std::find(couplingElemStencil.begin(), couplingElemStencil.end(), eIdxGlobal);
576 const auto idxInContext = std::distance(couplingElemStencil.begin(), it);
577 assert(it != couplingElemStencil.end());
579 auto& volVars = bulkContext_.lowDimVolVars[idxInContext];
580 const auto& fvGeom = bulkContext_.lowDimFvGeometries[idxInContext];
584 const auto elemGeom = element.geometry();
592 fvGeom.scv(eIdxGlobal) );
601 template<
class BulkLocalAssembler >
603 const BulkLocalAssembler& bulkLocalAssembler,
605 GridIndexType<bulkId> dofIdxGlobalJ,
606 const PrimaryVariables<bulkId>& priVarsJ,
610 ParentType::updateCouplingContext(domainI, bulkLocalAssembler, domainJ, dofIdxGlobalJ, priVarsJ, pvIdxJ);
618 template<
class LowDimLocalAssembler >
620 const LowDimLocalAssembler& lowDimLocalAssembler,
622 GridIndexType<bulkId> dofIdxGlobalJ,
623 const PrimaryVariables<bulkId>& priVarsJ,
627 ParentType::updateCouplingContext(domainI, lowDimLocalAssembler, domainJ, dofIdxGlobalJ, priVarsJ, pvIdxJ);
632 if (!LowDimLocalAssembler::isImplicit())
636 if (lowDimContext_.isSet)
638 assert(lowDimContext_.elementIdx == this->problem(lowDimId).gridGeometry().elementMapper().index(lowDimLocalAssembler.element()));
641 const auto& bulkGridGeom = this->problem(bulkId).gridGeometry();
642 const auto elementJ = bulkGridGeom.element(dofIdxGlobalJ);
645 const auto& scv = lowDimContext_.bulkFvGeometry->scv(dofIdxGlobalJ);
646 const auto elemSol =
elementSolution(elementJ, this->curSol()[bulkId], bulkGridGeom);
647 (*lowDimContext_.bulkElemVolVars)[dofIdxGlobalJ].update(elemSol, this->problem(bulkId), elementJ, scv);
650 if (dofIdxGlobalJ == bulkContext_.elementIdx)
651 lowDimContext_.bulkElemFluxVarsCache->update( elementJ, *lowDimContext_.bulkFvGeometry, *lowDimContext_.bulkElemVolVars);
653 lowDimContext_.bulkElemFluxVarsCache->update( this->problem(bulkId).gridGeometry().element(bulkContext_.elementIdx),
654 *lowDimContext_.bulkFvGeometry,
655 *lowDimContext_.bulkElemVolVars );
665 template<
class LowDimLocalAssembler >
667 const LowDimLocalAssembler& lowDimLocalAssembler,
668 LowDimIdType domainJ,
669 GridIndexType<lowDimId> dofIdxGlobalJ,
670 const PrimaryVariables<lowDimId>& priVarsJ,
674 ParentType::updateCouplingContext(domainI, lowDimLocalAssembler, domainJ, dofIdxGlobalJ, priVarsJ, pvIdxJ);
679 if (!LowDimLocalAssembler::isImplicit())
683 if (lowDimContext_.isSet)
685 const auto& ldSol = this->curSol()[lowDimId];
686 const auto& ldProblem = this->problem(lowDimId);
687 const auto& ldGridGeometry = this->problem(lowDimId).gridGeometry();
689 assert(bulkContext_.isSet);
690 assert(lowDimContext_.elementIdx == ldGridGeometry.elementMapper().index(lowDimLocalAssembler.element()));
693 const auto& bulkMap = couplingMapperPtr_->couplingMap(bulkGridId, lowDimGridId);
694 const auto& couplingElementStencil = bulkMap.find(bulkContext_.elementIdx)->second.couplingElementStencil;
695 auto it = std::find(couplingElementStencil.begin(), couplingElementStencil.end(), lowDimContext_.elementIdx);
696 assert(it != couplingElementStencil.end());
697 const auto idxInContext = std::distance(couplingElementStencil.begin(), it);
699 auto& volVars = bulkContext_.lowDimVolVars[idxInContext];
700 const auto& fvGeom = bulkContext_.lowDimFvGeometries[idxInContext];
701 const auto& element = lowDimLocalAssembler.element();
705 const auto elemGeom = element.geometry();
713 fvGeom.scv(lowDimContext_.elementIdx) );
716 const auto contextElem = this->problem(bulkId).gridGeometry().element(bulkContext_.elementIdx);
717 lowDimContext_.bulkElemFluxVarsCache->update(contextElem, *lowDimContext_.bulkFvGeometry, *lowDimContext_.bulkElemVolVars);
722 using ParentType::updateCoupledVariables;
728 template<
class BulkLocalAssembler,
class UpdatableFluxVarCache >
730 const BulkLocalAssembler& bulkLocalAssembler,
731 ElementVolumeVariables<bulkId>& elemVolVars,
732 UpdatableFluxVarCache& fluxVarsCache)
735 if (BulkLocalAssembler::isImplicit())
736 fluxVarsCache.update(bulkLocalAssembler.element(),
737 bulkLocalAssembler.fvGeometry(),
738 bulkLocalAssembler.curElemVolVars());
745 template<
class BulkLocalAssembler,
class UpdatableFluxVarCache >
747 const BulkLocalAssembler& bulkLocalAssembler,
748 GridVolumeVariables<bulkId>& gridVolVars,
749 UpdatableFluxVarCache& fluxVarsCache)
752 if (BulkLocalAssembler::isImplicit())
754 auto elemVolVars =
localView(gridVolVars);
755 elemVolVars.bind(bulkLocalAssembler.element(), bulkLocalAssembler.fvGeometry(), this->curSol()[bulkId]);
756 fluxVarsCache.update(bulkLocalAssembler.element(), bulkLocalAssembler.fvGeometry(), elemVolVars);
761 template<std::
size_t id, std::enable_if_t<(
id == bulkId ||
id == lowDimId),
int> = 0>
762 const typename CouplingMapper::template Stencil<id>&
764 {
return std::get<(id == bulkId ? 0 : 1)>(emptyStencilTuple_); }
776 template<
class BulkScvfIndices>
778 const FVElementGeometry<bulkId>& fvGeometry,
779 const ElementVolumeVariables<bulkId>& elemVolVars,
780 const ElementFluxVariablesCache<bulkId>& elemFluxVarsCache,
781 const LocalResidual<bulkId>& localResidual,
782 const BulkScvfIndices& scvfIndices)
const
785 NumEqVector<bulkId> coupledFluxes(0.0);
786 for (
const auto& scvfIdx : scvfIndices)
787 coupledFluxes += localResidual.evalFlux(this->problem(bulkId),
792 fvGeometry.scvf(scvfIdx));
793 return coupledFluxes;
797 std::shared_ptr<CouplingMapper> couplingMapperPtr_;
801 std::vector<bool> bulkElemIsCoupled_;
802 std::vector<bool> bulkScvfIsCoupled_;
805 using BulkStencil =
typename CouplingMapper::template Stencil<bulkId>;
806 using LowDimStencil =
typename CouplingMapper::template Stencil<lowDimId>;
807 std::tuple<BulkStencil, LowDimStencil> emptyStencilTuple_;
810 BulkCouplingContext bulkContext_;
811 LowDimCouplingContext lowDimContext_;
Defines the index types used for grid and local indices.
Element solution classes and factory functions.
The available discretization methods in Dumux.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethod::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:115
void makeInterpolatedVolVars(VolumeVariables &volVars, const Problem &problem, const SolutionVector &sol, const FVGeometry &fvGeometry, const typename FVGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const typename FVGeometry::GridGeometry::GridView::template Codim< 0 >::Entity::Geometry &elemGeom, const typename FVGeometry::GridGeometry::GridView::template Codim< 0 >::Entity::Geometry::GlobalCoordinate &pos)
Free function that allows the creation of a volume variables object interpolated to a given position ...
Definition: multidomain/facet/couplingmanager.hh:52
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
Struture to define the index types used for grid and local indices.
Definition: indextraits.hh:38
A vector of primary variables.
Definition: common/properties.hh:60
A vector of size number equations that can be used for Neumann fluxes, sources, residuals,...
Definition: common/properties.hh:62
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:70
Definition: common/properties.hh:89
Definition: common/properties.hh:113
If disabled, the volume variables are not stored (reduces memory, but is slower)
Definition: common/properties.hh:122
specifies if data on flux vars should be saved (faster, but more memory consuming)
Definition: common/properties.hh:132
The grid variables object managing variable data on the grid (volvars/fluxvars cache)
Definition: common/properties.hh:134
Definition: multidomain/couplingmanager.hh:46
void updateCouplingContext(BulkIdType domainI, const BulkLocalAssembler &bulkLocalAssembler, LowDimIdType domainJ, GridIndexType< lowDimId > dofIdxGlobalJ, const PrimaryVariables< lowDimId > &priVarsJ, unsigned int pvIdxJ)
After deflecting the solution of the lower-dimensional domain, we have to update the element volume v...
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:522
typename CouplingMapper::template Stencil< CouplingMapper::template gridId< GridView< j >::dimension >() > CouplingStencilType
types used for coupling stencils
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:156
typename MDTraits::SolutionVector SolutionVector
the type of the solution vector
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:159
void updateCoupledVariables(BulkIdType domainI, const BulkLocalAssembler &bulkLocalAssembler, ElementVolumeVariables< bulkId > &elemVolVars, UpdatableFluxVarCache &fluxVarsCache)
Update the transmissibilities in the bulk domain after the coupling context changed.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:729
const LowDimCouplingContext & lowDimCouplingContext() const
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:769
const CouplingStencilType< lowDimId > & couplingStencil(LowDimIdType domainI, const Element< lowDimId > &element, BulkIdType domainJ) const
The coupling stencil of the lower-dimensional domain with the bulk domain.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:220
LocalResidual< lowDimId >::ElementResidualVector evalCouplingResidual(LowDimIdType, const LowDimLocalAssembler &lowDimLocalAssembler, BulkIdType)
Evaluates the coupling element residual of a lower-dimensional domain element with respect to a dof i...
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:353
void bindCouplingContext(BulkIdType, const Element< bulkId > &element, const Assembler &assembler)
For the assembly of the element residual of a bulk domain element we need to prepare all variables of...
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:415
void updateCouplingContext(BulkIdType domainI, const BulkLocalAssembler &bulkLocalAssembler, BulkIdType domainJ, GridIndexType< bulkId > dofIdxGlobalJ, const PrimaryVariables< bulkId > &priVarsJ, unsigned int pvIdxJ)
Update the coupling context for a derivative bulk -> bulk. Here, we simply have to update the solutio...
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:602
BulkCouplingContext & bulkCouplingContext()
Return references to the bulk coupling contexts.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:772
bool isOnInteriorBoundary(const Element< bulkId > &element, const SubControlVolumeFace< bulkId > &scvf) const
returns true if a bulk scvf coincides with a facet element.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:243
const VolumeVariables< lowDimId > & getLowDimVolVars(const Element< bulkId > &element, const SubControlVolumeFace< bulkId > &scvf) const
returns the vol vars of a lower-dimensional element coinciding with a bulk scvf.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:250
LowDimCouplingContext & lowDimCouplingContext()
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:773
NumEqVector< bulkId > evalBulkFluxes(const Element< bulkId > &elementI, const FVElementGeometry< bulkId > &fvGeometry, const ElementVolumeVariables< bulkId > &elemVolVars, const ElementFluxVariablesCache< bulkId > &elemFluxVarsCache, const LocalResidual< bulkId > &localResidual, const BulkScvfIndices &scvfIndices) const
evaluates the bulk-facet exchange fluxes for a given facet element
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:777
const CouplingStencilType< bulkId > & couplingStencil(BulkIdType domainI, const Element< bulkId > &element, LowDimIdType domainJ) const
The coupling stencil of a given bulk domain element.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:200
const BulkCouplingContext & bulkCouplingContext() const
Return const references to the bulk coupling contexts.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:768
void updateCoupledVariables(BulkIdType domainI, const BulkLocalAssembler &bulkLocalAssembler, GridVolumeVariables< bulkId > &gridVolVars, UpdatableFluxVarCache &fluxVarsCache)
Update the transmissibilities in the bulk domain after the coupling context changed.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:746
const Element< lowDimId > getLowDimElement(const Element< bulkId > &element, const SubControlVolumeFace< bulkId > &scvf) const
returns the lower-dimensional element coinciding with a bulk scvf.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:267
const CouplingMapper::template Stencil< id > & getEmptyStencil(Dune::index_constant< id >) const
Empty stencil to be returned for elements that aren't coupled.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:763
void updateCouplingContext(LowDimIdType domainI, const LowDimLocalAssembler &lowDimLocalAssembler, BulkIdType domainJ, GridIndexType< bulkId > dofIdxGlobalJ, const PrimaryVariables< bulkId > &priVarsJ, unsigned int pvIdxJ)
After deflecting the solution of the bulk domain, we have to update the element volume variables and ...
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:619
const GridIndexType< lowDimId > getLowDimElementIndex(const Element< bulkId > &element, const SubControlVolumeFace< bulkId > &scvf) const
returns the index of the lower-dimensional element coinciding with a bulk scvf.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:277
void init(std::shared_ptr< Problem< bulkId > > bulkProblem, std::shared_ptr< Problem< lowDimId > > lowDimProblem, std::shared_ptr< CouplingMapper > couplingMapper, const SolutionVector &curSol)
Initialize the coupling manager.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:169
bool isCoupled(const Element< bulkId > &element, const SubControlVolumeFace< bulkId > &scvf) const
returns true if a bulk scvf flux depends on data in the facet domain.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:235
LocalResidual< lowDimId >::ElementResidualVector evalCouplingResidual(LowDimIdType, const LowDimLocalAssembler &lowDimLocalAssembler, BulkIdType, GridIndexType< bulkId > dofIdxGlobalJ)
Evaluates the coupling element residual of a lower-dimensional domain element with respect to a dof i...
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:340
LocalResidual< bulkId >::ElementResidualVector evalCouplingResidual(BulkIdType, const BulkLocalAssembler &bulkLocalAssembler, LowDimIdType, GridIndexType< lowDimId > dofIdxGlobalJ)
Evaluates the coupling element residual of a bulk domain element with respect to a dof in the lower-d...
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:305
void updateCouplingContext(LowDimIdType domainI, const LowDimLocalAssembler &lowDimLocalAssembler, LowDimIdType domainJ, GridIndexType< lowDimId > dofIdxGlobalJ, const PrimaryVariables< lowDimId > &priVarsJ, unsigned int pvIdxJ)
After deflecting the solution of the lower-dimensional domain has been deflected during the assembly ...
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:666
NumEqVector< lowDimId > evalSourcesFromBulk(const Element< lowDimId > &element, const FVElementGeometry< lowDimId > &fvGeometry, const ElementVolumeVariables< lowDimId > &elemVolVars, const SubControlVolume< lowDimId > &scv)
Computes the sources in a lower-dimensional element stemming from the bulk domain.
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:378
void bindCouplingContext(LowDimIdType, const Element< lowDimId > &element, const Assembler &assembler)
For the assembly of the element residual of a bulk domain element we need to prepare the local views ...
Definition: multidomain/facet/cellcentered/tpfa/couplingmanager.hh:476
Implementation for the coupling manager between two domains of dimension d and (d-1) for models consi...
Definition: multidomain/facet/couplingmanager.hh:95
Declares all properties used in Dumux.
The interface of the coupling manager for multi domain problems.