29#ifndef DUMUX_FACETCOUPLING_GRID_MANAGER_HH
30#define DUMUX_FACETCOUPLING_GRID_MANAGER_HH
38#include <dune/common/exceptions.hh>
39#include <dune/common/hybridutilities.hh>
40#include <dune/grid/common/gridfactory.hh>
50namespace FCGridManagerChecks {
58 template<
bool checkDimWorld,
typename G1,
typename G2>
61 return checkDimWorld ? int(G1::dimensionworld) == int(G2::dimensionworld)
62 : int(G1::dimension) > int(G2::dimension);
67 template<
bool checkDimWorld,
typename G1,
typename... Gs>
70 using G2 =
typename std::tuple_element_t<0, std::tuple<Gs...>>;
71 static constexpr bool value = evalCondition<checkDimWorld, G1, G2>() &&
FulfillConditions<checkDimWorld, Gs...>::value;
73 template<
bool checkDimWorld,
typename G1,
typename G2>
74 struct FulfillConditions<checkDimWorld, G1, G2> {
static constexpr bool value = evalCondition<checkDimWorld, G1, G2>(); };
85template<
typename... Grids>
93 static constexpr std::size_t numGrids =
sizeof...(Grids);
95 template<std::
size_t id>
using Grid =
typename std::tuple_element_t<id, std::tuple<Grids...>>;
97 template<std::
size_t id>
using GridDataPtr = std::shared_ptr< GridData<Grid<id>> >;
99 template<std::
size_t id>
using Intersection =
typename Grid<id>::LeafGridView::Intersection;
106 template<std::
size_t id>
109 std::get<id>(gridDataPtrTuple_) = std::make_shared<GridData<id>>( std::move(gridData) );
113 template<std::
size_t id>
115 {
return std::get<id>(gridDataPtrTuple_); }
118 template<std::
size_t id>
120 {
return std::get<id>(gridDataPtrTuple_)->getElementDomainMarker(element); }
123 template<std::
size_t id>
125 {
return std::get<id>(gridDataPtrTuple_)->getBoundaryDomainMarker(is); }
128 template<std::
size_t id>
130 {
return std::get<id>(gridDataPtrTuple_)->getBoundaryDomainMarker(boundarySegmentIndex); }
133 template<std::
size_t id>
135 {
return std::get<id>(gridDataPtrTuple_)->wasInserted(intersection); }
140 GridDataPtrTuple gridDataPtrTuple_;
155template<
typename... Grids>
163 template<std::
size_t id>
using Grid =
typename std::tuple_element_t<id, std::tuple<Grids...>>;
165 template<std::
size_t id>
using GridFactory =
typename Dune::GridFactory<Grid<id>>;
170 using EmbedmentMap = std::unordered_map<GIType, std::vector<GIType>>;
174 template<std::
size_t id>
using GridView =
typename Grid<id>::LeafGridView;
177 static constexpr std::size_t
numGrids =
sizeof...(Grids);
181 static constexpr int bulkDim = Grid<bulkGridId>::dimension;
189 template<std::
size_t id>
191 {
return *std::get<id>(gridViewPtrTuple_); }
194 template<std::
size_t id,
class Entity>
196 {
return std::get<id>(gridFactoryPtrTuple_)->insertionIndex(entity); }
199 template<std::
size_t id>
200 typename std::unordered_map< GridIndexType, std::vector<GridIndexType> >::mapped_type
203 const auto& map = embeddedEntityMaps_[id];
204 auto it = map.find( std::get<id>(gridFactoryPtrTuple_)->
insertionIndex(element) );
205 if (it != map.end())
return it->second;
206 else return typename std::unordered_map< GridIndexType, std::vector<GridIndexType> >::mapped_type();
210 template<std::
size_t id>
211 typename std::unordered_map< GridIndexType, std::vector<GridIndexType> >::mapped_type
214 const auto& map = adjoinedEntityMaps_[id];
215 auto it = map.find( std::get<id>(gridFactoryPtrTuple_)->
insertionIndex(element) );
216 if (it != map.end())
return it->second;
217 else return typename std::unordered_map< GridIndexType, std::vector<GridIndexType> >::mapped_type();
221 const std::unordered_map< GridIndexType, std::vector<GridIndexType> >&
embeddedEntityMap(std::size_t
id)
const
222 { assert(
id <
numGrids);
return embeddedEntityMaps_[id]; }
225 std::unordered_map< GridIndexType, std::vector<GridIndexType> >&
embeddedEntityMap(std::size_t
id)
226 { assert(
id <
numGrids);
return embeddedEntityMaps_[id]; }
229 const std::unordered_map< GridIndexType, std::vector<GridIndexType> >&
adjoinedEntityMap(std::size_t
id)
const
230 { assert(
id <
numGrids);
return adjoinedEntityMaps_[id]; }
233 std::unordered_map< GridIndexType, std::vector<GridIndexType> >&
adjoinedEntityMap(std::size_t
id)
234 { assert(
id <
numGrids);
return adjoinedEntityMaps_[id]; }
238 { assert(
id <
numGrids);
return gridVertexIndices_[id]; }
242 {
return numVerticesInHierarchy_; }
253 template<std::
size_t id>
254 void setData( std::shared_ptr<Grid<id>> gridPtr,
255 std::shared_ptr<GridFactory<id>> gridFactoryPtr,
258 std::vector<GridIndexType>&& gridVertexIndices,
261 std::get<id>(gridViewPtrTuple_) = std::make_shared<GridView<id>>(gridPtr->leafGridView());
262 std::get<id>(gridFactoryPtrTuple_) = gridFactoryPtr;
265 gridVertexIndices_[id] = std::move(gridVertexIndices);
271 std::array<EmbedmentMap, numGrids> embeddedEntityMaps_;
272 std::array<EmbedmentMap, numGrids> adjoinedEntityMaps_;
275 std::size_t numVerticesInHierarchy_;
276 std::array<std::vector<GridIndexType>,
numGrids> gridVertexIndices_;
279 using Indices = std::make_index_sequence<numGrids>;
280 template<std::
size_t id>
using GridViewPtr = std::shared_ptr<GridView<id>>;
282 GridPtrTuple gridViewPtrTuple_;
285 template<std::
size_t id>
using GridFactoryPtr = std::shared_ptr< Dune::GridFactory<Grid<id>> >;
287 GridFactoryPtrTuple gridFactoryPtrTuple_;
299template<
typename... Grids>
310 template<std::
size_t id>
using Grid =
typename std::tuple_element_t<id, std::tuple<Grids...>>;
312 template<std::
size_t id>
using GridPtr =
typename std::shared_ptr< Grid<id> >;
315 static constexpr std::size_t
numGrids =
sizeof...(Grids);
325 template<std::
size_t id>
327 {
return *std::get<id>(gridPtrTuple_); }
332 if (!enableEntityMarkers_)
333 DUNE_THROW(Dune::IOError,
"No grid data available");
339 {
return embeddingsPtr_; }
342 void init(
const std::string& paramGroup =
"")
345 gridDataPtr_ = std::make_shared<GridData>();
346 embeddingsPtr_ = std::make_shared<Embeddings>();
349 const auto fileName = getParamFromGroup<std::string>(paramGroup,
"Grid.File");
350 const auto ext = getFileExtension(fileName);
353 const bool verbose = getParamFromGroup<bool>(paramGroup,
"Grid.Verbosity",
false);
354 const bool domainMarkers = getParamFromGroup<bool>(paramGroup,
"Grid.DomainMarkers",
false);
355 const bool boundarySegments = getParamFromGroup<bool>(paramGroup,
"Grid.BoundarySegments",
false);
360 const auto thresh = getParamFromGroup<std::size_t>(paramGroup,
"Grid.GmshPhysicalEntityThreshold", 0);
362 gmshReader.
read(fileName, (boundarySegments ? thresh : 0), verbose);
363 passDataFromReader(gmshReader, domainMarkers, boundarySegments);
366 DUNE_THROW(Dune::NotImplemented,
"Reader for grid files of type ." + ext);
369 enableEntityMarkers_ = domainMarkers || boundarySegments;
375 using namespace Dune::Hybrid;
376 forEach(integralRange(Dune::Hybrid::size(gridPtrTuple_)), [&](
const auto id)
378 std::get<id>(this->gridPtrTuple_)->loadBalance();
384 template<std::
size_t id>
386 {
return *std::get<id>(gridPtrTuple_); }
390 {
return embeddingsPtr_; }
394 static std::string getFileExtension(
const std::string& fileName)
396 const auto pos = fileName.rfind(
'.', fileName.length());
397 if (pos != std::string::npos)
398 return(fileName.substr(pos+1, fileName.length() - pos));
400 DUNE_THROW(Dune::IOError,
"Please provide an extension for your grid file ('"<< fileName <<
"')!");
404 template<
typename MeshFileReader>
405 void passDataFromReader(MeshFileReader& reader,
bool domainMarkers,
bool boundarySegments)
407 const auto& vertices = reader.gridVertices();
409 using namespace Dune::Hybrid;
410 forEach(integralRange(Dune::Hybrid::size(gridPtrTuple_)), [&](
const auto id)
412 using GridFactory = Dune::GridFactory<Grid<id>>;
413 auto factoryPtr = std::make_shared<GridFactory>();
416 for (
const auto idx : reader.vertexIndices(
id))
417 factoryPtr->insertVertex(vertices[idx]);
420 for (
const auto& e : reader.elementData(
id))
421 factoryPtr->insertElement(e.gt, e.cornerIndices);
424 if (boundarySegments)
425 for (
const auto& segment : reader.boundarySegmentData(
id))
426 factoryPtr->insertBoundarySegment(segment);
429 auto gridPtr = std::shared_ptr<Grid<id>>(factoryPtr->createGrid());
432 if (domainMarkers || boundarySegments)
434 typename GridDataWrapper::template GridData<id> gridData( gridPtr,
436 std::move(reader.elementMarkerMap(
id)),
437 std::move(reader.boundaryMarkerMap(
id)) );
438 gridDataPtr_->template setGridData<id>( std::move(gridData) );
442 embeddingsPtr_->template setData<id>( gridPtr,
444 std::move(reader.embeddedEntityMap(
id)),
445 std::move(reader.adjoinedEntityMap(
id)),
446 std::move(reader.vertexIndices(
id)),
450 std::get<id>(gridPtrTuple_) = gridPtr;
455 using Indices = std::make_index_sequence<numGrids>;
456 using GridPtrTuple =
typename makeFromIndexedType<std::tuple, GridPtr, Indices>::type;
457 GridPtrTuple gridPtrTuple_;
460 bool enableEntityMarkers_;
461 std::shared_ptr<GridData> gridDataPtr_;
464 std::shared_ptr<Embeddings> embeddingsPtr_;
Utilities for template meta programming.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Defines the index types used for grid and local indices.
Class for grid data attached to dgf or gmsh grid files.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
static constexpr bool evalCondition()
Definition: multidomain/facet/gridmanager.hh:59
Structure 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:62
void read(const std::string &fileName, bool verbose=false)
Definition: gmshreader.hh:85
Definition: multidomain/facet/gridmanager.hh:66
typename std::tuple_element_t< 0, std::tuple< Gs... > > G2
Definition: multidomain/facet/gridmanager.hh:70
Grid data object to store element and boundary segment markers for all grids of the hierarchy.
Definition: multidomain/facet/gridmanager.hh:87
int getBoundaryDomainMarker(int boundarySegmentIndex) const
Returns the boundary marker for a given boundary segment index.
Definition: multidomain/facet/gridmanager.hh:129
std::shared_ptr< const GridData< id > > getSubDomainGridData() const
return the grid data for a specific grid
Definition: multidomain/facet/gridmanager.hh:114
GridData< Grid< id > > GridData
export the i-th grid data type
Definition: multidomain/facet/gridmanager.hh:103
bool wasInserted(const Intersection< id > &intersection) const
Returns true if an intersection was inserted during grid creation.
Definition: multidomain/facet/gridmanager.hh:134
int getElementDomainMarker(const typename Grid< id >::template Codim< 0 >::Entity &element) const
Returns domain marker of an element.
Definition: multidomain/facet/gridmanager.hh:119
void setGridData(GridData< id > &&gridData)
set the grid data object for the i-th grid
Definition: multidomain/facet/gridmanager.hh:107
int getBoundaryDomainMarker(const typename Grid< id >::LeafGridView::Intersection &is) const
Returns the boundary marker of an intersection.
Definition: multidomain/facet/gridmanager.hh:124
Contains the embeddings between grids with codimension one among the grid hierarchy....
Definition: multidomain/facet/gridmanager.hh:157
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:225
std::size_t numVerticesInHierarchy() const
Returns the number of vertices contained in the entire grid hierarch.
Definition: multidomain/facet/gridmanager.hh:241
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:237
static constexpr int bulkDim
state the dimension of the highest-dimensional grid
Definition: multidomain/facet/gridmanager.hh:181
GridView< bulkGridId > BulkGridView
export the bulk grid type
Definition: multidomain/facet/gridmanager.hh:184
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:229
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:212
static constexpr int bulkGridId
export the grid id of the bulk grid (descending grid dim -> always zero!)
Definition: multidomain/facet/gridmanager.hh:179
GridIndexType insertionIndex(const Entity &entity) const
return the insertion index of an entity of the i-th grid
Definition: multidomain/facet/gridmanager.hh:195
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:201
typename Grid< id >::LeafGridView GridView
export the i-th grid view type
Definition: multidomain/facet/gridmanager.hh:174
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:254
GIType GridIndexType
export the type used for indices
Definition: multidomain/facet/gridmanager.hh:186
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:221
static constexpr std::size_t numGrids
export the number of created grids
Definition: multidomain/facet/gridmanager.hh:177
const GridView< id > & gridView() const
return reference to the i-th grid view
Definition: multidomain/facet/gridmanager.hh:190
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:233
Creates the grids in the context of hybrid-dimensional coupled models, where the (n-1)-dimensional do...
Definition: multidomain/facet/gridmanager.hh:301
typename std::tuple_element_t< id, std::tuple< Grids... > > Grid
export the i-th grid type
Definition: multidomain/facet/gridmanager.hh:310
typename std::shared_ptr< Grid< id > > GridPtr
export the i-th grid pointer type
Definition: multidomain/facet/gridmanager.hh:312
static constexpr std::size_t numGrids
export the number of created grids
Definition: multidomain/facet/gridmanager.hh:315
std::shared_ptr< const Embeddings > getEmbeddings() const
return a pointer to the object containing embeddings
Definition: multidomain/facet/gridmanager.hh:338
Grid< id > & grid_()
return non-const reference to i-th grid
Definition: multidomain/facet/gridmanager.hh:385
void loadBalance()
Distributes the grid on all processes of a parallel computation.
Definition: multidomain/facet/gridmanager.hh:373
std::shared_ptr< Embeddings > getEmbeddings_()
return non-const pointer to the object containing embeddings
Definition: multidomain/facet/gridmanager.hh:389
static constexpr int bulkGridId
export the grid id of the bulk grid (descending grid dim -> always zero!)
Definition: multidomain/facet/gridmanager.hh:317
const Grid< id > & grid() const
returns the i-th grid
Definition: multidomain/facet/gridmanager.hh:326
std::shared_ptr< const GridData > getGridData() const
return a pointer to the grid data object
Definition: multidomain/facet/gridmanager.hh:330
void init(const std::string ¶mGroup="")
creates the grids from a file given in parameter tree
Definition: multidomain/facet/gridmanager.hh:342