24#ifndef DUMUX_GMSH_GRID_DATA_HANDLE_HH
25#define DUMUX_GMSH_GRID_DATA_HANDLE_HH
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>
38#include <dune/grid/uggrid.hh>
47template<
class Gr
id,
class Gr
idFactory,
class Data>
48struct GmshGridDataHandle :
public Dune::CommDataHandleIF<GmshGridDataHandle<Grid, GridFactory, Data>, typename Data::value_type>
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)
59 const auto& indexSet = gridView_.indexSet();
61 for (
const auto& element : elements(gridView_, Dune::Partitions::interior))
62 std::swap(elementMarkers_[gridFactory.insertionIndex(element)], data_[idSet_.id(element)]);
64 for (
const auto& face : entities(gridView_, Dune::Codim<1>()))
65 std::swap(faceMarkers_[indexSet.index(face)], data_[idSet_.id(face)]);
70 const auto& indexSet = gridView_.indexSet();
72 elementMarkers_.resize(indexSet.size(0));
73 for (
const auto& element : elements(gridView_))
74 std::swap(elementMarkers_[indexSet.index(element)], data_[idSet_.id(element)]);
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)]);
80 boundaryMarkers_.resize(gridView_.grid().numBoundarySegments(), 0);
81 for (
const auto& element : elements(gridView_.grid().leafGridView()))
83 for (
const auto& intersection : intersections(gridView_.grid().leafGridView(), element))
85 if (intersection.boundary())
87 const auto marker = faceMarkers_[indexSet.index(element.template subEntity<1>(intersection.indexInInside()))];
88 boundaryMarkers_[intersection.boundarySegmentIndex()] = marker;
94 Dune::CommDataHandleIF<GmshGridDataHandle<Grid, GridFactory, Data>,
typename Data::value_type>&
interface()
98 {
return codim == 0 || codim == 1; }
104 template<
class EntityType>
105 std::size_t
size (
const EntityType& e)
const
108 template<
class MessageBufferImp,
class EntityType>
109 void gather (MessageBufferImp& buff,
const EntityType& e)
const
110 { buff.write(data_[idSet_.id(e)]); }
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)]); }
117 using IdSet =
typename Grid::LocalIdSet;
121 Data& elementMarkers_;
122 Data& boundaryMarkers_;
124 mutable std::map< typename IdSet::IdType, typename Data::value_type> data_;
133template<
class Gr
idFactory,
class Data,
int dimgr
id>
134struct GmshGridDataHandle<
Dune::UGGrid<dimgrid>, GridFactory, Data>
135:
public Dune::CommDataHandleIF<GmshGridDataHandle<Dune::UGGrid<dimgrid>, GridFactory, Data>, typename Data::value_type>
137 using Grid = Dune::UGGrid<dimgrid>;
138 using GridView =
typename Grid::LevelGridView;
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)
146 for (
const auto& element : elements(gridView_, Dune::Partitions::interior))
147 std::swap(elementMarkers_[gridFactory.insertionIndex(element)], data_[idSet_.id(element)]);
154 auto bmSizeMin = boundaryMarkers_.size();
155 Dune::MPIHelper::getCommunication().min(&bmSizeMin, 1);
158 auto bmSize = boundaryMarkers_.size();
159 Dune::MPIHelper::getCommunication().broadcast(&bmSize, 1, 0);
160 boundaryMarkers_.resize(bmSize);
161 Dune::MPIHelper::getCommunication().broadcast(&boundaryMarkers_.front(), bmSize, 0);
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)]);
173 Dune::CommDataHandleIF<GmshGridDataHandle<Grid, GridFactory, Data>,
typename Data::value_type>&
interface()
176 bool contains (
int dim,
int codim)
const
177 {
return codim == 0 || codim == 1; }
183 template<
class EntityType>
184 std::size_t
size (
const EntityType& e)
const
187 template<
class MessageBufferImp,
class EntityType>
188 void gather (MessageBufferImp& buff,
const EntityType& e)
const
189 { buff.write(data_[idSet_.id(e)]); }
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)]); }
196 using IdSet =
typename Grid::LocalIdSet;
200 Data& elementMarkers_;
201 Data& boundaryMarkers_;
202 mutable std::map< typename IdSet::IdType, typename Data::value_type> data_;
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