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