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