25#ifndef DUMUX_IO_GRID_MANAGER_BASE_HH
26#define DUMUX_IO_GRID_MANAGER_BASE_HH
33#include <dune/common/exceptions.hh>
34#include <dune/common/classname.hh>
35#include <dune/common/parallel/communication.hh>
36#include <dune/common/parallel/mpihelper.hh>
37#include <dune/grid/io/file/dgfparser/dgfparser.hh>
38#include <dune/grid/io/file/gmshreader.hh>
39#include <dune/grid/common/gridfactory.hh>
40#include <dune/grid/utility/structuredgridfactory.hh>
65template <
class Gr
idType>
75 void init(
const std::string& modelParamGroup =
"")
78 "The header with the GridManager specialization for your grid type is not included "
79 "or no specialization has been implemented!"
80 " In case of the latter, consider providing your own GridManager."
100 if (Dune::MPIHelper::getCommunication().size() > 1)
114 auto dh =
gridData_->createGmshDataHandle();
115 gridPtr()->loadBalance(dh.interface());
116 gridPtr()->communicate(dh.interface(), Dune::InteriorBorder_All_Interface, Dune::ForwardCommunication);
130 DUNE_THROW(Dune::IOError,
"No grid data available");
151 DUNE_THROW(Dune::InvalidStateException,
"You are using DGF. To get the grid pointer use method dgfGridPtr()!");
162 DUNE_THROW(Dune::InvalidStateException,
"The DGF grid pointer is only available if the grid was constructed with a DGF file!");
170 std::size_t i = fileName.rfind(
'.', fileName.length());
171 if (i != std::string::npos)
173 return(fileName.substr(i+1, fileName.length() - i));
177 DUNE_THROW(Dune::IOError,
"Please provide and extension for your grid file ('"<< fileName <<
"')!");
189 const std::string& modelParamGroup)
193 if (extension !=
"dgf" && extension !=
"msh" && extension !=
"vtu" && extension !=
"vtp")
194 DUNE_THROW(Dune::IOError,
"Grid type " << Dune::className<Grid>() <<
" doesn't support grid files with extension: *."<< extension);
197 if (extension ==
"dgf")
200 dgfGridPtr() = Dune::GridPtr<Grid>(fileName.c_str(), Dune::MPIHelper::getCommunicator());
205 else if (extension ==
"msh")
208 const bool verbose = getParamFromGroup<bool>(modelParamGroup,
"Grid.Verbosity",
false);
209 const bool boundarySegments = getParamFromGroup<bool>(modelParamGroup,
"Grid.BoundarySegments",
false);
210 const bool domainMarkers = getParamFromGroup<bool>(modelParamGroup,
"Grid.DomainMarkers",
false);
218 std::vector<int> boundaryMarkers, elementMarkers;
219 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
220 Dune::GmshReader<Grid>::read(*gridFactory, fileName, boundaryMarkers, elementMarkers, verbose, boundarySegments);
221 gridPtr() = std::shared_ptr<Grid>(gridFactory->createGrid());
222 gridData_ = std::make_shared<GridData>(
gridPtr_, std::move(gridFactory), std::move(elementMarkers), std::move(boundaryMarkers));
226 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
227 Dune::GmshReader<Grid>::read(*gridFactory, fileName, verbose, boundarySegments);
228 gridPtr() = std::shared_ptr<Grid>(gridFactory->createGrid());
233 else if (extension ==
"vtu" || extension ==
"vtp")
235 if (Dune::MPIHelper::getCommunication().size() > 1)
236 DUNE_THROW(Dune::NotImplemented,
"Reading grids in parallel from VTK file formats is currently not supported!");
240 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
241 const bool verbose = getParamFromGroup<bool>(modelParamGroup,
"Grid.Verbosity",
false);
242 gridPtr() = vtkReader.
readGrid(*gridFactory, cellData, pointData, verbose);
243 gridData_ = std::make_shared<GridData>(
gridPtr_, std::move(gridFactory), std::move(cellData), std::move(pointData));
254 if(extension !=
"dgf")
255 DUNE_THROW(Dune::IOError,
"Grid type " << Dune::className<Grid>() <<
" only supports DGF (*.dgf) but the specified filename has extension: *."<< extension);
258 dgfGridPtr() = Dune::GridPtr<Grid>(fileName.c_str(), Dune::MPIHelper::getCommunicator());
270 template <
int dim,
int dimworld>
272 const std::string& modelParamGroup)
274 using GlobalPosition = Dune::FieldVector<typename Grid::ctype, dimworld>;
275 const auto upperRight = getParamFromGroup<GlobalPosition>(modelParamGroup,
"Grid.UpperRight");
276 const auto lowerLeft = getParamFromGroup<GlobalPosition>(modelParamGroup,
"Grid.LowerLeft", GlobalPosition(0.0));
278 using CellArray = std::array<unsigned int, dim>;
279 CellArray cells; cells.fill(1);
280 cells = getParamFromGroup<CellArray>(modelParamGroup,
"Grid.Cells", cells);
283 if (cellType == CellType::Cube)
285 gridPtr() = Dune::StructuredGridFactory<Grid>::createCubeGrid(lowerLeft, upperRight, cells);
287 else if (cellType == CellType::Simplex)
289 gridPtr() = Dune::StructuredGridFactory<Grid>::createSimplexGrid(lowerLeft, upperRight, cells);
293 DUNE_THROW(Dune::GridError,
"Unknown cell type for making structured grid! Choose Cube or Simplex.");
303 grid().globalRefine(getParamFromGroup<int>(modelParamGroup,
"Grid.Refinement"));
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
A vtk file reader using tinyxml2 as xml backend.
bool hasParamInGroup(const std::string ¶mGroup, const std::string ¶m)
Check whether a key exists in the parameter tree with a model group prefix.
Definition: parameters.hh:177
Template which always yields a false value.
Definition: typetraits.hh:35
Class for grid data attached to dgf or gmsh grid files.
Definition: griddata.hh:66
The grid manager (this is the class used by the user) for all supported grid managers that constructs...
Definition: gridmanager_base.hh:324
The grid manager base interface (public) and methods common to most grid manager specializations (pro...
Definition: gridmanager_base.hh:67
std::shared_ptr< GridData > getGridData() const
Get an owning pointer to grid data associated with the grid.
Definition: gridmanager_base.hh:127
GridType Grid
Definition: gridmanager_base.hh:69
Dune::GridPtr< Grid > dgfGridPtr_
Definition: gridmanager_base.hh:318
void maybeRefineGrid(const std::string &modelParamGroup)
Refines a grid after construction if GridParameterGroup.Refinement is set in the input file.
Definition: gridmanager_base.hh:300
void makeGridFromFile(const std::string &fileName, const std::string &modelParamGroup)
Makes a grid from a file. We currently support.
Definition: gridmanager_base.hh:188
void makeStructuredGrid(CellType cellType, const std::string &modelParamGroup)
Makes a structured cube grid using the structured grid factory.
Definition: gridmanager_base.hh:271
Dune::GridPtr< Grid > & dgfGridPtr()
Returns a reference to the DGF grid pointer (Dune::GridPtr<Grid>).
Definition: gridmanager_base.hh:157
bool hasGridData() const
Check whether there is data associated with the grid.
Definition: gridmanager_base.hh:138
std::shared_ptr< Grid > gridPtr_
Definition: gridmanager_base.hh:317
Grid & grid()
Returns a reference to the grid.
Definition: gridmanager_base.hh:87
std::shared_ptr< GridData > gridData_
Definition: gridmanager_base.hh:320
void makeGridFromDgfFile(const std::string &fileName)
Makes a grid from a DGF file. This is used by grid managers that only support DGF.
Definition: gridmanager_base.hh:250
bool enableGmshDomainMarkers_
A state variable if domain markers have been read from a Gmsh file.
Definition: gridmanager_base.hh:315
void loadBalance()
Call loadBalance() function of the grid.
Definition: gridmanager_base.hh:98
CellType
The cell types for structured grids.
Definition: gridmanager_base.hh:265
@ Simplex
Definition: gridmanager_base.hh:265
@ Cube
Definition: gridmanager_base.hh:265
std::shared_ptr< Grid > & gridPtr()
Returns a reference to the grid pointer (std::shared_ptr<Grid>)
Definition: gridmanager_base.hh:146
void init(const std::string &modelParamGroup="")
Make the grid. Implement this method in the specialization of this class for a grid type.
Definition: gridmanager_base.hh:75
bool enableDgfGridPointer_
A state variable if the DGF Dune::GridPtr has been enabled. It is always enabled if a DGF grid file w...
Definition: gridmanager_base.hh:310
std::string getFileExtension(const std::string &fileName) const
Returns the filename extension of a given filename.
Definition: gridmanager_base.hh:168
A vtk file reader using tinyxml2 as xml backend.
Definition: vtkreader.hh:51
std::unique_ptr< Grid > readGrid(bool verbose=false) const
Read a grid from a vtk/vtu/vtp file, ignoring cell and point data.
Definition: vtkreader.hh:128
std::unordered_map< std::string, std::vector< double > > Data
the cell / point data type for point data read from a grid file
Definition: vtkreader.hh:59
Class for grid data attached to dgf or gmsh grid files.