13#ifndef DUMUX_IO_GRID_MANAGER_BASE_HH
14#define DUMUX_IO_GRID_MANAGER_BASE_HH
21#include <dune/common/exceptions.hh>
22#include <dune/common/classname.hh>
23#include <dune/common/parallel/communication.hh>
24#include <dune/common/parallel/mpihelper.hh>
25#include <dune/grid/io/file/dgfparser/dgfparser.hh>
26#include <dune/grid/io/file/gmshreader.hh>
27#include <dune/grid/common/gridfactory.hh>
28#include <dune/grid/utility/structuredgridfactory.hh>
53template <
class Gr
idType>
63 void init(
const std::string& modelParamGroup =
"")
66 "The header with the GridManager specialization for your grid type is not included "
67 "or no specialization has been implemented!"
68 " In case of the latter, consider providing your own GridManager."
99 if (Dune::MPIHelper::getCommunication().size() > 1)
113 auto dh =
gridData_->createGmshDataHandle();
114 gridPtr()->loadBalance(dh.interface());
115 gridPtr()->communicate(dh.interface(), Dune::InteriorBorder_All_Interface, Dune::ForwardCommunication);
129 DUNE_THROW(Dune::IOError,
"No grid data available");
150 DUNE_THROW(Dune::InvalidStateException,
"You are using DGF. To get the grid pointer use method dgfGridPtr()!");
161 DUNE_THROW(Dune::InvalidStateException,
"The DGF grid pointer is only available if the grid was constructed with a DGF file!");
169 std::size_t i = fileName.rfind(
'.', fileName.length());
170 if (i != std::string::npos)
172 return(fileName.substr(i+1, fileName.length() - i));
176 DUNE_THROW(Dune::IOError,
"Please provide and extension for your grid file ('"<< fileName <<
"')!");
188 const std::string& modelParamGroup)
192 if (extension !=
"dgf" && extension !=
"msh" && extension !=
"vtu" && extension !=
"vtp")
193 DUNE_THROW(Dune::IOError,
"Grid type " << Dune::className<Grid>() <<
" doesn't support grid files with extension: *."<< extension);
196 if (extension ==
"dgf")
199 dgfGridPtr() = Dune::GridPtr<Grid>(fileName.c_str(), Dune::MPIHelper::getCommunicator());
204 else if (extension ==
"msh")
207 const bool verbose = getParamFromGroup<bool>(modelParamGroup,
"Grid.Verbosity",
false);
208 const bool boundarySegments = getParamFromGroup<bool>(modelParamGroup,
"Grid.BoundarySegments",
false);
209 const bool domainMarkers = getParamFromGroup<bool>(modelParamGroup,
"Grid.DomainMarkers",
false);
217 std::vector<int> boundaryMarkers, elementMarkers;
218 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
219 Dune::GmshReader<Grid>::read(*gridFactory, fileName, boundaryMarkers, elementMarkers, verbose, boundarySegments);
220 gridPtr() = std::shared_ptr<Grid>(gridFactory->createGrid());
221 gridData_ = std::make_shared<GridData>(
gridPtr_, std::move(gridFactory), std::move(elementMarkers), std::move(boundaryMarkers));
225 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
226 Dune::GmshReader<Grid>::read(*gridFactory, fileName, verbose, boundarySegments);
227 gridPtr() = std::shared_ptr<Grid>(gridFactory->createGrid());
232 else if (extension ==
"vtu" || extension ==
"vtp")
234 if (Dune::MPIHelper::getCommunication().size() > 1)
235 DUNE_THROW(Dune::NotImplemented,
"Reading grids in parallel from VTK file formats is currently not supported!");
239 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
240 const bool verbose = getParamFromGroup<bool>(modelParamGroup,
"Grid.Verbosity",
false);
241 gridPtr() = vtkReader.
readGrid(*gridFactory, cellData, pointData, verbose);
242 gridData_ = std::make_shared<GridData>(
gridPtr_, std::move(gridFactory), std::move(cellData), std::move(pointData));
253 if(extension !=
"dgf")
254 DUNE_THROW(Dune::IOError,
"Grid type " << Dune::className<Grid>() <<
" only supports DGF (*.dgf) but the specified filename has extension: *."<< extension);
257 dgfGridPtr() = Dune::GridPtr<Grid>(fileName.c_str(), Dune::MPIHelper::getCommunicator());
269 template <
int dim,
int dimworld>
271 const std::string& modelParamGroup)
273 using GlobalPosition = Dune::FieldVector<typename Grid::ctype, dimworld>;
274 const auto upperRight = getParamFromGroup<GlobalPosition>(modelParamGroup,
"Grid.UpperRight");
275 const auto lowerLeft = getParamFromGroup<GlobalPosition>(modelParamGroup,
"Grid.LowerLeft", GlobalPosition(0.0));
277 using CellArray = std::array<unsigned int, dim>;
278 CellArray cells; cells.fill(1);
279 cells = getParamFromGroup<CellArray>(modelParamGroup,
"Grid.Cells", cells);
282 if (cellType == CellType::Cube)
284 gridPtr() = Dune::StructuredGridFactory<Grid>::createCubeGrid(lowerLeft, upperRight, cells);
286 else if (cellType == CellType::Simplex)
288 gridPtr() = Dune::StructuredGridFactory<Grid>::createSimplexGrid(lowerLeft, upperRight, cells);
292 DUNE_THROW(Dune::GridError,
"Unknown cell type for making structured grid! Choose Cube or Simplex.");
302 grid().globalRefine(getParamFromGroup<int>(modelParamGroup,
"Grid.Refinement"));
Class for grid data attached to dgf or gmsh grid files.
Definition: griddata.hh:54
The grid manager base interface (public) and methods common to most grid manager specializations (pro...
Definition: gridmanager_base.hh:55
std::shared_ptr< GridData > getGridData() const
Get an owning pointer to grid data associated with the grid.
Definition: gridmanager_base.hh:126
GridType Grid
Definition: gridmanager_base.hh:57
Dune::GridPtr< Grid > dgfGridPtr_
Definition: gridmanager_base.hh:317
const Grid & grid() const
Returns a const reference to the grid.
Definition: gridmanager_base.hh:86
void maybeRefineGrid(const std::string &modelParamGroup)
Refines a grid after construction if GridParameterGroup.Refinement is set in the input file.
Definition: gridmanager_base.hh:299
void makeGridFromFile(const std::string &fileName, const std::string &modelParamGroup)
Makes a grid from a file. We currently support.
Definition: gridmanager_base.hh:187
void makeStructuredGrid(CellType cellType, const std::string &modelParamGroup)
Makes a structured cube grid using the structured grid factory.
Definition: gridmanager_base.hh:270
Dune::GridPtr< Grid > & dgfGridPtr()
Returns a reference to the DGF grid pointer (Dune::GridPtr<Grid>).
Definition: gridmanager_base.hh:156
bool hasGridData() const
Check whether there is data associated with the grid.
Definition: gridmanager_base.hh:137
std::shared_ptr< Grid > gridPtr_
Definition: gridmanager_base.hh:316
Grid & grid()
Returns a reference to the grid.
Definition: gridmanager_base.hh:75
std::shared_ptr< GridData > gridData_
Definition: gridmanager_base.hh:319
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:249
bool enableGmshDomainMarkers_
A state variable if domain markers have been read from a Gmsh file.
Definition: gridmanager_base.hh:314
void loadBalance()
Call loadBalance() function of the grid.
Definition: gridmanager_base.hh:97
CellType
The cell types for structured grids.
Definition: gridmanager_base.hh:264
@ Simplex
Definition: gridmanager_base.hh:264
@ Cube
Definition: gridmanager_base.hh:264
std::shared_ptr< Grid > & gridPtr()
Returns a reference to the grid pointer (std::shared_ptr<Grid>)
Definition: gridmanager_base.hh:145
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:63
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:309
std::string getFileExtension(const std::string &fileName) const
Returns the filename extension of a given filename.
Definition: gridmanager_base.hh:167
The grid manager (this is the class used by the user) for all supported grid managers that constructs...
Definition: gridmanager_base.hh:323
A vtk file reader using tinyxml2 as xml backend.
Definition: vtkreader.hh:39
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:116
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:47
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:165
The available discretization methods in Dumux.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Class for grid data attached to dgf or gmsh grid files.
Template which always yields a false value.
Definition: common/typetraits/typetraits.hh:24
A vtk file reader using tinyxml2 as xml backend.