3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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 *****************************************************************************/
25#ifndef DUMUX_STOKES_DARCY_COUPLINGMAPPER_HH
26#define DUMUX_STOKES_DARCY_COUPLINGMAPPER_HH
27
28#include <iostream>
29#include <fstream>
30#include <string>
31#include <utility>
32#include <memory>
33#include <unordered_map>
34#include <vector>
35
36#include <dune/common/timer.hh>
37#include <dune/common/exceptions.hh>
40
41namespace Dumux {
42
47template<class MDTraits>
49{
50 using Scalar = typename MDTraits::Scalar;
51
52public:
53 static constexpr auto stokesCellCenterIdx = typename MDTraits::template SubDomain<0>::Index();
54 static constexpr auto stokesFaceIdx = typename MDTraits::template SubDomain<1>::Index();
55 static constexpr auto cellCenterIdx = typename MDTraits::template SubDomain<0>::Index();
56 static constexpr auto faceIdx = typename MDTraits::template SubDomain<1>::Index();
57 static constexpr auto stokesIdx = stokesCellCenterIdx;
58 static constexpr auto darcyIdx = typename MDTraits::template SubDomain<2>::Index();
59
60
61private:
62 // obtain the type tags of the sub problems
63 using StokesTypeTag = typename MDTraits::template SubDomain<0>::TypeTag;
64 using DarcyTypeTag = typename MDTraits::template SubDomain<2>::TypeTag;
65
66 struct ElementMapInfo
67 {
68 std::size_t eIdx;
69 std::size_t scvfIdx;
70 std::size_t flipScvfIdx;
71 };
72
73 // the sub domain type tags
74 template<std::size_t id>
75 using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
77
79 "The free flow domain must use the staggered discretization");
80
82 "The Darcy domain must use the CCTpfa discretization");
83public:
84
88 StokesDarcyCouplingMapper(const CouplingManager& couplingManager) : couplingManager_(couplingManager) {}
89
93 template<class Stencils>
94 void computeCouplingMapsAndStencils(Stencils& darcyToStokesCellCenterStencils,
95 Stencils& darcyToStokesFaceStencils,
96 Stencils& stokesCellCenterToDarcyStencils,
97 Stencils& stokesFaceToDarcyStencils)
98 {
99 const auto& stokesProblem = couplingManager_.problem(stokesIdx);
100 const auto& darcyProblem = couplingManager_.problem(darcyIdx);
101
102 const auto& stokesFvGridGeometry = stokesProblem.gridGeometry();
103 const auto& darcyFvGridGeometry = darcyProblem.gridGeometry();
104
105 isCoupledDarcyScvf_.resize(darcyFvGridGeometry.numScvf(), false);
106
107 auto darcyFvGeometry = localView(darcyFvGridGeometry);
108 auto stokesFvGeometry = localView(stokesFvGridGeometry);
109 const auto& stokesGridView = stokesFvGridGeometry.gridView();
110
111 for(const auto& stokesElement : elements(stokesGridView))
112 {
113 stokesFvGeometry.bindElement(stokesElement);
114
115 for(const auto& scvf : scvfs(stokesFvGeometry))
116 {
117 // skip the DOF if it is not on the boundary
118 if(!scvf.boundary())
119 continue;
120
121 // get element intersecting with the scvf center
122 // for robustness, add epsilon in unit outer normal direction
123 const auto eps = (scvf.center() - stokesElement.geometry().center()).two_norm()*1e-8;
124 auto globalPos = scvf.center(); globalPos.axpy(eps, scvf.unitOuterNormal());
125 const auto darcyElementIdx = intersectingEntities(globalPos, darcyFvGridGeometry.boundingBoxTree());
126
127 // skip if no intersection was found
128 if(darcyElementIdx.empty())
129 continue;
130
131 // sanity check
132 if(darcyElementIdx.size() > 1)
133 DUNE_THROW(Dune::InvalidStateException, "Stokes face dof should only intersect with one Darcy element");
134
135 const auto stokesElementIdx = stokesFvGridGeometry.elementMapper().index(stokesElement);
136
137 const auto darcyDofIdx = darcyElementIdx[0];
138
139 stokesFaceToDarcyStencils[scvf.dofIndex()].push_back(darcyDofIdx);
140 stokesCellCenterToDarcyStencils[stokesElementIdx].push_back(darcyDofIdx);
141
142 darcyToStokesFaceStencils[darcyElementIdx[0]].push_back(scvf.dofIndex());
143 darcyToStokesCellCenterStencils[darcyElementIdx[0]].push_back(stokesElementIdx);
144
145 const auto& darcyElement = darcyFvGridGeometry.element(darcyElementIdx[0]);
146 darcyFvGeometry.bindElement(darcyElement);
147
148 // find the corresponding Darcy sub control volume face
149 for(const auto& darcyScvf : scvfs(darcyFvGeometry))
150 {
151 const Scalar distance = (darcyScvf.center() - scvf.center()).two_norm();
152
153 if(distance < eps)
154 {
155 isCoupledDarcyScvf_[darcyScvf.index()] = true;
156 darcyElementToStokesElementMap_[darcyElementIdx[0]].push_back({stokesElementIdx, scvf.index(), darcyScvf.index()});
157 stokesElementToDarcyElementMap_[stokesElementIdx].push_back({darcyElementIdx[0], darcyScvf.index(), scvf.index()});
158 }
159 }
160 }
161 }
162 }
163
167 bool isCoupledDarcyScvf(std::size_t darcyScvfIdx) const
168 {
169 return isCoupledDarcyScvf_[darcyScvfIdx];
170 }
171
176 {
177 return darcyElementToStokesElementMap_;
178 }
179
184 {
185 return stokesElementToDarcyElementMap_;
186 }
187
188private:
189 std::unordered_map<std::size_t, std::vector<ElementMapInfo>> darcyElementToStokesElementMap_;
190 std::unordered_map<std::size_t, std::vector<ElementMapInfo>> stokesElementToDarcyElementMap_;
191
192 std::vector<bool> isCoupledDarcyScvf_;
193
194 const CouplingManager& couplingManager_;
195};
196
197} // end namespace Dumux
198
199#endif
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
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:45
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
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
constexpr double eps
epsilon for checking direction of scvf normals
Definition: test_tpfafvgeometry_nonconforming.cc:44
Definition: common/properties.hh:150
Coupling mapper for Stokes and Darcy domains with equal dimension.
Definition: boundary/stokesdarcy/couplingmapper.hh:49
static constexpr auto stokesIdx
Definition: boundary/stokesdarcy/couplingmapper.hh:57
static constexpr auto darcyIdx
Definition: boundary/stokesdarcy/couplingmapper.hh:58
bool isCoupledDarcyScvf(std::size_t darcyScvfIdx) const
Returns whether a Darcy scvf is coupled to the other domain.
Definition: boundary/stokesdarcy/couplingmapper.hh:167
StokesDarcyCouplingMapper(const CouplingManager &couplingManager)
Constructor.
Definition: boundary/stokesdarcy/couplingmapper.hh:88
static constexpr auto faceIdx
Definition: boundary/stokesdarcy/couplingmapper.hh:56
static constexpr auto cellCenterIdx
Definition: boundary/stokesdarcy/couplingmapper.hh:55
const auto & darcyElementToStokesElementMap() const
A map that returns all Stokes elements coupled to a Darcy element.
Definition: boundary/stokesdarcy/couplingmapper.hh:175
const auto & stokesElementToDarcyElementMap() const
A map that returns all Darcy elements coupled to a Stokes element.
Definition: boundary/stokesdarcy/couplingmapper.hh:183
static constexpr auto stokesFaceIdx
Definition: boundary/stokesdarcy/couplingmapper.hh:54
void computeCouplingMapsAndStencils(Stencils &darcyToStokesCellCenterStencils, Stencils &darcyToStokesFaceStencils, Stencils &stokesCellCenterToDarcyStencils, Stencils &stokesFaceToDarcyStencils)
Main update routine.
Definition: boundary/stokesdarcy/couplingmapper.hh:94
static constexpr auto stokesCellCenterIdx
Definition: boundary/stokesdarcy/couplingmapper.hh:53
Definition: multidomain/couplingmanager.hh:46
Declares all properties used in Dumux.