version 3.10-dev
griddata.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_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>
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
35namespace Dumux {
36
37namespace Detail {
38
39template<class Grid>
40struct isUG : public std::false_type {};
41
42#if HAVE_DUNE_UGGRID
43template<int dim>
44struct isUG<Dune::UGGrid<dim>> : public std::true_type {};
45#endif
46
47} // end namespace Details
48
53template <class Grid>
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;
61
62 enum DataSourceType { dgf, gmsh, vtk };
63
64public:
66 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 , dataSourceType_(DataSourceType::gmsh)
74 {}
75
77 GridData(Dune::GridPtr<Grid> grid)
78 : dgfGrid_(grid)
79 , dataSourceType_(DataSourceType::dgf)
80 {}
81
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
96 // \{
97
102 const std::vector<double>& parameters(const Vertex& vertex) const
103 {
104 if (dataSourceType_ == DataSourceType::dgf)
105 {
106 if (vertex.level() != 0)
107 DUNE_THROW(Dune::IOError, "You can only obtain parameters for level 0 vertices!");
108
109 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
118 const std::vector<double>& parameters(const Element& element) const
119 {
120 if (dataSourceType_ == DataSourceType::dgf)
121 {
122 if (element.hasFather())
123 {
124 auto level0Element = element;
125 while(level0Element.hasFather())
126 level0Element = level0Element.father();
127
128 return dgfGrid_.parameters(level0Element);
129 }
130 else
131 {
132 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
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
156 // \{
157
163 int getBoundaryDomainMarker(int boundarySegmentIndex) const
164 {
165 if (dataSourceType_ != DataSourceType::gmsh)
166 DUNE_THROW(Dune::InvalidStateException, "Domain markers are only available for gmsh grids.");
167 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 return boundaryMarkers_[boundarySegmentIndex];
171 }
172
178 int getBoundaryDomainMarker(const Intersection& intersection) const
179 { return getBoundaryDomainMarker(intersection.boundarySegmentIndex()); }
180
184 bool wasInserted(const Intersection& intersection) const
185 { return gridFactory_->wasInserted(intersection); }
186
192 int getElementDomainMarker(const Element& element) const
193 {
194 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 auto level0element = element;
199 while (level0element.hasFather())
200 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 if (gridPtr_->comm().size() > 1 && !Detail::isUG<Grid>::value)
205 return elementMarkers_[gridPtr_->levelGridView(0).indexSet().index(level0element)];
206 else
207 return elementMarkers_[gridFactory_->insertionIndex(level0element)];
208 }
209
214 template<bool ug = Detail::isUG<Grid>::value, typename std::enable_if_t<!ug, int> = 0>
216 {
217 return GmshDataHandle(*gridPtr_, *gridFactory_, elementMarkers_, boundaryMarkers_, faceMarkers_);
218 }
219
224 template<bool ug = Detail::isUG<Grid>::value, typename std::enable_if_t<ug, int> = 0>
226 {
227 return GmshDataHandle(*gridPtr_, *gridFactory_, elementMarkers_, boundaryMarkers_);
228 }
229
230 // \}
231
235 // \{
236
242 double getParameter(const Element& element, const std::string& fieldName) const
243 { return getArrayParameter<double, 1>(element, fieldName)[0]; }
244
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
291 double getParameter(const Vertex& vertex, const std::string& fieldName) const
292 { return getArrayParameter<double, 1>(vertex, fieldName)[0]; }
293
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
337 {
338 return VtkDataHandle(*gridPtr_, *gridFactory_, cellData_, pointData_);
339 }
340
341 // \}
342
343private:
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
352 std::vector<int> elementMarkers_;
353 std::vector<int> boundaryMarkers_;
354 std::vector<int> faceMarkers_;
355
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
Class for grid data attached to dgf or gmsh grid files.
Definition: griddata.hh:55
GridData(std::shared_ptr< Grid > grid, std::shared_ptr< Dune::GridFactory< Grid > > factory, std::vector< int > &&elementMarkers, std::vector< int > &&boundaryMarkers, std::vector< int > &&faceMarkers=std::vector< int >{})
constructor for gmsh grid data
Definition: griddata.hh:66
double getParameter(const Vertex &vertex, const std::string &fieldName) const
Call the parameters function of the DGF grid pointer if available for vertex data.
Definition: griddata.hh:291
int getBoundaryDomainMarker(const Intersection &intersection) const
Return the boundary domain marker (Gmsh physical entity number) of an intersection Only available whe...
Definition: griddata.hh:178
double getParameter(const Element &element, const std::string &fieldName) const
Get an element parameter.
Definition: griddata.hh:242
GridData(std::shared_ptr< Grid > grid, std::shared_ptr< Dune::GridFactory< Grid > > factory, VTKReader::Data &&cellData, VTKReader::Data &&pointData)
constructor for VTK grid data
Definition: griddata.hh:83
const Dune::DGFBoundaryParameter::type & parameters(const Dune::Intersection< GridImp, IntersectionImp > &intersection) const
Call the parameters function of the DGF grid pointer if available.
Definition: griddata.hh:143
GridData(Dune::GridPtr< Grid > grid)
constructor for dgf grid data
Definition: griddata.hh:77
VtkDataHandle createVtkDataHandle()
Create a data handle for communication of the data in parallel simulations.
Definition: griddata.hh:336
const std::vector< double > & parameters(const Vertex &vertex) const
Call the parameters function of the DGF grid pointer if available for vertex data.
Definition: griddata.hh:102
std::array< T, size > getArrayParameter(const Vertex &vertex, const std::string &fieldName) const
Call the parameters function of the DGF grid pointer if available for vertex data.
Definition: griddata.hh:303
int getBoundaryDomainMarker(int boundarySegmentIndex) const
Return the boundary domain marker (Gmsh physical entity number) of an intersection Only available whe...
Definition: griddata.hh:163
GmshDataHandle createGmshDataHandle()
Create a data handle for communication of the data in parallel simulations.
Definition: griddata.hh:215
std::array< T, size > getArrayParameter(const Element &element, const std::string &fieldName) const
Get an element parameter that is an array.
Definition: griddata.hh:253
int getElementDomainMarker(const Element &element) const
Return the element domain marker (Gmsh physical entity number) of an element. Only available when usi...
Definition: griddata.hh:192
const std::vector< double > & parameters(const Element &element) const
Call the parameters function of the DGF grid pointer if available for element data.
Definition: griddata.hh:118
bool wasInserted(const Intersection &intersection) const
Returns true if an intersection was inserted during grid creation.
Definition: griddata.hh:184
std::unordered_map< std::string, std::vector< double > > Data
the cell / point data type for point data read from a grid file
Definition: vtkreader.hh:48
A data handle for commucating grid data for gmsh grids.
Definition: adapt.hh:17
Definition: common/pdesolver.hh:24
Definition: griddata.hh:40
A data handle for commucating grid data for gmsh grids.
Definition: gmshgriddatahandle.hh:37
A data handle for communicating grid data for VTK grids.
Definition: vtkgriddatahandle.hh:34
A data handle for communicating grid data for VTK grids.
A vtk file reader using tinyxml2 as xml backend.