3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
gmshgriddatahandle.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 *****************************************************************************/
24#ifndef DUMUX_GMSH_GRID_DATA_HANDLE_HH
25#define DUMUX_GMSH_GRID_DATA_HANDLE_HH
26
27#include <memory>
28#include <algorithm>
29#include <map>
30
31#include <dune/common/parallel/communication.hh>
32#include <dune/geometry/dimension.hh>
33#include <dune/grid/common/partitionset.hh>
34#include <dune/grid/common/datahandleif.hh>
35
36// UGGrid specific includes
37#if HAVE_DUNE_UGGRID
38#include <dune/grid/uggrid.hh>
39#endif
40
41namespace Dumux {
42
47template<class Grid, class GridFactory, class Data>
48struct GmshGridDataHandle : public Dune::CommDataHandleIF<GmshGridDataHandle<Grid, GridFactory, Data>, typename Data::value_type>
49{
50 using GridView = typename Grid::LevelGridView;
51
52 GmshGridDataHandle(const Grid& grid, const GridFactory& gridFactory, Data& elementMarkers, Data& boundaryMarkers, Data& faceMarkers)
53 : gridView_(grid.levelGridView(0))
54 , idSet_(grid.localIdSet())
55 , elementMarkers_(elementMarkers)
56 , boundaryMarkers_(boundaryMarkers)
57 , faceMarkers_(faceMarkers)
58 {
59 const auto& indexSet = gridView_.indexSet();
60
61 for (const auto& element : elements(gridView_, Dune::Partitions::interior))
62 std::swap(elementMarkers_[gridFactory.insertionIndex(element)], data_[idSet_.id(element)]);
63
64 for (const auto& face : entities(gridView_, Dune::Codim<1>()))
65 std::swap(faceMarkers_[indexSet.index(face)], data_[idSet_.id(face)]);
66 }
67
69 {
70 const auto& indexSet = gridView_.indexSet();
71
72 elementMarkers_.resize(indexSet.size(0));
73 for (const auto& element : elements(gridView_))
74 std::swap(elementMarkers_[indexSet.index(element)], data_[idSet_.id(element)]);
75
76 faceMarkers_.resize(indexSet.size(1));
77 for (const auto& face : entities(gridView_, Dune::Codim<1>()))
78 std::swap(faceMarkers_[indexSet.index(face)], data_[idSet_.id(face)]);
79
80 boundaryMarkers_.resize(gridView_.grid().numBoundarySegments(), 0);
81 for (const auto& element : elements(gridView_.grid().leafGridView()))
82 {
83 for (const auto& intersection : intersections(gridView_.grid().leafGridView(), element))
84 {
85 if (intersection.boundary())
86 {
87 const auto marker = faceMarkers_[indexSet.index(element.template subEntity<1>(intersection.indexInInside()))];
88 boundaryMarkers_[intersection.boundarySegmentIndex()] = marker;
89 }
90 }
91 }
92 }
93
94 Dune::CommDataHandleIF<GmshGridDataHandle<Grid, GridFactory, Data>, typename Data::value_type>& interface()
95 { return *this; }
96
97 bool contains (int dim, int codim) const
98 { return codim == 0 || codim == 1; }
99
101 bool fixedSize(int dim, int codim) const
102 { return true; }
103
104 template<class EntityType>
105 std::size_t size (const EntityType& e) const
106 { return 1; }
107
108 template<class MessageBufferImp, class EntityType>
109 void gather (MessageBufferImp& buff, const EntityType& e) const
110 { buff.write(data_[idSet_.id(e)]); }
111
112 template<class MessageBufferImp, class EntityType>
113 void scatter (MessageBufferImp& buff, const EntityType& e, std::size_t n)
114 { buff.read(data_[idSet_.id(e)]); }
115
116private:
117 using IdSet = typename Grid::LocalIdSet;
118
119 const GridView gridView_;
120 const IdSet &idSet_;
121 Data& elementMarkers_;
122 Data& boundaryMarkers_;
123 Data& faceMarkers_;
124 mutable std::map< typename IdSet::IdType, typename Data::value_type> data_;
125};
126
127#if HAVE_DUNE_UGGRID
128
133template<class GridFactory, class Data, int dimgrid>
134struct GmshGridDataHandle<Dune::UGGrid<dimgrid>, GridFactory, Data>
135: public Dune::CommDataHandleIF<GmshGridDataHandle<Dune::UGGrid<dimgrid>, GridFactory, Data>, typename Data::value_type>
136{
137 using Grid = Dune::UGGrid<dimgrid>;
138 using GridView = typename Grid::LevelGridView;
139
140 GmshGridDataHandle(const Grid& grid, const GridFactory& gridFactory, Data& elementMarkers, Data& boundaryMarkers)
141 : gridView_(grid.levelGridView(0))
142 , idSet_(grid.localIdSet())
143 , elementMarkers_(elementMarkers)
144 , boundaryMarkers_(boundaryMarkers)
145 {
146 for (const auto& element : elements(gridView_, Dune::Partitions::interior))
147 std::swap(elementMarkers_[gridFactory.insertionIndex(element)], data_[idSet_.id(element)]);
148
149 // Depending on the Dune version, the boundary markers are present on
150 // all processes (<= 2.6) or on the root process only (>= 2.7). Try to
151 // handle this in a flexible way: determine if the minimum size over
152 // all processes of the boundary markers vector is zero. If yes, assume
153 // that the root process contains all markers and broadcast them.
154 auto bmSizeMin = boundaryMarkers_.size();
155 Dune::MPIHelper::getCollectiveCommunication().min(&bmSizeMin, 1);
156 if (bmSizeMin == 0)
157 {
158 auto bmSize = boundaryMarkers_.size();
159 Dune::MPIHelper::getCollectiveCommunication().broadcast(&bmSize, 1, 0);
160 boundaryMarkers_.resize(bmSize);
161 Dune::MPIHelper::getCollectiveCommunication().broadcast(&boundaryMarkers_.front(), bmSize, 0);
162 }
163 }
164
166 {
167 const auto& indexSet = gridView_.indexSet();
168 elementMarkers_.resize(indexSet.size(0));
169 for (const auto& element : elements(gridView_))
170 std::swap(elementMarkers_[indexSet.index(element)], data_[idSet_.id(element)]);
171 }
172
173 Dune::CommDataHandleIF<GmshGridDataHandle<Grid, GridFactory, Data>, typename Data::value_type>& interface()
174 { return *this; }
175
176 bool contains (int dim, int codim) const
177 { return codim == 0 || codim == 1; }
178
180 bool fixedSize(int dim, int codim) const
181 { return true; }
182
183 template<class EntityType>
184 std::size_t size (const EntityType& e) const
185 { return 1; }
186
187 template<class MessageBufferImp, class EntityType>
188 void gather (MessageBufferImp& buff, const EntityType& e) const
189 { buff.write(data_[idSet_.id(e)]); }
190
191 template<class MessageBufferImp, class EntityType>
192 void scatter (MessageBufferImp& buff, const EntityType& e, std::size_t n)
193 { buff.read(data_[idSet_.id(e)]); }
194
195private:
196 using IdSet = typename Grid::LocalIdSet;
197
198 const GridView gridView_;
199 const IdSet &idSet_;
200 Data& elementMarkers_;
201 Data& boundaryMarkers_;
202 mutable std::map< typename IdSet::IdType, typename Data::value_type> data_;
203};
204
205#endif // HAVE_DUNE_UGGRID
206
207} // namespace Dumux
208
209#endif
Definition: adapt.hh:29
Definition: common/pdesolver.hh:36
A data handle for commucating grid data for gmsh grids.
Definition: gmshgriddatahandle.hh:49
void gather(MessageBufferImp &buff, const EntityType &e) const
Definition: gmshgriddatahandle.hh:109
std::size_t size(const EntityType &e) const
Definition: gmshgriddatahandle.hh:105
void scatter(MessageBufferImp &buff, const EntityType &e, std::size_t n)
Definition: gmshgriddatahandle.hh:113
GmshGridDataHandle(const Grid &grid, const GridFactory &gridFactory, Data &elementMarkers, Data &boundaryMarkers, Data &faceMarkers)
Definition: gmshgriddatahandle.hh:52
bool contains(int dim, int codim) const
Definition: gmshgriddatahandle.hh:97
bool fixedSize(int dim, int codim) const
returns true if size per entity of given dim and codim is a constant
Definition: gmshgriddatahandle.hh:101
~GmshGridDataHandle()
Definition: gmshgriddatahandle.hh:68
Dune::CommDataHandleIF< GmshGridDataHandle< Grid, GridFactory, Data >, typename Data::value_type > & interface()
Definition: gmshgriddatahandle.hh:94
typename Grid::LevelGridView GridView
Definition: gmshgriddatahandle.hh:50