3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_IO_GRID_DATA_HH
25#define DUMUX_IO_GRID_DATA_HH
26
27#include <vector>
28#include <memory>
29#include <type_traits>
30
31#include <dune/common/exceptions.hh>
32
33#include <dune/common/version.hh>
34#if DUNE_VERSION_LT(DUNE_COMMON,2,7)
35#include <dune/common/parallel/collectivecommunication.hh>
36#else
37#include <dune/common/parallel/communication.hh>
38#endif
39
40#include <dune/common/parallel/mpihelper.hh>
41#include <dune/grid/common/gridfactory.hh>
42#include <dune/grid/io/file/dgfparser/parser.hh>
43#include <dune/grid/io/file/dgfparser/gridptr.hh>
45
46// UGGrid specific includes
47#if HAVE_UG
48#include <dune/grid/uggrid.hh>
49#endif
50
51#include "gmshgriddatahandle.hh"
52
53namespace Dumux {
54
55namespace Detail {
56
57template<class Grid>
58struct isUG : public std::false_type {};
59
60#if HAVE_UG
61template<int dim>
62struct isUG<Dune::UGGrid<dim>> : public std::true_type {};
63#endif
64
65} // end namespace Details
66
71template <class Grid>
73{
74 using Intersection = typename Grid::LeafIntersection;
75 using Element = typename Grid::template Codim<0>::Entity;
76 using Vertex = typename Grid::template Codim<Grid::dimension>::Entity;
78
79 enum DataSourceType { dgf, gmsh, vtk };
80
81public:
83 GridData(std::shared_ptr<Grid> grid, std::shared_ptr<Dune::GridFactory<Grid>> factory,
84 std::vector<int>&& elementMarkers, std::vector<int>&& boundaryMarkers, std::vector<int>&& faceMarkers = std::vector<int>{})
85 : gridPtr_(grid)
86 , gridFactory_(factory)
87 , elementMarkers_(elementMarkers)
88 , boundaryMarkers_(boundaryMarkers)
89 , faceMarkers_(faceMarkers)
90 , dataSourceType_(DataSourceType::gmsh)
91 {}
92
94 GridData(Dune::GridPtr<Grid> grid)
95 : dgfGrid_(grid)
96 , dataSourceType_(DataSourceType::dgf)
97 {}
98
100 GridData(std::shared_ptr<Grid> grid, std::shared_ptr<Dune::GridFactory<Grid>> factory,
101 VTKReader::Data&& cellData, VTKReader::Data&& pointData)
102 : gridPtr_(grid)
103 , gridFactory_(factory)
104 , cellData_(cellData)
105 , pointData_(pointData)
106 , dataSourceType_(DataSourceType::vtk)
107 {}
108
109
113 // \{
114
119 const std::vector<double>& parameters(const Vertex& vertex) const
120 {
121 if (dataSourceType_ == DataSourceType::dgf)
122 {
123 if (vertex.level() != 0)
124 DUNE_THROW(Dune::IOError, "You can only obtain parameters for level 0 vertices!");
125
126 return dgfGrid_.parameters(vertex);
127 }
128 else
129 DUNE_THROW(Dune::InvalidStateException, "The parameters method is only available if the grid was constructed with a DGF file.");
130 }
131
135 const std::vector<double>& parameters(const Element& element) const
136 {
137 if (dataSourceType_ == DataSourceType::dgf)
138 {
139 if (element.hasFather())
140 {
141 auto level0Element = element;
142 while(level0Element.hasFather())
143 level0Element = level0Element.father();
144
145 return dgfGrid_.parameters(level0Element);
146 }
147 else
148 {
149 return dgfGrid_.parameters(element);
150 }
151 }
152 else
153 DUNE_THROW(Dune::InvalidStateException, "The parameters method is only available if the grid was constructed with a DGF file.");
154 }
155
159 template <class GridImp, class IntersectionImp>
160 const Dune::DGFBoundaryParameter::type& parameters(const Dune::Intersection<GridImp, IntersectionImp>& intersection) const
161 {
162 if (dataSourceType_ == DataSourceType::dgf)
163 return dgfGrid_.parameters(intersection);
164 else
165 DUNE_THROW(Dune::InvalidStateException, "The parameters method is only available if the grid was constructed with a DGF file.");
166 }
167
168 // \}
169
173 // \{
174
180 int getBoundaryDomainMarker(int boundarySegmentIndex) const
181 {
182 if (dataSourceType_ != DataSourceType::gmsh)
183 DUNE_THROW(Dune::InvalidStateException, "Domain markers are only available for gmsh grids.");
184 if (boundarySegmentIndex >= boundaryMarkers_.size())
185 DUNE_THROW(Dune::RangeError, "Boundary segment index "<< boundarySegmentIndex << " bigger than number of boundary segments in grid.\n"
186 "Make sure to call this function only for boundaries that were defined as physical entities in gmsh.");
187 return boundaryMarkers_[boundarySegmentIndex];
188 }
189
195 int getBoundaryDomainMarker(const Intersection& intersection) const
196 { return getBoundaryDomainMarker(intersection.boundarySegmentIndex()); }
197
201 bool wasInserted(const Intersection& intersection) const
202 { return gridFactory_->wasInserted(intersection); }
203
209 int getElementDomainMarker(const Element& element) const
210 {
211 if (dataSourceType_ != DataSourceType::gmsh)
212 DUNE_THROW(Dune::InvalidStateException, "Domain markers are only available for gmsh grids.");
213
214 // parameters are only given for level 0 elements
215 auto level0element = element;
216 while (level0element.hasFather())
217 level0element = level0element.father();
218
219 // in the parallel case the data is load balanced and then accessed with indices of the index set
220 // for UGGrid element data is read on all processes since UGGrid can't communicate element data (yet)
221 if (gridPtr_->comm().size() > 1 && !Detail::isUG<Grid>::value)
222 return elementMarkers_[gridPtr_->levelGridView(0).indexSet().index(level0element)];
223 else
224 return elementMarkers_[gridFactory_->insertionIndex(level0element)];
225 }
226
231 template<bool ug = Detail::isUG<Grid>::value, typename std::enable_if_t<!ug, int> = 0>
233 {
234 return DataHandle(*gridPtr_, *gridFactory_, elementMarkers_, boundaryMarkers_, faceMarkers_);
235 }
236
241 template<bool ug = Detail::isUG<Grid>::value, typename std::enable_if_t<ug, int> = 0>
243 {
244 return DataHandle(*gridPtr_, *gridFactory_, elementMarkers_, boundaryMarkers_);
245 }
246
247 // \}
248
252 // \{
253
259 double getParameter(const Element& element, const std::string& fieldName) const
260 {
261 if (dataSourceType_ != DataSourceType::vtk)
262 DUNE_THROW(Dune::InvalidStateException, "This access function is only available for data from VTK files.");
263
264 if (cellData_.count(fieldName) == 0)
265 DUNE_THROW(Dune::IOError, "No field with the name " << fieldName << " found in cell data");
266
267 // parameters are only given for level 0 elements
268 auto level0element = element;
269 while (level0element.hasFather())
270 level0element = level0element.father();
271
272 return cellData_.at(fieldName)[gridFactory_->insertionIndex(level0element)];
273 }
274
281 double getParameter(const Vertex& vertex, const std::string& fieldName) const
282 {
283 if (dataSourceType_ != DataSourceType::vtk)
284 DUNE_THROW(Dune::InvalidStateException, "This access function is only available for data from VTK files.");
285
286 if (vertex.level() != 0)
287 DUNE_THROW(Dune::IOError, "You can only obtain parameters for level 0 vertices!");
288
289 if (pointData_.count(fieldName) == 0)
290 DUNE_THROW(Dune::IOError, "No field with the name " << fieldName << " found in point data");
291
292 return pointData_.at(fieldName)[gridFactory_->insertionIndex(vertex)];
293 }
294
295 // \}
296
297private:
298 // grid and grid factor for gmsh grid data / vtk grid data
299 std::shared_ptr<Grid> gridPtr_;
300 std::shared_ptr<Dune::GridFactory<Grid>> gridFactory_;
301
306 std::vector<int> elementMarkers_;
307 std::vector<int> boundaryMarkers_;
308 std::vector<int> faceMarkers_;
309
313 VTKReader::Data cellData_, pointData_;
314
315 // dgf grid data
316 Dune::GridPtr<Grid> dgfGrid_;
317
318 // specify which type of data we have
319 // TODO unfortunately all grid readers provide different data types, should be streamlined (changes in Dune)
320 DataSourceType dataSourceType_;
321};
322
323} // namespace Dumux
324
325#endif
A data handle for commucating grid data for gmsh grids.
A vtk file reader using tinyxml2 as xml backend.
Definition: adapt.hh:29
Definition: common/pdesolver.hh:35
A data handle for commucating grid data for gmsh grids.
Definition: gmshgriddatahandle.hh:55
Definition: griddata.hh:58
Class for grid data attached to dgf or gmsh grid files.
Definition: griddata.hh:73
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:83
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:281
int getBoundaryDomainMarker(const Intersection &intersection) const
Return the boundary domain marker (Gmsh physical entity number) of an intersection Only available whe...
Definition: griddata.hh:195
double getParameter(const Element &element, const std::string &fieldName) const
Get a element parameter.
Definition: griddata.hh:259
GridData(std::shared_ptr< Grid > grid, std::shared_ptr< Dune::GridFactory< Grid > > factory, VTKReader::Data &&cellData, VTKReader::Data &&pointData)
constructor for gmsh grid data
Definition: griddata.hh:100
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:160
GridData(Dune::GridPtr< Grid > grid)
constructor for dgf grid data
Definition: griddata.hh:94
DataHandle createGmshDataHandle()
Create a data handle for communication of the data in parallel simulations.
Definition: griddata.hh:232
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:119
int getBoundaryDomainMarker(int boundarySegmentIndex) const
Return the boundary domain marker (Gmsh physical entity number) of an intersection Only available whe...
Definition: griddata.hh:180
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:209
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:135
bool wasInserted(const Intersection &intersection) const
Returns true if an intersection was inserted during grid creation.
Definition: griddata.hh:201
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:57