27#ifndef DUMUX_FACETCOUPLING_GRID_MANAGER_HH
28#define DUMUX_FACETCOUPLING_GRID_MANAGER_HH
36#include <dune/common/exceptions.hh>
37#include <dune/common/hybridutilities.hh>
38#include <dune/grid/common/gridfactory.hh>
48namespace FCGridManagerChecks {
56 template<
bool checkDimWorld,
typename G1,
typename G2>
59 return checkDimWorld ? int(G1::dimensionworld) == int(G2::dimensionworld)
60 : int(G1::dimension) > int(G2::dimension);
65 template<
bool checkDimWorld,
typename G1,
typename... Gs>
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;
71 template<
bool checkDimWorld,
typename G1,
typename G2>
72 struct FulfillConditions<checkDimWorld, G1, G2> {
static constexpr bool value = evalCondition<checkDimWorld, G1, G2>(); };
83template<
typename... Grids>
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;
104 template<std::
size_t id>
107 std::get<id>(gridDataPtrTuple_) = std::make_shared<GridData<id>>( std::move(gridData) );
111 template<std::
size_t id>
113 {
return std::get<id>(gridDataPtrTuple_); }
116 template<std::
size_t id>
118 {
return std::get<id>(gridDataPtrTuple_)->getElementDomainMarker(element); }
121 template<std::
size_t id>
123 {
return std::get<id>(gridDataPtrTuple_)->getBoundaryDomainMarker(is); }
126 template<std::
size_t id>
128 {
return std::get<id>(gridDataPtrTuple_)->getBoundaryDomainMarker(boundarySegmentIndex); }
131 template<std::
size_t id>
133 {
return std::get<id>(gridDataPtrTuple_)->wasInserted(intersection); }
138 GridDataPtrTuple gridDataPtrTuple_;
153template<
typename... Grids>
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>>;
168 using EmbedmentMap = std::unordered_map<GIType, std::vector<GIType>>;
172 template<std::
size_t id>
using GridView =
typename Grid<id>::LeafGridView;
175 static constexpr std::size_t
numGrids =
sizeof...(Grids);
179 static constexpr int bulkDim = Grid<bulkGridId>::dimension;
187 template<std::
size_t id>
189 {
return *std::get<id>(gridViewPtrTuple_); }
192 template<std::
size_t id,
class Entity>
194 {
return std::get<id>(gridFactoryPtrTuple_)->insertionIndex(entity); }
197 template<std::
size_t id>
198 typename std::unordered_map< GridIndexType, std::vector<GridIndexType> >::mapped_type
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();
208 template<std::
size_t id>
209 typename std::unordered_map< GridIndexType, std::vector<GridIndexType> >::mapped_type
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();
219 const std::unordered_map< GridIndexType, std::vector<GridIndexType> >&
embeddedEntityMap(std::size_t
id)
const
220 { assert(
id <
numGrids);
return embeddedEntityMaps_[id]; }
223 std::unordered_map< GridIndexType, std::vector<GridIndexType> >&
embeddedEntityMap(std::size_t
id)
224 { assert(
id <
numGrids);
return embeddedEntityMaps_[id]; }
227 const std::unordered_map< GridIndexType, std::vector<GridIndexType> >&
adjoinedEntityMap(std::size_t
id)
const
228 { assert(
id <
numGrids);
return adjoinedEntityMaps_[id]; }
231 std::unordered_map< GridIndexType, std::vector<GridIndexType> >&
adjoinedEntityMap(std::size_t
id)
232 { assert(
id <
numGrids);
return adjoinedEntityMaps_[id]; }
236 { assert(
id <
numGrids);
return gridVertexIndices_[id]; }
240 {
return numVerticesInHierarchy_; }
251 template<std::
size_t id>
252 void setData( std::shared_ptr<Grid<id>> gridPtr,
253 std::shared_ptr<GridFactory<id>> gridFactoryPtr,
256 std::vector<GridIndexType>&& gridVertexIndices,
259 std::get<id>(gridViewPtrTuple_) = std::make_shared<GridView<id>>(gridPtr->leafGridView());
260 std::get<id>(gridFactoryPtrTuple_) = gridFactoryPtr;
263 gridVertexIndices_[id] = std::move(gridVertexIndices);
269 std::array<EmbedmentMap, numGrids> embeddedEntityMaps_;
270 std::array<EmbedmentMap, numGrids> adjoinedEntityMaps_;
273 std::size_t numVerticesInHierarchy_;
274 std::array<std::vector<GridIndexType>,
numGrids> gridVertexIndices_;
277 using Indices = std::make_index_sequence<numGrids>;
278 template<std::
size_t id>
using GridViewPtr = std::shared_ptr<GridView<id>>;
280 GridPtrTuple gridViewPtrTuple_;
283 template<std::
size_t id>
using GridFactoryPtr = std::shared_ptr< Dune::GridFactory<Grid<id>> >;
285 GridFactoryPtrTuple gridFactoryPtrTuple_;
297template<
typename... Grids>
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> >;
313 static constexpr std::size_t
numGrids =
sizeof...(Grids);
323 template<std::
size_t id>
325 {
return *std::get<id>(gridPtrTuple_); }
330 if (!enableEntityMarkers_)
331 DUNE_THROW(Dune::IOError,
"No grid data available");
337 {
return embeddingsPtr_; }
340 void init(
const std::string& paramGroup =
"")
343 gridDataPtr_ = std::make_shared<GridData>();
344 embeddingsPtr_ = std::make_shared<Embeddings>();
347 const auto fileName = getParamFromGroup<std::string>(paramGroup,
"Grid.File");
348 const auto ext = getFileExtension(fileName);
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);
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);
364 DUNE_THROW(Dune::NotImplemented,
"Reader for grid files of type ." + ext);
367 enableEntityMarkers_ = domainMarkers || boundarySegments;
373 using namespace Dune::Hybrid;
374 forEach(integralRange(Dune::Hybrid::size(gridPtrTuple_)), [&](
const auto id)
376 std::get<id>(this->gridPtrTuple_)->loadBalance();
382 template<std::
size_t id>
384 {
return *std::get<id>(gridPtrTuple_); }
388 {
return embeddingsPtr_; }
392 static std::string getFileExtension(
const std::string& fileName)
394 const auto pos = fileName.rfind(
'.', fileName.length());
395 if (pos != std::string::npos)
396 return(fileName.substr(pos+1, fileName.length() - pos));
398 DUNE_THROW(Dune::IOError,
"Please provide an extension for your grid file ('"<< fileName <<
"')!");
402 template<
typename MeshFileReader>
403 void passDataFromReader(MeshFileReader& reader,
bool domainMarkers,
bool boundarySegments)
405 const auto& vertices = reader.gridVertices();
407 using namespace Dune::Hybrid;
408 forEach(integralRange(Dune::Hybrid::size(gridPtrTuple_)), [&](
const auto id)
410 using GridFactory = Dune::GridFactory<Grid<id>>;
411 auto factoryPtr = std::make_shared<GridFactory>();
414 for (
const auto idx : reader.vertexIndices(
id))
415 factoryPtr->insertVertex(vertices[idx]);
418 for (
const auto& e : reader.elementData(
id))
419 factoryPtr->insertElement(e.gt, e.cornerIndices);
422 if (boundarySegments)
423 for (
const auto& segment : reader.boundarySegmentData(
id))
424 factoryPtr->insertBoundarySegment(segment);
427 auto gridPtr = std::shared_ptr<Grid<id>>(factoryPtr->createGrid());
430 if (domainMarkers || boundarySegments)
432 typename GridDataWrapper::template GridData<id> gridData( gridPtr,
434 std::move(reader.elementMarkerMap(
id)),
435 std::move(reader.boundaryMarkerMap(
id)) );
436 gridDataPtr_->template setGridData<id>( std::move(gridData) );
440 embeddingsPtr_->template setData<id>( gridPtr,
442 std::move(reader.embeddedEntityMap(
id)),
443 std::move(reader.adjoinedEntityMap(
id)),
444 std::move(reader.vertexIndices(
id)),
448 std::get<id>(gridPtrTuple_) = gridPtr;
453 using Indices = std::make_index_sequence<numGrids>;
454 using GridPtrTuple =
typename makeFromIndexedType<std::tuple, GridPtr, Indices>::type;
455 GridPtrTuple gridPtrTuple_;
458 bool enableEntityMarkers_;
459 std::shared_ptr<GridData> gridDataPtr_;
462 std::shared_ptr<Embeddings> embeddingsPtr_;
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.
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 ¶mGroup="")
creates the grids from a file given in parameter tree
Definition: multidomain/facet/gridmanager.hh:340