Line | Branch | Exec | Source |
---|---|---|---|
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 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup InputOutput | ||
10 | * \brief Grid manager specialization for UGGrid | ||
11 | */ | ||
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 | ||
21 | #include <dumux/io/grid/gridmanager_base.hh> | ||
22 | #endif | ||
23 | |||
24 | #include <dumux/common/gridcapabilities.hh> | ||
25 | |||
26 | namespace Dumux { | ||
27 | |||
28 | #if HAVE_DUNE_UGGRID | ||
29 | |||
30 | /*! | ||
31 | * \ingroup InputOutput | ||
32 | * \brief Provides a grid manager for UGGrids | ||
33 | * from information in the input file | ||
34 | * | ||
35 | * All keys are expected to be in group GridParameterGroup. | ||
36 | |||
37 | * The following keys are recognized: | ||
38 | * - File : A DGF or gmsh file to load from, type detection by file extension | ||
39 | * - LowerLeft : lowerLeft corner of a structured grid | ||
40 | * - UpperRight : upperright corner of a structured grid | ||
41 | * - Cells : number of elements in a structured grid | ||
42 | * - CellType : "Cube" or "Simplex" to be used for structured grids | ||
43 | * - Refinement : the number of global refines to perform | ||
44 | * - Verbosity : whether the grid construction should output to standard out | ||
45 | * - BoundarySegments : whether to insert boundary segments into the grid | ||
46 | * | ||
47 | */ | ||
48 | template<int dim> | ||
49 |
0/4✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
156 | class GridManager<Dune::UGGrid<dim>> |
50 | : public GridManagerBase<Dune::UGGrid<dim>> | ||
51 | { | ||
52 | public: | ||
53 | using Grid = typename Dune::UGGrid<dim>; | ||
54 | using ParentType = GridManagerBase<Grid>; | ||
55 | using Element = typename Grid::template Codim<0>::Entity; | ||
56 | |||
57 | /*! | ||
58 | * \brief Make the UGGrid. | ||
59 | */ | ||
60 | 78 | void init(const std::string& modelParamGroup = "") | |
61 | { | ||
62 | |||
63 | // try to create it from a DGF or msh file in GridParameterGroup.File | ||
64 |
8/14✓ Branch 1 taken 78 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 78 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 78 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 78 times.
✓ Branch 11 taken 34 times.
✓ Branch 12 taken 44 times.
✓ Branch 13 taken 34 times.
✓ Branch 14 taken 44 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
|
234 | if (hasParamInGroup(modelParamGroup, "Grid.File")) |
65 | { | ||
66 | 34 | preProcessing_(modelParamGroup); | |
67 |
3/6✓ Branch 2 taken 34 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 32 times.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
34 | ParentType::makeGridFromFile(getParamFromGroup<std::string>(modelParamGroup, "Grid.File"), modelParamGroup); |
68 | 34 | postProcessing_(modelParamGroup); | |
69 | 34 | return; | |
70 | } | ||
71 | |||
72 | // Then look for the necessary keys to construct from the input file | ||
73 |
6/14✓ Branch 1 taken 44 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 44 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 44 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 44 times.
✓ Branch 11 taken 44 times.
✗ Branch 12 not taken.
✓ Branch 13 taken 44 times.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
|
132 | else if (hasParamInGroup(modelParamGroup, "Grid.UpperRight")) |
74 | { | ||
75 | 44 | preProcessing_(modelParamGroup); | |
76 | // make the grid | ||
77 |
1/4✗ Branch 1 not taken.
✓ Branch 2 taken 44 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
88 | const auto cellType = getParamFromGroup<std::string>(modelParamGroup, "Grid.CellType", "Cube"); |
78 |
2/2✓ Branch 1 taken 42 times.
✓ Branch 2 taken 2 times.
|
44 | if (cellType == "Cube") |
79 |
1/2✓ Branch 1 taken 42 times.
✗ Branch 2 not taken.
|
42 | ParentType::template makeStructuredGrid<dim, dim>(ParentType::CellType::Cube, modelParamGroup); |
80 |
1/2✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
|
2 | else if (cellType == "Simplex") |
81 |
1/2✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
|
2 | 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 |
1/2✓ Branch 1 taken 44 times.
✗ Branch 2 not taken.
|
44 | 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 | |||
98 | /*! | ||
99 | * \brief Call loadBalance() function of the grid. | ||
100 | * \note This overload is necessary because UGGrid cannot load balance element data yet. | ||
101 | * For dgf grids using element data in parallel will (hopefully) through an error | ||
102 | * For gmsh grids the parameters are read on every process so we use too much memory but | ||
103 | * the parameters are available via the insertionIndex of the level 0 element. | ||
104 | */ | ||
105 | 78 | void loadBalance() | |
106 | { | ||
107 |
4/4✓ Branch 1 taken 14 times.
✓ Branch 2 taken 64 times.
✓ Branch 3 taken 14 times.
✓ Branch 4 taken 64 times.
|
78 | if (Dune::MPIHelper::getCommunication().size() > 1) |
108 | { | ||
109 | // if we may have dgf parameters use load balancing of the dgf pointer | ||
110 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 8 times.
|
14 | if(ParentType::enableDgfGridPointer_) |
111 | { | ||
112 | 6 | ParentType::dgfGridPtr().loadBalance(); | |
113 | // update the grid data | ||
114 |
2/4✓ Branch 2 taken 6 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 6 times.
|
12 | 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 |
2/2✓ Branch 0 taken 4 times.
✓ Branch 1 taken 4 times.
|
8 | else if (ParentType::enableGmshDomainMarkers_) |
119 | { | ||
120 | // element and face markers are communicated during load balance | ||
121 | 12 | auto dh = ParentType::gridData_->createGmshDataHandle(); | |
122 |
3/6✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 4 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 4 times.
✗ Branch 8 not taken.
|
4 | 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 | |||
127 | // if we have VTK parameters we have to manually load balance the data | ||
128 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 4 times.
|
4 | else if (ParentType::enableVtkData_) |
129 | { | ||
130 | // cell and point data is communicated during load balance | ||
131 | ✗ | auto dh = ParentType::gridData_->createVtkDataHandle(); | |
132 | ✗ | ParentType::gridPtr()->loadBalance(dh.interface()); | |
133 | // Right now, UGGrid cannot communicate element data. If this gets implemented, communicate the data here: | ||
134 | // gridPtr()->communicate(dh.interface(), Dune::InteriorBorder_All_Interface, Dune::ForwardCommunication); | ||
135 | } | ||
136 | |||
137 | else | ||
138 | 4 | ParentType::gridPtr()->loadBalance(); | |
139 | } | ||
140 | 78 | } | |
141 | |||
142 | private: | ||
143 | /*! | ||
144 | * \brief Do some operations before making the grid | ||
145 | */ | ||
146 | ✗ | void preProcessing_(const std::string& modelParamGroup) | |
147 | ✗ | {} | |
148 | |||
149 | /*! | ||
150 | * \brief Do some operatrion after making the grid, like global refinement | ||
151 | */ | ||
152 | 78 | void postProcessing_(const std::string& modelParamGroup) | |
153 | { | ||
154 | // Set refinement type | ||
155 |
0/2✗ Branch 1 not taken.
✗ Branch 2 not taken.
|
78 | const auto refType = getParamFromGroup<std::string>(modelParamGroup, "Grid.RefinementType", "Local"); |
156 |
1/2✓ Branch 1 taken 78 times.
✗ Branch 2 not taken.
|
78 | if (refType == "Local") |
157 |
1/2✓ Branch 1 taken 78 times.
✗ Branch 2 not taken.
|
78 | ParentType::grid().setRefinementType(Dune::UGGrid<dim>::RefinementType::LOCAL); |
158 | ✗ | else if (refType == "Copy") | |
159 | ✗ | ParentType::grid().setRefinementType(Dune::UGGrid<dim>::RefinementType::COPY); | |
160 | else | ||
161 | ✗ | DUNE_THROW(Dune::IOError, "UGGrid only supports 'Local' or 'Copy' as refinement type. Not '"<< refType<<"'!"); | |
162 | |||
163 | // Set closure type | ||
164 |
3/8✓ Branch 1 taken 78 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 78 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 78 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
|
156 | const auto closureType = getParamFromGroup<std::string>(modelParamGroup, "Grid.ClosureType", "Green"); |
165 |
2/2✓ Branch 1 taken 77 times.
✓ Branch 2 taken 1 times.
|
78 | if (closureType == "Green") |
166 |
1/2✓ Branch 1 taken 77 times.
✗ Branch 2 not taken.
|
77 | ParentType::grid().setClosureType(Dune::UGGrid<dim>::ClosureType::GREEN); |
167 |
1/2✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
|
1 | else if (closureType == "None") |
168 |
1/2✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
|
1 | ParentType::grid().setClosureType(Dune::UGGrid<dim>::ClosureType::NONE); |
169 | else | ||
170 | ✗ | DUNE_THROW(Dune::IOError, "UGGrid only supports 'Green' or 'None' as closure type. Not '"<< closureType<<"'!"); | |
171 | |||
172 | // Check if should refine the grid | ||
173 |
1/2✓ Branch 1 taken 78 times.
✗ Branch 2 not taken.
|
78 | ParentType::maybeRefineGrid(modelParamGroup); |
174 | // do load balancing | ||
175 |
1/2✓ Branch 1 taken 78 times.
✗ Branch 2 not taken.
|
78 | loadBalance(); |
176 | 78 | } | |
177 | }; | ||
178 | |||
179 | namespace Grid::Capabilities { | ||
180 | |||
181 | // To the best of our knowledge UGGrid is view thread-safe for sequential runs | ||
182 | // TODO / Deprecation notice: This specialization may be removed after we depend on Dune release 2.10 if is guaranteed by UGGrid itself by then | ||
183 | template<int dim> | ||
184 | struct MultithreadingSupported<Dune::UGGrid<dim>> | ||
185 | { | ||
186 | template<class GV> | ||
187 | static bool eval(const GV& gv) // default is independent of the grid view | ||
188 |
3/4✗ Branch 0 not taken.
✓ Branch 1 taken 58 times.
✓ Branch 2 taken 52 times.
✓ Branch 3 taken 6 times.
|
58 | { return gv.comm().size() <= 1; } |
189 | }; | ||
190 | |||
191 | } // end namespace Grid::Capabilities | ||
192 | |||
193 | #endif // HAVE_DUNE_UGGRID | ||
194 | |||
195 | } // end namespace Dumux | ||
196 | |||
197 | #endif | ||
198 |