version 3.11-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-FileCopyrightText: 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#include "gridinput_.hh"
37
38namespace Dumux {
39
46template <class Grid>
47class GridManager;
48
54template <class GridType>
56{
57public:
58 using Grid = GridType;
60
64 void init(const std::string& modelParamGroup = "")
65 {
66 static_assert(AlwaysFalse<GridType>::value,
67 "The header with the GridManager specialization for your grid type is not included "
68 "or no specialization has been implemented!"
69 " In case of the latter, consider providing your own GridManager."
70 );
71 }
72
77 {
79 return *dgfGridPtr();
80 else
81 return *gridPtr();
82 }
83
87 const Grid& grid() const
88 {
90 return *dgfGridPtr();
91 else
92 return *gridPtr();
93 }
94
99 {
100 if (Dune::MPIHelper::instance().size() > 1)
101 {
102 // if we may have dgf parameters use load balancing of the dgf pointer
104 {
105 dgfGridPtr().loadBalance();
106 // update the grid data
107 gridData_ = std::make_shared<GridData>(dgfGridPtr());
108 }
109
110 // if we have gmsh parameters we have to manually load balance the data
112 {
113 // element and face markers are communicated during load balance
114 auto dh = gridData_->createGmshDataHandle();
115 gridPtr()->loadBalance(dh.interface());
116 gridPtr()->communicate(dh.interface(), Dune::InteriorBorder_All_Interface, Dune::ForwardCommunication);
117 }
118
119 // if we have VTK parameters we have to manually load balance the data
120 else if (enableVtkData_)
121 {
122 if constexpr (Detail::GridData::CartesianGrid<Grid>)
123 {
124 gridData_->communicateStructuredVtkData();
125 }
126 else
127 {
128 // cell and point data is communicated during load balance
129 auto dh = gridData_->createVtkDataHandle();
130 gridPtr()->loadBalance(dh.interface());
131 gridPtr()->communicate(dh.interface(), Dune::InteriorBorder_All_Interface, Dune::ForwardCommunication);
132 }
133 }
134
135 else
136 gridPtr()->loadBalance();
137 }
138 }
139
144 std::shared_ptr<GridData> getGridData() const
145 {
146 if (!gridData_)
147 DUNE_THROW(Dune::IOError, "No grid data available");
148
149 return gridData_;
150 }
151
155 bool hasGridData() const
156 { return static_cast<bool>(gridData_); }
157
158protected:
159
163 std::shared_ptr<Grid>& gridPtr()
164 {
166 return gridPtr_;
167 else
168 DUNE_THROW(Dune::InvalidStateException, "You are using DGF. To get the grid pointer use method dgfGridPtr()!");
169 }
170
174 Dune::GridPtr<Grid>& dgfGridPtr()
175 {
177 return dgfGridPtr_;
178 else
179 DUNE_THROW(Dune::InvalidStateException, "The DGF grid pointer is only available if the grid was constructed with a DGF file!");
180 }
181
185 std::string getFileExtension(const std::string& fileName) const
186 {
187 std::size_t i = fileName.rfind('.', fileName.length());
188 if (i != std::string::npos)
189 {
190 return(fileName.substr(i+1, fileName.length() - i));
191 }
192 else
193 {
194 DUNE_THROW(Dune::IOError, "Please provide and extension for your grid file ('"<< fileName << "')!");
195 }
196 return "";
197 }
198
205 void makeGridFromFile(const std::string& fileName,
206 const std::string& modelParamGroup)
207 {
208 // We found a file in the input file...does it have a supported extension?
209 const std::string extension = getFileExtension(fileName);
210 if (extension != "dgf" && extension != "msh" && extension != "vtu" && extension != "vtp")
211 DUNE_THROW(Dune::IOError, "Grid type " << Dune::className<Grid>() << " doesn't support grid files with extension: *."<< extension);
212
213 // Dune Grid Format (DGF) files
214 if (extension == "dgf")
215 {
217 dgfGridPtr() = Dune::GridPtr<Grid>(fileName.c_str(), Dune::MPIHelper::getCommunicator());
218 gridData_ = std::make_shared<GridData>(dgfGridPtr_);
219 }
220
221 // Gmsh mesh format
222 else if (extension == "msh")
223 {
224 // get some optional parameters
225 const bool verbose = getParamFromGroup<bool>(modelParamGroup, "Grid.Verbosity", false);
226 const bool boundarySegments = getParamFromGroup<bool>(modelParamGroup, "Grid.BoundarySegments", false);
227 const bool domainMarkers = getParamFromGroup<bool>(modelParamGroup, "Grid.DomainMarkers", false);
228
229 if (domainMarkers)
230 {
232 std::vector<int> boundaryMarkers, elementMarkers;
233 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
234 // the gmsh reader reads the grid on rank 0
235 Dune::GmshReader<Grid>::read(*gridFactory, fileName, boundaryMarkers, elementMarkers, verbose, boundarySegments);
236 gridPtr() = std::shared_ptr<Grid>(gridFactory->createGrid());
237 gridData_ = std::make_shared<GridData>(gridPtr_, std::move(gridFactory), std::move(elementMarkers), std::move(boundaryMarkers));
238 }
239 else
240 {
241 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
242 Dune::GmshReader<Grid>::read(*gridFactory, fileName, verbose, boundarySegments);
243 gridPtr() = std::shared_ptr<Grid>(gridFactory->createGrid());
244 }
245 }
246
247 // VTK file formats for unstructured grids
248 else if (extension == "vtu" || extension == "vtp")
249 {
250 if (Dune::MPIHelper::getCommunication().size() > 1)
251 DUNE_THROW(Dune::NotImplemented, "Reading grids in parallel from VTK file formats is currently not supported!");
252
253 VTKReader vtkReader(fileName);
254 VTKReader::Data cellData, pointData;
255 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
256 const bool verbose = getParamFromGroup<bool>(modelParamGroup, "Grid.Verbosity", false);
257 gridPtr() = vtkReader.readGrid(*gridFactory, cellData, pointData, verbose);
258 gridData_ = std::make_shared<GridData>(gridPtr_, std::move(gridFactory), std::move(cellData), std::move(pointData));
259 }
260 }
261
265 void makeGridFromDgfFile(const std::string& fileName)
266 {
267 // We found a file in the input file...does it have a supported extension?
268 const std::string extension = getFileExtension(fileName);
269 if(extension != "dgf")
270 DUNE_THROW(Dune::IOError, "Grid type " << Dune::className<Grid>() << " only supports DGF (*.dgf) but the specified filename has extension: *."<< extension);
271
273 dgfGridPtr() = Dune::GridPtr<Grid>(fileName.c_str(), Dune::MPIHelper::getCommunicator());
274 gridData_ = std::make_shared<GridData>(dgfGridPtr_);
275 }
276
281
285 template <int dim, int dimworld>
287 const std::string& modelParamGroup)
288 {
289 using GlobalPosition = Dune::FieldVector<typename Grid::ctype, dimworld>;
290 const auto upperRight = getParamFromGroup<GlobalPosition>(modelParamGroup, "Grid.UpperRight");
291 const auto lowerLeft = getParamFromGroup<GlobalPosition>(modelParamGroup, "Grid.LowerLeft", GlobalPosition(0.0));
292
293 using CellArray = std::array<unsigned int, dim>;
294 CellArray cells; cells.fill(1);
295 cells = getParamFromGroup<CellArray>(modelParamGroup, "Grid.Cells", cells);
296
297 // make the grid
298 if (cellType == CellType::Cube)
299 {
300 gridPtr() = Dune::StructuredGridFactory<Grid>::createCubeGrid(lowerLeft, upperRight, cells);
301 }
302 else if (cellType == CellType::Simplex)
303 {
304 gridPtr() = Dune::StructuredGridFactory<Grid>::createSimplexGrid(lowerLeft, upperRight, cells);
305 }
306 else
307 {
308 DUNE_THROW(Dune::GridError, "Unknown cell type for making structured grid! Choose Cube or Simplex.");
309 }
310 }
311
315 void maybeRefineGrid(const std::string& modelParamGroup)
316 {
317 if (hasParamInGroup(modelParamGroup, "Grid.Refinement"))
318 grid().globalRefine(getParamFromGroup<int>(modelParamGroup, "Grid.Refinement"));
319 }
320
326
331
335 bool enableVtkData_ = false;
336
337 std::shared_ptr<Grid> gridPtr_;
338 Dune::GridPtr<Grid> dgfGridPtr_;
339
340 std::shared_ptr<GridData> gridData_;
341};
342
343template <class Grid>
344class GridManager : public GridManagerBase<Grid> {};
345
346} // end namespace Dumux
347
348#endif
Class for grid data attached to dgf or gmsh grid files.
Definition: griddata.hh:57
The grid manager base interface (public) and methods common to most grid manager specializations (pro...
Definition: gridmanager_base.hh:56
std::shared_ptr< GridData > getGridData() const
Get an owning pointer to grid data associated with the grid.
Definition: gridmanager_base.hh:144
GridType Grid
Definition: gridmanager_base.hh:58
Dune::GridPtr< Grid > dgfGridPtr_
Definition: gridmanager_base.hh:338
const Grid & grid() const
Returns a const reference to the grid.
Definition: gridmanager_base.hh:87
void maybeRefineGrid(const std::string &modelParamGroup)
Refines a grid after construction if GridParameterGroup.Refinement is set in the input file.
Definition: gridmanager_base.hh:315
bool enableVtkData_
A state variable if cell or point data have been read from a VTK file.
Definition: gridmanager_base.hh:335
void makeGridFromFile(const std::string &fileName, const std::string &modelParamGroup)
Makes a grid from a file. We currently support.
Definition: gridmanager_base.hh:205
void makeStructuredGrid(CellType cellType, const std::string &modelParamGroup)
Makes a structured cube grid using the structured grid factory.
Definition: gridmanager_base.hh:286
Dune::GridPtr< Grid > & dgfGridPtr()
Returns a reference to the DGF grid pointer (Dune::GridPtr<Grid>).
Definition: gridmanager_base.hh:174
bool hasGridData() const
Check whether there is data associated with the grid.
Definition: gridmanager_base.hh:155
std::shared_ptr< Grid > gridPtr_
Definition: gridmanager_base.hh:337
Grid & grid()
Returns a reference to the grid.
Definition: gridmanager_base.hh:76
std::shared_ptr< GridData > gridData_
Definition: gridmanager_base.hh:340
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:265
bool enableGmshDomainMarkers_
A state variable if domain markers have been read from a Gmsh file.
Definition: gridmanager_base.hh:330
void loadBalance()
Call loadBalance() function of the grid.
Definition: gridmanager_base.hh:98
CellType
The cell types for structured grids.
Definition: gridmanager_base.hh:280
@ Simplex
Definition: gridmanager_base.hh:280
@ Cube
Definition: gridmanager_base.hh:280
std::shared_ptr< Grid > & gridPtr()
Returns a reference to the grid pointer (std::shared_ptr<Grid>)
Definition: gridmanager_base.hh:163
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:64
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:325
std::string getFileExtension(const std::string &fileName) const
Returns the filename extension of a given filename.
Definition: gridmanager_base.hh:185
The grid manager (this is the class used by the user) for all supported grid managers that constructs...
Definition: gridmanager_base.hh:344
A vtk file reader using tinyxml2 as xml backend.
Definition: vtkreader.hh:342
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:423
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:350
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.