24#ifndef DUMUX_GMSH_GRID_DATA_HANDLE_HH
25#define DUMUX_GMSH_GRID_DATA_HANDLE_HH
31#include <dune/common/version.hh>
32#if DUNE_VERSION_LT(DUNE_COMMON,2,7)
33#include <dune/common/parallel/collectivecommunication.hh>
35#include <dune/common/parallel/communication.hh>
38#include <dune/geometry/dimension.hh>
39#include <dune/grid/common/partitionset.hh>
40#include <dune/grid/common/datahandleif.hh>
44#include <dune/grid/uggrid.hh>
53template<
class Gr
id,
class Gr
idFactory,
class Data>
54struct GmshGridDataHandle :
public Dune::CommDataHandleIF<GmshGridDataHandle<Grid, GridFactory, Data>, typename Data::value_type>
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)
65 const auto& indexSet = gridView_.indexSet();
67 for (
const auto& element : elements(gridView_, Dune::Partitions::interior))
68 std::swap(elementMarkers_[gridFactory.insertionIndex(element)], data_[idSet_.id(element)]);
70 for (
const auto& face : entities(gridView_, Dune::Codim<1>()))
71 std::swap(faceMarkers_[indexSet.index(face)], data_[idSet_.id(face)]);
76 const auto& indexSet = gridView_.indexSet();
78 elementMarkers_.resize(indexSet.size(0));
79 for (
const auto& element : elements(gridView_))
80 std::swap(elementMarkers_[indexSet.index(element)], data_[idSet_.id(element)]);
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)]);
86 boundaryMarkers_.resize(gridView_.grid().numBoundarySegments(), 0);
87 for (
const auto& element : elements(gridView_.grid().leafGridView()))
89 for (
const auto& intersection : intersections(gridView_.grid().leafGridView(), element))
91 if (intersection.boundary())
93 const auto marker = faceMarkers_[indexSet.index(element.template subEntity<1>(intersection.indexInInside()))];
94 boundaryMarkers_[intersection.boundarySegmentIndex()] = marker;
100 Dune::CommDataHandleIF<GmshGridDataHandle<Grid, GridFactory, Data>,
typename Data::value_type>&
interface()
104 {
return codim == 0 || codim == 1; }
106#if DUNE_VERSION_GTE(DUNE_GRID,2,7)
108 bool fixedSize(
int dim,
int codim)
const
116 template<
class EntityType>
117 std::size_t
size (
const EntityType& e)
const
120 template<
class MessageBufferImp,
class EntityType>
121 void gather (MessageBufferImp& buff,
const EntityType& e)
const
122 { buff.write(data_[idSet_.id(e)]); }
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)]); }
129 using IdSet =
typename Grid::LocalIdSet;
133 Data& elementMarkers_;
134 Data& boundaryMarkers_;
136 mutable std::map< typename IdSet::IdType, typename Data::value_type> data_;
145template<
class Gr
idFactory,
class Data,
int dimgr
id>
146struct GmshGridDataHandle<
Dune::UGGrid<dimgrid>, GridFactory, Data>
147:
public Dune::CommDataHandleIF<GmshGridDataHandle<Dune::UGGrid<dimgrid>, GridFactory, Data>, typename Data::value_type>
149 using Grid = Dune::UGGrid<dimgrid>;
150 using GridView =
typename Grid::LevelGridView;
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)
158 for (
const auto& element : elements(gridView_, Dune::Partitions::interior))
159 std::swap(elementMarkers_[gridFactory.insertionIndex(element)], data_[idSet_.id(element)]);
166 auto bmSizeMin = boundaryMarkers_.size();
167 Dune::MPIHelper::getCollectiveCommunication().min(&bmSizeMin, 1);
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);
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)]);
185 Dune::CommDataHandleIF<GmshGridDataHandle<Grid, GridFactory, Data>,
typename Data::value_type>&
interface()
188 bool contains (
int dim,
int codim)
const
189 {
return codim == 0 || codim == 1; }
191#if DUNE_VERSION_GTE(DUNE_GRID,2,7)
193 bool fixedSize(
int dim,
int codim)
const
201 template<
class EntityType>
202 std::size_t
size (
const EntityType& e)
const
205 template<
class MessageBufferImp,
class EntityType>
206 void gather (MessageBufferImp& buff,
const EntityType& e)
const
207 { buff.write(data_[idSet_.id(e)]); }
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)]); }
214 using IdSet =
typename Grid::LocalIdSet;
218 Data& elementMarkers_;
219 Data& boundaryMarkers_;
220 mutable std::map< typename IdSet::IdType, typename Data::value_type> data_;
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