version 3.8
gmshreader.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_FACETCOUPLING_GMSH_READER_HH
13#define DUMUX_FACETCOUPLING_GMSH_READER_HH
14
15#include <algorithm>
16#include <cassert>
17#include <fstream>
18#include <iostream>
19#include <sstream>
20#include <typeinfo>
21#include <unordered_map>
22
23#include <dune/common/timer.hh>
24#include <dune/common/fvector.hh>
25#include <dune/geometry/type.hh>
26
28
29namespace Dumux {
30
48template <class BulkGrid, int numGrids>
50{
51 // extract some necessary info from bulk grid
52 static constexpr int bulkDim = BulkGrid::dimension;
53 static constexpr int bulkDimWorld = BulkGrid::dimensionworld;
54 using ctype = typename BulkGrid::ctype;
56 using GlobalPosition = Dune::FieldVector<ctype, bulkDimWorld>;
57
58 // determine minimum dimension for which a grid is created
59 static constexpr int minGridDim = BulkGrid::dimension - numGrids + 1;
60 static_assert(minGridDim >= 1, "Grids with dim < 1 cannot be read!");
61
62 // structure to store data on an element
63 using VertexIndexSet = std::vector<GridIndexType>;
64 struct ElementData
65 {
66 Dune::GeometryType gt;
67 VertexIndexSet cornerIndices;
68 };
69
70public:
73 void read(const std::string& fileName, bool verbose = false)
74 {
75 read(fileName, 0, verbose);
76 }
77
79 void read(const std::string& fileName, std::size_t boundarySegThresh, bool verbose = false)
80 {
81 Dune::Timer watch;
82 if (verbose) std::cout << "Opening " << fileName << std::endl;
83 std::ifstream gridFile(fileName);
84 if (gridFile.fail())
85 DUNE_THROW(Dune::InvalidStateException, "Could not open the given .msh file. Make sure it exists");
86
87 // currently we only support version 2 file format
88 std::string line;
89 std::getline(gridFile, line);
90 if (line.find("$MeshFormat") == std::string::npos)
91 DUNE_THROW(Dune::InvalidStateException, "Expected $MeshFormat in the first line of the grid file!");
92 std::getline(gridFile, line);
93 if (line.find_first_of("2") != 0)
94 DUNE_THROW(Dune::InvalidStateException, "Currently only Gmsh mesh file format version 2 is supported!");
95
96 // read file until we get to the list of nodes
97 while (line.find("$Nodes") == std::string::npos)
98 std::getline(gridFile, line);
99
100 // read all vertices
101 std::getline(gridFile, line);
102 const auto numVertices = convertString<std::size_t>(line);
103 gridVertices_.resize(numVertices);
104
105 std::getline(gridFile, line);
106 std::size_t vertexCount = 0;
107 while (line.find("$EndNodes") == std::string::npos)
108 {
109 // drop first entry in line (vertex index) and read coordinates
110 std::istringstream stream(line);
111 std::string buf; stream >> buf;
112 GlobalPosition v;
113 for (auto& coord : v)
114 {
115 stream >> coord;
116 if (stream.fail()) DUNE_THROW(Dune::IOError, "Could not read vertex coordinate");
117 }
118
119 // insert vertex into container and move to next line
120 gridVertices_[vertexCount++] = v;
121 std::getline(gridFile, line);
122 }
123
124 // we should always find as many vertices as the mesh file states
125 if (vertexCount != numVertices)
126 DUNE_THROW(Dune::InvalidStateException, "Couldn't find as many vertices as stated in the .msh file");
127
128 // read file until we get to the list of elements
129 while(line.find("$Elements") == std::string::npos)
130 std::getline(gridFile, line);
131
132 // read elements
133 std::getline(gridFile, line);
134 const auto numElements = convertString<std::size_t>(line);
135
136 // keep track of number of elements
137 std::array<std::size_t, numGrids> elementCount;
138 std::fill(elementCount.begin(), elementCount.end(), 0);
139
140 // Construct maps that map bulk grid vertex
141 // indices to lowDim vertex indices. -1 indicates non-initialized status
142 std::size_t elemCount = 0;
143 std::array<std::size_t, numGrids> gridVertexCount;
144 std::array<std::vector<GridIndexType>, numGrids> gridVertexMap;
145 std::array<std::vector<bool>, numGrids> idxIsAssigned;
146 std::fill(gridVertexCount.begin(), gridVertexCount.end(), 0);
147 std::fill(gridVertexMap.begin(), gridVertexMap.end(), std::vector<GridIndexType>(vertexCount));
148 std::fill(idxIsAssigned.begin(), idxIsAssigned.end(), std::vector<bool>(vertexCount, false));
149 std::getline(gridFile, line);
150 while (line.find("$EndElements") == std::string::npos)
151 {
152 // pass all indices into vector
153 std::istringstream stream(line);
154 std::string buf;
155 std::vector<std::size_t> lineData;
156 while (stream >> buf) lineData.push_back(convertString<std::size_t>(buf));
157 assert(lineData.size() >= 4 && "Grid format erroneous or unsupported");
158
159 // obtain geometry type
160 const auto gt = obtainGeometryType( lineData[1] );
161 const std::size_t physicalIndex = lineData[3];
162 const auto geoDim = gt.dim();
163 const bool isBoundarySeg = geoDim != bulkDim && physicalIndex < boundarySegThresh;
164 if (geoDim >= minGridDim-1)
165 {
166 // insert boundary segment
167 if ((isBoundarySeg || geoDim == minGridDim-1))
168 {
169 const unsigned int nextLevelGridIdx = bulkDim-geoDim-1;
170
171 VertexIndexSet corners;
172 auto it = lineData.begin()+2+lineData[2]+1;
173 for (; it != lineData.end(); ++it)
174 {
175 *it -= 1; // gmsh indices start from 1
176
177 // insert map if vertex is not inserted yet
178 if (!idxIsAssigned[nextLevelGridIdx][*it])
179 {
180 gridVertexMap[nextLevelGridIdx][*it] = gridVertexCount[nextLevelGridIdx]++;
181 idxIsAssigned[nextLevelGridIdx][*it] = true;
182 vertexIndices_[nextLevelGridIdx].push_back(*it);
183 }
184
185 corners.push_back(gridVertexMap[nextLevelGridIdx][*it]);
186 }
187
188 // marker = physical entity index
189 boundaryMarkerMaps_[nextLevelGridIdx].push_back(physicalIndex);
190 boundarySegments_[nextLevelGridIdx].push_back(corners);
191 }
192
193 // insert element
194 else
195 {
196 const unsigned int gridIdx = bulkDim-geoDim;
197
198 VertexIndexSet corners;
199 auto it = lineData.begin()+2+lineData[2]+1;
200 for (; it != lineData.end(); ++it)
201 {
202 *it -= 1; // gmsh indices start from 1
203
204 // insert map if vertex is not inserted yet
205 if (!idxIsAssigned[gridIdx][*it])
206 {
207 gridVertexMap[gridIdx][*it] = gridVertexCount[gridIdx]++;
208 idxIsAssigned[gridIdx][*it] = true;
209 vertexIndices_[gridIdx].push_back(*it);
210 }
211
212 corners.push_back(gridVertexMap[gridIdx][*it]);
213 }
214
215 // add data to embedments/embeddings
216 if (geoDim > minGridDim)
217 {
218 const auto gridElemCount = elementData_[gridIdx].size();
219 const auto& embeddedVIndices = vertexIndices_[gridIdx+1];
220 const auto& embeddedIndicesAssigned = idxIsAssigned[gridIdx+1];
221
222 VertexIndexSet cornerIndicesGlobal(corners.size());
223 for (unsigned int i = 0; i < corners.size(); ++i)
224 cornerIndicesGlobal[i] = vertexIndices_[gridIdx][corners[i]];
225 addEmbeddings(cornerIndicesGlobal, gridIdx, gridElemCount, embeddedVIndices, embeddedIndicesAssigned);
226 }
227
228 // ensure dune-specific corner ordering
229 reorder(gt, corners);
230
231 // insert element data to grid's container
232 elementMarkerMaps_[gridIdx].push_back(physicalIndex);
233 elementData_[gridIdx].emplace_back(ElementData({gt, corners}));
234 }
235 }
236
237 // get next line
238 std::getline(gridFile, line);
239 elemCount++;
240 }
241
242 // make sure we read all elements
243 if (elemCount != numElements)
244 DUNE_THROW(Dune::InvalidStateException, "Didn't read as many elements as stated in the .msh file");
245
246 if (verbose)
247 {
248 std::cout << "Finished reading gmsh file" << std::endl;
249 for (std::size_t id = 0; id < numGrids; ++id)
250 {
251 std::cout << elementData_[id].size() << " "
252 << bulkDim-id << "-dimensional elements comprising of "
253 << gridVertexCount[id] << " vertices";
254 if (id < numGrids-1) std::cout << "," << std::endl;
255 }
256 std::cout << " have been read in " << watch.elapsed() << " seconds." << std::endl;
257 }
258 }
259
261 const std::vector<GlobalPosition>& gridVertices() const
262 { return gridVertices_; }
263
265 VertexIndexSet& vertexIndices(std::size_t id)
266 {
267 assert(id < numGrids && "Index exceeds number of grids provided");
268 return vertexIndices_[id];
269 }
270
272 const std::vector<ElementData>& elementData(std::size_t id) const
273 {
274 assert(id < numGrids && "Index exceeds number of grids provided");
275 return elementData_[id];
276 }
277
279 const std::vector<VertexIndexSet>& boundarySegmentData(std::size_t id) const
280 {
281 assert(id < numGrids && "Index exceeds number of grids provided");
282 return boundarySegments_[id];
283 }
284
286 std::vector<int>& elementMarkerMap(std::size_t id)
287 {
288 assert(id < numGrids && "Index exceeds number of grids provided");
289 return elementMarkerMaps_[id];
290 }
291
293 std::vector<int>& boundaryMarkerMap(std::size_t id)
294 {
295 assert(id < numGrids && "Index exceeds number of grids provided");
296 return boundaryMarkerMaps_[id];
297 }
298
300 std::unordered_map< GridIndexType, std::vector<GridIndexType> >& embeddedEntityMap(std::size_t id)
301 {
302 assert(id < numGrids && "Index exceeds number of grids provided");
303 return embeddedEntityMaps_[id];
304 }
305
307 std::unordered_map< GridIndexType, std::vector<GridIndexType> >& adjoinedEntityMap(std::size_t id)
308 {
309 assert(id < numGrids && "Index exceeds number of grids provided");
310 return adjoinedEntityMaps_[id];
311 }
312
313private:
315 template<class T>
316 T convertString(const std::string& string) const
317 {
318 T value;
319 std::istringstream stream(string);
320 stream >> value;
321 if (stream.fail())
322 DUNE_THROW(Dune::InvalidStateException, "Couldn't convert string: " << string << "to type: " << typeid(T).name());
323 return value;
324 }
325
327 Dune::GeometryType obtainGeometryType(std::size_t gmshElemType) const
328 {
329 switch (gmshElemType)
330 {
331 case 15: return Dune::GeometryTypes::vertex; // points
332 case 1: return Dune::GeometryTypes::line; // lines
333 case 2: return Dune::GeometryTypes::triangle; // triangle
334 case 3: return Dune::GeometryTypes::quadrilateral; // quadrilateral
335 case 4: return Dune::GeometryTypes::tetrahedron; // tetrahedron
336 case 5: return Dune::GeometryTypes::hexahedron; // hexahedron
337 default:
338 DUNE_THROW(Dune::NotImplemented, "FacetCoupling gmsh reader for gmsh element type " << gmshElemType);
339 }
340 }
341
343 void reorder(const Dune::GeometryType gt, VertexIndexSet& cornerIndices) const
344 {
345 // triangles, lines & tetrahedra need no reordering
346 if (gt == Dune::GeometryTypes::hexahedron)
347 {
348 using std::swap;
349 assert(cornerIndices.size() == 8);
350 swap(cornerIndices[2], cornerIndices[3]);
351 swap(cornerIndices[6], cornerIndices[7]);
352 }
353 else if (gt == Dune::GeometryTypes::quadrilateral)
354 {
355 using std::swap;
356 assert(cornerIndices.size() == 4);
357 swap(cornerIndices[2], cornerIndices[3]);
358 }
359 }
360
362 void addEmbeddings(const VertexIndexSet& globalCornerIndices,
363 unsigned int gridIdx,
364 std::size_t curElemIdx,
365 const std::vector<GridIndexType>& embeddedVIndices,
366 const std::vector<bool>& embededIndexIsAssigned)
367 {
368 const unsigned int embeddedGridIdx = gridIdx+1;
369
370 // check for embedments only if any of the corners exist in the embedded grid
371 for (auto globalCornerIdx : globalCornerIndices)
372 {
373 if (embededIndexIsAssigned[globalCornerIdx])
374 {
375 for (std::size_t i = 0; i < elementData_[embeddedGridIdx].size(); ++i)
376 {
377 const auto& e = elementData_[embeddedGridIdx][i];
378
379 auto vertIsContained = [&embeddedVIndices, &globalCornerIndices] (auto eCornerIdx)
380 {
381 return std::find(globalCornerIndices.begin(),
382 globalCornerIndices.end(),
383 embeddedVIndices[eCornerIdx]) != globalCornerIndices.end();
384 };
385
386 // if all corners are contained within this element, it is embedded
387 if ( std::all_of(e.cornerIndices.begin(), e.cornerIndices.end(), vertIsContained) )
388 {
389 embeddedEntityMaps_[gridIdx][curElemIdx].push_back(i);
390 adjoinedEntityMaps_[embeddedGridIdx][i].push_back(curElemIdx);
391 }
392 }
393
394 return;
395 }
396 }
397 }
398
400 std::vector<GlobalPosition> gridVertices_;
401 std::array<VertexIndexSet, numGrids> vertexIndices_;
402 std::array<std::vector<ElementData>, numGrids> elementData_;
403 std::array<std::vector<VertexIndexSet>, numGrids> boundarySegments_;
404
406 std::array< std::unordered_map< GridIndexType, std::vector<GridIndexType> >, numGrids > embeddedEntityMaps_;
407 std::array< std::unordered_map< GridIndexType, std::vector<GridIndexType> >, numGrids > adjoinedEntityMaps_;
408
410 std::array< std::vector<int>, numGrids > elementMarkerMaps_;
411 std::array< std::vector<int>, numGrids > boundaryMarkerMaps_;
412};
413
414} // end namespace Dumux
415
416#endif
Reads gmsh files where (n-1)-dimensional grids are defined on the faces or edges of n-dimensional gri...
Definition: gmshreader.hh:50
const std::vector< VertexIndexSet > & boundarySegmentData(std::size_t id) const
Returns the vector of read elements for a grid.
Definition: gmshreader.hh:279
std::unordered_map< GridIndexType, std::vector< GridIndexType > > & embeddedEntityMap(std::size_t id)
Returns the maps of the embedded entities.
Definition: gmshreader.hh:300
void read(const std::string &fileName, bool verbose=false)
Definition: gmshreader.hh:73
std::vector< int > & elementMarkerMap(std::size_t id)
Returns the maps of element markers.
Definition: gmshreader.hh:286
void read(const std::string &fileName, std::size_t boundarySegThresh, bool verbose=false)
Reads the data from a given mesh file.
Definition: gmshreader.hh:79
std::unordered_map< GridIndexType, std::vector< GridIndexType > > & adjoinedEntityMap(std::size_t id)
Returns the maps of the embedments.
Definition: gmshreader.hh:307
const std::vector< ElementData > & elementData(std::size_t id) const
Returns the vector of read elements for a grid.
Definition: gmshreader.hh:272
const std::vector< GlobalPosition > & gridVertices() const
Returns the vector with all grid vertices (entire hierarchy)
Definition: gmshreader.hh:261
std::vector< int > & boundaryMarkerMap(std::size_t id)
Returns the maps of domain markers.
Definition: gmshreader.hh:293
VertexIndexSet & vertexIndices(std::size_t id)
Returns a grid's vertex indices.
Definition: gmshreader.hh:265
Defines the index types used for grid and local indices.
constexpr Line line
Definition: couplingmanager1d3d_line.hh:31
Definition: adapt.hh:17
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:27