version 3.10-dev
boundary/stokesdarcy/couplingmapper.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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_STOKES_DARCY_COUPLINGMAPPER_HH
14#define DUMUX_STOKES_DARCY_COUPLINGMAPPER_HH
15
16#include <type_traits>
17#include <unordered_map>
18#include <vector>
19
20#include <dune/common/exceptions.hh>
22
23namespace Dumux {
24
30{
31 struct ElementMapInfo
32 {
33 std::size_t eIdx;
34 std::size_t scvfIdx;
35 std::size_t flipScvfIdx;
36 };
37
38public:
39
43 template<class CouplingManager, class Stencils>
45 Stencils& darcyToStokesCellCenterStencils,
46 Stencils& darcyToStokesFaceStencils,
47 Stencils& stokesCellCenterToDarcyStencils,
48 Stencils& stokesFaceToDarcyStencils)
49 {
50 const auto& stokesFvGridGeometry = couplingManager.problem(CouplingManager::stokesIdx).gridGeometry();
51 const auto& darcyFvGridGeometry = couplingManager.problem(CouplingManager::darcyIdx).gridGeometry();
52
53 static_assert(std::decay_t<decltype(stokesFvGridGeometry)>::discMethod == DiscretizationMethods::staggered,
54 "The free flow domain must use the staggered discretization");
55
56 static_assert(std::decay_t<decltype(darcyFvGridGeometry)>::discMethod == DiscretizationMethods::cctpfa,
57 "The Darcy domain must use the CCTpfa discretization");
58
59 isCoupledDarcyScvf_.resize(darcyFvGridGeometry.numScvf(), false);
60
61 const auto& stokesGridView = stokesFvGridGeometry.gridView();
62 auto stokesFvGeometry = localView(stokesFvGridGeometry);
63 for (const auto& stokesElement : elements(stokesGridView))
64 {
65 stokesFvGeometry.bindElement(stokesElement);
66
67 for (const auto& scvf : scvfs(stokesFvGeometry))
68 {
69 // skip the DOF if it is not on the boundary
70 if (!scvf.boundary())
71 continue;
72
73 // get element intersecting with the scvf center
74 // for robustness, add epsilon in unit outer normal direction
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());
78
79 // skip if no intersection was found
80 if (darcyElementIdx.empty())
81 continue;
82
83 // sanity check
84 if (darcyElementIdx.size() > 1)
85 DUNE_THROW(Dune::InvalidStateException, "Stokes face dof should only intersect with one Darcy element");
86
87 const auto stokesElementIdx = stokesFvGridGeometry.elementMapper().index(stokesElement);
88
89 const auto darcyDofIdx = darcyElementIdx[0];
90
91 stokesFaceToDarcyStencils[scvf.dofIndex()].push_back(darcyDofIdx);
92 stokesCellCenterToDarcyStencils[stokesElementIdx].push_back(darcyDofIdx);
93
94 darcyToStokesFaceStencils[darcyElementIdx[0]].push_back(scvf.dofIndex());
95 darcyToStokesCellCenterStencils[darcyElementIdx[0]].push_back(stokesElementIdx);
96
97 const auto& darcyElement = darcyFvGridGeometry.element(darcyElementIdx[0]);
98 const auto darcyFvGeometry = localView(darcyFvGridGeometry).bindElement(darcyElement);
99
100 // find the corresponding Darcy sub control volume face
101 for (const auto& darcyScvf : scvfs(darcyFvGeometry))
102 {
103 const auto distance = (darcyScvf.center() - scvf.center()).two_norm();
104
105 if (distance < eps)
106 {
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()});
110 }
111 }
112 }
113 }
114 }
115
119 bool isCoupledDarcyScvf(std::size_t darcyScvfIdx) const
120 {
121 return isCoupledDarcyScvf_[darcyScvfIdx];
122 }
123
128 {
129 return darcyElementToStokesElementMap_;
130 }
131
136 {
137 return stokesElementToDarcyElementMap_;
138 }
139
140private:
141 std::unordered_map<std::size_t, std::vector<ElementMapInfo>> darcyElementToStokesElementMap_;
142 std::unordered_map<std::size_t, std::vector<ElementMapInfo>> stokesElementToDarcyElementMap_;
143
144 std::vector<bool> isCoupledDarcyScvf_;
145};
146
147} // end namespace Dumux
148
149#endif
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:37
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:298
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
Definition: adapt.hh:17