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