version 3.10-dev
gridmanager_base.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_IO_GRID_MANAGER_BASE_HH
14#define DUMUX_IO_GRID_MANAGER_BASE_HH
15
16#include <array>
17#include <bitset>
18#include <memory>
19#include <sstream>
20
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>
29
34
35#include "griddata.hh"
36
37namespace Dumux {
38
45template <class Grid>
46class GridManager;
47
53template <class GridType>
55{
56public:
57 using Grid = GridType;
59
63 void init(const std::string& modelParamGroup = "")
64 {
65 static_assert(AlwaysFalse<GridType>::value,
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."
69 );
70 }
71
76 {
78 return *dgfGridPtr();
79 else
80 return *gridPtr();
81 }
82
86 const Grid& grid() const
87 {
89 return *dgfGridPtr();
90 else
91 return *gridPtr();
92 }
93
98 {
99 if (Dune::MPIHelper::instance().size() > 1)
100 {
101 // if we may have dgf parameters use load balancing of the dgf pointer
103 {
104 dgfGridPtr().loadBalance();
105 // update the grid data
106 gridData_ = std::make_shared<GridData>(dgfGridPtr());
107 }
108
109 // if we have gmsh parameters we have to manually load balance the data
111 {
112 // element and face markers are communicated during load balance
113 auto dh = gridData_->createGmshDataHandle();
114 gridPtr()->loadBalance(dh.interface());
115 gridPtr()->communicate(dh.interface(), Dune::InteriorBorder_All_Interface, Dune::ForwardCommunication);
116 }
117
118 // if we have VTK parameters we have to manually load balance the data
119 else if (enableVtkData_)
120 {
121 // cell and point data is communicated during load balance
122 auto dh = gridData_->createVtkDataHandle();
123 gridPtr()->loadBalance(dh.interface());
124 gridPtr()->communicate(dh.interface(), Dune::InteriorBorder_All_Interface, Dune::ForwardCommunication);
125 }
126
127 else
128 gridPtr()->loadBalance();
129 }
130 }
131
136 std::shared_ptr<GridData> getGridData() const
137 {
138 if (!gridData_)
139 DUNE_THROW(Dune::IOError, "No grid data available");
140
141 return gridData_;
142 }
143
147 bool hasGridData() const
148 { return static_cast<bool>(gridData_); }
149
150protected:
151
155 std::shared_ptr<Grid>& gridPtr()
156 {
158 return gridPtr_;
159 else
160 DUNE_THROW(Dune::InvalidStateException, "You are using DGF. To get the grid pointer use method dgfGridPtr()!");
161 }
162
166 Dune::GridPtr<Grid>& dgfGridPtr()
167 {
169 return dgfGridPtr_;
170 else
171 DUNE_THROW(Dune::InvalidStateException, "The DGF grid pointer is only available if the grid was constructed with a DGF file!");
172 }
173
177 std::string getFileExtension(const std::string& fileName) const
178 {
179 std::size_t i = fileName.rfind('.', fileName.length());
180 if (i != std::string::npos)
181 {
182 return(fileName.substr(i+1, fileName.length() - i));
183 }
184 else
185 {
186 DUNE_THROW(Dune::IOError, "Please provide and extension for your grid file ('"<< fileName << "')!");
187 }
188 return "";
189 }
190
197 void makeGridFromFile(const std::string& fileName,
198 const std::string& modelParamGroup)
199 {
200 // We found a file in the input file...does it have a supported extension?
201 const std::string extension = getFileExtension(fileName);
202 if (extension != "dgf" && extension != "msh" && extension != "vtu" && extension != "vtp")
203 DUNE_THROW(Dune::IOError, "Grid type " << Dune::className<Grid>() << " doesn't support grid files with extension: *."<< extension);
204
205 // Dune Grid Format (DGF) files
206 if (extension == "dgf")
207 {
209 dgfGridPtr() = Dune::GridPtr<Grid>(fileName.c_str(), Dune::MPIHelper::getCommunicator());
210 gridData_ = std::make_shared<GridData>(dgfGridPtr_);
211 }
212
213 // Gmsh mesh format
214 else if (extension == "msh")
215 {
216 // get some optional parameters
217 const bool verbose = getParamFromGroup<bool>(modelParamGroup, "Grid.Verbosity", false);
218 const bool boundarySegments = getParamFromGroup<bool>(modelParamGroup, "Grid.BoundarySegments", false);
219 const bool domainMarkers = getParamFromGroup<bool>(modelParamGroup, "Grid.DomainMarkers", false);
220
221 if (domainMarkers)
222 {
224 std::vector<int> boundaryMarkers, elementMarkers;
225 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
226 // the gmsh reader reads the grid on rank 0
227 Dune::GmshReader<Grid>::read(*gridFactory, fileName, boundaryMarkers, elementMarkers, verbose, boundarySegments);
228 gridPtr() = std::shared_ptr<Grid>(gridFactory->createGrid());
229 gridData_ = std::make_shared<GridData>(gridPtr_, std::move(gridFactory), std::move(elementMarkers), std::move(boundaryMarkers));
230 }
231 else
232 {
233 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
234 Dune::GmshReader<Grid>::read(*gridFactory, fileName, verbose, boundarySegments);
235 gridPtr() = std::shared_ptr<Grid>(gridFactory->createGrid());
236 }
237 }
238
239 // VTK file formats for unstructured grids
240 else if (extension == "vtu" || extension == "vtp")
241 {
242 if (Dune::MPIHelper::getCommunication().size() > 1)
243 DUNE_THROW(Dune::NotImplemented, "Reading grids in parallel from VTK file formats is currently not supported!");
244
245 VTKReader vtkReader(fileName);
246 VTKReader::Data cellData, pointData;
247 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
248 const bool verbose = getParamFromGroup<bool>(modelParamGroup, "Grid.Verbosity", false);
249 gridPtr() = vtkReader.readGrid(*gridFactory, cellData, pointData, verbose);
250 gridData_ = std::make_shared<GridData>(gridPtr_, std::move(gridFactory), std::move(cellData), std::move(pointData));
251 }
252 }
253
257 void makeGridFromDgfFile(const std::string& fileName)
258 {
259 // We found a file in the input file...does it have a supported extension?
260 const std::string extension = getFileExtension(fileName);
261 if(extension != "dgf")
262 DUNE_THROW(Dune::IOError, "Grid type " << Dune::className<Grid>() << " only supports DGF (*.dgf) but the specified filename has extension: *."<< extension);
263
265 dgfGridPtr() = Dune::GridPtr<Grid>(fileName.c_str(), Dune::MPIHelper::getCommunicator());
266 gridData_ = std::make_shared<GridData>(dgfGridPtr_);
267 }
268
273
277 template <int dim, int dimworld>
279 const std::string& modelParamGroup)
280 {
281 using GlobalPosition = Dune::FieldVector<typename Grid::ctype, dimworld>;
282 const auto upperRight = getParamFromGroup<GlobalPosition>(modelParamGroup, "Grid.UpperRight");
283 const auto lowerLeft = getParamFromGroup<GlobalPosition>(modelParamGroup, "Grid.LowerLeft", GlobalPosition(0.0));
284
285 using CellArray = std::array<unsigned int, dim>;
286 CellArray cells; cells.fill(1);
287 cells = getParamFromGroup<CellArray>(modelParamGroup, "Grid.Cells", cells);
288
289 // make the grid
290 if (cellType == CellType::Cube)
291 {
292 gridPtr() = Dune::StructuredGridFactory<Grid>::createCubeGrid(lowerLeft, upperRight, cells);
293 }
294 else if (cellType == CellType::Simplex)
295 {
296 gridPtr() = Dune::StructuredGridFactory<Grid>::createSimplexGrid(lowerLeft, upperRight, cells);
297 }
298 else
299 {
300 DUNE_THROW(Dune::GridError, "Unknown cell type for making structured grid! Choose Cube or Simplex.");
301 }
302 }
303
307 void maybeRefineGrid(const std::string& modelParamGroup)
308 {
309 if (hasParamInGroup(modelParamGroup, "Grid.Refinement"))
310 grid().globalRefine(getParamFromGroup<int>(modelParamGroup, "Grid.Refinement"));
311 }
312
318
323
327 bool enableVtkData_ = false;
328
329 std::shared_ptr<Grid> gridPtr_;
330 Dune::GridPtr<Grid> dgfGridPtr_;
331
332 std::shared_ptr<GridData> gridData_;
333};
334
335template <class Grid>
336class GridManager : public GridManagerBase<Grid> {};
337
338} // end namespace Dumux
339
340#endif
Class for grid data attached to dgf or gmsh grid files.
Definition: griddata.hh:55
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:136
GridType Grid
Definition: gridmanager_base.hh:57
Dune::GridPtr< Grid > dgfGridPtr_
Definition: gridmanager_base.hh:330
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:307
bool enableVtkData_
A state variable if cell or point data have been read from a VTK file.
Definition: gridmanager_base.hh:327
void makeGridFromFile(const std::string &fileName, const std::string &modelParamGroup)
Makes a grid from a file. We currently support.
Definition: gridmanager_base.hh:197
void makeStructuredGrid(CellType cellType, const std::string &modelParamGroup)
Makes a structured cube grid using the structured grid factory.
Definition: gridmanager_base.hh:278
Dune::GridPtr< Grid > & dgfGridPtr()
Returns a reference to the DGF grid pointer (Dune::GridPtr<Grid>).
Definition: gridmanager_base.hh:166
bool hasGridData() const
Check whether there is data associated with the grid.
Definition: gridmanager_base.hh:147
std::shared_ptr< Grid > gridPtr_
Definition: gridmanager_base.hh:329
Grid & grid()
Returns a reference to the grid.
Definition: gridmanager_base.hh:75
std::shared_ptr< GridData > gridData_
Definition: gridmanager_base.hh:332
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:257
bool enableGmshDomainMarkers_
A state variable if domain markers have been read from a Gmsh file.
Definition: gridmanager_base.hh:322
void loadBalance()
Call loadBalance() function of the grid.
Definition: gridmanager_base.hh:97
CellType
The cell types for structured grids.
Definition: gridmanager_base.hh:272
@ Simplex
Definition: gridmanager_base.hh:272
@ Cube
Definition: gridmanager_base.hh:272
std::shared_ptr< Grid > & gridPtr()
Returns a reference to the grid pointer (std::shared_ptr<Grid>)
Definition: gridmanager_base.hh:155
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:317
std::string getFileExtension(const std::string &fileName) const
Returns the filename extension of a given filename.
Definition: gridmanager_base.hh:177
The grid manager (this is the class used by the user) for all supported grid managers that constructs...
Definition: gridmanager_base.hh:336
A vtk file reader using tinyxml2 as xml backend.
Definition: vtkreader.hh:40
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:121
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:48
Type traits.
bool hasParamInGroup(const std::string &paramGroup, const std::string &param)
Check whether a key exists in the parameter tree with a model group prefix.
Definition: parameters.hh:165
The available discretization methods in Dumux.
Definition: adapt.hh:17
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.