version 3.8
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-FileCopyrightInfo: 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")
127 DUNE_THROW(Dune::IOError, "Grid type " << Dune::className<Grid>() << " only supports DGF (*.dgf) and Gmsh (*.msh) grid files but the specified filename has extension: *."<< extension);
128
129 // make the grid
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 if (extension == "msh")
137 {
138 // get some optional parameters
139 const bool verbose = getParamFromGroup<bool>(modelParamGroup, "Grid.Verbosity", false);
140 const bool boundarySegments = getParamFromGroup<bool>(modelParamGroup, "Grid.BoundarySegments", false);
141 const bool domainMarkers = getParamFromGroup<bool>(modelParamGroup, "Grid.DomainMarkers", false);
142
143 if (domainMarkers)
144 ParentType::enableGmshDomainMarkers_ = true;
145
146 // only fill the factory for rank 0
147 if (domainMarkers)
148 {
149 std::vector<int> boundaryMarkersInsertionIndex, boundaryMarkers, faceMarkers, elementMarkers;
150 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
151
152 Dune::GmshReader<Grid>::read(*gridFactory, fileName, boundaryMarkersInsertionIndex, elementMarkers, verbose, boundarySegments);
153 ParentType::gridPtr() = std::shared_ptr<Grid>(gridFactory->createGrid());
154
155 // reorder boundary markers according to boundarySegmentIndex
156 boundaryMarkers.resize(ParentType::gridPtr()->numBoundarySegments(), 0);
157 faceMarkers.resize(ParentType::gridPtr()->leafGridView().size(1), 0);
158 const auto& indexSet = ParentType::gridPtr()->leafGridView().indexSet();
159 for (const auto& element : elements(ParentType::gridPtr()->leafGridView()))
160 {
161 for (const auto& intersection : intersections(ParentType::gridPtr()->leafGridView(), element))
162 {
163 if (intersection.boundary() && gridFactory->wasInserted(intersection))
164 {
165 auto marker = boundaryMarkersInsertionIndex[gridFactory->insertionIndex(intersection)];
166 boundaryMarkers[intersection.boundarySegmentIndex()] = marker;
167 faceMarkers[indexSet.index(element.template subEntity<1>(intersection.indexInInside()))] = marker;
168 }
169 }
170 }
171
172 ParentType::gridData_ = std::make_shared<typename ParentType::GridData>(ParentType::gridPtr(), std::move(gridFactory),
173 std::move(elementMarkers), std::move(boundaryMarkers), std::move(faceMarkers));
174 }
175 else
176 {
177 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
178
179 if (Dune::MPIHelper::getCommunication().rank() == 0)
180 Dune::GmshReader<Grid>::read(*gridFactory, fileName, verbose, boundarySegments);
181
182 ParentType::gridPtr() = std::shared_ptr<Grid>(gridFactory->createGrid());
183 }
184 }
185 }
186
190 template <int dimension, int dimensionworld, std::enable_if_t<dimension != dimensionworld, int> = 0>
191 void makeStructuredGrid(typename ParentType::CellType cellType,
192 const std::string& modelParamGroup)
193 {
194 DUNE_THROW(Dune::IOError, "ALUGrid currently only supports the creation of structured grids with dimension == dimensionworld. Consider reading in a grid file instead.");
195 }
196
200 template <int dimension, int dimensionworld, std::enable_if_t<dimension == dimensionworld, int> = 0>
201 void makeStructuredGrid(typename ParentType::CellType cellType,
202 const std::string& modelParamGroup)
203 {
204 // make a structured grid
205 if (elType == Dune::cube)
206 ParentType::template makeStructuredGrid<dimension, dimensionworld>(ParentType::CellType::Cube, modelParamGroup);
207 else if (elType == Dune::simplex)
208 ParentType::template makeStructuredGrid<dimension, dimensionworld>(ParentType::CellType::Simplex, modelParamGroup);
209 else
210 DUNE_THROW(Dune::IOError, "ALUGrid only supports Dune::cube or Dune::simplex as cell type!");
211 }
212};
213
219template<int dim, int dimworld, Dune::ALUGridElementType elType, Dune::ALUGridRefinementType refinementType>
220class BoundaryFlag<Dune::ALUGrid<dim, dimworld, elType, refinementType>>
221{
222public:
223 BoundaryFlag() : flag_(-1) {}
224
225 template<class Intersection>
226 BoundaryFlag(const Intersection& i) : flag_(-1)
227 {
228 if (i.boundary())
229 flag_ = i.impl().boundaryId();
230 }
231
232 using value_type = int;
233
234 value_type get() const { return flag_; }
235
236private:
237 int flag_;
238};
239
240namespace Grid::Capabilities {
241
242// To the best of our knowledge ALUGrid is view thread-safe
243// This specialization can be removed after we depend on Dune release 2.9 in which this is guaranteed by ALUGrid itself
244template<int dim, int dimworld, Dune::ALUGridElementType elType, Dune::ALUGridRefinementType refinementType>
245struct MultithreadingSupported<Dune::ALUGrid<dim, dimworld, elType, refinementType>>
246{
247 template<class GV>
248 static bool eval(const GV&) // default is independent of the grid view
249 { return true; }
250};
251
252} // end namespace Grid::Capabilities
253
254#endif // HAVE_DUNE_ALUGRID
255
256} // end namespace Dumux
257
258#endif
Boundary flag to store e.g. in sub control volume faces.
std::size_t value_type
Definition: boundaryflag.hh:39
value_type get() const
Definition: boundaryflag.hh:41
Grid Grid
Definition: gridmanager_base.hh:57
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
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
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:69