version 3.11-dev
Loading...
Searching...
No Matches
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-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_UG_HH
14#define DUMUX_IO_GRID_MANAGER_UG_HH
15
16#if HAVE_DUNE_UGGRID
17#include <dune/grid/uggrid.hh>
18#include <dune/grid/io/file/dgfparser/dgfug.hh>
19#endif
20
21#ifndef DUMUX_IO_GRID_MANAGER_BASE_HH
23#endif
24
26
27namespace Dumux {
28
29#if HAVE_DUNE_UGGRID
30
51template<int dim>
52class GridManager<Dune::UGGrid<dim>>
53: public GridManagerBase<Dune::UGGrid<dim>>
54{
55public:
56 using Grid = typename Dune::UGGrid<dim>;
58 using Element = typename Grid::template Codim<0>::Entity;
59
63 void init(const std::string& modelParamGroup = "")
64 {
65
66 // First, try to create it from a file in GridParameterGroup.File
67 if (hasParamInGroup(modelParamGroup, "Grid.File"))
68 {
69 preProcessing_(modelParamGroup);
70 ParentType::makeGridFromFile(getParamFromGroup<std::string>(modelParamGroup, "Grid.File"), modelParamGroup);
71 postProcessing_(modelParamGroup);
72 return;
73 }
74
75 // Then look for the necessary keys to construct from the input file
76 else if (hasParamInGroup(modelParamGroup, "Grid.UpperRight"))
77 {
78 preProcessing_(modelParamGroup);
79 // make the grid
80 const auto cellType = getParamFromGroup<std::string>(modelParamGroup, "Grid.CellType", "Cube");
81 if (cellType == "Cube")
82 ParentType::template makeStructuredGrid<dim, dim>(ParentType::CellType::Cube, modelParamGroup);
83 else if (cellType == "Simplex")
84 ParentType::template makeStructuredGrid<dim, dim>(ParentType::CellType::Simplex, modelParamGroup);
85 else
86 DUNE_THROW(Dune::IOError, "UGGrid only supports 'Cube' or 'Simplex' as cell type. Not '"<< cellType<<"'!");
87 postProcessing_(modelParamGroup);
88 }
89
90 // Didn't find a way to construct the grid
91 else
92 {
93 const auto prefix = modelParamGroup.empty() ? modelParamGroup : modelParamGroup + ".";
94 DUNE_THROW(ParameterException, "Please supply one of the parameters "
95 << prefix + "Grid.UpperRight"
96 << ", or a grid file in " << prefix + "Grid.File");
97
98 }
99 }
100
109 {
110 if (Dune::MPIHelper::getCommunication().size() > 1)
111 {
112 // if we may have dgf parameters use load balancing of the dgf pointer
114 {
115 ParentType::dgfGridPtr().loadBalance();
116 // update the grid data
117 ParentType::gridData_ = std::make_shared<typename ParentType::GridData>(ParentType::dgfGridPtr());
118 }
119
120 // if we have gmsh parameters we have to manually load balance the data
122 {
123 // element and face markers are communicated during load balance
124 auto dh = ParentType::gridData_->createGmshDataHandle();
125 ParentType::gridPtr()->loadBalance(dh.interface());
126 // Right now, UGGrid cannot communicate element data. If this gets implemented, communicate the data here:
127 // ParentType::gridPtr()->communicate(dh.interface(), Dune::InteriorBorder_All_Interface, Dune::ForwardCommunication);
128 }
129
130 // if we have VTK parameters we have to manually load balance the data
132 {
133 // cell and point data is communicated during load balance
134 auto dh = ParentType::gridData_->createVtkDataHandle();
135 ParentType::gridPtr()->loadBalance(dh.interface());
136 // Right now, UGGrid cannot communicate element data. If this gets implemented, communicate the data here:
137 // gridPtr()->communicate(dh.interface(), Dune::InteriorBorder_All_Interface, Dune::ForwardCommunication);
138 }
139
140 else
141 ParentType::gridPtr()->loadBalance();
142 }
143 }
144
145private:
149 void preProcessing_(const std::string& modelParamGroup)
150 {}
151
155 void postProcessing_(const std::string& modelParamGroup)
156 {
157 // Set refinement type
158 const auto refType = getParamFromGroup<std::string>(modelParamGroup, "Grid.RefinementType", "Local");
159 if (refType == "Local")
160 ParentType::grid().setRefinementType(Dune::UGGrid<dim>::RefinementType::LOCAL);
161 else if (refType == "Copy")
162 ParentType::grid().setRefinementType(Dune::UGGrid<dim>::RefinementType::COPY);
163 else
164 DUNE_THROW(Dune::IOError, "UGGrid only supports 'Local' or 'Copy' as refinement type. Not '"<< refType<<"'!");
165
166 // Set closure type
167 const auto closureType = getParamFromGroup<std::string>(modelParamGroup, "Grid.ClosureType", "Green");
168 if (closureType == "Green")
169 ParentType::grid().setClosureType(Dune::UGGrid<dim>::ClosureType::GREEN);
170 else if (closureType == "None")
171 ParentType::grid().setClosureType(Dune::UGGrid<dim>::ClosureType::NONE);
172 else
173 DUNE_THROW(Dune::IOError, "UGGrid only supports 'Green' or 'None' as closure type. Not '"<< closureType<<"'!");
174
175 // Check if should refine the grid
176 ParentType::maybeRefineGrid(modelParamGroup);
177 // do load balancing
178 loadBalance();
179 }
180};
181
182namespace Grid::Capabilities {
183
184// To the best of our knowledge UGGrid is view thread-safe for sequential runs
185// TODO / Deprecation notice: This specialization may be removed after we depend on Dune release 2.10 if is guaranteed by UGGrid itself by then
186template<int dim>
187struct MultithreadingSupported<Dune::UGGrid<dim>>
188{
189 template<class GV>
190 static bool eval(const GV& gv) // default is independent of the grid view
191 { return gv.comm().size() <= 1; }
192};
193
194} // end namespace Grid::Capabilities
195
196#endif // HAVE_DUNE_UGGRID
197
198} // end namespace Dumux
199
200#endif
typename Grid::template Codim< 0 >::Entity Element
Definition gridmanager_ug.hh:58
GridManagerBase< Grid > ParentType
Definition gridmanager_ug.hh:57
void init(const std::string &modelParamGroup="")
Make the UGGrid.
Definition gridmanager_ug.hh:63
typename Dune::UGGrid< dim > Grid
Definition gridmanager_ug.hh:56
void loadBalance()
Call loadBalance() function of the grid.
Definition gridmanager_ug.hh:108
The grid manager base interface (public) and methods common to most grid manager specializations (pro...
Definition gridmanager_base.hh:56
bool enableVtkData_
Definition gridmanager_base.hh:335
void makeGridFromFile(const std::string &fileName, const std::string &modelParamGroup)
Definition gridmanager_base.hh:205
void makeStructuredGrid(CellType cellType, const std::string &modelParamGroup)
Definition gridmanager_base.hh:286
Dune::GridPtr< Grid > & dgfGridPtr()
Definition gridmanager_base.hh:174
std::shared_ptr< GridData > gridData_
Definition gridmanager_base.hh:340
bool enableGmshDomainMarkers_
Definition gridmanager_base.hh:330
void loadBalance()
Definition gridmanager_base.hh:98
@ Simplex
Definition gridmanager_base.hh:280
@ Cube
Definition gridmanager_base.hh:280
std::shared_ptr< Grid > & gridPtr()
Definition gridmanager_base.hh:163
bool enableDgfGridPointer_
Definition gridmanager_base.hh:325
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
dune-grid capabilities compatibility layer
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 gridcapabilities.hh:17
Definition adapt.hh:17
Definition common/pdesolver.hh:24
static bool eval(const GV &gv)
Definition gridmanager_ug.hh:190