12#ifndef DUMUX_IO_GRID_MANAGER_SP_HH
13#define DUMUX_IO_GRID_MANAGER_SP_HH
17#include <dune/grid/spgrid.hh>
18#include <dune/grid/spgrid/dgfparser.hh>
21#ifndef DUMUX_IO_GRID_MANAGER_BASE_HH
37template<
class ct,
int dim,
template<
int >
class Ref,
class Comm>
38class GridManager<
Dune::SPGrid<ct, dim, Ref, Comm>>
39:
public GridManagerBase<Dune::SPGrid<ct, dim, Ref, Comm>>
43 using ParentType = GridManagerBase<Grid>;
48 void init(
const std::string& paramGroup =
"")
50 const auto overlap = getParamFromGroup<int>(paramGroup,
"Grid.Overlap", 1);
52 DUNE_THROW(Dune::NotImplemented,
"dune-spgrid does currently not support zero overlap!");
57 ParentType::makeGridFromDgfFile(getParamFromGroup<std::string>(paramGroup,
"Grid.File"));
58 ParentType::maybeRefineGrid(paramGroup);
59 ParentType::loadBalance();
65 using GlobalPosition = Dune::FieldVector<ct, dim>;
66 const auto lowerLeft = getParamFromGroup<GlobalPosition>(paramGroup,
"Grid.LowerLeft", GlobalPosition(0.0));
67 const auto upperRight = getParamFromGroup<GlobalPosition>(paramGroup,
"Grid.UpperRight");
69 using IntArray = std::array<int, dim>;
70 IntArray cells; cells.fill(1);
71 cells = getParamFromGroup<IntArray>(paramGroup,
"Grid.Cells", cells);
73 const auto periodic = getParamFromGroup<std::bitset<dim>>(paramGroup,
"Grid.Periodic", std::bitset<dim>{});
74 init(lowerLeft, upperRight, cells, paramGroup, overlap, periodic);
78 const auto prefix = paramGroup.empty() ? paramGroup : paramGroup +
".";
79 DUNE_THROW(ParameterException,
"Please supply a grid file in " << prefix <<
"Grid.File or " << prefix <<
"Grid.UpperRight/Cells.");
83 void init(
const Dune::FieldVector<ct, dim>& lowerLeft,
84 const Dune::FieldVector<ct, dim>& upperRight,
85 const std::array<int, dim>& cells,
86 const std::string& paramGroup =
"",
87 const int overlap = 1,
88 const std::bitset<dim> periodic = std::bitset<dim>{})
91 DUNE_THROW(Dune::NotImplemented,
"dune-spgrid does currently not support zero overlap!");
92 using IntArray = std::array<int, dim>;
93 IntArray spOverlap; spOverlap.fill(overlap);
94 using Domain =
typename Grid::Domain;
95 std::vector< typename Domain::Cube > cubes;
96 cubes.push_back(
typename Domain::Cube( lowerLeft, upperRight ) );
97 Domain domain( cubes,
typename Domain::Topology(
static_cast<unsigned int>(periodic.to_ulong()) ) );
98 ParentType::gridPtr() = std::make_shared<Grid>( domain, cells, spOverlap );
99 ParentType::maybeRefineGrid(paramGroup);
100 ParentType::loadBalance();
Grid Grid
Definition: gridmanager_base.hh:57
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
Definition: momentum/velocitygradients.hh:26
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