24#ifndef DUMUX_IO_GRID_MANAGER_UG_HH
25#define DUMUX_IO_GRID_MANAGER_UG_HH
28#include <dune/grid/uggrid.hh>
29#include <dune/grid/io/file/dgfparser/dgfug.hh>
32#ifndef DUMUX_IO_GRID_MANAGER_BASE_HH
61class GridManager<
Dune::UGGrid<dim>>
62:
public GridManagerBase<Dune::UGGrid<dim>>
65 using Grid =
typename Dune::UGGrid<dim>;
66 using ParentType = GridManagerBase<Grid>;
67 using Element =
typename Grid::template Codim<0>::Entity;
72 void init(
const std::string& modelParamGroup =
"")
78 preProcessing_(modelParamGroup);
79 ParentType::makeGridFromFile(getParamFromGroup<std::string>(modelParamGroup,
"Grid.File"), modelParamGroup);
80 postProcessing_(modelParamGroup);
87 preProcessing_(modelParamGroup);
89 const auto cellType = getParamFromGroup<std::string>(modelParamGroup,
"Grid.CellType",
"Cube");
90 if (cellType ==
"Cube")
91 ParentType::template makeStructuredGrid<dim, dim>(ParentType::CellType::Cube, modelParamGroup);
92 else if (cellType ==
"Simplex")
93 ParentType::template makeStructuredGrid<dim, dim>(ParentType::CellType::Simplex, modelParamGroup);
95 DUNE_THROW(Dune::IOError,
"UGGrid only supports 'Cube' or 'Simplex' as cell type. Not '"<< cellType<<
"'!");
96 postProcessing_(modelParamGroup);
102 const auto prefix = modelParamGroup.empty() ? modelParamGroup : modelParamGroup +
".";
103 DUNE_THROW(ParameterException,
"Please supply one of the parameters "
104 << prefix +
"Grid.UpperRight"
105 <<
", or a grid file in " << prefix +
"Grid.File");
119 if (Dune::MPIHelper::getCommunication().size() > 1)
122 if(ParentType::enableDgfGridPointer_)
124 ParentType::dgfGridPtr().loadBalance();
126 ParentType::gridData_ = std::make_shared<typename ParentType::GridData>(ParentType::dgfGridPtr());
130 else if (ParentType::enableGmshDomainMarkers_)
133 auto dh = ParentType::gridData_->createGmshDataHandle();
134 ParentType::gridPtr()->loadBalance(dh.interface());
139 ParentType::gridPtr()->loadBalance();
147 void preProcessing_(
const std::string& modelParamGroup)
153 void postProcessing_(
const std::string& modelParamGroup)
156 const auto refType = getParamFromGroup<std::string>(modelParamGroup,
"Grid.RefinementType",
"Local");
157 if (refType ==
"Local")
158 ParentType::grid().setRefinementType(Dune::UGGrid<dim>::RefinementType::LOCAL);
159 else if (refType ==
"Copy")
160 ParentType::grid().setRefinementType(Dune::UGGrid<dim>::RefinementType::COPY);
162 DUNE_THROW(Dune::IOError,
"UGGrid only supports 'Local' or 'Copy' as refinement type. Not '"<< refType<<
"'!");
165 const auto closureType = getParamFromGroup<std::string>(modelParamGroup,
"Grid.ClosureType",
"Green");
166 if (closureType ==
"Green")
167 ParentType::grid().setClosureType(Dune::UGGrid<dim>::ClosureType::GREEN);
168 else if (closureType ==
"None")
169 ParentType::grid().setClosureType(Dune::UGGrid<dim>::ClosureType::NONE);
171 DUNE_THROW(Dune::IOError,
"UGGrid only supports 'Green' or 'None' as closure type. Not '"<< closureType<<
"'!");
174 ParentType::maybeRefineGrid(modelParamGroup);
180namespace Grid::Capabilities {
185struct MultithreadingSupported<
Dune::UGGrid<dim>>
188 static bool eval(
const GV& gv)
189 {
return gv.comm().size() <= 1; }
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:177
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
Definition: deprecated.hh:149
static bool eval(const GV &)
Definition: gridcapabilities.hh:81
Grid Grid
Definition: gridmanager_base.hh:69
void loadBalance()
Call loadBalance() function of the grid.
Definition: gridmanager_base.hh:98
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