13#ifndef DUMUX_STOKES_DARCY_COUPLINGMAPPER_HH
14#define DUMUX_STOKES_DARCY_COUPLINGMAPPER_HH
17#include <unordered_map>
20#include <dune/common/exceptions.hh>
35 std::size_t flipScvfIdx;
43 template<
class CouplingManager,
class Stencils>
45 Stencils& darcyToStokesCellCenterStencils,
46 Stencils& darcyToStokesFaceStencils,
47 Stencils& stokesCellCenterToDarcyStencils,
48 Stencils& stokesFaceToDarcyStencils)
50 const auto& stokesFvGridGeometry = couplingManager.
problem(CouplingManager::stokesIdx).gridGeometry();
51 const auto& darcyFvGridGeometry = couplingManager.
problem(CouplingManager::darcyIdx).gridGeometry();
54 "The free flow domain must use the staggered discretization");
57 "The Darcy domain must use the CCTpfa discretization");
59 isCoupledDarcyScvf_.resize(darcyFvGridGeometry.numScvf(),
false);
61 const auto& stokesGridView = stokesFvGridGeometry.gridView();
62 auto stokesFvGeometry =
localView(stokesFvGridGeometry);
63 for (
const auto& stokesElement : elements(stokesGridView))
65 stokesFvGeometry.bindElement(stokesElement);
67 for (
const auto& scvf : scvfs(stokesFvGeometry))
75 const auto eps = (scvf.center() - stokesElement.geometry().center()).two_norm()*1e-8;
76 auto globalPos = scvf.center(); globalPos.axpy(eps, scvf.unitOuterNormal());
77 const auto darcyElementIdx =
intersectingEntities(globalPos, darcyFvGridGeometry.boundingBoxTree());
80 if (darcyElementIdx.empty())
84 if (darcyElementIdx.size() > 1)
85 DUNE_THROW(Dune::InvalidStateException,
"Stokes face dof should only intersect with one Darcy element");
87 const auto stokesElementIdx = stokesFvGridGeometry.elementMapper().index(stokesElement);
89 const auto darcyDofIdx = darcyElementIdx[0];
91 stokesFaceToDarcyStencils[scvf.dofIndex()].push_back(darcyDofIdx);
92 stokesCellCenterToDarcyStencils[stokesElementIdx].push_back(darcyDofIdx);
94 darcyToStokesFaceStencils[darcyElementIdx[0]].push_back(scvf.dofIndex());
95 darcyToStokesCellCenterStencils[darcyElementIdx[0]].push_back(stokesElementIdx);
97 const auto& darcyElement = darcyFvGridGeometry.element(darcyElementIdx[0]);
98 const auto darcyFvGeometry =
localView(darcyFvGridGeometry).bindElement(darcyElement);
101 for (
const auto& darcyScvf : scvfs(darcyFvGeometry))
103 const auto distance = (darcyScvf.center() - scvf.center()).two_norm();
107 isCoupledDarcyScvf_[darcyScvf.index()] =
true;
108 darcyElementToStokesElementMap_[darcyElementIdx[0]].push_back({stokesElementIdx, scvf.index(), darcyScvf.index()});
109 stokesElementToDarcyElementMap_[stokesElementIdx].push_back({darcyElementIdx[0], darcyScvf.index(), scvf.index()});
121 return isCoupledDarcyScvf_[darcyScvfIdx];
129 return darcyElementToStokesElementMap_;
137 return stokesElementToDarcyElementMap_;
141 std::unordered_map<std::size_t, std::vector<ElementMapInfo>> darcyElementToStokesElementMap_;
142 std::unordered_map<std::size_t, std::vector<ElementMapInfo>> stokesElementToDarcyElementMap_;
144 std::vector<bool> isCoupledDarcyScvf_;
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:48
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:309
Coupling mapper for Stokes and Darcy domains with equal dimension.
Definition: boundary/stokesdarcy/couplingmapper.hh:30
bool isCoupledDarcyScvf(std::size_t darcyScvfIdx) const
Returns whether a Darcy scvf is coupled to the other domain.
Definition: boundary/stokesdarcy/couplingmapper.hh:119
void computeCouplingMapsAndStencils(const CouplingManager &couplingManager, Stencils &darcyToStokesCellCenterStencils, Stencils &darcyToStokesFaceStencils, Stencils &stokesCellCenterToDarcyStencils, Stencils &stokesFaceToDarcyStencils)
Main update routine.
Definition: boundary/stokesdarcy/couplingmapper.hh:44
const auto & darcyElementToStokesElementMap() const
A map that returns all Stokes elements coupled to a Darcy element.
Definition: boundary/stokesdarcy/couplingmapper.hh:127
const auto & stokesElementToDarcyElementMap() const
A map that returns all Darcy elements coupled to a Stokes element.
Definition: boundary/stokesdarcy/couplingmapper.hh:135
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
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:282
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: intersectingentities.hh:102
The available discretization methods in Dumux.
constexpr CCTpfa cctpfa
Definition: method.hh:145
constexpr Staggered staggered
Definition: method.hh:149