3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
multidomain/facet/gridmanager.hh
Go to the documentation of this file.
1/*****************************************************************************
2 * See the file COPYING for full copying permissions. *
3 * *
4 * This program is free software: you can redistribute it and/or modify *
5 * it under the terms of the GNU General Public License as published by *
6 * the Free Software Foundation, either version 3 of the License, or *
7 * (at your option) any later version. *
8 * *
9 * This program is distributed in the hope that it will be useful, *
10 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
11 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
12 * GNU General Public License for more details. *
13 * *
14 * You should have received a copy of the GNU General Public License *
15 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
16 *****************************************************************************/
27#ifndef DUMUX_FACETCOUPLING_GRID_MANAGER_HH
28#define DUMUX_FACETCOUPLING_GRID_MANAGER_HH
29
30#include <cassert>
31#include <map>
32#include <tuple>
33#include <type_traits>
34#include <vector>
35
36#include <dune/common/exceptions.hh>
37#include <dune/common/hybridutilities.hh>
38#include <dune/grid/common/gridfactory.hh>
39
44
45#include "gmshreader.hh"
46
47namespace Dumux {
48namespace FCGridManagerChecks {
49
50 // The grid creator and grid data classes provided below
51 // require that the grids passed as template arguments are
52 // ordered in descending grid dimension and they must have
53 // the same world dimension. The following helper structs
54 // verify these conditions and are reused below.
55 // evaluates if dim is descending or dimworld is equal for two grids
56 template<bool checkDimWorld, typename G1, typename G2>
57 static constexpr bool evalCondition()
58 {
59 return checkDimWorld ? int(G1::dimensionworld) == int(G2::dimensionworld)
60 : int(G1::dimension) > int(G2::dimension);
61 }
62
63 // helper structs to evaluate conditions on the grid
64 template<bool checkDimWorld, typename... Gs> struct FulfillConditions;
65 template<bool checkDimWorld, typename G1, typename... Gs>
66 struct FulfillConditions<checkDimWorld, G1, Gs...>
67 {
68 using G2 = typename std::tuple_element_t<0, std::tuple<Gs...>>;
69 static constexpr bool value = evalCondition<checkDimWorld, G1, G2>() && FulfillConditions<checkDimWorld, Gs...>::value;
70 };
71 template<bool checkDimWorld, typename G1, typename G2>
72 struct FulfillConditions<checkDimWorld, G1, G2> { static constexpr bool value = evalCondition<checkDimWorld, G1, G2>(); };
73}
74
83template<typename... Grids>
85{
86 // make sure all grids have the same world dimension and are ordered in descending dimension
87 static_assert(FCGridManagerChecks::FulfillConditions<false, Grids...>::value, "All grids must have the same world dimension!");
88 static_assert(FCGridManagerChecks::FulfillConditions<true, Grids...>::value, "Grids must be ordered w.r.t the dimension in descending order!");
89
91 static constexpr std::size_t numGrids = sizeof...(Grids);
93 template<std::size_t id> using Grid = typename std::tuple_element_t<id, std::tuple<Grids...>>;
95 template<std::size_t id> using GridDataPtr = std::shared_ptr< GridData<Grid<id>> >;
97 template<std::size_t id> using Intersection = typename Grid<id>::LeafGridView::Intersection;
98
99public:
101 template<std::size_t id> using GridData = GridData<Grid<id>>;
102
104 template<std::size_t id>
105 void setGridData(GridData<id>&& gridData)
106 {
107 std::get<id>(gridDataPtrTuple_) = std::make_shared<GridData<id>>( std::move(gridData) );
108 }
109
111 template<std::size_t id>
112 std::shared_ptr<const GridData<id>> getSubDomainGridData() const
113 { return std::get<id>(gridDataPtrTuple_); }
114
116 template<std::size_t id>
117 int getElementDomainMarker(const typename Grid<id>::template Codim<0>::Entity& element) const
118 { return std::get<id>(gridDataPtrTuple_)->getElementDomainMarker(element); }
119
121 template<std::size_t id>
122 int getBoundaryDomainMarker(const typename Grid<id>::LeafGridView::Intersection& is) const
123 { return std::get<id>(gridDataPtrTuple_)->getBoundaryDomainMarker(is); }
124
126 template<std::size_t id>
127 int getBoundaryDomainMarker(int boundarySegmentIndex) const
128 { return std::get<id>(gridDataPtrTuple_)->getBoundaryDomainMarker(boundarySegmentIndex); }
129
131 template<std::size_t id>
132 bool wasInserted(const Intersection<id>& intersection) const
133 { return std::get<id>(gridDataPtrTuple_)->wasInserted(intersection); }
134
135private:
138 GridDataPtrTuple gridDataPtrTuple_;
139};
140
153template<typename... Grids>
155{
156 // make sure all grids have the same world dimension and are ordered in descending dimension
157 static_assert(FCGridManagerChecks::FulfillConditions<false, Grids...>::value, "All grids must have the same world dimension!");
158 static_assert(FCGridManagerChecks::FulfillConditions<true, Grids...>::value, "Grids must be ordered w.r.t the dimension in descending order!");
159
161 template<std::size_t id> using Grid = typename std::tuple_element_t<id, std::tuple<Grids...>>;
163 template<std::size_t id> using GridFactory = typename Dune::GridFactory<Grid<id>>;
164
166 using GIType = typename IndexTraits< typename Grid<0>::LeafGridView >::GridIndex;
168 using EmbedmentMap = std::unordered_map<GIType, std::vector<GIType>>;
169
170public:
172 template<std::size_t id> using GridView = typename Grid<id>::LeafGridView;
173
175 static constexpr std::size_t numGrids = sizeof...(Grids);
177 static constexpr int bulkGridId = 0;
179 static constexpr int bulkDim = Grid<bulkGridId>::dimension;
180
184 using GridIndexType = GIType;
185
187 template<std::size_t id>
188 const GridView<id>& gridView() const
189 { return *std::get<id>(gridViewPtrTuple_); }
190
192 template<std::size_t id, class Entity>
193 GridIndexType insertionIndex(const Entity& entity) const
194 { return std::get<id>(gridFactoryPtrTuple_)->insertionIndex(entity); }
195
197 template<std::size_t id>
198 typename std::unordered_map< GridIndexType, std::vector<GridIndexType> >::mapped_type
199 embeddedEntityIndices(const typename Grid<id>::template Codim<0>::Entity& element) const
200 {
201 const auto& map = embeddedEntityMaps_[id];
202 auto it = map.find( std::get<id>(gridFactoryPtrTuple_)->insertionIndex(element) );
203 if (it != map.end()) return it->second;
204 else return typename std::unordered_map< GridIndexType, std::vector<GridIndexType> >::mapped_type();
205 }
206
208 template<std::size_t id>
209 typename std::unordered_map< GridIndexType, std::vector<GridIndexType> >::mapped_type
210 adjoinedEntityIndices(const typename Grid<id>::template Codim<0>::Entity& element) const
211 {
212 const auto& map = adjoinedEntityMaps_[id];
213 auto it = map.find( std::get<id>(gridFactoryPtrTuple_)->insertionIndex(element) );
214 if (it != map.end()) return it->second;
215 else return typename std::unordered_map< GridIndexType, std::vector<GridIndexType> >::mapped_type();
216 }
217
219 const std::unordered_map< GridIndexType, std::vector<GridIndexType> >& embeddedEntityMap(std::size_t id) const
220 { assert(id < numGrids); return embeddedEntityMaps_[id]; }
221
223 std::unordered_map< GridIndexType, std::vector<GridIndexType> >& embeddedEntityMap(std::size_t id)
224 { assert(id < numGrids); return embeddedEntityMaps_[id]; }
225
227 const std::unordered_map< GridIndexType, std::vector<GridIndexType> >& adjoinedEntityMap(std::size_t id) const
228 { assert(id < numGrids); return adjoinedEntityMaps_[id]; }
229
231 std::unordered_map< GridIndexType, std::vector<GridIndexType> >& adjoinedEntityMap(std::size_t id)
232 { assert(id < numGrids); return adjoinedEntityMaps_[id]; }
233
235 const std::vector<GridIndexType>& gridHierarchyIndices(std::size_t id) const
236 { assert(id < numGrids); return gridVertexIndices_[id]; }
237
239 std::size_t numVerticesInHierarchy() const
240 { return numVerticesInHierarchy_; }
241
251 template<std::size_t id>
252 void setData( std::shared_ptr<Grid<id>> gridPtr,
253 std::shared_ptr<GridFactory<id>> gridFactoryPtr,
254 EmbedmentMap&& embeddedEntityMap,
255 EmbedmentMap&& adjoinedEntityMap,
256 std::vector<GridIndexType>&& gridVertexIndices,
257 std::size_t numVerticesInHierarchy )
258 {
259 std::get<id>(gridViewPtrTuple_) = std::make_shared<GridView<id>>(gridPtr->leafGridView());
260 std::get<id>(gridFactoryPtrTuple_) = gridFactoryPtr;
261 embeddedEntityMaps_[id] = std::move(embeddedEntityMap);
262 adjoinedEntityMaps_[id] = std::move(adjoinedEntityMap);
263 gridVertexIndices_[id] = std::move(gridVertexIndices);
264 numVerticesInHierarchy_ = numVerticesInHierarchy;
265 }
266
267private:
269 std::array<EmbedmentMap, numGrids> embeddedEntityMaps_;
270 std::array<EmbedmentMap, numGrids> adjoinedEntityMaps_;
271
273 std::size_t numVerticesInHierarchy_;
274 std::array<std::vector<GridIndexType>, numGrids> gridVertexIndices_;
275
277 using Indices = std::make_index_sequence<numGrids>;
278 template<std::size_t id> using GridViewPtr = std::shared_ptr<GridView<id>>;
280 GridPtrTuple gridViewPtrTuple_;
281
283 template<std::size_t id> using GridFactoryPtr = std::shared_ptr< Dune::GridFactory<Grid<id>> >;
284 using GridFactoryPtrTuple = typename makeFromIndexedType<std::tuple, GridFactoryPtr, Indices>::type;
285 GridFactoryPtrTuple gridFactoryPtrTuple_;
286};
287
297template<typename... Grids>
299{
300 // make sure all grids have the same world dimension and are ordered in descending dimension
301 static_assert(FCGridManagerChecks::FulfillConditions<false, Grids...>::value, "All grids must have the same world dimension!");
302 static_assert(FCGridManagerChecks::FulfillConditions<true, Grids...>::value, "Grids must be ordered w.r.t the dimension in descending order!");
303
304 // we use a wrapper class for the grid data containing the data on all grids
306public:
308 template<std::size_t id> using Grid = typename std::tuple_element_t<id, std::tuple<Grids...>>;
310 template<std::size_t id> using GridPtr = typename std::shared_ptr< Grid<id> >;
311
313 static constexpr std::size_t numGrids = sizeof...(Grids);
315 static constexpr int bulkGridId = 0;
316
321
323 template<std::size_t id>
324 const Grid<id>& grid() const
325 { return *std::get<id>(gridPtrTuple_); }
326
328 std::shared_ptr<const GridData> getGridData() const
329 {
330 if (!enableEntityMarkers_)
331 DUNE_THROW(Dune::IOError, "No grid data available");
332 return gridDataPtr_;
333 }
334
336 std::shared_ptr<const Embeddings> getEmbeddings() const
337 { return embeddingsPtr_; }
338
340 void init(const std::string& paramGroup = "")
341 {
342 // reset the grid & embedding data
343 gridDataPtr_ = std::make_shared<GridData>();
344 embeddingsPtr_ = std::make_shared<Embeddings>();
345
346 // get filename and determine grid file extension
347 const auto fileName = getParamFromGroup<std::string>(paramGroup, "Grid.File");
348 const auto ext = getFileExtension(fileName);
349
350 // get some parameters
351 const bool verbose = getParamFromGroup<bool>(paramGroup, "Grid.Verbosity", false);
352 const bool domainMarkers = getParamFromGroup<bool>(paramGroup, "Grid.DomainMarkers", false);
353 const bool boundarySegments = getParamFromGroup<bool>(paramGroup, "Grid.BoundarySegments", false);
354
355 // forward to the corresponding reader
356 if (ext == "msh")
357 {
358 const auto thresh = getParamFromGroup<std::size_t>(paramGroup, "Grid.GmshPhysicalEntityThreshold", 0);
360 gmshReader.read(fileName, (boundarySegments ? thresh : 0), verbose);
361 passDataFromReader(gmshReader, domainMarkers, boundarySegments);
362 }
363 else
364 DUNE_THROW(Dune::NotImplemented, "Reader for grid files of type ." + ext);
365
366 // find out if entity markers are active
367 enableEntityMarkers_ = domainMarkers || boundarySegments;
368 }
369
372 {
373 using namespace Dune::Hybrid;
374 forEach(integralRange(Dune::Hybrid::size(gridPtrTuple_)), [&](const auto id)
375 {
376 std::get<id>(this->gridPtrTuple_)->loadBalance();
377 });
378 }
379
380protected:
382 template<std::size_t id>
384 { return *std::get<id>(gridPtrTuple_); }
385
387 std::shared_ptr<Embeddings> getEmbeddings_()
388 { return embeddingsPtr_; }
389
390private:
392 static std::string getFileExtension(const std::string& fileName)
393 {
394 const auto pos = fileName.rfind('.', fileName.length());
395 if (pos != std::string::npos)
396 return(fileName.substr(pos+1, fileName.length() - pos));
397 else
398 DUNE_THROW(Dune::IOError, "Please provide an extension for your grid file ('"<< fileName << "')!");
399 }
400
402 template<typename MeshFileReader>
403 void passDataFromReader(MeshFileReader& reader, bool domainMarkers, bool boundarySegments)
404 {
405 const auto& vertices = reader.gridVertices();
406
407 using namespace Dune::Hybrid;
408 forEach(integralRange(Dune::Hybrid::size(gridPtrTuple_)), [&](const auto id)
409 {
410 using GridFactory = Dune::GridFactory<Grid<id>>;
411 auto factoryPtr = std::make_shared<GridFactory>();
412
413 // insert grid vertices
414 for (const auto idx : reader.vertexIndices(id))
415 factoryPtr->insertVertex(vertices[idx]);
416
417 // insert elements
418 for (const auto& e : reader.elementData(id))
419 factoryPtr->insertElement(e.gt, e.cornerIndices);
420
421 // insert boundary segments
422 if (boundarySegments)
423 for (const auto& segment : reader.boundarySegmentData(id))
424 factoryPtr->insertBoundarySegment(segment);
425
426 // make grid
427 auto gridPtr = std::shared_ptr<Grid<id>>(factoryPtr->createGrid());
428
429 // maybe create and set grid data object
430 if (domainMarkers || boundarySegments)
431 {
432 typename GridDataWrapper::template GridData<id> gridData( gridPtr,
433 factoryPtr,
434 std::move(reader.elementMarkerMap(id)),
435 std::move(reader.boundaryMarkerMap(id)) );
436 gridDataPtr_->template setGridData<id>( std::move(gridData) );
437 }
438
439 // copy the embeddings
440 embeddingsPtr_->template setData<id>( gridPtr,
441 factoryPtr,
442 std::move(reader.embeddedEntityMap(id)),
443 std::move(reader.adjoinedEntityMap(id)),
444 std::move(reader.vertexIndices(id)),
445 vertices.size() );
446
447 // set the grid pointer
448 std::get<id>(gridPtrTuple_) = gridPtr;
449 });
450 }
451
453 using Indices = std::make_index_sequence<numGrids>;
454 using GridPtrTuple = typename makeFromIndexedType<std::tuple, GridPtr, Indices>::type;
455 GridPtrTuple gridPtrTuple_;
456
458 bool enableEntityMarkers_;
459 std::shared_ptr<GridData> gridDataPtr_;
460
462 std::shared_ptr<Embeddings> embeddingsPtr_;
463};
464
465} // end namespace Dumux
466
467#endif
Utilities for template meta programming.
Defines the index types used for grid and local indices.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Class for grid data attached to dgf or gmsh grid files.
Definition: adapt.hh:29
static constexpr bool evalCondition()
Definition: multidomain/facet/gridmanager.hh:57
Struture to define the index types used for grid and local indices.
Definition: indextraits.hh:38
Definition: utility.hh:40
Reads gmsh files where (n-1)-dimensional grids are defined on the faces or edges of n-dimensional gri...
Definition: gmshreader.hh:61
void read(const std::string &fileName, bool verbose=false)
Definition: gmshreader.hh:84
Definition: multidomain/facet/gridmanager.hh:64
typename std::tuple_element_t< 0, std::tuple< Gs... > > G2
Definition: multidomain/facet/gridmanager.hh:68
Grid data object to store element and boundary segment markers for all grids of the hierarchy.
Definition: multidomain/facet/gridmanager.hh:85
int getBoundaryDomainMarker(int boundarySegmentIndex) const
Returns the boundary marker for a given bounday segment index.
Definition: multidomain/facet/gridmanager.hh:127
std::shared_ptr< const GridData< id > > getSubDomainGridData() const
return the grid data for a specific grid
Definition: multidomain/facet/gridmanager.hh:112
GridData< Grid< id > > GridData
export the i-th grid data type
Definition: multidomain/facet/gridmanager.hh:101
bool wasInserted(const Intersection< id > &intersection) const
Returns true if an intersection was inserted during grid creation.
Definition: multidomain/facet/gridmanager.hh:132
int getElementDomainMarker(const typename Grid< id >::template Codim< 0 >::Entity &element) const
Returns domain marker of an element.
Definition: multidomain/facet/gridmanager.hh:117
void setGridData(GridData< id > &&gridData)
set the grid data object for the i-th grid
Definition: multidomain/facet/gridmanager.hh:105
int getBoundaryDomainMarker(const typename Grid< id >::LeafGridView::Intersection &is) const
Returns the boundary marker of an intersection.
Definition: multidomain/facet/gridmanager.hh:122
Contains the embeddings between grids with codimension one among the grid hierarchy....
Definition: multidomain/facet/gridmanager.hh:155
std::unordered_map< GridIndexType, std::vector< GridIndexType > > & embeddedEntityMap(std::size_t id)
Returns non-const reference to maps of the embedded entities.
Definition: multidomain/facet/gridmanager.hh:223
std::size_t numVerticesInHierarchy() const
Returns the number of vertices contained in the entire grid hierarch.
Definition: multidomain/facet/gridmanager.hh:239
const std::vector< GridIndexType > & gridHierarchyIndices(std::size_t id) const
Returns the hierachy's insertion indices that make up the grid for the given id.
Definition: multidomain/facet/gridmanager.hh:235
static constexpr int bulkDim
state the dimension of the highest-dimensional grid
Definition: multidomain/facet/gridmanager.hh:179
GridView< bulkGridId > BulkGridView
export the bulk grid type
Definition: multidomain/facet/gridmanager.hh:182
const std::unordered_map< GridIndexType, std::vector< GridIndexType > > & adjoinedEntityMap(std::size_t id) const
Returns const reference to the maps of the adjoined entities of dimension d+1.
Definition: multidomain/facet/gridmanager.hh:227
std::unordered_map< GridIndexType, std::vector< GridIndexType > >::mapped_type adjoinedEntityIndices(const typename Grid< id >::template Codim< 0 >::Entity &element) const
Returns the insertion indices of the entities in which the element is embedded.
Definition: multidomain/facet/gridmanager.hh:210
static constexpr int bulkGridId
export the grid id of the bulk grid (descending grid dim -> always zero!)
Definition: multidomain/facet/gridmanager.hh:177
GridIndexType insertionIndex(const Entity &entity) const
return the insertion index of an entity of the i-th grid
Definition: multidomain/facet/gridmanager.hh:193
std::unordered_map< GridIndexType, std::vector< GridIndexType > >::mapped_type embeddedEntityIndices(const typename Grid< id >::template Codim< 0 >::Entity &element) const
Returns the insertion indices of the entities embedded in given element.
Definition: multidomain/facet/gridmanager.hh:199
typename Grid< id >::LeafGridView GridView
export the i-th grid view type
Definition: multidomain/facet/gridmanager.hh:172
void setData(std::shared_ptr< Grid< id > > gridPtr, std::shared_ptr< GridFactory< id > > gridFactoryPtr, EmbedmentMap &&embeddedEntityMap, EmbedmentMap &&adjoinedEntityMap, std::vector< GridIndexType > &&gridVertexIndices, std::size_t numVerticesInHierarchy)
Sets the required data for a specific grid on the hierarchy.
Definition: multidomain/facet/gridmanager.hh:252
GIType GridIndexType
export the type used for indices
Definition: multidomain/facet/gridmanager.hh:184
const std::unordered_map< GridIndexType, std::vector< GridIndexType > > & embeddedEntityMap(std::size_t id) const
Returns const reference to maps of the embedded entities.
Definition: multidomain/facet/gridmanager.hh:219
static constexpr std::size_t numGrids
export the number of created grids
Definition: multidomain/facet/gridmanager.hh:175
const GridView< id > & gridView() const
return reference to the i-th grid view
Definition: multidomain/facet/gridmanager.hh:188
std::unordered_map< GridIndexType, std::vector< GridIndexType > > & adjoinedEntityMap(std::size_t id)
Returns non-const reference to the maps of the adjoined entities of dimension d+1.
Definition: multidomain/facet/gridmanager.hh:231
Creates the grids in the context of hybrid-dimensional coupled models, where the (n-1)-dimensional do...
Definition: multidomain/facet/gridmanager.hh:299
typename std::tuple_element_t< id, std::tuple< Grids... > > Grid
export the i-th grid type
Definition: multidomain/facet/gridmanager.hh:308
typename std::shared_ptr< Grid< id > > GridPtr
export the i-th grid pointer type
Definition: multidomain/facet/gridmanager.hh:310
static constexpr std::size_t numGrids
export the number of created grids
Definition: multidomain/facet/gridmanager.hh:313
std::shared_ptr< const Embeddings > getEmbeddings() const
return a pointer to the object containing embeddings
Definition: multidomain/facet/gridmanager.hh:336
Grid< id > & grid_()
return non-const reference to i-th grid
Definition: multidomain/facet/gridmanager.hh:383
void loadBalance()
Distributes the grid on all processes of a parallel computation.
Definition: multidomain/facet/gridmanager.hh:371
std::shared_ptr< Embeddings > getEmbeddings_()
return non-const pointer to the object containing embeddings
Definition: multidomain/facet/gridmanager.hh:387
static constexpr int bulkGridId
export the grid id of the bulk grid (descending grid dim -> always zero!)
Definition: multidomain/facet/gridmanager.hh:315
const Grid< id > & grid() const
returns the i-th grid
Definition: multidomain/facet/gridmanager.hh:324
std::shared_ptr< const GridData > getGridData() const
return a pointer to the grid data object
Definition: multidomain/facet/gridmanager.hh:328
void init(const std::string &paramGroup="")
creates the grids from a file given in parameter tree
Definition: multidomain/facet/gridmanager.hh:340