12#ifndef DUMUX_IO_GRID_MANAGER_UG_HH
13#define DUMUX_IO_GRID_MANAGER_UG_HH
16#include <dune/grid/uggrid.hh>
17#include <dune/grid/io/file/dgfparser/dgfug.hh>
20#ifndef DUMUX_IO_GRID_MANAGER_BASE_HH
49class GridManager<
Dune::UGGrid<dim>>
50:
public GridManagerBase<Dune::UGGrid<dim>>
53 using Grid =
typename Dune::UGGrid<dim>;
54 using ParentType = GridManagerBase<Grid>;
55 using Element =
typename Grid::template Codim<0>::Entity;
60 void init(
const std::string& modelParamGroup =
"")
66 preProcessing_(modelParamGroup);
67 ParentType::makeGridFromFile(getParamFromGroup<std::string>(modelParamGroup,
"Grid.File"), modelParamGroup);
68 postProcessing_(modelParamGroup);
75 preProcessing_(modelParamGroup);
77 const auto cellType = getParamFromGroup<std::string>(modelParamGroup,
"Grid.CellType",
"Cube");
78 if (cellType ==
"Cube")
79 ParentType::template makeStructuredGrid<dim, dim>(ParentType::CellType::Cube, modelParamGroup);
80 else if (cellType ==
"Simplex")
81 ParentType::template makeStructuredGrid<dim, dim>(ParentType::CellType::Simplex, modelParamGroup);
83 DUNE_THROW(Dune::IOError,
"UGGrid only supports 'Cube' or 'Simplex' as cell type. Not '"<< cellType<<
"'!");
84 postProcessing_(modelParamGroup);
90 const auto prefix = modelParamGroup.empty() ? modelParamGroup : modelParamGroup +
".";
91 DUNE_THROW(ParameterException,
"Please supply one of the parameters "
92 << prefix +
"Grid.UpperRight"
93 <<
", or a grid file in " << prefix +
"Grid.File");
107 if (Dune::MPIHelper::getCommunication().size() > 1)
110 if(ParentType::enableDgfGridPointer_)
112 ParentType::dgfGridPtr().loadBalance();
114 ParentType::gridData_ = std::make_shared<typename ParentType::GridData>(ParentType::dgfGridPtr());
118 else if (ParentType::enableGmshDomainMarkers_)
121 auto dh = ParentType::gridData_->createGmshDataHandle();
122 ParentType::gridPtr()->loadBalance(dh.interface());
128 else if (ParentType::enableVtkData_)
131 auto dh = ParentType::gridData_->createVtkDataHandle();
132 ParentType::gridPtr()->loadBalance(dh.interface());
138 ParentType::gridPtr()->loadBalance();
146 void preProcessing_(
const std::string& modelParamGroup)
152 void postProcessing_(
const std::string& modelParamGroup)
155 const auto refType = getParamFromGroup<std::string>(modelParamGroup,
"Grid.RefinementType",
"Local");
156 if (refType ==
"Local")
157 ParentType::grid().setRefinementType(Dune::UGGrid<dim>::RefinementType::LOCAL);
158 else if (refType ==
"Copy")
159 ParentType::grid().setRefinementType(Dune::UGGrid<dim>::RefinementType::COPY);
161 DUNE_THROW(Dune::IOError,
"UGGrid only supports 'Local' or 'Copy' as refinement type. Not '"<< refType<<
"'!");
164 const auto closureType = getParamFromGroup<std::string>(modelParamGroup,
"Grid.ClosureType",
"Green");
165 if (closureType ==
"Green")
166 ParentType::grid().setClosureType(Dune::UGGrid<dim>::ClosureType::GREEN);
167 else if (closureType ==
"None")
168 ParentType::grid().setClosureType(Dune::UGGrid<dim>::ClosureType::NONE);
170 DUNE_THROW(Dune::IOError,
"UGGrid only supports 'Green' or 'None' as closure type. Not '"<< closureType<<
"'!");
173 ParentType::maybeRefineGrid(modelParamGroup);
179namespace Grid::Capabilities {
184struct MultithreadingSupported<
Dune::UGGrid<dim>>
187 static bool eval(
const GV& gv)
188 {
return gv.comm().size() <= 1; }
Grid Grid
Definition: gridmanager_base.hh:57
void loadBalance()
Call loadBalance() function of the grid.
Definition: gridmanager_base.hh:97
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
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
Definition: common/pdesolver.hh:24
static bool eval(const GV &)
Definition: gridcapabilities.hh:69