24#ifndef DUMUX_GMSH_GRID_DATA_HANDLE_HH
25#define DUMUX_GMSH_GRID_DATA_HANDLE_HH
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>
40#include <dune/grid/uggrid.hh>
49template<
class Gr
id,
class Gr
idFactory,
class Data>
50struct GmshGridDataHandle :
public Dune::CommDataHandleIF<GmshGridDataHandle<Grid, GridFactory, Data>, typename Data::value_type>
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)
61 const auto& indexSet = gridView_.indexSet();
63 for (
const auto& element : elements(gridView_, Dune::Partitions::interior))
64 std::swap(elementMarkers_[gridFactory.insertionIndex(element)], data_[idSet_.id(element)]);
66 for (
const auto& face : entities(gridView_, Dune::Codim<1>()))
67 std::swap(faceMarkers_[indexSet.index(face)], data_[idSet_.id(face)]);
72 const auto& indexSet = gridView_.indexSet();
74 elementMarkers_.resize(indexSet.size(0));
75 for (
const auto& element : elements(gridView_))
76 std::swap(elementMarkers_[indexSet.index(element)], data_[idSet_.id(element)]);
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)]);
82 boundaryMarkers_.resize(gridView_.grid().numBoundarySegments(), 0);
83 for (
const auto& element : elements(gridView_.grid().leafGridView()))
85 for (
const auto& intersection :
intersections(gridView_.grid().leafGridView(), element))
87 if (intersection.boundary())
89 const auto marker = faceMarkers_[indexSet.index(element.template subEntity<1>(intersection.indexInInside()))];
90 boundaryMarkers_[intersection.boundarySegmentIndex()] = marker;
96 Dune::CommDataHandleIF<GmshGridDataHandle<Grid, GridFactory, Data>,
typename Data::value_type>&
interface()
100 {
return codim == 0 || codim == 1; }
105 template<
class EntityType>
106 std::size_t
size (
const EntityType& e)
const
109 template<
class MessageBufferImp,
class EntityType>
110 void gather (MessageBufferImp& buff,
const EntityType& e)
const
111 { buff.write(data_[idSet_.id(e)]); }
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)]); }
118 using IdSet =
typename Grid::LocalIdSet;
122 Data& elementMarkers_;
123 Data& boundaryMarkers_;
125 mutable std::map< typename IdSet::IdType, typename Data::value_type> data_;
134template<
class Gr
idFactory,
class Data,
int dimgr
id>
135struct GmshGridDataHandle<
Dune::UGGrid<dimgrid>, GridFactory, Data>
136:
public Dune::CommDataHandleIF<GmshGridDataHandle<Dune::UGGrid<dimgrid>, GridFactory, Data>, typename Data::value_type>
138 using Grid = Dune::UGGrid<dimgrid>;
139 using GridView =
typename Grid::LevelGridView;
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)
147 for (
const auto& element : elements(gridView_, Dune::Partitions::interior))
148 std::swap(elementMarkers_[gridFactory.insertionIndex(element)], data_[idSet_.id(element)]);
155 auto bmSizeMin = boundaryMarkers_.size();
156 Dune::MPIHelper::getCollectiveCommunication().min(&bmSizeMin, 1);
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);
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)]);
174 Dune::CommDataHandleIF<GmshGridDataHandle<Grid, GridFactory, Data>,
typename Data::value_type>&
interface()
177 bool contains (
int dim,
int codim)
const
178 {
return codim == 0 || codim == 1; }
180 bool fixedSize (
int dim,
int codim)
const
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_;
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