version 3.11-dev
gridmanager_alu.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//
12#ifndef DUMUX_IO_GRID_MANAGER_ALU_HH
13#define DUMUX_IO_GRID_MANAGER_ALU_HH
14
15// ALUGrid specific includes
16#if HAVE_DUNE_ALUGRID
17#include <dune/alugrid/grid.hh>
18#include <dune/alugrid/dgf.hh>
19#endif
20
21#ifndef DUMUX_IO_GRID_MANAGER_BASE_HH
23#endif
24
27
28namespace Dumux {
29
30#if HAVE_DUNE_ALUGRID
31
49template<int dim, int dimworld, Dune::ALUGridElementType elType, Dune::ALUGridRefinementType refinementType>
50class GridManager<Dune::ALUGrid<dim, dimworld, elType, refinementType>>
51: public GridManagerBase<Dune::ALUGrid<dim, dimworld, elType, refinementType>>
52{
53public:
54 using Grid = Dune::ALUGrid<dim, dimworld, elType, refinementType>;
55 using ParentType = GridManagerBase<Grid>;
56
60 void init(const std::string& modelParamGroup = "", bool adaptiveRestart = false)
61 {
62 // restarting an adaptive grid using Dune's BackupRestoreFacility
63 // TODO: the part after first || is backward compatibility with old sequential models remove once sequential adaptive restart is replaced
64 if (adaptiveRestart || hasParam("Restart") || hasParam("TimeManager.Restart"))
65 {
66 auto restartTime = getParamFromGroup<double>(modelParamGroup, "TimeLoop.Restart", 0.0);
67 // TODO: backward compatibility with old sequential models remove once sequential adaptive restart is replaced
68 if (hasParam("Restart") || hasParam("TimeManager.Restart"))
69 {
70 restartTime = getParamFromGroup<double>("TimeManager", "Restart");
71 std::cerr << "Warning: You are using a deprecated restart mechanism. The usage will change in the future.\n";
72 }
73
74 const int rank = Dune::MPIHelper::getCommunication().rank();
75 const std::string name = getParamFromGroup<std::string>(modelParamGroup, "Problem.Name");
76 std::ostringstream oss;
77 oss << name << "_time=" << restartTime << "_rank=" << rank << ".grs";
78 std::cout << "Restoring an ALUGrid from " << oss.str() << std::endl;
79 ParentType::gridPtr() = std::shared_ptr<Grid>(Dune::BackupRestoreFacility<Grid>::restore(oss.str()));
80 ParentType::loadBalance();
81 return;
82 }
83
84 // try to create it from a DGF or msh file in GridParameterGroup.File
85 else if (hasParamInGroup(modelParamGroup, "Grid.File"))
86 {
87 makeGridFromFile(getParamFromGroup<std::string>(modelParamGroup, "Grid.File"), modelParamGroup);
88 ParentType::maybeRefineGrid(modelParamGroup);
89 ParentType::loadBalance();
90 return;
91 }
92
93 // Then look for the necessary keys to construct from the input file
94 else if (hasParamInGroup(modelParamGroup, "Grid.UpperRight"))
95 {
96 if (elType == Dune::cube)
97 makeStructuredGrid<dim, dimworld>(ParentType::CellType::Cube, modelParamGroup);
98 else if (elType == Dune::simplex)
99 makeStructuredGrid<dim, dimworld>(ParentType::CellType::Simplex, modelParamGroup);
100 else
101 DUNE_THROW(Dune::IOError, "ALUGrid only supports Dune::cube or Dune::simplex as cell type!");
102
103 ParentType::maybeRefineGrid(modelParamGroup);
104 ParentType::loadBalance();
105 }
106
107 // Didn't find a way to construct the grid
108 else
109 {
110 const auto prefix = modelParamGroup.empty() ? modelParamGroup : modelParamGroup + ".";
111 DUNE_THROW(ParameterException, "Please supply one of the parameters "
112 << prefix + "Grid.UpperRight"
113 << ", or a grid file in " << prefix + "Grid.File");
114
115 }
116 }
117
121 void makeGridFromFile(const std::string& fileName,
122 const std::string& modelParamGroup)
123 {
124 // We found a file in the input file...does it have a supported extension?
125 const std::string extension = ParentType::getFileExtension(fileName);
126 if (extension != "dgf" && extension != "msh" && extension != "vtu" && extension != "vti")
127 DUNE_THROW(Dune::IOError, "Grid type " << Dune::className<Grid>() << " doesn't support grid files with extension: *."<< extension);
128
129 // Dune Grid Format (DGF) files
130 if (extension == "dgf")
131 {
132 ParentType::enableDgfGridPointer_ = true;
133 ParentType::dgfGridPtr() = Dune::GridPtr<Grid>(fileName.c_str(), Dune::MPIHelper::getCommunicator());
134 ParentType::gridData_ = std::make_shared<typename ParentType::GridData>(ParentType::dgfGridPtr());
135 }
136
137 // Gmsh mesh format
138 else if (extension == "msh")
139 {
140 // get some optional parameters
141 const bool verbose = getParamFromGroup<bool>(modelParamGroup, "Grid.Verbosity", false);
142 const bool boundarySegments = getParamFromGroup<bool>(modelParamGroup, "Grid.BoundarySegments", false);
143 const bool domainMarkers = getParamFromGroup<bool>(modelParamGroup, "Grid.DomainMarkers", false);
144
145 if (domainMarkers)
146 {
147 ParentType::enableGmshDomainMarkers_ = true;
148 std::vector<int> boundaryMarkersInsertionIndex, boundaryMarkers, faceMarkers, elementMarkers;
149 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
150 Dune::GmshReader<Grid>::read(*gridFactory, fileName, boundaryMarkersInsertionIndex, elementMarkers, verbose, boundarySegments);
151 ParentType::gridPtr() = std::shared_ptr<Grid>(gridFactory->createGrid());
152
153 // reorder boundary markers according to boundarySegmentIndex
154 boundaryMarkers.resize(ParentType::gridPtr()->numBoundarySegments(), 0);
155 faceMarkers.resize(ParentType::gridPtr()->leafGridView().size(1), 0);
156 const auto& indexSet = ParentType::gridPtr()->leafGridView().indexSet();
157 for (const auto& element : elements(ParentType::gridPtr()->leafGridView()))
158 {
159 for (const auto& intersection : intersections(ParentType::gridPtr()->leafGridView(), element))
160 {
161 if (intersection.boundary() && gridFactory->wasInserted(intersection))
162 {
163 auto marker = boundaryMarkersInsertionIndex[gridFactory->insertionIndex(intersection)];
164 boundaryMarkers[intersection.boundarySegmentIndex()] = marker;
165 faceMarkers[indexSet.index(element.template subEntity<1>(intersection.indexInInside()))] = marker;
166 }
167 }
168 }
169
170 ParentType::gridData_ = std::make_shared<typename ParentType::GridData>(ParentType::gridPtr(), std::move(gridFactory),
171 std::move(elementMarkers), std::move(boundaryMarkers), std::move(faceMarkers));
172 }
173 else
174 {
175 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
176 Dune::GmshReader<Grid>::read(*gridFactory, fileName, verbose, boundarySegments);
177 ParentType::gridPtr() = std::shared_ptr<Grid>(gridFactory->createGrid());
178 }
179 }
180
181 // VTK file formats for unstructured grids
182 // (can be constructed from both unstructured or structured grid data file formats)
183 else if (extension == "vtu" || extension == "vti")
184 {
185 VTKReader vtkReader(fileName);
186 VTKReader::Data cellData, pointData;
187 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
188 const bool verbose = getParamFromGroup<bool>(modelParamGroup, "Grid.Verbosity", false);
189 ParentType::gridPtr() = vtkReader.readGrid(*gridFactory, cellData, pointData, verbose);
190 ParentType::gridData_ = std::make_shared<typename ParentType::GridData>(ParentType::gridPtr(), std::move(gridFactory), std::move(cellData), std::move(pointData));
191 ParentType::enableVtkData_ = true;
192 }
193 }
194
198 template <int dimension, int dimensionworld, std::enable_if_t<dimension != dimensionworld, int> = 0>
199 void makeStructuredGrid(typename ParentType::CellType cellType,
200 const std::string& modelParamGroup)
201 {
202 DUNE_THROW(Dune::IOError, "ALUGrid currently only supports the creation of structured grids with dimension == dimensionworld. Consider reading in a grid file instead.");
203 }
204
208 template <int dimension, int dimensionworld, std::enable_if_t<dimension == dimensionworld, int> = 0>
209 void makeStructuredGrid(typename ParentType::CellType cellType,
210 const std::string& modelParamGroup)
211 {
212 // make a structured grid
213 if (elType == Dune::cube)
214 ParentType::template makeStructuredGrid<dimension, dimensionworld>(ParentType::CellType::Cube, modelParamGroup);
215 else if (elType == Dune::simplex)
216 ParentType::template makeStructuredGrid<dimension, dimensionworld>(ParentType::CellType::Simplex, modelParamGroup);
217 else
218 DUNE_THROW(Dune::IOError, "ALUGrid only supports Dune::cube or Dune::simplex as cell type!");
219 }
220};
221
227template<int dim, int dimworld, Dune::ALUGridElementType elType, Dune::ALUGridRefinementType refinementType>
228class BoundaryFlag<Dune::ALUGrid<dim, dimworld, elType, refinementType>>
229{
230public:
231 BoundaryFlag() : flag_(invalidFlag_) {}
232
233 template<class Intersection>
234 BoundaryFlag(const Intersection& i) : flag_(invalidFlag_)
235 {
236 if (i.boundary())
237 flag_ = i.impl().boundaryId();
238 }
239
240 using value_type = int;
241
242 value_type get() const { return flag_; }
243
244 operator bool() const { return flag_ != invalidFlag_; }
245
246private:
247 static constexpr value_type invalidFlag_ = -1;
248 value_type flag_;
249};
250
251namespace Grid::Capabilities {
252
253// To the best of our knowledge ALUGrid is view thread-safe
254// This specialization can be removed after we depend on Dune release 2.9 in which this is guaranteed by ALUGrid itself
255template<int dim, int dimworld, Dune::ALUGridElementType elType, Dune::ALUGridRefinementType refinementType>
256struct MultithreadingSupported<Dune::ALUGrid<dim, dimworld, elType, refinementType>>
257{
258 template<class GV>
259 static bool eval(const GV&) // default is independent of the grid view
260 { return true; }
261};
262
263} // end namespace Grid::Capabilities
264
265#endif // HAVE_DUNE_ALUGRID
266
267} // end namespace Dumux
268
269#endif
Boundary flag to store e.g. in sub control volume faces.
std::size_t value_type
Definition: boundaryflag.hh:28
value_type get() const
Definition: boundaryflag.hh:41
Grid Grid
Definition: gridmanager_base.hh:58
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
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
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
dune-grid capabilities compatibility layer
Provides a grid manager for all supported grid managers with input file interfaces....
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
bool hasParam(const std::string &param)
Check whether a key exists in the parameter tree.
Definition: parameters.hh:157
Definition: adapt.hh:17
Definition: common/pdesolver.hh:24
static bool eval(const GV &)
Definition: gridcapabilities.hh:29