Dune Core Modules (2.6.0)

gmshreader.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3
4#ifndef DUNE_GMSHREADER_HH
5#define DUNE_GMSHREADER_HH
6
7#include <cstdarg>
8#include <cstdio>
9#include <cstring>
10#include <fstream>
11#include <iostream>
12#include <map>
13#include <memory>
14#include <string>
15#include <vector>
16
19
20#include <dune/geometry/type.hh>
21
24
25namespace Dune
26{
27
35 {
41 };
42 };
43
44 namespace {
45
46 // arbitrary dimension, implementation is in specialization
47 template< int dimension, int dimWorld = dimension >
48 class GmshReaderQuadraticBoundarySegment
49 {};
50
51 // quadratic boundary segments in 1d
52 /*
53 Note the points
54
55 (0) (alpha) (1)
56
57 are mapped to the points in global coordinates
58
59 p0 p2 p1
60
61 alpha is determined automatically from the given points.
62 */
63 template< int dimWorld >
64 struct GmshReaderQuadraticBoundarySegment< 2, dimWorld >
65 : public Dune::BoundarySegment< 2, dimWorld >
66 {
67 typedef Dune::FieldVector< double, dimWorld > GlobalVector;
68
69 GmshReaderQuadraticBoundarySegment ( const GlobalVector &p0_, const GlobalVector &p1_, const GlobalVector &p2_)
70 : p0(p0_), p1(p1_), p2(p2_)
71 {
72 GlobalVector d1 = p1;
73 d1 -= p0;
74 GlobalVector d2 = p2;
75 d2 -= p1;
76
77 alpha=d1.two_norm()/(d1.two_norm()+d2.two_norm());
78 if (alpha<1E-6 || alpha>1-1E-6)
79 DUNE_THROW(Dune::IOError, "ration in quadratic boundary segment bad");
80 }
81
82 virtual GlobalVector operator() ( const Dune::FieldVector<double,1> &local ) const
83 {
84 GlobalVector y;
85 y = 0.0;
86 y.axpy((local[0]-alpha)*(local[0]-1.0)/alpha,p0);
87 y.axpy(local[0]*(local[0]-1.0)/(alpha*(alpha-1.0)),p1);
88 y.axpy(local[0]*(local[0]-alpha)/(1.0-alpha),p2);
89 return y;
90 }
91
92 private:
93 GlobalVector p0,p1,p2;
94 double alpha;
95 };
96
97
98 // quadratic boundary segments in 2d
99 /* numbering of points corresponding to gmsh:
100
101 2
102
103 5 4
104
105 0 3 1
106
107 Note: The vertices 3, 4, 5 are not necessarily at the edge midpoints but can
108 be placed with parameters alpha, beta , gamma at the following positions
109 in local coordinates:
110
111
112 2 = (0,1)
113
114 5 = (0,beta) 4 = (1-gamma/sqrt(2),gamma/sqrt(2))
115
116 0 = (0,0) 3 = (alpha,0) 1 = (1,0)
117
118 The parameters alpha, beta, gamma are determined from the given vertices in
119 global coordinates.
120 */
121 template<>
122 class GmshReaderQuadraticBoundarySegment< 3, 3 >
123 : public Dune::BoundarySegment< 3 >
124 {
125 public:
126 GmshReaderQuadraticBoundarySegment (Dune::FieldVector<double,3> p0_, Dune::FieldVector<double,3> p1_,
129 : p0(p0_), p1(p1_), p2(p2_), p3(p3_), p4(p4_), p5(p5_)
130 {
131 sqrt2 = sqrt(2.0);
133
134 d1 = p3; d1 -= p0;
135 d2 = p1; d2 -= p3;
136 alpha=d1.two_norm()/(d1.two_norm()+d2.two_norm());
137 if (alpha<1E-6 || alpha>1-1E-6)
138 DUNE_THROW(Dune::IOError, "alpha in quadratic boundary segment bad");
139
140 d1 = p5; d1 -= p0;
141 d2 = p2; d2 -= p5;
142 beta=d1.two_norm()/(d1.two_norm()+d2.two_norm());
143 if (beta<1E-6 || beta>1-1E-6)
144 DUNE_THROW(Dune::IOError, "beta in quadratic boundary segment bad");
145
146 d1 = p4; d1 -= p1;
147 d2 = p2; d2 -= p4;
148 gamma=sqrt2*(d1.two_norm()/(d1.two_norm()+d2.two_norm()));
149 if (gamma<1E-6 || gamma>1-1E-6)
150 DUNE_THROW(Dune::IOError, "gamma in quadratic boundary segment bad");
151 }
152
153 virtual Dune::FieldVector<double,3> operator() (const Dune::FieldVector<double,2>& local) const
154 {
156 y = 0.0;
157 y.axpy(phi0(local),p0);
158 y.axpy(phi1(local),p1);
159 y.axpy(phi2(local),p2);
160 y.axpy(phi3(local),p3);
161 y.axpy(phi4(local),p4);
162 y.axpy(phi5(local),p5);
163 return y;
164 }
165
166 private:
167 // The six Lagrange basis function on the reference element
168 // for the points given above
169
170 double phi0 (const Dune::FieldVector<double,2>& local) const
171 {
172 return (alpha*beta-beta*local[0]-alpha*local[1])*(1-local[0]-local[1])/(alpha*beta);
173 }
174 double phi3 (const Dune::FieldVector<double,2>& local) const
175 {
176 return local[0]*(1-local[0]-local[1])/(alpha*(1-alpha));
177 }
178 double phi1 (const Dune::FieldVector<double,2>& local) const
179 {
180 return local[0]*(gamma*local[0]-(sqrt2-gamma-sqrt2*alpha)*local[1]-alpha*gamma)/(gamma*(1-alpha));
181 }
182 double phi5 (const Dune::FieldVector<double,2>& local) const
183 {
184 return local[1]*(1-local[0]-local[1])/(beta*(1-beta));
185 }
186 double phi4 (const Dune::FieldVector<double,2>& local) const
187 {
188 return local[0]*local[1]/((1-gamma/sqrt2)*gamma/sqrt2);
189 }
190 double phi2 (const Dune::FieldVector<double,2>& local) const
191 {
192 return local[1]*(beta*(1-gamma/sqrt2)-local[0]*(beta-gamma/sqrt2)-local[1]*(1-gamma/sqrt2))/((1-gamma/sqrt2)*(beta-1));
193 }
194
195 Dune::FieldVector<double,3> p0,p1,p2,p3,p4,p5;
196 double alpha,beta,gamma,sqrt2;
197 };
198
199 } // end empty namespace
200
202 template<typename GridType>
204 {
205 protected:
206 // private data
208 bool verbose;
209 bool insert_boundary_segments;
210 unsigned int number_of_real_vertices;
211 int boundary_element_count;
212 int element_count;
213 // read buffer
214 char buf[512];
215 std::string fileName;
216 // exported data
217 std::vector<int> boundary_id_to_physical_entity;
218 std::vector<int> element_index_to_physical_entity;
219
220 // static data
221 static const int dim = GridType::dimension;
222 static const int dimWorld = GridType::dimensionworld;
223 static_assert( (dimWorld <= 3), "GmshReader requires dimWorld <= 3." );
224
225 // typedefs
227
228 // don't use something like
229 // readfile(file, 1, "%s\n", buf);
230 // to skip the rest of of the line -- that will only skip the next
231 // whitespace-separated word! Use skipline() instead.
232 void readfile(FILE * file, int cnt, const char * format,
233 void* t1, void* t2 = 0, void* t3 = 0, void* t4 = 0,
234 void* t5 = 0, void* t6 = 0, void* t7 = 0, void* t8 = 0,
235 void* t9 = 0, void* t10 = 0)
236 {
237 off_t pos = ftello(file);
238 int c = fscanf(file, format, t1, t2, t3, t4, t5, t6, t7, t8, t9, t10);
239 if (c != cnt)
240 DUNE_THROW(Dune::IOError, "Error parsing " << fileName << " "
241 "file pos " << pos
242 << ": Expected '" << format << "', only read " << c << " entries instead of " << cnt << ".");
243 }
244
245 // skip over the rest of the line, including the terminating newline
246 void skipline(FILE * file)
247 {
248 int c;
249 do {
250 c = std::fgetc(file);
251 } while(c != '\n' && c != EOF);
252 }
253
254 public:
255
256 GmshReaderParser(Dune::GridFactory<GridType>& _factory, bool v, bool i) :
257 factory(_factory), verbose(v), insert_boundary_segments(i) {}
258
259 std::vector<int> & boundaryIdMap()
260 {
261 return boundary_id_to_physical_entity;
262 }
263
264 std::vector<int> & elementIndexMap()
265 {
266 return element_index_to_physical_entity;
267 }
268
269 void read (const std::string& f)
270 {
271 if (verbose) std::cout << "Reading " << dim << "d Gmsh grid..." << std::endl;
272
273 // open file name, we use C I/O
274 fileName = f;
275 FILE* file = fopen(fileName.c_str(),"r");
276 if (file==0)
277 DUNE_THROW(Dune::IOError, "Could not open " << fileName);
278
279 //=========================================
280 // Header: Read vertices into vector
281 // Check vertices that are needed
282 //=========================================
283
284 number_of_real_vertices = 0;
285 boundary_element_count = 0;
286 element_count = 0;
287
288 // process header
289 double version_number;
290 int file_type, data_size;
291
292 readfile(file,1,"%s\n",buf);
293 if (strcmp(buf,"$MeshFormat")!=0)
294 DUNE_THROW(Dune::IOError, "expected $MeshFormat in first line");
295 readfile(file,3,"%lg %d %d\n",&version_number,&file_type,&data_size);
296 if( (version_number < 2.0) || (version_number > 2.2) )
297 DUNE_THROW(Dune::IOError, "can only read Gmsh version 2 files");
298 if (verbose) std::cout << "version " << version_number << " Gmsh file detected" << std::endl;
299 readfile(file,1,"%s\n",buf);
300 if (strcmp(buf,"$EndMeshFormat")!=0)
301 DUNE_THROW(Dune::IOError, "expected $EndMeshFormat");
302
303 // node section
304 int number_of_nodes;
305
306 readfile(file,1,"%s\n",buf);
307 if (strcmp(buf,"$Nodes")!=0)
308 DUNE_THROW(Dune::IOError, "expected $Nodes");
309 readfile(file,1,"%d\n",&number_of_nodes);
310 if (verbose) std::cout << "file contains " << number_of_nodes << " nodes" << std::endl;
311
312 // read nodes
313 // The '+1' is due to the fact that gmsh numbers node starting from 1 rather than from 0
314 std::vector< GlobalVector > nodes( number_of_nodes+1 );
315 {
316 int id;
317 double x[ 3 ];
318 for( int i = 1; i <= number_of_nodes; ++i )
319 {
320 readfile(file,4, "%d %lg %lg %lg\n", &id, &x[ 0 ], &x[ 1 ], &x[ 2 ] );
321
322 if (id > number_of_nodes) {
324 "Only dense sequences of node indices are currently supported (node index "
325 << id << " is invalid).");
326 }
327
328 // just store node position
329 for( int j = 0; j < dimWorld; ++j )
330 nodes[ id ][ j ] = x[ j ];
331 }
332 readfile(file,1,"%s\n",buf);
333 if (strcmp(buf,"$EndNodes")!=0)
334 DUNE_THROW(Dune::IOError, "expected $EndNodes");
335 }
336
337 // element section
338 readfile(file,1,"%s\n",buf);
339 if (strcmp(buf,"$Elements")!=0)
340 DUNE_THROW(Dune::IOError, "expected $Elements");
341 int number_of_elements;
342 readfile(file,1,"%d\n",&number_of_elements);
343 if (verbose) std::cout << "file contains " << number_of_elements << " elements" << std::endl;
344
345 //=========================================
346 // Pass 1: Select and insert those vertices in the file that
347 // actually occur as corners of an element.
348 //=========================================
349
350 long section_element_offset = ftell(file);
351 std::map<int,unsigned int> renumber;
352 for (int i=1; i<=number_of_elements; i++)
353 {
354 int id, elm_type, number_of_tags;
355 readfile(file,3,"%d %d %d ",&id,&elm_type,&number_of_tags);
356 for (int k=1; k<=number_of_tags; k++)
357 {
358 int blub;
359 readfile(file,1,"%d ",&blub);
360 // k == 1: physical entity (not used here)
361 // k == 2: elementary entity (not used here either)
362 // if version_number < 2.2:
363 // k == 3: mesh partition 0
364 // else
365 // k == 3: number of mesh partitions
366 // k => 4: mesh partition k-4
367 }
368 pass1HandleElement(file, elm_type, renumber, nodes);
369 }
370 if (verbose) std::cout << "number of real vertices = " << number_of_real_vertices << std::endl;
371 if (verbose) std::cout << "number of boundary elements = " << boundary_element_count << std::endl;
372 if (verbose) std::cout << "number of elements = " << element_count << std::endl;
373 readfile(file,1,"%s\n",buf);
374 if (strcmp(buf,"$EndElements")!=0)
375 DUNE_THROW(Dune::IOError, "expected $EndElements");
376 boundary_id_to_physical_entity.resize(boundary_element_count);
377 element_index_to_physical_entity.resize(element_count);
378
379 //==============================================
380 // Pass 2: Insert boundary segments and elements
381 //==============================================
382
383 fseek(file, section_element_offset, SEEK_SET);
384 boundary_element_count = 0;
385 element_count = 0;
386 for (int i=1; i<=number_of_elements; i++)
387 {
388 int id, elm_type, number_of_tags;
389 readfile(file,3,"%d %d %d ",&id,&elm_type,&number_of_tags);
390 int physical_entity = -1;
391
392 for (int k=1; k<=number_of_tags; k++)
393 {
394 int blub;
395 readfile(file,1,"%d ",&blub);
396 if (k==1) physical_entity = blub;
397 }
398 pass2HandleElement(file, elm_type, renumber, nodes, physical_entity);
399 }
400 readfile(file,1,"%s\n",buf);
401 if (strcmp(buf,"$EndElements")!=0)
402 DUNE_THROW(Dune::IOError, "expected $EndElements");
403
404 fclose(file);
405 }
406
412 void pass1HandleElement(FILE* file, const int elm_type,
413 std::map<int,unsigned int> & renumber,
414 const std::vector< GlobalVector > & nodes)
415 {
416 // some data about gmsh elements
417 const int nDofs[12] = {-1, 2, 3, 4, 4, 8, 6, 5, 3, 6, -1, 10};
418 const int nVertices[12] = {-1, 2, 3, 4, 4, 8, 6, 5, 2, 3, -1, 4};
419 const int elementDim[12] = {-1, 1, 2, 2, 3, 3, 3, 3, 1, 2, -1, 3};
420
421 // test whether we support the element type
422 if ( not (elm_type >= 0 && elm_type < 12 // index in suitable range?
423 && (elementDim[elm_type] == dim || elementDim[elm_type] == (dim-1) ) ) ) // real element or boundary element?
424 {
425 skipline(file); // skip rest of line if element is unknown
426 return;
427 }
428
429 // The format string for parsing is n times '%d' in a row
430 std::string formatString = "%d";
431 for (int i=1; i<nDofs[elm_type]; i++)
432 formatString += " %d";
433 formatString += "\n";
434
435 // '10' is the largest number of dofs we may encounter in a .msh file
436 std::vector<int> elementDofs(10);
437
438 readfile(file,nDofs[elm_type], formatString.c_str(),
439 &(elementDofs[0]),&(elementDofs[1]),&(elementDofs[2]),
440 &(elementDofs[3]),&(elementDofs[4]),&(elementDofs[5]),
441 &(elementDofs[6]),&(elementDofs[7]),&(elementDofs[8]),
442 &(elementDofs[9]));
443
444 // insert each vertex if it hasn't been inserted already
445 for (int i=0; i<nVertices[elm_type]; i++)
446 if (renumber.find(elementDofs[i])==renumber.end())
447 {
448 renumber[elementDofs[i]] = number_of_real_vertices++;
449 factory.insertVertex(nodes[elementDofs[i]]);
450 }
451
452 // count elements and boundary elements
453 if (elementDim[elm_type] == dim)
454 element_count++;
455 else
456 boundary_element_count++;
457
458 }
459
460
461
462 // generic-case: This is not supposed to be used at runtime.
463 template <class E, class V, class V2>
464 void boundarysegment_insert(
465 const V& nodes,
466 const E& elementDofs,
467 const V2& vertices
468 )
469 {
470 DUNE_THROW(Dune::IOError, "tried to create a 3D boundary segment in a non-3D Grid");
471 }
472
473 // 3d-case:
474 template <class E, class V>
475 void boundarysegment_insert(
476 const std::vector<FieldVector<double, 3> >& nodes,
477 const E& elementDofs,
478 const V& vertices
479 )
480 {
481 std::array<FieldVector<double,dimWorld>, 6> v;
482 for (int i=0; i<6; i++)
483 for (int j=0; j<dimWorld; j++)
484 v[i][j] = nodes[elementDofs[i]][j];
485
486 BoundarySegment<dim,dimWorld>* newBoundarySegment
487 = (BoundarySegment<dim,dimWorld>*) new GmshReaderQuadraticBoundarySegment< 3, 3 >( v[0], v[1], v[2],
488 v[3], v[4], v[5] );
489
490 factory.insertBoundarySegment( vertices,
491 std::shared_ptr<BoundarySegment<dim,dimWorld> >(newBoundarySegment) );
492 }
493
494
495
500 virtual void pass2HandleElement(FILE* file, const int elm_type,
501 std::map<int,unsigned int> & renumber,
502 const std::vector< GlobalVector > & nodes,
503 const int physical_entity)
504 {
505 // some data about gmsh elements
506 const int nDofs[12] = {-1, 2, 3, 4, 4, 8, 6, 5, 3, 6, -1, 10};
507 const int nVertices[12] = {-1, 2, 3, 4, 4, 8, 6, 5, 2, 3, -1, 4};
508 const int elementDim[12] = {-1, 1, 2, 2, 3, 3, 3, 3, 1, 2, -1, 3};
509
510 // test whether we support the element type
511 if ( not (elm_type >= 0 && elm_type < 12 // index in suitable range?
512 && (elementDim[elm_type] == dim || elementDim[elm_type] == (dim-1) ) ) ) // real element or boundary element?
513 {
514 skipline(file); // skip rest of line if element is unknown
515 return;
516 }
517
518 // The format string for parsing is n times '%d' in a row
519 std::string formatString = "%d";
520 for (int i=1; i<nDofs[elm_type]; i++)
521 formatString += " %d";
522 formatString += "\n";
523
524 // '10' is the largest number of dofs we may encounter in a .msh file
525 std::vector<int> elementDofs(10);
526
527 readfile(file,nDofs[elm_type], formatString.c_str(),
528 &(elementDofs[0]),&(elementDofs[1]),&(elementDofs[2]),
529 &(elementDofs[3]),&(elementDofs[4]),&(elementDofs[5]),
530 &(elementDofs[6]),&(elementDofs[7]),&(elementDofs[8]),
531 &(elementDofs[9]));
532
533 // correct differences between gmsh and Dune in the local vertex numbering
534 switch (elm_type)
535 {
536 case 3 : // 4-node quadrilateral
537 std::swap(elementDofs[2],elementDofs[3]);
538 break;
539 case 5 : // 8-node hexahedron
540 std::swap(elementDofs[2],elementDofs[3]);
541 std::swap(elementDofs[6],elementDofs[7]);
542 break;
543 case 7 : // 5-node pyramid
544 std::swap(elementDofs[2],elementDofs[3]);
545 break;
546 }
547
548 // renumber corners to account for the explicitly given vertex
549 // numbering in the file
550 std::vector<unsigned int> vertices(nVertices[elm_type]);
551
552 for (int i=0; i<nVertices[elm_type]; i++)
553 vertices[i] = renumber[elementDofs[i]];
554
555 // If it is an element, insert it as such
556 if (elementDim[elm_type] == dim) {
557
558 switch (elm_type)
559 {
560 case 1 : // 2-node line
562 break;
563 case 2 : // 3-node triangle
565 break;
566 case 3 : // 4-node quadrilateral
568 break;
569 case 4 : // 4-node tetrahedron
571 break;
572 case 5 : // 8-node hexahedron
574 break;
575 case 6 : // 6-node prism
577 break;
578 case 7 : // 5-node pyramid
580 break;
581 case 9 : // 6-node triangle
583 break;
584 case 11 : // 10-node tetrahedron
586 break;
587 }
588
589 } else {
590 // it must be a boundary segment then
591 if (insert_boundary_segments) {
592
593 switch (elm_type)
594 {
595 case 1 : // 2-node line
596 factory.insertBoundarySegment(vertices);
597 break;
598
599 case 2 : // 3-node triangle
600 factory.insertBoundarySegment(vertices);
601 break;
602
603 case 8 : { // 3-node line
604 std::array<FieldVector<double,dimWorld>, 3> v;
605 for (int i=0; i<dimWorld; i++) {
606 v[0][i] = nodes[elementDofs[0]][i];
607 v[1][i] = nodes[elementDofs[2]][i]; // yes, the renumbering is intended!
608 v[2][i] = nodes[elementDofs[1]][i];
609 }
610 BoundarySegment<dim,dimWorld>* newBoundarySegment
611 = (BoundarySegment<dim,dimWorld>*) new GmshReaderQuadraticBoundarySegment< 2, dimWorld >(v[0], v[1], v[2]);
612 factory.insertBoundarySegment(vertices,
613 std::shared_ptr<BoundarySegment<dim,dimWorld> >(newBoundarySegment));
614 break;
615 }
616 case 9 : { // 6-node triangle
617 boundarysegment_insert(nodes, elementDofs, vertices);
618 break;
619 }
620
621 }
622
623 }
624 }
625
626 // count elements and boundary elements
627 if (elementDim[elm_type] == dim) {
628 element_index_to_physical_entity[element_count] = physical_entity;
629 element_count++;
630 } else {
631 boundary_id_to_physical_entity[boundary_element_count] = physical_entity;
632 boundary_element_count++;
633 }
634
635 }
636
637 };
638
655 template<typename GridType>
657 {
658 public:
659 typedef GridType Grid;
660
662 static Grid* read (const std::string& fileName, bool verbose = true, bool insertBoundarySegments=true)
663 {
664 // make a grid factory
666
667 // create parse object
668 GmshReaderParser<Grid> parser(factory,verbose,insertBoundarySegments);
669 parser.read(fileName);
670
671 return factory.createGrid();
672 }
673
675 static Grid* read (const std::string& fileName,
676 std::vector<int>& boundarySegmentToPhysicalEntity,
677 std::vector<int>& elementToPhysicalEntity,
678 bool verbose = true, bool insertBoundarySegments=true)
679 {
680 // make a grid factory
682
683 // create parse object
684 GmshReaderParser<Grid> parser(factory,verbose,insertBoundarySegments);
685 parser.read(fileName);
686
687 boundarySegmentToPhysicalEntity.swap(parser.boundaryIdMap());
688 elementToPhysicalEntity.swap(parser.elementIndexMap());
689
690 return factory.createGrid();
691 }
692
694 static void read (Dune::GridFactory<Grid>& factory, const std::string& fileName,
695 bool verbose = true, bool insertBoundarySegments=true)
696 {
697 // create parse object
698 GmshReaderParser<Grid> parser(factory,verbose,insertBoundarySegments);
699 parser.read(fileName);
700 }
701
703 static void read (Dune::GridFactory<Grid>& factory,
704 const std::string& fileName,
705 std::vector<int>& boundarySegmentToPhysicalEntity,
706 std::vector<int>& elementToPhysicalEntity,
707 bool verbose = true, bool insertBoundarySegments=true)
708 {
709 // create parse object
710 GmshReaderParser<Grid> parser(factory,verbose,insertBoundarySegments);
711 parser.read(fileName);
712
713 boundarySegmentToPhysicalEntity.swap(parser.boundaryIdMap());
714 elementToPhysicalEntity.swap(parser.elementIndexMap());
715 }
716 };
717
720} // namespace Dune
721
722#endif
Base class for grid boundary segments of arbitrary geometry.
FieldTraits< value_type >::real_type two_norm() const
two norm sqrt(sum over squared values of entries)
Definition: densevector.hh:590
derived_type & axpy(const field_type &a, const DenseVector< Other > &y)
vector space axpy operation ( *this += a y )
Definition: densevector.hh:524
vector space out of a tensor product of fields.
Definition: fvector.hh:93
dimension independent parts for GmshReaderParser
Definition: gmshreader.hh:204
void pass1HandleElement(FILE *file, const int elm_type, std::map< int, unsigned int > &renumber, const std::vector< GlobalVector > &nodes)
Process one element during the first pass through the list of all elements.
Definition: gmshreader.hh:412
virtual void pass2HandleElement(FILE *file, const int elm_type, std::map< int, unsigned int > &renumber, const std::vector< GlobalVector > &nodes, const int physical_entity)
Process one element during the second pass through the list of all elements.
Definition: gmshreader.hh:500
Read Gmsh mesh file.
Definition: gmshreader.hh:657
static Grid * read(const std::string &fileName, std::vector< int > &boundarySegmentToPhysicalEntity, std::vector< int > &elementToPhysicalEntity, bool verbose=true, bool insertBoundarySegments=true)
Definition: gmshreader.hh:675
static Grid * read(const std::string &fileName, bool verbose=true, bool insertBoundarySegments=true)
Definition: gmshreader.hh:662
static void read(Dune::GridFactory< Grid > &factory, const std::string &fileName, bool verbose=true, bool insertBoundarySegments=true)
Definition: gmshreader.hh:694
static void read(Dune::GridFactory< Grid > &factory, const std::string &fileName, std::vector< int > &boundarySegmentToPhysicalEntity, std::vector< int > &elementToPhysicalEntity, bool verbose=true, bool insertBoundarySegments=true)
Definition: gmshreader.hh:703
Provide a generic factory class for unstructured grids.
Definition: gridfactory.hh:263
virtual void insertElement(const GeometryType &type, const std::vector< unsigned int > &vertices)
Insert an element into the coarse grid.
Definition: gridfactory.hh:290
virtual GridType * createGrid()
Finalize grid creation and hand over the grid.
Definition: gridfactory.hh:316
virtual void insertVertex(const FieldVector< ctype, dimworld > &pos)
Insert a vertex into the coarse grid.
Definition: gridfactory.hh:279
virtual void insertBoundarySegment(const std::vector< unsigned int > &vertices)
insert a boundary segment
Definition: gridfactory.hh:308
Default exception class for I/O errors.
Definition: exceptions.hh:229
Provide a generic factory class for unstructured grids.
A few common exception classes.
Implements a vector constructed from a given type representing a field and a compile-time given size.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
constexpr GeometryType line
GeometryType representing a line.
Definition: type.hh:733
constexpr GeometryType prism
GeometryType representing a 3D prism.
Definition: type.hh:763
constexpr GeometryType triangle
GeometryType representing a triangle.
Definition: type.hh:739
constexpr GeometryType quadrilateral
GeometryType representing a quadrilateral (a square).
Definition: type.hh:745
constexpr GeometryType hexahedron
GeometryType representing a hexahedron.
Definition: type.hh:769
constexpr GeometryType pyramid
GeometryType representing a 3D pyramid.
Definition: type.hh:757
constexpr GeometryType tetrahedron
GeometryType representing a tetrahedron.
Definition: type.hh:751
static std::string formatString(const std::string &s, const T &... args)
Format values according to printf format string.
Definition: stringutility.hh:72
Dune namespace.
Definition: alignedallocator.hh:10
Base class for classes implementing geometries of boundary segments.
Definition: boundarysegment.hh:30
Options for read operation.
Definition: gmshreader.hh:35
GeometryOrder
Definition: gmshreader.hh:36
@ firstOrder
edges are straight lines.
Definition: gmshreader.hh:38
@ secondOrder
quadratic boundary approximation.
Definition: gmshreader.hh:40
A unique label for each type of element that can occur in a grid.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)