version 3.11-dev
gridmanager_sp.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_IO_GRID_MANAGER_SP_HH
13#define DUMUX_IO_GRID_MANAGER_SP_HH
14
15#include <type_traits>
16
17// SPGrid specific includes
18#if HAVE_DUNE_SPGRID
19#include <dune/grid/spgrid.hh>
20#include <dune/grid/spgrid/dgfparser.hh>
21#endif
22
23#ifndef DUMUX_IO_GRID_MANAGER_BASE_HH
25#endif
27
28namespace Dumux {
29
30#if HAVE_DUNE_SPGRID
31
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>>
43{
44public:
45 using Grid = Dune::SPGrid<ct, dim, Ref, Comm>;
46 using ParentType = GridManagerBase<Grid>;
47
51 void init(const std::string& paramGroup = "")
52 {
53 const auto overlap = getParamFromGroup<int>(paramGroup, "Grid.Overlap", 1);
54 if (overlap == 0)
55 DUNE_THROW(Dune::NotImplemented, "dune-spgrid does currently not support zero overlap!");
56
57 // try to create it from file
58 if (hasParamInGroup(paramGroup, "Grid.File"))
59 {
60 ParentType::makeGridFromDgfFile(getParamFromGroup<std::string>(paramGroup, "Grid.File"));
61 ParentType::maybeRefineGrid(paramGroup);
62 ParentType::loadBalance();
63 return;
64 }
65 // Didn't find a way to construct the grid
66 else if (hasParamInGroup(paramGroup, "Grid.UpperRight"))
67 {
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");
71
72 using IntArray = std::array<int, dim>;
73 IntArray cells; cells.fill(1);
74 cells = getParamFromGroup<IntArray>(paramGroup, "Grid.Cells", cells);
75
76 const auto periodic = getParamFromGroup<std::bitset<dim>>(paramGroup, "Grid.Periodic", std::bitset<dim>{});
77 init(lowerLeft, upperRight, cells, paramGroup, overlap, periodic);
78 }
79 else
80 {
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.");
83 }
84 }
85
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>{})
92 {
93 if (overlap == 0)
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();
104 }
105};
106
107template<class ct, int dim, template< int > class Ref, class Comm>
108struct PeriodicGridTraits<Dune::SPGrid<ct, dim, Ref, Comm>>
109{
110private:
111 using Grid = Dune::SPGrid<ct, dim, Ref, Comm>;
112public:
113 struct SupportsPeriodicity : public std::true_type {};
114
115 PeriodicGridTraits(const Grid& grid) {};
116
117 bool isPeriodic (const typename Grid::LeafIntersection& intersection) const
118 {
119 return intersection.neighbor() && intersection.boundary();
120 }
121};
122
123#endif // HAVE_DUNE_SPGRID
124
125} // end namespace Dumux
126
127#endif
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 &paramGroup, const std::string &param)
Check whether a key exists in the parameter tree with a model group prefix.
Definition: parameters.hh:165
Definition: adapt.hh:17
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