GCC Code Coverage Report


Directory: ../../../builds/dumux-repositories/
File: /builds/dumux-repositories/dumux/dumux/io/grid/griddata.hh
Date: 2024-09-21 20:52:54
Exec Total Coverage
Lines: 34 83 41.0%
Functions: 26 82 31.7%
Branches: 59 454 13.0%

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 Class for grid data attached to dgf or gmsh grid files
11 */
12 #ifndef DUMUX_IO_GRID_DATA_HH
13 #define DUMUX_IO_GRID_DATA_HH
14
15 #include <vector>
16 #include <memory>
17 #include <type_traits>
18
19 #include <dune/common/exceptions.hh>
20 #include <dune/common/parallel/communication.hh>
21 #include <dune/common/parallel/mpihelper.hh>
22 #include <dune/grid/common/gridfactory.hh>
23 #include <dune/grid/io/file/dgfparser/parser.hh>
24 #include <dune/grid/io/file/dgfparser/gridptr.hh>
25 #include <dumux/io/vtk/vtkreader.hh>
26
27 // UGGrid specific includes
28 #if HAVE_DUNE_UGGRID
29 #include <dune/grid/uggrid.hh>
30 #endif
31
32 #include "gmshgriddatahandle.hh"
33 #include "vtkgriddatahandle.hh"
34
35 namespace Dumux {
36
37 namespace Detail {
38
39 template<class Grid>
40 struct isUG : public std::false_type {};
41
42 #if HAVE_DUNE_UGGRID
43 template<int dim>
44 struct isUG<Dune::UGGrid<dim>> : public std::true_type {};
45 #endif
46
47 } // end namespace Details
48
49 /*!
50 * \ingroup InputOutput
51 * \brief Class for grid data attached to dgf or gmsh grid files
52 */
53 template <class Grid>
54 class GridData
55 {
56 using Intersection = typename Grid::LeafIntersection;
57 using Element = typename Grid::template Codim<0>::Entity;
58 using Vertex = typename Grid::template Codim<Grid::dimension>::Entity;
59 using GmshDataHandle = GmshGridDataHandle<Grid, Dune::GridFactory<Grid>, std::vector<int>>;
60 using VtkDataHandle = VtkGridDataHandle<Grid, Dune::GridFactory<Grid>, std::vector<double>>;
61
62 enum DataSourceType { dgf, gmsh, vtk };
63
64 public:
65 //! constructor for gmsh grid data
66 65 GridData(std::shared_ptr<Grid> grid, std::shared_ptr<Dune::GridFactory<Grid>> factory,
67 std::vector<int>&& elementMarkers, std::vector<int>&& boundaryMarkers, std::vector<int>&& faceMarkers = std::vector<int>{})
68 : gridPtr_(grid)
69 , gridFactory_(factory)
70 , elementMarkers_(elementMarkers)
71 , boundaryMarkers_(boundaryMarkers)
72 , faceMarkers_(faceMarkers)
73
6/22
✓ Branch 3 taken 41 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 41 times.
✗ Branch 7 not taken.
✓ Branch 9 taken 41 times.
✗ Branch 10 not taken.
✓ Branch 12 taken 41 times.
✗ Branch 13 not taken.
✓ Branch 15 taken 41 times.
✗ Branch 16 not taken.
✓ Branch 18 taken 41 times.
✗ Branch 19 not taken.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✗ Branch 31 not taken.
65 , dataSourceType_(DataSourceType::gmsh)
74 65 {}
75
76 //! constructor for dgf grid data
77 81 GridData(Dune::GridPtr<Grid> grid)
78 : dgfGrid_(grid)
79
8/26
✓ Branch 1 taken 75 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 75 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 75 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 75 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 75 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 75 times.
✗ Branch 17 not taken.
✓ Branch 19 taken 75 times.
✗ Branch 20 not taken.
✓ Branch 22 taken 75 times.
✗ Branch 23 not taken.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✗ Branch 31 not taken.
✗ Branch 32 not taken.
✗ Branch 33 not taken.
✗ Branch 34 not taken.
✗ Branch 35 not taken.
648 , dataSourceType_(DataSourceType::dgf)
80 81 {}
81
82 //! constructor for VTK grid data
83 GridData(std::shared_ptr<Grid> grid, std::shared_ptr<Dune::GridFactory<Grid>> factory,
84 VTKReader::Data&& cellData, VTKReader::Data&& pointData)
85 : gridPtr_(grid)
86 , gridFactory_(factory)
87 , cellData_(cellData)
88 , pointData_(pointData)
89 , dataSourceType_(DataSourceType::vtk)
90 {}
91
92
93 /*!
94 * \name DGF interface functions
95 */
96 // \{
97
98 /*!
99 * \brief Call the parameters function of the DGF grid pointer if available for vertex data
100 * \note You can only pass vertices that exist on level 0!
101 */
102 15686 const std::vector<double>& parameters(const Vertex& vertex) const
103 {
104
1/2
✓ Branch 0 taken 15686 times.
✗ Branch 1 not taken.
15686 if (dataSourceType_ == DataSourceType::dgf)
105 {
106
2/4
✗ Branch 0 not taken.
✓ Branch 1 taken 15686 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 15686 times.
31372 if (vertex.level() != 0)
107 DUNE_THROW(Dune::IOError, "You can only obtain parameters for level 0 vertices!");
108
109 15686 return dgfGrid_.parameters(vertex);
110 }
111 else
112 DUNE_THROW(Dune::InvalidStateException, "The parameters method is only available if the grid was constructed with a DGF file.");
113 }
114
115 /*!
116 * \brief Call the parameters function of the DGF grid pointer if available for element data
117 */
118 134784 const std::vector<double>& parameters(const Element& element) const
119 {
120
1/2
✓ Branch 0 taken 134784 times.
✗ Branch 1 not taken.
134784 if (dataSourceType_ == DataSourceType::dgf)
121 {
122
4/4
✓ Branch 0 taken 29519 times.
✓ Branch 1 taken 53135 times.
✓ Branch 2 taken 81649 times.
✓ Branch 3 taken 23125 times.
187428 if (element.hasFather())
123 {
124 89539 auto level0Element = element;
125
5/5
✓ Branch 0 taken 29519 times.
✓ Branch 1 taken 90887 times.
✓ Branch 2 taken 29519 times.
✓ Branch 3 taken 60877 times.
✓ Branch 4 taken 30010 times.
179444 while(level0Element.hasFather())
126
1/2
✓ Branch 1 taken 31358 times.
✗ Branch 2 not taken.
60877 level0Element = level0Element.father();
127
128
1/2
✓ Branch 1 taken 30010 times.
✗ Branch 2 not taken.
59529 return dgfGrid_.parameters(level0Element);
129 }
130 else
131 {
132 75255 return dgfGrid_.parameters(element);
133 }
134 }
135 else
136 DUNE_THROW(Dune::InvalidStateException, "The parameters method is only available if the grid was constructed with a DGF file.");
137 }
138
139 /*!
140 * \brief Call the parameters function of the DGF grid pointer if available
141 */
142 template <class GridImp, class IntersectionImp>
143 const Dune::DGFBoundaryParameter::type& parameters(const Dune::Intersection<GridImp, IntersectionImp>& intersection) const
144 {
145 if (dataSourceType_ == DataSourceType::dgf)
146 return dgfGrid_.parameters(intersection);
147 else
148 DUNE_THROW(Dune::InvalidStateException, "The parameters method is only available if the grid was constructed with a DGF file.");
149 }
150
151 // \}
152
153 /*!
154 * \name Gmsh interface functions
155 */
156 // \{
157
158 /*!
159 * \brief Return the boundary domain marker (Gmsh physical entity number) of an intersection
160 Only available when using Gmsh with GridParameterGroup.DomainMarkers = 1.
161 * \param boundarySegmentIndex The boundary segment index of the intersection (intersection.boundarySegmentIndex()
162 */
163 759198 int getBoundaryDomainMarker(int boundarySegmentIndex) const
164 {
165
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 759198 times.
759198 if (dataSourceType_ != DataSourceType::gmsh)
166 DUNE_THROW(Dune::InvalidStateException, "Domain markers are only available for gmsh grids.");
167
2/4
✗ Branch 0 not taken.
✓ Branch 1 taken 759198 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 759198 times.
1518396 if (boundarySegmentIndex >= boundaryMarkers_.size())
168 DUNE_THROW(Dune::RangeError, "Boundary segment index "<< boundarySegmentIndex << " bigger than number of boundary segments in grid.\n"
169 "Make sure to call this function only for boundaries that were defined as physical entities in gmsh.");
170 1518396 return boundaryMarkers_[boundarySegmentIndex];
171 }
172
173 /*!
174 * \brief Return the boundary domain marker (Gmsh physical entity number) of an intersection
175 Only available when using Gmsh with GridParameterGroup.DomainMarkers = 1.
176 * \param intersection The intersection to be evaluated
177 */
178 int getBoundaryDomainMarker(const Intersection& intersection) const
179
6/12
✓ Branch 1 taken 3781 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 3781 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 17921 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 17921 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 84 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 84 times.
✗ Branch 17 not taken.
21786 { return getBoundaryDomainMarker(intersection.boundarySegmentIndex()); }
180
181 /*!
182 * \brief Returns true if an intersection was inserted during grid creation
183 */
184 bool wasInserted(const Intersection& intersection) const
185
2/4
✓ Branch 1 taken 480 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 480 times.
✗ Branch 5 not taken.
960 { return gridFactory_->wasInserted(intersection); }
186
187 /*!
188 * \brief Return the element domain marker (Gmsh physical entity number) of an element.
189 Only available when using Gmsh with GridParameterGroup.DomainMarkers = 1.
190 * \param element The element to be evaluated
191 */
192 143962 int getElementDomainMarker(const Element& element) const
193 {
194
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 142912 times.
143962 if (dataSourceType_ != DataSourceType::gmsh)
195 DUNE_THROW(Dune::InvalidStateException, "Domain markers are only available for gmsh grids.");
196
197 // parameters are only given for level 0 elements
198 207363 auto level0element = element;
199
5/5
✓ Branch 0 taken 63791 times.
✓ Branch 1 taken 189024 times.
✓ Branch 2 taken 63791 times.
✓ Branch 3 taken 126114 times.
✓ Branch 4 taken 62910 times.
398217 while (level0element.hasFather())
200
1/2
✓ Branch 1 taken 46112 times.
✗ Branch 2 not taken.
109903 level0element = level0element.father();
201
202 // in the parallel case the data is load balanced and then accessed with indices of the index set
203 // for UGGrid element data is read on all processes since UGGrid can't communicate element data (yet)
204
6/8
✗ Branch 0 not taken.
✓ Branch 1 taken 62910 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 62910 times.
✓ Branch 4 taken 32092 times.
✓ Branch 5 taken 30818 times.
✓ Branch 6 taken 32092 times.
✓ Branch 7 taken 30818 times.
368417 if (gridPtr_->comm().size() > 1 && !Detail::isUG<Grid>::value)
205
4/8
✓ Branch 1 taken 32092 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 32092 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 32092 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 32092 times.
✗ Branch 11 not taken.
96276 return elementMarkers_[gridPtr_->levelGridView(0).indexSet().index(level0element)];
206 else
207
2/4
✓ Branch 1 taken 30818 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 30818 times.
✗ Branch 5 not taken.
223740 return elementMarkers_[gridFactory_->insertionIndex(level0element)];
208 }
209
210 /*!
211 * \brief Create a data handle for communication of the data in parallel simulations
212 * \note this data handle is the default
213 */
214 template<bool ug = Detail::isUG<Grid>::value, typename std::enable_if_t<!ug, int> = 0>
215 GmshDataHandle createGmshDataHandle()
216 {
217 12 return GmshDataHandle(*gridPtr_, *gridFactory_, elementMarkers_, boundaryMarkers_, faceMarkers_);
218 }
219
220 /*!
221 * \brief Create a data handle for communication of the data in parallel simulations
222 * \note this data handle is the specialized for UGGrid since UGGrid can't communicate element data (yet)
223 */
224 template<bool ug = Detail::isUG<Grid>::value, typename std::enable_if_t<ug, int> = 0>
225 GmshDataHandle createGmshDataHandle()
226 {
227 12 return GmshDataHandle(*gridPtr_, *gridFactory_, elementMarkers_, boundaryMarkers_);
228 }
229
230 // \}
231
232 /*!
233 * \name VTK interface functions
234 */
235 // \{
236
237 /*!
238 * \brief Get an element parameter
239 * \param element the element
240 * \param fieldName the name of the field to read from the vtk data
241 */
242 double getParameter(const Element& element, const std::string& fieldName) const
243 { return getArrayParameter<double, 1>(element, fieldName)[0]; }
244
245 /*!
246 * \brief Get an element parameter that is an array
247 * \param element the element
248 * \param fieldName the name of the field to read from the vtk data
249 * \tparam T the parameter scalar type
250 * \tparam size the number of parameters in the array
251 */
252 template<class T, std::size_t size>
253 std::array<T, size> getArrayParameter(const Element& element, const std::string& fieldName) const
254 {
255 if (dataSourceType_ != DataSourceType::vtk)
256 DUNE_THROW(Dune::InvalidStateException, "This access function is only available for data from VTK files.");
257
258 if (cellData_.count(fieldName) == 0)
259 DUNE_THROW(Dune::IOError, "No field with the name " << fieldName << " found in cell data.");
260
261 // parameters are only given for level 0 elements
262 auto level0element = element;
263 while (level0element.hasFather())
264 level0element = level0element.father();
265
266 // different index depending on whether we have communicated user data (parallel setting) or not (sequential setting)
267 const auto index = [&]{
268 if (gridPtr_->comm().size() > 1)
269 return gridPtr_->levelGridView(0).indexSet().index(level0element);
270 else
271 return gridFactory_->insertionIndex(level0element);
272 }();
273
274 std::array<T, size> param;
275 const auto& data = cellData_.at(fieldName);
276
277 if (const std::size_t nc = data.size()/gridPtr_->levelGridView(0).size(0); size != nc)
278 DUNE_THROW(Dune::IOError, "Array size " << size << " does not match number of components of field " << nc);
279
280 for (std::size_t i = 0; i < size; ++i)
281 param[i] = static_cast<T>(data[i + size*index]);
282 return param;
283 }
284
285 /*!
286 * \brief Call the parameters function of the DGF grid pointer if available for vertex data
287 * \param vertex the vertex
288 * \param fieldName the name of the field to read from the vtk data
289 * \note You can only pass vertices that exist on level 0!
290 */
291 double getParameter(const Vertex& vertex, const std::string& fieldName) const
292 { return getArrayParameter<double, 1>(vertex, fieldName)[0]; }
293
294 /*!
295 * \brief Call the parameters function of the DGF grid pointer if available for vertex data
296 * \param vertex the vertex
297 * \param fieldName the name of the field to read from the vtk data
298 * \tparam T the parameter scalar type
299 * \tparam size the number of parameters in the array
300 * \note You can only pass vertices that exist on level 0 (otherwise: undefined behaviour).
301 */
302 template<class T, std::size_t size>
303 std::array<T, size> getArrayParameter(const Vertex& vertex, const std::string& fieldName) const
304 {
305 if (dataSourceType_ != DataSourceType::vtk)
306 DUNE_THROW(Dune::InvalidStateException, "This access function is only available for data from VTK files.");
307
308 if (vertex.level() != 0)
309 DUNE_THROW(Dune::IOError, "You can only obtain parameters for level 0 vertices!");
310
311 if (pointData_.count(fieldName) == 0)
312 DUNE_THROW(Dune::IOError, "No field with the name " << fieldName << " found in point data");
313
314 // different index depending on whether we have communicated user data (parallel setting) or not (sequential setting)
315 const auto index = [&]{
316 if (gridPtr_->comm().size() > 1)
317 return gridPtr_->levelGridView(0).indexSet().index(vertex);
318 else
319 return gridFactory_->insertionIndex(vertex);
320 }();
321
322 std::array<T, size> param;
323 const auto& data = pointData_.at(fieldName);
324
325 if (const std::size_t nc = data.size()/gridPtr_->levelGridView(0).size(Grid::dimension); size != nc)
326 DUNE_THROW(Dune::IOError, "Array size " << size << " does not match number of components of field " << nc);
327
328 for (std::size_t i = 0; i < size; ++i)
329 param[i] = static_cast<T>(data[i + size*index]);
330 return param;
331 }
332
333 /*!
334 * \brief Create a data handle for communication of the data in parallel simulations
335 */
336 VtkDataHandle createVtkDataHandle()
337 {
338 return VtkDataHandle(*gridPtr_, *gridFactory_, cellData_, pointData_);
339 }
340
341 // \}
342
343 private:
344 // grid and grid factor for gmsh grid data / vtk grid data
345 std::shared_ptr<Grid> gridPtr_;
346 std::shared_ptr<Dune::GridFactory<Grid>> gridFactory_;
347
348 /*!
349 * \brief Element and boundary domain markers obtained from Gmsh physical entities
350 * They map from element indices / boundary ids to the physical entity number
351 */
352 std::vector<int> elementMarkers_;
353 std::vector<int> boundaryMarkers_;
354 std::vector<int> faceMarkers_;
355
356 /*!
357 * \brief Cell and vertex data obtained from VTK files
358 */
359 VTKReader::Data cellData_, pointData_;
360
361 // dgf grid data
362 Dune::GridPtr<Grid> dgfGrid_;
363
364 // specify which type of data we have
365 // TODO unfortunately all grid readers provide different data types, should be streamlined (changes in Dune)
366 DataSourceType dataSourceType_;
367 };
368
369 } // namespace Dumux
370
371 #endif
372