12#ifndef DUMUX_IO_GRID_MANAGER_SP_HH
13#define DUMUX_IO_GRID_MANAGER_SP_HH
19#include <dune/grid/spgrid.hh>
20#include <dune/grid/spgrid/dgfparser.hh>
23#ifndef DUMUX_IO_GRID_MANAGER_BASE_HH
40template<
class ct,
int dim,
template<
int >
class Ref,
class Comm>
41class GridManager<
Dune::SPGrid<ct, dim, Ref, Comm>>
42:
public GridManagerBase<Dune::SPGrid<ct, dim, Ref, Comm>>
45 using Grid = Dune::SPGrid<ct, dim, Ref, Comm>;
46 using ParentType = GridManagerBase<Grid>;
51 void init(
const std::string& paramGroup =
"")
53 const auto overlap = getParamFromGroup<int>(paramGroup,
"Grid.Overlap", 1);
55 DUNE_THROW(Dune::NotImplemented,
"dune-spgrid does currently not support zero overlap!");
60 ParentType::makeGridFromDgfFile(getParamFromGroup<std::string>(paramGroup,
"Grid.File"));
61 ParentType::maybeRefineGrid(paramGroup);
62 ParentType::loadBalance();
68 using GlobalPosition = Dune::FieldVector<ct, dim>;
69 const auto lowerLeft = getParamFromGroup<GlobalPosition>(paramGroup,
"Grid.LowerLeft", GlobalPosition(0.0));
70 const auto upperRight = getParamFromGroup<GlobalPosition>(paramGroup,
"Grid.UpperRight");
72 using IntArray = std::array<int, dim>;
73 IntArray cells; cells.fill(1);
74 cells = getParamFromGroup<IntArray>(paramGroup,
"Grid.Cells", cells);
76 const auto periodic = getParamFromGroup<std::bitset<dim>>(paramGroup,
"Grid.Periodic", std::bitset<dim>{});
77 init(lowerLeft, upperRight, cells, paramGroup, overlap, periodic);
81 const auto prefix = paramGroup.empty() ? paramGroup : paramGroup +
".";
82 DUNE_THROW(ParameterException,
"Please supply a grid file in " << prefix <<
"Grid.File or " << prefix <<
"Grid.UpperRight/Cells.");
86 void init(
const Dune::FieldVector<ct, dim>& lowerLeft,
87 const Dune::FieldVector<ct, dim>& upperRight,
88 const std::array<int, dim>& cells,
89 const std::string& paramGroup =
"",
90 const int overlap = 1,
91 const std::bitset<dim> periodic = std::bitset<dim>{})
94 DUNE_THROW(Dune::NotImplemented,
"dune-spgrid does currently not support zero overlap!");
95 using IntArray = std::array<int, dim>;
96 IntArray spOverlap; spOverlap.fill(overlap);
97 using Domain =
typename Grid::Domain;
98 std::vector< typename Domain::Cube > cubes;
99 cubes.push_back(
typename Domain::Cube( lowerLeft, upperRight ) );
100 Domain domain( cubes,
typename Domain::Topology(
static_cast<unsigned int>(periodic.to_ulong()) ) );
101 ParentType::gridPtr() = std::make_shared<Grid>( domain, cells, spOverlap );
102 ParentType::maybeRefineGrid(paramGroup);
103 ParentType::loadBalance();
107template<
class ct,
int dim,
template<
int >
class Ref,
class Comm>
108struct PeriodicGridTraits<
Dune::SPGrid<ct, dim, Ref, Comm>>
111 using Grid = Dune::SPGrid<ct, dim, Ref, Comm>;
113 struct SupportsPeriodicity :
public std::true_type {};
117 bool isPeriodic (
const typename Grid::LeafIntersection& intersection)
const
119 return intersection.neighbor() && intersection.boundary();
Grid Grid
Definition: gridmanager_base.hh:58
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:64
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
Grid properties related to periodicity.
bool isPeriodic(const typename Grid::LeafIntersection &intersection) const
Definition: periodicgridtraits.hh:28
PeriodicGridTraits(const Grid &grid)
Definition: periodicgridtraits.hh:26