12#ifndef DUMUX_IO_GRID_MANAGER_ALU_HH
13#define DUMUX_IO_GRID_MANAGER_ALU_HH
17#include <dune/alugrid/grid.hh>
18#include <dune/alugrid/dgf.hh>
21#ifndef DUMUX_IO_GRID_MANAGER_BASE_HH
49template<
int dim,
int dimworld, Dune::ALUGr
idElementType elType, Dune::ALUGr
idRefinementType refinementType>
50class GridManager<
Dune::ALUGrid<dim, dimworld, elType, refinementType>>
51:
public GridManagerBase<Dune::ALUGrid<dim, dimworld, elType, refinementType>>
54 using Grid = Dune::ALUGrid<dim, dimworld, elType, refinementType>;
55 using ParentType = GridManagerBase<Grid>;
60 void init(
const std::string& modelParamGroup =
"",
bool adaptiveRestart =
false)
64 if (adaptiveRestart ||
hasParam(
"Restart") ||
hasParam(
"TimeManager.Restart"))
66 auto restartTime = getParamFromGroup<double>(modelParamGroup,
"TimeLoop.Restart", 0.0);
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";
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();
87 makeGridFromFile(getParamFromGroup<std::string>(modelParamGroup,
"Grid.File"), modelParamGroup);
88 ParentType::maybeRefineGrid(modelParamGroup);
89 ParentType::loadBalance();
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);
101 DUNE_THROW(Dune::IOError,
"ALUGrid only supports Dune::cube or Dune::simplex as cell type!");
103 ParentType::maybeRefineGrid(modelParamGroup);
104 ParentType::loadBalance();
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");
122 const std::string& modelParamGroup)
125 const std::string extension = ParentType::getFileExtension(fileName);
126 if (extension !=
"dgf" && extension !=
"msh" && extension !=
"vtu")
127 DUNE_THROW(Dune::IOError,
"Grid type " << Dune::className<Grid>() <<
" doesn't support grid files with extension: *."<< extension);
130 if (extension ==
"dgf")
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());
138 else if (extension ==
"msh")
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);
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());
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()))
159 for (
const auto& intersection : intersections(ParentType::gridPtr()->leafGridView(), element))
161 if (intersection.boundary() && gridFactory->wasInserted(intersection))
163 auto marker = boundaryMarkersInsertionIndex[gridFactory->insertionIndex(intersection)];
164 boundaryMarkers[intersection.boundarySegmentIndex()] = marker;
165 faceMarkers[indexSet.index(
element.template subEntity<1>(intersection.indexInInside()))] = marker;
170 ParentType::gridData_ = std::make_shared<typename ParentType::GridData>(ParentType::gridPtr(), std::move(gridFactory),
171 std::move(elementMarkers), std::move(boundaryMarkers), std::move(faceMarkers));
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());
182 else if (extension ==
"vtu")
184 VTKReader vtkReader(fileName);
186 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
187 const bool verbose = getParamFromGroup<bool>(modelParamGroup,
"Grid.Verbosity",
false);
188 ParentType::gridPtr() = vtkReader.readGrid(*gridFactory, cellData, pointData, verbose);
189 ParentType::gridData_ = std::make_shared<typename ParentType::GridData>(ParentType::gridPtr(), std::move(gridFactory), std::move(cellData), std::move(pointData));
190 ParentType::enableVtkData_ =
true;
197 template <
int dimension,
int dimensionworld, std::enable_if_t<dimension != dimensionworld,
int> = 0>
199 const std::string& modelParamGroup)
201 DUNE_THROW(Dune::IOError,
"ALUGrid currently only supports the creation of structured grids with dimension == dimensionworld. Consider reading in a grid file instead.");
207 template <
int dimension,
int dimensionworld, std::enable_if_t<dimension == dimensionworld,
int> = 0>
209 const std::string& modelParamGroup)
212 if (elType == Dune::cube)
213 ParentType::template makeStructuredGrid<dimension, dimensionworld>(ParentType::CellType::Cube, modelParamGroup);
214 else if (elType == Dune::simplex)
215 ParentType::template makeStructuredGrid<dimension, dimensionworld>(ParentType::CellType::Simplex, modelParamGroup);
217 DUNE_THROW(Dune::IOError,
"ALUGrid only supports Dune::cube or Dune::simplex as cell type!");
226template<
int dim,
int dimworld, Dune::ALUGr
idElementType elType, Dune::ALUGr
idRefinementType refinementType>
227class BoundaryFlag<
Dune::ALUGrid<dim, dimworld, elType, refinementType>>
230 BoundaryFlag() : flag_(-1) {}
232 template<
class Intersection>
233 BoundaryFlag(
const Intersection& i) : flag_(-1)
236 flag_ = i.impl().boundaryId();
247namespace Grid::Capabilities {
251template<
int dim,
int dimworld, Dune::ALUGr
idElementType elType, Dune::ALUGr
idRefinementType refinementType>
252struct MultithreadingSupported<
Dune::ALUGrid<dim, dimworld, elType, refinementType>>
255 static bool eval(
const GV&)
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:197
void makeStructuredGrid(CellType cellType, const std::string &modelParamGroup)
Makes a structured cube grid using the structured grid factory.
Definition: gridmanager_base.hh:278
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
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:48
dune-grid capabilities compatibility layer
Provides a grid manager for all supported grid managers with input file interfaces....
bool hasParamInGroup(const std::string ¶mGroup, const std::string ¶m)
Check whether a key exists in the parameter tree with a model group prefix.
Definition: parameters.hh:165
bool hasParam(const std::string ¶m)
Check whether a key exists in the parameter tree.
Definition: parameters.hh:157
Definition: common/pdesolver.hh:24
static bool eval(const GV &)
Definition: gridcapabilities.hh:69