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
60template<
int dim,
int dimworld, Dune::ALUGr
idElementType elType, Dune::ALUGr
idRefinementType refinementType>
61class GridManager<
Dune::ALUGrid<dim, dimworld, elType, refinementType>>
62:
public GridManagerBase<Dune::ALUGrid<dim, dimworld, elType, refinementType>>
65 using Grid = Dune::ALUGrid<dim, dimworld, elType, refinementType>;
66 using ParentType = GridManagerBase<Grid>;
71 void init(
const std::string& modelParamGroup =
"",
bool adaptiveRestart =
false)
75 if (adaptiveRestart ||
hasParam(
"Restart") ||
hasParam(
"TimeManager.Restart"))
77 auto restartTime = getParamFromGroup<double>(modelParamGroup,
"TimeLoop.Restart", 0.0);
81 restartTime = getParamFromGroup<double>(
"TimeManager",
"Restart");
82 std::cerr <<
"Warning: You are using a deprecated restart mechanism. The usage will change in the future.\n";
85 const int rank = Dune::MPIHelper::getCollectiveCommunication().rank();
86 const std::string name = getParamFromGroup<std::string>(modelParamGroup,
"Problem.Name");
87 std::ostringstream oss;
88 oss << name <<
"_time=" << restartTime <<
"_rank=" << rank <<
".grs";
89 std::cout <<
"Restoring an ALUGrid from " << oss.str() << std::endl;
90 ParentType::gridPtr() = std::shared_ptr<Grid>(Dune::BackupRestoreFacility<Grid>::restore(oss.str()));
91 ParentType::loadBalance();
98 makeGridFromFile(getParamFromGroup<std::string>(modelParamGroup,
"Grid.File"), modelParamGroup);
99 ParentType::maybeRefineGrid(modelParamGroup);
100 ParentType::loadBalance();
107 if (elType == Dune::cube)
108 makeStructuredGrid<dim, dimworld>(ParentType::CellType::Cube, modelParamGroup);
109 else if (elType == Dune::simplex)
110 makeStructuredGrid<dim, dimworld>(ParentType::CellType::Simplex, modelParamGroup);
112 DUNE_THROW(Dune::IOError,
"ALUGrid only supports Dune::cube or Dune::simplex as cell type!");
114 ParentType::maybeRefineGrid(modelParamGroup);
115 ParentType::loadBalance();
121 const auto prefix = modelParamGroup.empty() ? modelParamGroup : modelParamGroup +
".";
122 DUNE_THROW(ParameterException,
"Please supply one of the parameters "
123 << prefix +
"Grid.UpperRight"
124 <<
", or a grid file in " << prefix +
"Grid.File");
133 const std::string& modelParamGroup)
136 const std::string extension = ParentType::getFileExtension(fileName);
137 if(extension !=
"dgf" && extension !=
"msh")
138 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 if (extension ==
"dgf")
143 ParentType::enableDgfGridPointer_ =
true;
144 ParentType::dgfGridPtr() = Dune::GridPtr<Grid>(fileName.c_str(), Dune::MPIHelper::getCommunicator());
145 ParentType::gridData_ = std::make_shared<typename ParentType::GridData>(ParentType::dgfGridPtr());
147 if (extension ==
"msh")
150 const bool verbose = getParamFromGroup<bool>(modelParamGroup,
"Grid.Verbosity",
false);
151 const bool boundarySegments = getParamFromGroup<bool>(modelParamGroup,
"Grid.BoundarySegments",
false);
152 const bool domainMarkers = getParamFromGroup<bool>(modelParamGroup,
"Grid.DomainMarkers",
false);
155 ParentType::enableGmshDomainMarkers_ =
true;
160 std::vector<int> boundaryMarkersInsertionIndex, boundaryMarkers, faceMarkers, elementMarkers;
161 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
163 Dune::GmshReader<Grid>::read(*gridFactory, fileName, boundaryMarkersInsertionIndex, elementMarkers, verbose, boundarySegments);
164 ParentType::gridPtr() = std::shared_ptr<Grid>(gridFactory->createGrid());
167 boundaryMarkers.resize(ParentType::gridPtr()->numBoundarySegments(), 0);
168 faceMarkers.resize(ParentType::gridPtr()->leafGridView().size(1), 0);
169 const auto& indexSet = ParentType::gridPtr()->leafGridView().indexSet();
170 for (
const auto& element : elements(ParentType::gridPtr()->leafGridView()))
172 for (
const auto& intersection : intersections(ParentType::gridPtr()->leafGridView(), element))
174 if (intersection.boundary() && gridFactory->wasInserted(intersection))
176 auto marker = boundaryMarkersInsertionIndex[gridFactory->insertionIndex(intersection)];
177 boundaryMarkers[intersection.boundarySegmentIndex()] = marker;
178 faceMarkers[indexSet.index(
element.template subEntity<1>(intersection.indexInInside()))] = marker;
183 ParentType::gridData_ = std::make_shared<typename ParentType::GridData>(ParentType::gridPtr(), std::move(gridFactory),
184 std::move(elementMarkers), std::move(boundaryMarkers), std::move(faceMarkers));
188 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
190 if (Dune::MPIHelper::getCollectiveCommunication().rank() == 0)
191 Dune::GmshReader<Grid>::read(*gridFactory, fileName, verbose, boundarySegments);
193 ParentType::gridPtr() = std::shared_ptr<Grid>(gridFactory->createGrid());
201 template <
int dimension,
int dimensionworld, std::enable_if_t<dimension != dimensionworld,
int> = 0>
203 const std::string& modelParamGroup)
205 DUNE_THROW(Dune::IOError,
"ALUGrid currently only supports the creation of structured grids with dimension == dimensionworld. Consider reading in a grid file instead.");
211 template <
int dimension,
int dimensionworld, std::enable_if_t<dimension == dimensionworld,
int> = 0>
213 const std::string& modelParamGroup)
216 if (elType == Dune::cube)
217 ParentType::template makeStructuredGrid<dimension, dimensionworld>(ParentType::CellType::Cube, modelParamGroup);
218 else if (elType == Dune::simplex)
219 ParentType::template makeStructuredGrid<dimension, dimensionworld>(ParentType::CellType::Simplex, modelParamGroup);
221 DUNE_THROW(Dune::IOError,
"ALUGrid only supports Dune::cube or Dune::simplex as cell type!");
230template<
int dim,
int dimworld, Dune::ALUGr
idElementType elType, Dune::ALUGr
idRefinementType refinementType>
231class BoundaryFlag<
Dune::ALUGrid<dim, dimworld, elType, refinementType>>
234 BoundaryFlag() : flag_(-1) {}
236 template<
class Intersection>
237 BoundaryFlag(
const Intersection& i) : flag_(-1)
240 flag_ = i.impl().boundaryId();
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:374
bool hasParam(const std::string ¶m)
Check whether a key exists in the parameter tree.
Definition: parameters.hh:366
Definition: common/pdesolver.hh:36
std::size_t value_type
Definition: boundaryflag.hh:51
value_type get() const
Definition: boundaryflag.hh:53
Grid Grid
Definition: gridmanager_base.hh:68
void makeGridFromFile(const std::string &fileName, const std::string &modelParamGroup)
Makes a grid from a file. We currently support.
Definition: gridmanager_base.hh:186
void makeStructuredGrid(CellType cellType, const std::string &modelParamGroup)
Makes a structured cube grid using the structured grid factory.
Definition: gridmanager_base.hh:269
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:74