version 3.11-dev
Loading...
Searching...
No Matches
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//
13#ifndef DUMUX_IO_GRID_MANAGER_SP_HH
14#define DUMUX_IO_GRID_MANAGER_SP_HH
15
16#include <type_traits>
17
18// SPGrid specific includes
19#if HAVE_DUNE_SPGRID
20#include <dune/grid/spgrid.hh>
21#include <dune/grid/spgrid/dgfparser.hh>
22#endif
23
24#ifndef DUMUX_IO_GRID_MANAGER_BASE_HH
26#endif
28
29namespace Dumux {
30
31#if HAVE_DUNE_SPGRID
32
46template<class ct, int dim, template< int > class Ref, class Comm>
47class GridManager<Dune::SPGrid<ct, dim, Ref, Comm>>
48: public GridManagerBase<Dune::SPGrid<ct, dim, Ref, Comm>>
49{
50public:
51 using Grid = Dune::SPGrid<ct, dim, Ref, Comm>;
53
57 void init(const std::string& paramGroup = "")
58 {
59 const auto overlap = getParamFromGroup<int>(paramGroup, "Grid.Overlap", 1);
60 if (overlap == 0)
61 DUNE_THROW(Dune::NotImplemented, "dune-spgrid does currently not support zero overlap!");
62
63 // First, try to create it from file
64 if (hasParamInGroup(paramGroup, "Grid.File"))
65 {
69 return;
70 }
71
72 // Then look for the necessary keys to construct a structured grid from the input file
73 else if (hasParamInGroup(paramGroup, "Grid.UpperRight"))
74 {
75 using GlobalPosition = Dune::FieldVector<ct, dim>;
76 const auto lowerLeft = getParamFromGroup<GlobalPosition>(paramGroup, "Grid.LowerLeft", GlobalPosition(0.0));
77 const auto upperRight = getParamFromGroup<GlobalPosition>(paramGroup, "Grid.UpperRight");
78
79 using IntArray = std::array<int, dim>;
80 IntArray cells; cells.fill(1);
81 cells = getParamFromGroup<IntArray>(paramGroup, "Grid.Cells", cells);
82
83 const auto periodic = getParamFromGroup<std::bitset<dim>>(paramGroup, "Grid.Periodic", std::bitset<dim>{});
84 init(lowerLeft, upperRight, cells, paramGroup, overlap, periodic);
85 }
86
87 // Didn't find a way to construct the grid
88 else
89 {
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.");
92 }
93 }
94
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>{})
101 {
102 if (overlap == 0)
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();
113 }
114};
115
116template<class ct, int dim, template< int > class Ref, class Comm>
117struct PeriodicGridTraits<Dune::SPGrid<ct, dim, Ref, Comm>>
118{
119private:
120 using Grid = Dune::SPGrid<ct, dim, Ref, Comm>;
121public:
122 struct SupportsPeriodicity : public std::true_type {};
123
124 PeriodicGridTraits(const Grid& grid) {};
125
126 bool isPeriodic (const typename Grid::LeafIntersection& intersection) const
127 {
128 return intersection.neighbor() && intersection.boundary();
129 }
130};
131
132#endif // HAVE_DUNE_SPGRID
133
134} // end namespace Dumux
135
136#endif
void init(const Dune::FieldVector< ct, dim > &lowerLeft, const Dune::FieldVector< ct, dim > &upperRight, const std::array< int, dim > &cells, const std::string &paramGroup="", 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 &paramGroup="")
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 &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 gridmanager_sp.hh:126
PeriodicGridTraits(const Grid &grid)
Definition gridmanager_sp.hh:124