version 3.8
gridmanager_ug.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-FileCopyrightInfo: 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_UG_HH
13#define DUMUX_IO_GRID_MANAGER_UG_HH
14
15#if HAVE_DUNE_UGGRID
16#include <dune/grid/uggrid.hh>
17#include <dune/grid/io/file/dgfparser/dgfug.hh>
18#endif
19
20#ifndef DUMUX_IO_GRID_MANAGER_BASE_HH
22#endif
23
25
26namespace Dumux {
27
28#if HAVE_DUNE_UGGRID
29
48template<int dim>
49class GridManager<Dune::UGGrid<dim>>
50: public GridManagerBase<Dune::UGGrid<dim>>
51{
52public:
53 using Grid = typename Dune::UGGrid<dim>;
54 using ParentType = GridManagerBase<Grid>;
55 using Element = typename Grid::template Codim<0>::Entity;
56
60 void init(const std::string& modelParamGroup = "")
61 {
62
63 // try to create it from a DGF or msh file in GridParameterGroup.File
64 if (hasParamInGroup(modelParamGroup, "Grid.File"))
65 {
66 preProcessing_(modelParamGroup);
67 ParentType::makeGridFromFile(getParamFromGroup<std::string>(modelParamGroup, "Grid.File"), modelParamGroup);
68 postProcessing_(modelParamGroup);
69 return;
70 }
71
72 // Then look for the necessary keys to construct from the input file
73 else if (hasParamInGroup(modelParamGroup, "Grid.UpperRight"))
74 {
75 preProcessing_(modelParamGroup);
76 // make the grid
77 const auto cellType = getParamFromGroup<std::string>(modelParamGroup, "Grid.CellType", "Cube");
78 if (cellType == "Cube")
79 ParentType::template makeStructuredGrid<dim, dim>(ParentType::CellType::Cube, modelParamGroup);
80 else if (cellType == "Simplex")
81 ParentType::template makeStructuredGrid<dim, dim>(ParentType::CellType::Simplex, modelParamGroup);
82 else
83 DUNE_THROW(Dune::IOError, "UGGrid only supports 'Cube' or 'Simplex' as cell type. Not '"<< cellType<<"'!");
84 postProcessing_(modelParamGroup);
85 }
86
87 // Didn't find a way to construct the grid
88 else
89 {
90 const auto prefix = modelParamGroup.empty() ? modelParamGroup : modelParamGroup + ".";
91 DUNE_THROW(ParameterException, "Please supply one of the parameters "
92 << prefix + "Grid.UpperRight"
93 << ", or a grid file in " << prefix + "Grid.File");
94
95 }
96 }
97
105 void loadBalance()
106 {
107 if (Dune::MPIHelper::getCommunication().size() > 1)
108 {
109 // if we may have dgf parameters use load balancing of the dgf pointer
110 if(ParentType::enableDgfGridPointer_)
111 {
112 ParentType::dgfGridPtr().loadBalance();
113 // update the grid data
114 ParentType::gridData_ = std::make_shared<typename ParentType::GridData>(ParentType::dgfGridPtr());
115 }
116
117 // if we have gmsh parameters we have to manually load balance the data
118 else if (ParentType::enableGmshDomainMarkers_)
119 {
120 // element and face markers are communicated during load balance
121 auto dh = ParentType::gridData_->createGmshDataHandle();
122 ParentType::gridPtr()->loadBalance(dh.interface());
123 // Right now, UGGrid cannot communicate element data. If this gets implemented, communicate the data here:
124 // ParentType::gridPtr()->communicate(dh.interface(), Dune::InteriorBorder_All_Interface, Dune::ForwardCommunication);
125 }
126 else
127 ParentType::gridPtr()->loadBalance();
128 }
129 }
130
131private:
135 void preProcessing_(const std::string& modelParamGroup)
136 {}
137
141 void postProcessing_(const std::string& modelParamGroup)
142 {
143 // Set refinement type
144 const auto refType = getParamFromGroup<std::string>(modelParamGroup, "Grid.RefinementType", "Local");
145 if (refType == "Local")
146 ParentType::grid().setRefinementType(Dune::UGGrid<dim>::RefinementType::LOCAL);
147 else if (refType == "Copy")
148 ParentType::grid().setRefinementType(Dune::UGGrid<dim>::RefinementType::COPY);
149 else
150 DUNE_THROW(Dune::IOError, "UGGrid only supports 'Local' or 'Copy' as refinement type. Not '"<< refType<<"'!");
151
152 // Set closure type
153 const auto closureType = getParamFromGroup<std::string>(modelParamGroup, "Grid.ClosureType", "Green");
154 if (closureType == "Green")
155 ParentType::grid().setClosureType(Dune::UGGrid<dim>::ClosureType::GREEN);
156 else if (closureType == "None")
157 ParentType::grid().setClosureType(Dune::UGGrid<dim>::ClosureType::NONE);
158 else
159 DUNE_THROW(Dune::IOError, "UGGrid only supports 'Green' or 'None' as closure type. Not '"<< closureType<<"'!");
160
161 // Check if should refine the grid
162 ParentType::maybeRefineGrid(modelParamGroup);
163 // do load balancing
164 loadBalance();
165 }
166};
167
168namespace Grid::Capabilities {
169
170// To the best of our knowledge UGGrid is view thread-safe for sequential runs
171// This specialization maybe be removed after we depend on Dune release 2.9 if is guaranteed by UGGrid itself by then
172template<int dim>
173struct MultithreadingSupported<Dune::UGGrid<dim>>
174{
175 template<class GV>
176 static bool eval(const GV& gv) // default is independent of the grid view
177 { return gv.comm().size() <= 1; }
178};
179
180} // end namespace Grid::Capabilities
181
182#endif // HAVE_DUNE_UGGRID
183
184} // end namespace Dumux
185
186#endif
Grid Grid
Definition: gridmanager_base.hh:57
void loadBalance()
Call loadBalance() function of the grid.
Definition: gridmanager_base.hh:97
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:63
dune-grid capabilities compatibility layer
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
static bool eval(const GV &)
Definition: gridcapabilities.hh:69