3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_IO_GRID_MANAGER_ALU_HH
25#define DUMUX_IO_GRID_MANAGER_ALU_HH
26
27#include <dune/common/version.hh>
28
29// ALUGrid specific includes
30#if HAVE_DUNE_ALUGRID
31#include <dune/alugrid/grid.hh>
32#include <dune/alugrid/dgf.hh>
33#endif
34
35#ifndef DUMUX_IO_GRID_MANAGER_BASE_HH
37#endif
38
40
41namespace Dumux {
42
43#if HAVE_DUNE_ALUGRID
44
62template<int dim, int dimworld, Dune::ALUGridElementType elType, Dune::ALUGridRefinementType refinementType>
63class GridManager<Dune::ALUGrid<dim, dimworld, elType, refinementType>>
64: public GridManagerBase<Dune::ALUGrid<dim, dimworld, elType, refinementType>>
65{
66public:
67 using Grid = Dune::ALUGrid<dim, dimworld, elType, refinementType>;
68 using ParentType = GridManagerBase<Grid>;
69
73 void init(const std::string& modelParamGroup = "", bool adaptiveRestart = false)
74 {
75 // restarting an adaptive grid using Dune's BackupRestoreFacility
76 // TODO: the part after first || is backward compatibilty with old sequential models remove once sequential adpative restart is replaced
77 if (adaptiveRestart || hasParam("Restart") || hasParam("TimeManager.Restart"))
78 {
79 auto restartTime = getParamFromGroup<double>(modelParamGroup, "TimeLoop.Restart", 0.0);
80 // TODO: backward compatibilty with old sequential models remove once sequential adpative restart is replaced
81 if (hasParam("Restart") || hasParam("TimeManager.Restart"))
82 {
83 restartTime = getParamFromGroup<double>("TimeManager", "Restart");
84 std::cerr << "Warning: You are using a deprecated restart mechanism. The usage will change in the future.\n";
85 }
86
87 const int rank = Dune::MPIHelper::getCollectiveCommunication().rank();
88 const std::string name = getParamFromGroup<std::string>(modelParamGroup, "Problem.Name");
89 std::ostringstream oss;
90 oss << name << "_time=" << restartTime << "_rank=" << rank << ".grs";
91 std::cout << "Restoring an ALUGrid from " << oss.str() << std::endl;
92 ParentType::gridPtr() = std::shared_ptr<Grid>(Dune::BackupRestoreFacility<Grid>::restore(oss.str()));
93 ParentType::loadBalance();
94 return;
95 }
96
97 // try to create it from a DGF or msh file in GridParameterGroup.File
98 else if (hasParamInGroup(modelParamGroup, "Grid.File"))
99 {
100 makeGridFromFile(getParamFromGroup<std::string>(modelParamGroup, "Grid.File"), modelParamGroup);
101 ParentType::maybeRefineGrid(modelParamGroup);
102 ParentType::loadBalance();
103 return;
104 }
105
106 // Then look for the necessary keys to construct from the input file
107 else if (hasParamInGroup(modelParamGroup, "Grid.UpperRight"))
108 {
109 if (elType == Dune::cube)
110 makeStructuredGrid<dim, dimworld>(ParentType::CellType::Cube, modelParamGroup);
111 else if (elType == Dune::simplex)
112 makeStructuredGrid<dim, dimworld>(ParentType::CellType::Simplex, modelParamGroup);
113 else
114 DUNE_THROW(Dune::IOError, "ALUGrid only supports Dune::cube or Dune::simplex as cell type!");
115
116 ParentType::maybeRefineGrid(modelParamGroup);
117 ParentType::loadBalance();
118 }
119
120 // Didn't find a way to construct the grid
121 else
122 {
123 const auto prefix = modelParamGroup.empty() ? modelParamGroup : modelParamGroup + ".";
124 DUNE_THROW(ParameterException, "Please supply one of the parameters "
125 << prefix + "Grid.UpperRight"
126 << ", or a grid file in " << prefix + "Grid.File");
127
128 }
129 }
130
134 void makeGridFromFile(const std::string& fileName,
135 const std::string& modelParamGroup)
136 {
137 // We found a file in the input file...does it have a supported extension?
138 const std::string extension = ParentType::getFileExtension(fileName);
139 if(extension != "dgf" && extension != "msh")
140 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);
141
142 // make the grid
143 if (extension == "dgf")
144 {
145 ParentType::enableDgfGridPointer_ = true;
146 ParentType::dgfGridPtr() = Dune::GridPtr<Grid>(fileName.c_str(), Dune::MPIHelper::getCommunicator());
147 ParentType::gridData_ = std::make_shared<typename ParentType::GridData>(ParentType::dgfGridPtr());
148 }
149 if (extension == "msh")
150 {
151 // get some optional parameters
152 const bool verbose = getParamFromGroup<bool>(modelParamGroup, "Grid.Verbosity", false);
153 const bool boundarySegments = getParamFromGroup<bool>(modelParamGroup, "Grid.BoundarySegments", false);
154 const bool domainMarkers = getParamFromGroup<bool>(modelParamGroup, "Grid.DomainMarkers", false);
155
156 if (domainMarkers)
157 ParentType::enableGmshDomainMarkers_ = true;
158
159 // only filll the factory for rank 0
160 if (domainMarkers)
161 {
162 std::vector<int> boundaryMarkersInsertionIndex, boundaryMarkers, faceMarkers, elementMarkers;
163 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
164
165// Older versions of Dune(-Alugrid) require that the Gmsh file is read only on process 0
166#if DUNE_VERSION_LT(DUNE_GRID, 2, 7)
167 if (Dune::MPIHelper::getCollectiveCommunication().rank() == 0)
168#endif
169 Dune::GmshReader<Grid>::read(*gridFactory, fileName, boundaryMarkersInsertionIndex, elementMarkers, verbose, boundarySegments);
170
171 ParentType::gridPtr() = std::shared_ptr<Grid>(gridFactory->createGrid());
172
173 // reorder boundary markers according to boundarySegmentIndex
174 boundaryMarkers.resize(ParentType::gridPtr()->numBoundarySegments(), 0);
175 faceMarkers.resize(ParentType::gridPtr()->leafGridView().size(1), 0);
176 const auto& indexSet = ParentType::gridPtr()->leafGridView().indexSet();
177 for (const auto& element : elements(ParentType::gridPtr()->leafGridView()))
178 {
179 for (const auto& intersection : intersections(ParentType::gridPtr()->leafGridView(), element))
180 {
181 if (intersection.boundary() && gridFactory->wasInserted(intersection))
182 {
183 auto marker = boundaryMarkersInsertionIndex[gridFactory->insertionIndex(intersection)];
184 boundaryMarkers[intersection.boundarySegmentIndex()] = marker;
185 faceMarkers[indexSet.index(element.template subEntity<1>(intersection.indexInInside()))] = marker;
186 }
187 }
188 }
189
190 ParentType::gridData_ = std::make_shared<typename ParentType::GridData>(ParentType::gridPtr(), std::move(gridFactory),
191 std::move(elementMarkers), std::move(boundaryMarkers), std::move(faceMarkers));
192 }
193 else
194 {
195 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
196
197// Older versions of Dune(-Alugrid) require that the Gmsh file is read only on process 0
198#if DUNE_VERSION_LT(DUNE_GRID, 2, 7)
199 if (Dune::MPIHelper::getCollectiveCommunication().rank() == 0)
200#endif
201 Dune::GmshReader<Grid>::read(*gridFactory, fileName, verbose, boundarySegments);
202
203 ParentType::gridPtr() = std::shared_ptr<Grid>(gridFactory->createGrid());
204 }
205 }
206 }
207
211 template <int dimension, int dimensionworld, std::enable_if_t<dimension != dimensionworld, int> = 0>
212 void makeStructuredGrid(typename ParentType::CellType cellType,
213 const std::string& modelParamGroup)
214 {
215 DUNE_THROW(Dune::IOError, "ALUGrid currently only supports the creation of structured grids with dimension == dimensionworld. Consider reading in a grid file instead.");
216 }
217
221 template <int dimension, int dimensionworld, std::enable_if_t<dimension == dimensionworld, int> = 0>
222 void makeStructuredGrid(typename ParentType::CellType cellType,
223 const std::string& modelParamGroup)
224 {
225 // make a structured grid
226 if (elType == Dune::cube)
227 ParentType::template makeStructuredGrid<dimension, dimensionworld>(ParentType::CellType::Cube, modelParamGroup);
228 else if (elType == Dune::simplex)
229 ParentType::template makeStructuredGrid<dimension, dimensionworld>(ParentType::CellType::Simplex, modelParamGroup);
230 else
231 DUNE_THROW(Dune::IOError, "ALUGrid only supports Dune::cube or Dune::simplex as cell type!");
232 }
233};
234
240template<int dim, int dimworld, Dune::ALUGridElementType elType, Dune::ALUGridRefinementType refinementType>
241class BoundaryFlag<Dune::ALUGrid<dim, dimworld, elType, refinementType>>
242{
243public:
244 BoundaryFlag() : flag_(-1) {}
245
246 template<class Intersection>
247 BoundaryFlag(const Intersection& i) : flag_(-1)
248 {
249 if (i.boundary())
250 flag_ = i.impl().boundaryId();
251 }
252
253 using value_type = int;
254
255 value_type get() const { return flag_; }
256
257private:
258 int flag_;
259};
260
261#endif // HAVE_DUNE_ALUGRID
262
263} // end namespace Dumux
264
265#endif
Boundary flag to store e.g. in sub control volume faces.
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:391
bool hasParam(const std::string &param)
Check whether a key exists in the parameter tree.
Definition: parameters.hh:383
Definition: adapt.hh:29
Definition: common/pdesolver.hh:35
std::size_t value_type
Definition: boundaryflag.hh:51
value_type get() const
Definition: boundaryflag.hh:53
Grid Grid
Definition: gridmanager_base.hh:75
void makeGridFromFile(const std::string &fileName, const std::string &modelParamGroup)
Makes a grid from a file. We currently support.
Definition: gridmanager_base.hh:184
void makeStructuredGrid(CellType cellType, const std::string &modelParamGroup)
Makes a structured cube grid using the structured grid factory.
Definition: gridmanager_base.hh:267
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:81