3.1-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
61template<int dim, int dimworld, Dune::ALUGridElementType elType, Dune::ALUGridRefinementType refinementType>
62class GridManager<Dune::ALUGrid<dim, dimworld, elType, refinementType>>
63: public GridManagerBase<Dune::ALUGrid<dim, dimworld, elType, refinementType>>
64{
65public:
66 using Grid = Dune::ALUGrid<dim, dimworld, elType, refinementType>;
67 using ParentType = GridManagerBase<Grid>;
68
72 void init(const std::string& modelParamGroup = "", bool adaptiveRestart = false)
73 {
74 // restarting an adaptive grid using Dune's BackupRestoreFacility
75 // TODO: the part after first || is backward compatibilty with old sequential models remove once sequential adpative restart is replaced
76 if (adaptiveRestart || hasParam("Restart") || hasParam("TimeManager.Restart"))
77 {
78 auto restartTime = getParamFromGroup<double>(modelParamGroup, "TimeLoop.Restart", 0.0);
79 // TODO: backward compatibilty with old sequential models remove once sequential adpative restart is replaced
80 if (hasParam("Restart") || hasParam("TimeManager.Restart"))
81 {
82 restartTime = getParamFromGroup<double>("TimeManager", "Restart");
83 std::cerr << "Warning: You are using a deprecated restart mechanism. The usage will change in the future.\n";
84 }
85
86 const int rank = Dune::MPIHelper::getCollectiveCommunication().rank();
87 const std::string name = getParamFromGroup<std::string>(modelParamGroup, "Problem.Name");
88 std::ostringstream oss;
89 oss << name << "_time=" << restartTime << "_rank=" << rank << ".grs";
90 std::cout << "Restoring an ALUGrid from " << oss.str() << std::endl;
91 ParentType::gridPtr() = std::shared_ptr<Grid>(Dune::BackupRestoreFacility<Grid>::restore(oss.str()));
92 ParentType::loadBalance();
93 return;
94 }
95
96 // try to create it from a DGF or msh file in GridParameterGroup.File
97 else if (hasParamInGroup(modelParamGroup, "Grid.File"))
98 {
99 makeGridFromFile(getParamFromGroup<std::string>(modelParamGroup, "Grid.File"), modelParamGroup);
100 ParentType::maybeRefineGrid(modelParamGroup);
101 ParentType::loadBalance();
102 return;
103 }
104
105 // Then look for the necessary keys to construct from the input file
106 else if (hasParamInGroup(modelParamGroup, "Grid.UpperRight"))
107 {
108 if (elType == Dune::cube)
109 makeStructuredGrid<dim, dimworld>(ParentType::CellType::Cube, modelParamGroup);
110 else if (elType == Dune::simplex)
111 makeStructuredGrid<dim, dimworld>(ParentType::CellType::Simplex, modelParamGroup);
112 else
113 DUNE_THROW(Dune::IOError, "ALUGrid only supports Dune::cube or Dune::simplex as cell type!");
114
115 ParentType::maybeRefineGrid(modelParamGroup);
116 ParentType::loadBalance();
117 }
118
119 // Didn't find a way to construct the grid
120 else
121 {
122 const auto prefix = modelParamGroup == "" ? modelParamGroup : modelParamGroup + ".";
123 DUNE_THROW(ParameterException, "Please supply one of the parameters "
124 << prefix + "Grid.UpperRight"
125 << ", or a grid file in " << prefix + "Grid.File");
126
127 }
128 }
129
133 void makeGridFromFile(const std::string& fileName,
134 const std::string& modelParamGroup)
135 {
136 // We found a file in the input file...does it have a supported extension?
137 const std::string extension = ParentType::getFileExtension(fileName);
138 if(extension != "dgf" && extension != "msh")
139 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);
140
141 // make the grid
142 if (extension == "dgf")
143 {
144 ParentType::enableDgfGridPointer_ = true;
145 ParentType::dgfGridPtr() = Dune::GridPtr<Grid>(fileName.c_str(), Dune::MPIHelper::getCommunicator());
146 ParentType::gridData_ = std::make_shared<typename ParentType::GridData>(ParentType::dgfGridPtr());
147 }
148 if (extension == "msh")
149 {
150 // get some optional parameters
151 const bool verbose = getParamFromGroup<bool>(modelParamGroup, "Grid.Verbosity", false);
152 const bool boundarySegments = getParamFromGroup<bool>(modelParamGroup, "Grid.BoundarySegments", false);
153 const bool domainMarkers = getParamFromGroup<bool>(modelParamGroup, "Grid.DomainMarkers", false);
154
155 if (domainMarkers)
156 ParentType::enableGmshDomainMarkers_ = true;
157
158 // only filll the factory for rank 0
159 if (domainMarkers)
160 {
161 std::vector<int> boundaryMarkersInsertionIndex, boundaryMarkers, faceMarkers, elementMarkers;
162 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
163
164// Older versions of Dune(-Alugrid) require that the Gmsh file is read only on process 0
165#if DUNE_VERSION_LT(DUNE_GRID, 2, 7)
166 if (Dune::MPIHelper::getCollectiveCommunication().rank() == 0)
167#endif
168 Dune::GmshReader<Grid>::read(*gridFactory, fileName, boundaryMarkersInsertionIndex, elementMarkers, verbose, boundarySegments);
169
170 ParentType::gridPtr() = std::shared_ptr<Grid>(gridFactory->createGrid());
171
172 // reorder boundary markers according to boundarySegmentIndex
173 boundaryMarkers.resize(ParentType::gridPtr()->numBoundarySegments(), 0);
174 faceMarkers.resize(ParentType::gridPtr()->leafGridView().size(1), 0);
175 const auto& indexSet = ParentType::gridPtr()->leafGridView().indexSet();
176 for (const auto& element : elements(ParentType::gridPtr()->leafGridView()))
177 {
178 for (const auto& intersection : intersections(ParentType::gridPtr()->leafGridView(), element))
179 {
180 if (intersection.boundary() && gridFactory->wasInserted(intersection))
181 {
182 auto marker = boundaryMarkersInsertionIndex[gridFactory->insertionIndex(intersection)];
183 boundaryMarkers[intersection.boundarySegmentIndex()] = marker;
184 faceMarkers[indexSet.index(element.template subEntity<1>(intersection.indexInInside()))] = marker;
185 }
186 }
187 }
188
189 ParentType::gridData_ = std::make_shared<typename ParentType::GridData>(ParentType::gridPtr(), std::move(gridFactory),
190 std::move(elementMarkers), std::move(boundaryMarkers), std::move(faceMarkers));
191 }
192 else
193 {
194 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
195
196// Older versions of Dune(-Alugrid) require that the Gmsh file is read only on process 0
197#if DUNE_VERSION_LT(DUNE_GRID, 2, 7)
198 if (Dune::MPIHelper::getCollectiveCommunication().rank() == 0)
199#endif
200 Dune::GmshReader<Grid>::read(*gridFactory, fileName, verbose, boundarySegments);
201
202 ParentType::gridPtr() = std::shared_ptr<Grid>(gridFactory->createGrid());
203 }
204 }
205 }
206
210 template <int dimension, int dimensionworld, std::enable_if_t<dimension != dimensionworld, int> = 0>
211 void makeStructuredGrid(typename ParentType::CellType cellType,
212 const std::string& modelParamGroup)
213 {
214 DUNE_THROW(Dune::IOError, "ALUGrid currently only supports the creation of structured grids with dimension == dimensionworld. Consider reading in a grid file instead.");
215 }
216
220 template <int dimension, int dimensionworld, std::enable_if_t<dimension == dimensionworld, int> = 0>
221 void makeStructuredGrid(typename ParentType::CellType cellType,
222 const std::string& modelParamGroup)
223 {
224 // make a structured grid
225 if (elType == Dune::cube)
226 ParentType::template makeStructuredGrid<dimension, dimensionworld>(ParentType::CellType::Cube, modelParamGroup);
227 else if (elType == Dune::simplex)
228 ParentType::template makeStructuredGrid<dimension, dimensionworld>(ParentType::CellType::Simplex, modelParamGroup);
229 else
230 DUNE_THROW(Dune::IOError, "ALUGrid only supports Dune::cube or Dune::simplex as cell type!");
231 }
232};
233
234#if DUNE_GRID_EXPERIMENTAL_GRID_EXTENSIONS
235
237template<int dim, int dimworld, Dune::ALUGridElementType elType, Dune::ALUGridRefinementType refinementType>
238class BoundaryFlag<Dune::ALUGrid<dim, dimworld, elType, refinementType>>
239{
240public:
241 BoundaryFlag() : flag_(-1) {}
242
243 template<class Intersection>
244 BoundaryFlag(const Intersection& i) : flag_(-1)
245 {
246 if (i.boundary())
247 flag_ = i.impl().boundaryId();
248 }
249
250 using value_type = int;
251
252 value_type get() const { return flag_; }
253
254private:
255 int flag_;
256};
257
258#endif // DUNE_GRID_EXPERIMENTAL_GRID_EXTENSIONS
259#endif // HAVE_DUNE_ALUGRID
260
261} // end namespace Dumux
262
263#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:454
bool hasParam(const std::string &param)
Check whether a key exists in the parameter tree.
Definition: parameters.hh:446
Dune::IteratorRange< typename MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper >::Intersections::const_iterator > intersections(const MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper > &glue)
Range generator to iterate with range-based for loops over all intersections as follows: for (const a...
Definition: glue.hh:62
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Definition: common/properties/model.hh:34
std::size_t value_type
Definition: boundaryflag.hh:51
value_type get() const
Definition: boundaryflag.hh:53
Grid Grid
Definition: gridmanager_base.hh:67
void makeGridFromFile(const std::string &fileName, const std::string &modelParamGroup)
Makes a grid from a file. We currently support.
Definition: gridmanager_base.hh:176
void makeStructuredGrid(CellType cellType, const std::string &modelParamGroup)
Makes a structured cube grid using the structured grid factory.
Definition: gridmanager_base.hh:259
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:73