Dune Core Modules (2.4.2)

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