24#ifndef DUMUX_IO_GRID_MANAGER_ALU_HH
25#define DUMUX_IO_GRID_MANAGER_ALU_HH
29#include <dune/alugrid/grid.hh>
30#include <dune/alugrid/dgf.hh>
33#ifndef DUMUX_IO_GRID_MANAGER_BASE_HH
61template<
int dim,
int dimworld, Dune::ALUGr
idElementType elType, Dune::ALUGr
idRefinementType refinementType>
62class GridManager<
Dune::ALUGrid<dim, dimworld, elType, refinementType>>
63:
public GridManagerBase<Dune::ALUGrid<dim, dimworld, elType, refinementType>>
66 using Grid = Dune::ALUGrid<dim, dimworld, elType, refinementType>;
67 using ParentType = GridManagerBase<Grid>;
72 void init(
const std::string& modelParamGroup =
"",
bool adaptiveRestart =
false)
76 if (adaptiveRestart ||
hasParam(
"Restart") ||
hasParam(
"TimeManager.Restart"))
78 auto restartTime = getParamFromGroup<double>(modelParamGroup,
"TimeLoop.Restart", 0.0);
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";
86 const int rank = Dune::MPIHelper::getCommunication().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();
99 makeGridFromFile(getParamFromGroup<std::string>(modelParamGroup,
"Grid.File"), modelParamGroup);
100 ParentType::maybeRefineGrid(modelParamGroup);
101 ParentType::loadBalance();
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);
113 DUNE_THROW(Dune::IOError,
"ALUGrid only supports Dune::cube or Dune::simplex as cell type!");
115 ParentType::maybeRefineGrid(modelParamGroup);
116 ParentType::loadBalance();
122 const auto prefix = modelParamGroup.empty() ? 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");
134 const std::string& modelParamGroup)
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);
142 if (extension ==
"dgf")
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());
148 if (extension ==
"msh")
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);
156 ParentType::enableGmshDomainMarkers_ =
true;
161 std::vector<int> boundaryMarkersInsertionIndex, boundaryMarkers, faceMarkers, elementMarkers;
162 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
164 Dune::GmshReader<Grid>::read(*gridFactory, fileName, boundaryMarkersInsertionIndex, elementMarkers, verbose, boundarySegments);
165 ParentType::gridPtr() = std::shared_ptr<Grid>(gridFactory->createGrid());
168 boundaryMarkers.resize(ParentType::gridPtr()->numBoundarySegments(), 0);
169 faceMarkers.resize(ParentType::gridPtr()->leafGridView().size(1), 0);
170 const auto& indexSet = ParentType::gridPtr()->leafGridView().indexSet();
171 for (
const auto& element : elements(ParentType::gridPtr()->leafGridView()))
173 for (
const auto& intersection : intersections(ParentType::gridPtr()->leafGridView(), element))
175 if (intersection.boundary() && gridFactory->wasInserted(intersection))
177 auto marker = boundaryMarkersInsertionIndex[gridFactory->insertionIndex(intersection)];
178 boundaryMarkers[intersection.boundarySegmentIndex()] = marker;
179 faceMarkers[indexSet.index(
element.template subEntity<1>(intersection.indexInInside()))] = marker;
184 ParentType::gridData_ = std::make_shared<typename ParentType::GridData>(ParentType::gridPtr(), std::move(gridFactory),
185 std::move(elementMarkers), std::move(boundaryMarkers), std::move(faceMarkers));
189 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
191 if (Dune::MPIHelper::getCommunication().rank() == 0)
192 Dune::GmshReader<Grid>::read(*gridFactory, fileName, verbose, boundarySegments);
194 ParentType::gridPtr() = std::shared_ptr<Grid>(gridFactory->createGrid());
202 template <
int dimension,
int dimensionworld, std::enable_if_t<dimension != dimensionworld,
int> = 0>
204 const std::string& modelParamGroup)
206 DUNE_THROW(Dune::IOError,
"ALUGrid currently only supports the creation of structured grids with dimension == dimensionworld. Consider reading in a grid file instead.");
212 template <
int dimension,
int dimensionworld, std::enable_if_t<dimension == dimensionworld,
int> = 0>
214 const std::string& modelParamGroup)
217 if (elType == Dune::cube)
218 ParentType::template makeStructuredGrid<dimension, dimensionworld>(ParentType::CellType::Cube, modelParamGroup);
219 else if (elType == Dune::simplex)
220 ParentType::template makeStructuredGrid<dimension, dimensionworld>(ParentType::CellType::Simplex, modelParamGroup);
222 DUNE_THROW(Dune::IOError,
"ALUGrid only supports Dune::cube or Dune::simplex as cell type!");
231template<
int dim,
int dimworld, Dune::ALUGr
idElementType elType, Dune::ALUGr
idRefinementType refinementType>
232class BoundaryFlag<
Dune::ALUGrid<dim, dimworld, elType, refinementType>>
235 BoundaryFlag() : flag_(-1) {}
237 template<
class Intersection>
238 BoundaryFlag(
const Intersection& i) : flag_(-1)
241 flag_ = i.impl().boundaryId();
252namespace Grid::Capabilities {
256template<
int dim,
int dimworld, Dune::ALUGr
idElementType elType, Dune::ALUGr
idRefinementType refinementType>
257struct MultithreadingSupported<
Dune::ALUGrid<dim, dimworld, elType, refinementType>>
260 static bool eval(
const GV&)
dune-grid capabilities compatibility layer
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 ¶mGroup, const std::string ¶m)
Check whether a key exists in the parameter tree with a model group prefix.
Definition: parameters.hh:177
bool hasParam(const std::string ¶m)
Check whether a key exists in the parameter tree.
Definition: parameters.hh:169
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
Definition: deprecated.hh:149
std::size_t value_type
Definition: boundaryflag.hh:51
value_type get() const
Definition: boundaryflag.hh:53
static bool eval(const GV &)
Definition: gridcapabilities.hh:81
Grid Grid
Definition: gridmanager_base.hh:69
void makeGridFromFile(const std::string &fileName, const std::string &modelParamGroup)
Makes a grid from a file. We currently support.
Definition: gridmanager_base.hh:188
void makeStructuredGrid(CellType cellType, const std::string &modelParamGroup)
Makes a structured cube grid using the structured grid factory.
Definition: gridmanager_base.hh:271
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:75