3.1-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/parallel/collectivecommunication.hh>
32#include <dune/geometry/dimension.hh>
33#include <dune/grid/common/partitionset.hh>
34#include <dune/grid/common/datahandleif.hh>
35
36
37
38// UGGrid specific includes
39#if HAVE_UG
40#include <dune/grid/uggrid.hh>
41#endif
42
43namespace Dumux {
44
49template<class Grid, class GridFactory, class Data>
50struct GmshGridDataHandle : public Dune::CommDataHandleIF<GmshGridDataHandle<Grid, GridFactory, Data>, typename Data::value_type>
51{
52 using GridView = typename Grid::LevelGridView;
53
54 GmshGridDataHandle(const Grid& grid, const GridFactory& gridFactory, Data& elementMarkers, Data& boundaryMarkers, Data& faceMarkers)
55 : gridView_(grid.levelGridView(0))
56 , idSet_(grid.localIdSet())
57 , elementMarkers_(elementMarkers)
58 , boundaryMarkers_(boundaryMarkers)
59 , faceMarkers_(faceMarkers)
60 {
61 const auto& indexSet = gridView_.indexSet();
62
63 for (const auto& element : elements(gridView_, Dune::Partitions::interior))
64 std::swap(elementMarkers_[gridFactory.insertionIndex(element)], data_[idSet_.id(element)]);
65
66 for (const auto& face : entities(gridView_, Dune::Codim<1>()))
67 std::swap(faceMarkers_[indexSet.index(face)], data_[idSet_.id(face)]);
68 }
69
71 {
72 const auto& indexSet = gridView_.indexSet();
73
74 elementMarkers_.resize(indexSet.size(0));
75 for (const auto& element : elements(gridView_))
76 std::swap(elementMarkers_[indexSet.index(element)], data_[idSet_.id(element)]);
77
78 faceMarkers_.resize(indexSet.size(1));
79 for (const auto& face : entities(gridView_, Dune::Codim<1>()))
80 std::swap(faceMarkers_[indexSet.index(face)], data_[idSet_.id(face)]);
81
82 boundaryMarkers_.resize(gridView_.grid().numBoundarySegments(), 0);
83 for (const auto& element : elements(gridView_.grid().leafGridView()))
84 {
85 for (const auto& intersection : intersections(gridView_.grid().leafGridView(), element))
86 {
87 if (intersection.boundary())
88 {
89 const auto marker = faceMarkers_[indexSet.index(element.template subEntity<1>(intersection.indexInInside()))];
90 boundaryMarkers_[intersection.boundarySegmentIndex()] = marker;
91 }
92 }
93 }
94 }
95
96 Dune::CommDataHandleIF<GmshGridDataHandle<Grid, GridFactory, Data>, typename Data::value_type>& interface()
97 { return *this; }
98
99 bool contains (int dim, int codim) const
100 { return codim == 0 || codim == 1; }
101
102 bool fixedSize (int dim, int codim) const
103 { return true; }
104
105 template<class EntityType>
106 std::size_t size (const EntityType& e) const
107 { return 1; }
108
109 template<class MessageBufferImp, class EntityType>
110 void gather (MessageBufferImp& buff, const EntityType& e) const
111 { buff.write(data_[idSet_.id(e)]); }
112
113 template<class MessageBufferImp, class EntityType>
114 void scatter (MessageBufferImp& buff, const EntityType& e, std::size_t n)
115 { buff.read(data_[idSet_.id(e)]); }
116
117private:
118 using IdSet = typename Grid::LocalIdSet;
119
120 const GridView gridView_;
121 const IdSet &idSet_;
122 Data& elementMarkers_;
123 Data& boundaryMarkers_;
124 Data& faceMarkers_;
125 mutable std::map< typename IdSet::IdType, typename Data::value_type> data_;
126};
127
128#if HAVE_UG
129
134template<class GridFactory, class Data, int dimgrid>
135struct GmshGridDataHandle<Dune::UGGrid<dimgrid>, GridFactory, Data>
136: public Dune::CommDataHandleIF<GmshGridDataHandle<Dune::UGGrid<dimgrid>, GridFactory, Data>, typename Data::value_type>
137{
138 using Grid = Dune::UGGrid<dimgrid>;
139 using GridView = typename Grid::LevelGridView;
140
141 GmshGridDataHandle(const Grid& grid, const GridFactory& gridFactory, Data& elementMarkers, Data& boundaryMarkers)
142 : gridView_(grid.levelGridView(0))
143 , idSet_(grid.localIdSet())
144 , elementMarkers_(elementMarkers)
145 , boundaryMarkers_(boundaryMarkers)
146 {
147 for (const auto& element : elements(gridView_, Dune::Partitions::interior))
148 std::swap(elementMarkers_[gridFactory.insertionIndex(element)], data_[idSet_.id(element)]);
149
150 // Depending on the Dune version, the boundary markers are present on
151 // all processes (<= 2.6) or on the root process only (>= 2.7). Try to
152 // handle this in a flexible way: determine if the minimum size over
153 // all processes of the boundary markers vector is zero. If yes, assume
154 // that the root process contains all markers and broadcast them.
155 auto bmSizeMin = boundaryMarkers_.size();
156 Dune::MPIHelper::getCollectiveCommunication().min(&bmSizeMin, 1);
157 if (bmSizeMin == 0)
158 {
159 auto bmSize = boundaryMarkers_.size();
160 Dune::MPIHelper::getCollectiveCommunication().broadcast(&bmSize, 1, 0);
161 boundaryMarkers_.resize(bmSize);
162 Dune::MPIHelper::getCollectiveCommunication().broadcast(&boundaryMarkers_.front(), bmSize, 0);
163 }
164 }
165
167 {
168 const auto& indexSet = gridView_.indexSet();
169 elementMarkers_.resize(indexSet.size(0));
170 for (const auto& element : elements(gridView_))
171 std::swap(elementMarkers_[indexSet.index(element)], data_[idSet_.id(element)]);
172 }
173
174 Dune::CommDataHandleIF<GmshGridDataHandle<Grid, GridFactory, Data>, typename Data::value_type>& interface()
175 { return *this; }
176
177 bool contains (int dim, int codim) const
178 { return codim == 0 || codim == 1; }
179
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_UG
206
207} // namespace Dumux
208
209#endif
Dune::IteratorRange< typename MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper >::Intersections::const_iterator > intersections(const MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper > &glue)
Range generator to iterate with range-based for loops over all intersections as follows: for (const a...
Definition: glue.hh:62
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Definition: common/properties/model.hh:34
A data handle for commucating grid data for gmsh grids.
Definition: gmshgriddatahandle.hh:51
void gather(MessageBufferImp &buff, const EntityType &e) const
Definition: gmshgriddatahandle.hh:110
std::size_t size(const EntityType &e) const
Definition: gmshgriddatahandle.hh:106
void scatter(MessageBufferImp &buff, const EntityType &e, std::size_t n)
Definition: gmshgriddatahandle.hh:114
GmshGridDataHandle(const Grid &grid, const GridFactory &gridFactory, Data &elementMarkers, Data &boundaryMarkers, Data &faceMarkers)
Definition: gmshgriddatahandle.hh:54
bool contains(int dim, int codim) const
Definition: gmshgriddatahandle.hh:99
bool fixedSize(int dim, int codim) const
Definition: gmshgriddatahandle.hh:102
~GmshGridDataHandle()
Definition: gmshgriddatahandle.hh:70
Dune::CommDataHandleIF< GmshGridDataHandle< Grid, GridFactory, Data >, typename Data::value_type > & interface()
Definition: gmshgriddatahandle.hh:96
typename Grid::LevelGridView GridView
Definition: gmshgriddatahandle.hh:52