24#ifndef DUMUX_IO_GRID_MANAGER_ALU_HH
25#define DUMUX_IO_GRID_MANAGER_ALU_HH
27#include <dune/common/version.hh>
31#include <dune/alugrid/grid.hh>
32#include <dune/alugrid/dgf.hh>
35#ifndef DUMUX_IO_GRID_MANAGER_BASE_HH
62template<
int dim,
int dimworld, Dune::ALUGr
idElementType elType, Dune::ALUGr
idRefinementType refinementType>
63class GridManager<
Dune::ALUGrid<dim, dimworld, elType, refinementType>>
64:
public GridManagerBase<Dune::ALUGrid<dim, dimworld, elType, refinementType>>
67 using Grid = Dune::ALUGrid<dim, dimworld, elType, refinementType>;
68 using ParentType = GridManagerBase<Grid>;
73 void init(
const std::string& modelParamGroup =
"",
bool adaptiveRestart =
false)
77 if (adaptiveRestart ||
hasParam(
"Restart") ||
hasParam(
"TimeManager.Restart"))
79 auto restartTime = getParamFromGroup<double>(modelParamGroup,
"TimeLoop.Restart", 0.0);
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";
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();
100 makeGridFromFile(getParamFromGroup<std::string>(modelParamGroup,
"Grid.File"), modelParamGroup);
101 ParentType::maybeRefineGrid(modelParamGroup);
102 ParentType::loadBalance();
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);
114 DUNE_THROW(Dune::IOError,
"ALUGrid only supports Dune::cube or Dune::simplex as cell type!");
116 ParentType::maybeRefineGrid(modelParamGroup);
117 ParentType::loadBalance();
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");
135 const std::string& modelParamGroup)
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);
143 if (extension ==
"dgf")
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());
149 if (extension ==
"msh")
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);
157 ParentType::enableGmshDomainMarkers_ =
true;
162 std::vector<int> boundaryMarkersInsertionIndex, boundaryMarkers, faceMarkers, elementMarkers;
163 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
166#if DUNE_VERSION_LT(DUNE_GRID, 2, 7)
167 if (Dune::MPIHelper::getCollectiveCommunication().rank() == 0)
169 Dune::GmshReader<Grid>::read(*gridFactory, fileName, boundaryMarkersInsertionIndex, elementMarkers, verbose, boundarySegments);
171 ParentType::gridPtr() = std::shared_ptr<Grid>(gridFactory->createGrid());
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()))
179 for (
const auto& intersection : intersections(ParentType::gridPtr()->leafGridView(), element))
181 if (intersection.boundary() && gridFactory->wasInserted(intersection))
183 auto marker = boundaryMarkersInsertionIndex[gridFactory->insertionIndex(intersection)];
184 boundaryMarkers[intersection.boundarySegmentIndex()] = marker;
185 faceMarkers[indexSet.index(element.template subEntity<1>(intersection.indexInInside()))] = marker;
190 ParentType::gridData_ = std::make_shared<typename ParentType::GridData>(ParentType::gridPtr(), std::move(gridFactory),
191 std::move(elementMarkers), std::move(boundaryMarkers), std::move(faceMarkers));
195 auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>();
198#if DUNE_VERSION_LT(DUNE_GRID, 2, 7)
199 if (Dune::MPIHelper::getCollectiveCommunication().rank() == 0)
201 Dune::GmshReader<Grid>::read(*gridFactory, fileName, verbose, boundarySegments);
203 ParentType::gridPtr() = std::shared_ptr<Grid>(gridFactory->createGrid());
211 template <
int dimension,
int dimensionworld, std::enable_if_t<dimension != dimensionworld,
int> = 0>
213 const std::string& modelParamGroup)
215 DUNE_THROW(Dune::IOError,
"ALUGrid currently only supports the creation of structured grids with dimension == dimensionworld. Consider reading in a grid file instead.");
221 template <
int dimension,
int dimensionworld, std::enable_if_t<dimension == dimensionworld,
int> = 0>
223 const std::string& modelParamGroup)
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);
231 DUNE_THROW(Dune::IOError,
"ALUGrid only supports Dune::cube or Dune::simplex as cell type!");
240template<
int dim,
int dimworld, Dune::ALUGr
idElementType elType, Dune::ALUGr
idRefinementType refinementType>
241class BoundaryFlag<
Dune::ALUGrid<dim, dimworld, elType, refinementType>>
244 BoundaryFlag() : flag_(-1) {}
246 template<
class Intersection>
247 BoundaryFlag(
const Intersection& i) : flag_(-1)
250 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:391
bool hasParam(const std::string ¶m)
Check whether a key exists in the parameter tree.
Definition: parameters.hh:383
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