3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_IO_GRID_MANAGER_SP_HH
25#define DUMUX_IO_GRID_MANAGER_SP_HH
26
27// SPGrid specific includes
28#if HAVE_DUNE_SPGRID
29#include <dune/grid/spgrid.hh>
30#include <dune/grid/spgrid/dgfparser.hh>
31#endif
32
33#ifndef DUMUX_IO_GRID_MANAGER_BASE_HH
35#endif
36
37namespace Dumux {
38
39#if HAVE_DUNE_SPGRID
40
48template<class ct, int dim, template< int > class Ref, class Comm>
49class GridManager<Dune::SPGrid<ct, dim, Ref, Comm>>
50: public GridManagerBase<Dune::SPGrid<ct, dim, Ref, Comm>>
51{
52public:
53 using Grid = Dune::SPGrid<ct, dim, Ref, Comm>;
54 using ParentType = GridManagerBase<Grid>;
55
59 void init(const std::string& paramGroup = "")
60 {
61 const auto overlap = getParamFromGroup<int>(paramGroup, "Grid.Overlap", 1);
62 if (overlap == 0)
63 DUNE_THROW(Dune::NotImplemented, "dune-spgrid does currently not support zero overlap!");
64
65 // try to create it from file
66 if (hasParamInGroup(paramGroup, "Grid.File"))
67 {
68 ParentType::makeGridFromDgfFile(getParamFromGroup<std::string>(paramGroup, "Grid.File"));
69 ParentType::maybeRefineGrid(paramGroup);
70 ParentType::loadBalance();
71 return;
72 }
73 // Didn't find a way to construct the grid
74 else if (hasParamInGroup(paramGroup, "Grid.UpperRight"))
75 {
76 using GlobalPosition = Dune::FieldVector<ct, dim>;
77 const auto lowerLeft = getParamFromGroup<GlobalPosition>(paramGroup, "Grid.LowerLeft", GlobalPosition(0.0));
78 const auto upperRight = getParamFromGroup<GlobalPosition>(paramGroup, "Grid.UpperRight");
79
80 using IntArray = std::array<int, dim>;
81 IntArray cells; cells.fill(1);
82 cells = getParamFromGroup<IntArray>(paramGroup, "Grid.Cells", cells);
83
84 const auto periodic = getParamFromGroup<std::bitset<dim>>(paramGroup, "Grid.Periodic", std::bitset<dim>{});
85 IntArray spOverlap; spOverlap.fill(overlap);
86
87 using Domain = typename Grid::Domain;
88 std::vector< typename Domain::Cube > cubes;
89 cubes.push_back( typename Domain::Cube( lowerLeft, upperRight ) );
90 Domain domain( cubes, typename Domain::Topology( static_cast<unsigned int>(periodic.to_ulong()) ) );
91 ParentType::gridPtr() = std::make_shared<Grid>( domain, cells, spOverlap );
92 ParentType::maybeRefineGrid(paramGroup);
93 ParentType::loadBalance();
94 }
95 else
96 {
97 const auto prefix = paramGroup == "" ? paramGroup : paramGroup + ".";
98 DUNE_THROW(ParameterException, "Please supply a grid file in " << prefix << "Grid.File or " << prefix << "Grid.UpperRight/Cells.");
99 }
100 }
101};
102
103#endif // HAVE_DUNE_SPGRID
104
105} // end namespace Dumux
106
107#endif
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:454
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Definition: common/properties/model.hh:34
Grid Grid
Definition: gridmanager_base.hh:67
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:73