13#ifndef DUMUX_IO_GRID_MANAGER_SP_HH
14#define DUMUX_IO_GRID_MANAGER_SP_HH
20#include <dune/grid/spgrid.hh>
21#include <dune/grid/spgrid/dgfparser.hh>
24#ifndef DUMUX_IO_GRID_MANAGER_BASE_HH
46template<
class ct,
int dim,
template<
int >
class Ref,
class Comm>
51 using Grid = Dune::SPGrid<ct, dim, Ref, Comm>;
57 void init(
const std::string& paramGroup =
"")
61 DUNE_THROW(Dune::NotImplemented,
"dune-spgrid does currently not support zero overlap!");
75 using GlobalPosition = Dune::FieldVector<ct, dim>;
79 using IntArray = std::array<int, dim>;
80 IntArray cells; cells.fill(1);
84 init(lowerLeft, upperRight, cells, paramGroup, overlap, periodic);
90 const auto prefix = paramGroup.empty() ? paramGroup : paramGroup +
".";
91 DUNE_THROW(
ParameterException,
"Please supply a grid file in " << prefix <<
"Grid.File or " << prefix <<
"Grid.UpperRight/Cells.");
95 void init(
const Dune::FieldVector<ct, dim>& lowerLeft,
96 const Dune::FieldVector<ct, dim>& upperRight,
97 const std::array<int, dim>& cells,
98 const std::string& paramGroup =
"",
99 const int overlap = 1,
100 const std::bitset<dim> periodic = std::bitset<dim>{})
103 DUNE_THROW(Dune::NotImplemented,
"dune-spgrid does currently not support zero overlap!");
104 using IntArray = std::array<int, dim>;
105 IntArray spOverlap; spOverlap.fill(overlap);
106 using Domain =
typename Grid::Domain;
107 std::vector< typename Domain::Cube > cubes;
108 cubes.push_back(
typename Domain::Cube( lowerLeft, upperRight ) );
109 Domain domain( cubes,
typename Domain::Topology(
static_cast<unsigned int>(periodic.to_ulong()) ) );
110 ParentType::gridPtr() = std::make_shared<Grid>( domain, cells, spOverlap );
111 ParentType::maybeRefineGrid(paramGroup);
112 ParentType::loadBalance();
116template<
class ct,
int dim,
template<
int >
class Ref,
class Comm>
120 using Grid = Dune::SPGrid<ct, dim, Ref, Comm>;
126 bool isPeriodic (
const typename Grid::LeafIntersection& intersection)
const
128 return intersection.neighbor() && intersection.boundary();
void init(const Dune::FieldVector< ct, dim > &lowerLeft, const Dune::FieldVector< ct, dim > &upperRight, const std::array< int, dim > &cells, const std::string ¶mGroup="", const int overlap=1, const std::bitset< dim > periodic=std::bitset< dim >{})
Definition gridmanager_sp.hh:95
Dune::SPGrid< ct, dim, Ref, Comm > Grid
Definition gridmanager_sp.hh:51
GridManagerBase< Grid > ParentType
Definition gridmanager_sp.hh:52
void init(const std::string ¶mGroup="")
Make the grid. This is implemented by specializations of this method.
Definition gridmanager_sp.hh:57
The grid manager base interface (public) and methods common to most grid manager specializations (pro...
Definition gridmanager_base.hh:56
void maybeRefineGrid(const std::string &modelParamGroup)
Definition gridmanager_base.hh:315
void makeGridFromDgfFile(const std::string &fileName)
Definition gridmanager_base.hh:265
void loadBalance()
Definition gridmanager_base.hh:98
void init(const std::string &modelParamGroup="")
Definition gridmanager_base.hh:64
The grid manager (this is the class used by the user) for all supported grid managers that constructs...
Definition gridmanager_base.hh:344
Exception thrown if a run-time parameter is not specified correctly.
Definition exceptions.hh:48
Provides a grid manager for all supported grid managers with input file interfaces....
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition parameters.hh:149
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.
Definition gridmanager_sp.hh:122
bool isPeriodic(const typename Grid::LeafIntersection &intersection) const
Definition gridmanager_sp.hh:126
PeriodicGridTraits(const Grid &grid)
Definition gridmanager_sp.hh:124