25#ifndef DUMUX_STOKES_DARCY_COUPLINGMAPPER_HH
26#define DUMUX_STOKES_DARCY_COUPLINGMAPPER_HH
29#include <unordered_map>
32#include <dune/common/exceptions.hh>
47 std::size_t flipScvfIdx;
55 template<
class CouplingManager,
class Stencils>
57 Stencils& darcyToStokesCellCenterStencils,
58 Stencils& darcyToStokesFaceStencils,
59 Stencils& stokesCellCenterToDarcyStencils,
60 Stencils& stokesFaceToDarcyStencils)
62 const auto& stokesFvGridGeometry = couplingManager.
problem(CouplingManager::stokesIdx).gridGeometry();
63 const auto& darcyFvGridGeometry = couplingManager.
problem(CouplingManager::darcyIdx).gridGeometry();
66 "The free flow domain must use the staggered discretization");
69 "The Darcy domain must use the CCTpfa discretization");
71 isCoupledDarcyScvf_.resize(darcyFvGridGeometry.numScvf(),
false);
73 auto darcyFvGeometry =
localView(darcyFvGridGeometry);
74 auto stokesFvGeometry =
localView(stokesFvGridGeometry);
75 const auto& stokesGridView = stokesFvGridGeometry.gridView();
77 for (
const auto& stokesElement : elements(stokesGridView))
79 stokesFvGeometry.bindElement(stokesElement);
81 for (
const auto& scvf : scvfs(stokesFvGeometry))
89 const auto eps = (scvf.center() - stokesElement.geometry().center()).two_norm()*1e-8;
90 auto globalPos = scvf.center(); globalPos.axpy(eps, scvf.unitOuterNormal());
91 const auto darcyElementIdx =
intersectingEntities(globalPos, darcyFvGridGeometry.boundingBoxTree());
94 if (darcyElementIdx.empty())
98 if (darcyElementIdx.size() > 1)
99 DUNE_THROW(Dune::InvalidStateException,
"Stokes face dof should only intersect with one Darcy element");
101 const auto stokesElementIdx = stokesFvGridGeometry.elementMapper().index(stokesElement);
103 const auto darcyDofIdx = darcyElementIdx[0];
105 stokesFaceToDarcyStencils[scvf.dofIndex()].push_back(darcyDofIdx);
106 stokesCellCenterToDarcyStencils[stokesElementIdx].push_back(darcyDofIdx);
108 darcyToStokesFaceStencils[darcyElementIdx[0]].push_back(scvf.dofIndex());
109 darcyToStokesCellCenterStencils[darcyElementIdx[0]].push_back(stokesElementIdx);
111 const auto& darcyElement = darcyFvGridGeometry.element(darcyElementIdx[0]);
112 darcyFvGeometry.bindElement(darcyElement);
115 for (
const auto& darcyScvf : scvfs(darcyFvGeometry))
117 const auto distance = (darcyScvf.center() - scvf.center()).two_norm();
121 isCoupledDarcyScvf_[darcyScvf.index()] =
true;
122 darcyElementToStokesElementMap_[darcyElementIdx[0]].push_back({stokesElementIdx, scvf.index(), darcyScvf.index()});
123 stokesElementToDarcyElementMap_[stokesElementIdx].push_back({darcyElementIdx[0], darcyScvf.index(), scvf.index()});
135 return isCoupledDarcyScvf_[darcyScvfIdx];
143 return darcyElementToStokesElementMap_;
151 return stokesElementToDarcyElementMap_;
155 std::unordered_map<std::size_t, std::vector<ElementMapInfo>> darcyElementToStokesElementMap_;
156 std::unordered_map<std::size_t, std::vector<ElementMapInfo>> stokesElementToDarcyElementMap_;
158 std::vector<bool> isCoupledDarcyScvf_;
The available discretization methods in Dumux.
std::vector< std::size_t > intersectingEntities(const Dune::FieldVector< ctype, dimworld > &point, const BoundingBoxTree< EntitySet > &tree, bool isCartesianGrid=false)
Compute all intersections between entities and a point.
Definition: geometry/intersectingentities.hh:100
ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: geometry/distance.hh:138
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
Coupling mapper for Stokes and Darcy domains with equal dimension.
Definition: boundary/stokesdarcy/couplingmapper.hh:42
bool isCoupledDarcyScvf(std::size_t darcyScvfIdx) const
Returns whether a Darcy scvf is coupled to the other domain.
Definition: boundary/stokesdarcy/couplingmapper.hh:133
void computeCouplingMapsAndStencils(const CouplingManager &couplingManager, Stencils &darcyToStokesCellCenterStencils, Stencils &darcyToStokesFaceStencils, Stencils &stokesCellCenterToDarcyStencils, Stencils &stokesFaceToDarcyStencils)
Main update routine.
Definition: boundary/stokesdarcy/couplingmapper.hh:56
const auto & darcyElementToStokesElementMap() const
A map that returns all Stokes elements coupled to a Darcy element.
Definition: boundary/stokesdarcy/couplingmapper.hh:141
const auto & stokesElementToDarcyElementMap() const
A map that returns all Darcy elements coupled to a Stokes element.
Definition: boundary/stokesdarcy/couplingmapper.hh:149
Definition: multidomain/couplingmanager.hh:46
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:264