Dune Core Modules (2.3.1)

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 dune_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 if( id != i )
322 DUNE_THROW( Dune::IOError, "Expected id " << i << "(got id " << id << "." );
323
324 // just store node position
325 for( int j = 0; j < dimWorld; ++j )
326 nodes[ i ][ j ] = x[ j ];
327 }
328 readfile(file,1,"%s\n",buf);
329 if (strcmp(buf,"$EndNodes")!=0)
330 DUNE_THROW(Dune::IOError, "expected $EndNodes");
331 }
332
333 // element section
334 readfile(file,1,"%s\n",buf);
335 if (strcmp(buf,"$Elements")!=0)
336 DUNE_THROW(Dune::IOError, "expected $Elements");
337 int number_of_elements;
338 readfile(file,1,"%d\n",&number_of_elements);
339 if (verbose) std::cout << "file contains " << number_of_elements << " elements" << std::endl;
340
341 //=========================================
342 // Pass 1: Renumber needed vertices
343 //=========================================
344
345 long section_element_offset = ftell(file);
346 std::map<int,unsigned int> renumber;
347 for (int i=1; i<=number_of_elements; i++)
348 {
349 int id, elm_type, number_of_tags;
350 readfile(file,3,"%d %d %d ",&id,&elm_type,&number_of_tags);
351 for (int k=1; k<=number_of_tags; k++)
352 {
353 int blub;
354 readfile(file,1,"%d ",&blub);
355 // k == 1: physical entity (not used here)
356 // k == 2: elementary entity (not used here either)
357 // if version_number < 2.2:
358 // k == 3: mesh partition 0
359 // else
360 // k == 3: number of mesh partitions
361 // k => 4: mesh partition k-4
362 }
363 pass1HandleElement(file, elm_type, renumber, nodes);
364 }
365 if (verbose) std::cout << "number of real vertices = " << number_of_real_vertices << std::endl;
366 if (verbose) std::cout << "number of boundary elements = " << boundary_element_count << std::endl;
367 if (verbose) std::cout << "number of elements = " << element_count << std::endl;
368 readfile(file,1,"%s\n",buf);
369 if (strcmp(buf,"$EndElements")!=0)
370 DUNE_THROW(Dune::IOError, "expected $EndElements");
371 boundary_id_to_physical_entity.resize(boundary_element_count);
372 element_index_to_physical_entity.resize(element_count);
373
374 //==============================================
375 // Pass 2: Insert boundary segments and elements
376 //==============================================
377
378 fseek(file, section_element_offset, SEEK_SET);
379 boundary_element_count = 0;
380 element_count = 0;
381 for (int i=1; i<=number_of_elements; i++)
382 {
383 int id, elm_type, number_of_tags;
384 readfile(file,3,"%d %d %d ",&id,&elm_type,&number_of_tags);
385 int physical_entity = -1;
386 std::vector<int> mesh_partitions;
387 if ( version_number < 2.2 )
388 {
389 mesh_partitions.resize(1);
390 }
391 for (int k=1; k<=number_of_tags; k++)
392 {
393 int blub;
394 readfile(file,1,"%d ",&blub);
395 if (k==1) physical_entity = blub;
396 // k == 2: elementary entity (not used here)
397 if ( version_number < 2.2 )
398 {
399 if (k==3) mesh_partitions[0] = blub;
400 }
401 else
402 {
403 if (k > 3)
404 mesh_partitions[k-4] = blub;
405 else
406 mesh_partitions.resize(blub);
407 }
408 }
409 pass2HandleElement(file, elm_type, renumber, nodes, physical_entity);
410 }
411 readfile(file,1,"%s\n",buf);
412 if (strcmp(buf,"$EndElements")!=0)
413 DUNE_THROW(Dune::IOError, "expected $EndElements");
414
415 fclose(file);
416 }
417
418 // dimension dependent routines
419 void pass1HandleElement(FILE* file, const int elm_type,
420 std::map<int,unsigned int> & renumber,
421 const std::vector< GlobalVector > & nodes)
422 {
423 // some data about gmsh elements
424 const int nDofs[12] = {-1, 2, 3, 4, 4, 8, 6, 5, 3, 6, -1, 10};
425 const int nVertices[12] = {-1, 2, 3, 4, 4, 8, 6, 5, 2, 3, -1, 4};
426 const int elementDim[12] = {-1, 1, 2, 2, 3, 3, 3, 3, 1, 2, -1, 3};
427
428 // test whether we support the element type
429 if ( not (elm_type >= 0 && elm_type < 12 // index in suitable range?
430 && (elementDim[elm_type] == dim || elementDim[elm_type] == (dim-1) ) ) ) // real element or boundary element?
431 {
432 skipline(file); // skip rest of line if element is unknown
433 return;
434 }
435
436 // The format string for parsing is n times '%d' in a row
437 std::string formatString = "%d";
438 for (int i=1; i<nDofs[elm_type]; i++)
439 formatString += " %d";
440 formatString += "\n";
441
442 // '10' is the largest number of dofs we may encounter in a .msh file
443 std::vector<int> elementDofs(10);
444
445 readfile(file,nDofs[elm_type], formatString.c_str(),
446 &(elementDofs[0]),&(elementDofs[1]),&(elementDofs[2]),
447 &(elementDofs[3]),&(elementDofs[4]),&(elementDofs[5]),
448 &(elementDofs[6]),&(elementDofs[7]),&(elementDofs[8]),
449 &(elementDofs[9]));
450
451 // insert each vertex if it hasn't been inserted already
452 for (int i=0; i<nVertices[elm_type]; i++)
453 if (renumber.find(elementDofs[i])==renumber.end())
454 {
455 renumber[elementDofs[i]] = number_of_real_vertices++;
456 factory.insertVertex(nodes[elementDofs[i]]);
457 }
458
459 // count elements and boundary elements
460 if (elementDim[elm_type] == dim)
461 element_count++;
462 else
463 boundary_element_count++;
464
465 }
466
467
468
469 // generic-case: This is not supposed to be used at runtime.
470 template <class E, class V, class V2>
471 void boundarysegment_insert(
472 const V& nodes,
473 const E& elementDofs,
474 const V2& vertices
475 )
476 {
477 DUNE_THROW(Dune::IOError, "tried to create a 3D boundary segment in a non-3D Grid");
478 }
479
480 // 3d-case:
481 template <class E, class V>
482 void boundarysegment_insert(
483 const std::vector<FieldVector<double, 3> >& nodes,
484 const E& elementDofs,
485 const V& vertices
486 )
487 {
489 for (int i=0; i<6; i++)
490 for (int j=0; j<dimWorld; j++)
491 v[i][j] = nodes[elementDofs[i]][j];
492
493 BoundarySegment<dim,dimWorld>* newBoundarySegment
494 = (BoundarySegment<dim,dimWorld>*) new GmshReaderQuadraticBoundarySegment< 3, 3 >( v[0], v[1], v[2],
495 v[3], v[4], v[5] );
496
497 factory.insertBoundarySegment( vertices,
498 shared_ptr<BoundarySegment<dim,dimWorld> >(newBoundarySegment) );
499 }
500
501
502
503
504
505 virtual void pass2HandleElement(FILE* file, const int elm_type,
506 std::map<int,unsigned int> & renumber,
507 const std::vector< GlobalVector > & nodes,
508 const int physical_entity)
509 {
510 // some data about gmsh elements
511 const int nDofs[12] = {-1, 2, 3, 4, 4, 8, 6, 5, 3, 6, -1, 10};
512 const int nVertices[12] = {-1, 2, 3, 4, 4, 8, 6, 5, 2, 3, -1, 4};
513 const int elementDim[12] = {-1, 1, 2, 2, 3, 3, 3, 3, 1, 2, -1, 3};
514
515 // test whether we support the element type
516 if ( not (elm_type >= 0 && elm_type < 12 // index in suitable range?
517 && (elementDim[elm_type] == dim || elementDim[elm_type] == (dim-1) ) ) ) // real element or boundary element?
518 {
519 skipline(file); // skip rest of line if element is unknown
520 return;
521 }
522
523 // The format string for parsing is n times '%d' in a row
524 std::string formatString = "%d";
525 for (int i=1; i<nDofs[elm_type]; i++)
526 formatString += " %d";
527 formatString += "\n";
528
529 // '10' is the largest number of dofs we may encounter in a .msh file
530 std::vector<int> elementDofs(10);
531
532 readfile(file,nDofs[elm_type], formatString.c_str(),
533 &(elementDofs[0]),&(elementDofs[1]),&(elementDofs[2]),
534 &(elementDofs[3]),&(elementDofs[4]),&(elementDofs[5]),
535 &(elementDofs[6]),&(elementDofs[7]),&(elementDofs[8]),
536 &(elementDofs[9]));
537
538 // correct differences between gmsh and Dune in the local vertex numbering
539 switch (elm_type)
540 {
541 case 3 : // 4-node quadrilateral
542 std::swap(elementDofs[2],elementDofs[3]);
543 break;
544 case 5 : // 8-node hexahedron
545 std::swap(elementDofs[2],elementDofs[3]);
546 std::swap(elementDofs[6],elementDofs[7]);
547 break;
548 case 7 : // 5-node pyramid
549 std::swap(elementDofs[2],elementDofs[3]);
550 break;
551 }
552
553 // renumber corners to account for the explicitly given vertex
554 // numbering in the file
555 std::vector<unsigned int> vertices(nVertices[elm_type]);
556
557 for (int i=0; i<nVertices[elm_type]; i++)
558 vertices[i] = renumber[elementDofs[i]];
559
560 // If it is an element, insert it as such
561 if (elementDim[elm_type] == dim) {
562
563 switch (elm_type)
564 {
565 case 1 : // 2-node line
567 break;
568 case 2 : // 3-node triangle
570 break;
571 case 3 : // 4-node quadrilateral
573 break;
574 case 4 : // 4-node tetrahedron
576 break;
577 case 5 : // 8-node hexahedron
579 break;
580 case 6 : // 6-node prism
582 break;
583 case 7 : // 5-node pyramid
585 break;
586 case 9 : // 6-node triangle
588 break;
589 case 11 : // 10-node tetrahedron
591 break;
592 }
593
594 } else {
595 // it must be a boundary segment then
596 if (insert_boundary_segments) {
597
598 switch (elm_type)
599 {
600 case 1 : // 2-node line
601 factory.insertBoundarySegment(vertices);
602 break;
603
604 case 2 : // 3-node triangle
605 factory.insertBoundarySegment(vertices);
606 break;
607
608 case 8 : { // 3-node line
610 for (int i=0; i<dimWorld; i++) {
611 v[0][i] = nodes[elementDofs[0]][i];
612 v[1][i] = nodes[elementDofs[2]][i]; // yes, the renumbering is intended!
613 v[2][i] = nodes[elementDofs[1]][i];
614 }
615 BoundarySegment<dim,dimWorld>* newBoundarySegment
616 = (BoundarySegment<dim,dimWorld>*) new GmshReaderQuadraticBoundarySegment< 2, dimWorld >(v[0], v[1], v[2]);
617 factory.insertBoundarySegment(vertices,
618 shared_ptr<BoundarySegment<dim,dimWorld> >(newBoundarySegment));
619 break;
620 }
621 case 9 : { // 6-node triangle
622 boundarysegment_insert(nodes, elementDofs, vertices);
623 break;
624 }
625
626 }
627
628 }
629 }
630
631 // count elements and boundary elements
632 if (elementDim[elm_type] == dim) {
633 element_index_to_physical_entity[element_count] = physical_entity;
634 element_count++;
635 } else {
636 boundary_id_to_physical_entity[boundary_element_count] = physical_entity;
637 boundary_element_count++;
638 }
639
640 }
641
642 };
643
660 template<typename GridType>
662 {
663 public:
664 typedef GridType Grid;
665
667 static Grid* read (const std::string& fileName, bool verbose = true, bool insert_boundary_segments=true)
668 {
669 // make a grid factory
671
672 // create parse object
673 GmshReaderParser<Grid> parser(factory,verbose,insert_boundary_segments);
674 parser.read(fileName);
675
676 return factory.createGrid();
677 }
678
680 static Grid* read (const std::string& fileName,
681 std::vector<int>& boundary_id_to_physical_entity,
682 std::vector<int>& element_index_to_physical_entity,
683 bool verbose = true, bool insert_boundary_segments=true)
684 {
685 // make a grid factory
687
688 // create parse object
689 GmshReaderParser<Grid> parser(factory,verbose,insert_boundary_segments);
690 parser.read(fileName);
691
692 boundary_id_to_physical_entity.swap(parser.boundaryIdMap());
693 element_index_to_physical_entity.swap(parser.elementIndexMap());
694
695 return factory.createGrid();
696 }
697
699 static void read (Dune::GridFactory<Grid>& factory, const std::string& fileName,
700 bool verbose = true, bool insert_boundary_segments=true)
701 {
702 // create parse object
703 GmshReaderParser<Grid> parser(factory,verbose,insert_boundary_segments);
704 parser.read(fileName);
705 }
706
708 static void read (Dune::GridFactory<Grid>& factory,
709 const std::string& fileName,
710 std::vector<int>& boundary_id_to_physical_entity,
711 std::vector<int>& element_index_to_physical_entity,
712 bool verbose = true, bool insert_boundary_segments=true)
713 {
714 // create parse object
715 GmshReaderParser<Grid> parser(factory,verbose,insert_boundary_segments);
716 parser.read(fileName);
717
718 boundary_id_to_physical_entity.swap(parser.boundaryIdMap());
719 element_index_to_physical_entity.swap(parser.elementIndexMap());
720 }
721 };
722
725} // namespace Dune
726
727#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:521
derived_type & axpy(const value_type &a, const DenseVector< Other > &y)
vector space axpy operation ( *this += a y )
Definition: densevector.hh:456
vector space out of a tensor product of fields.
Definition: fvector.hh:92
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:662
static Grid * read(const std::string &fileName, bool verbose=true, bool insert_boundary_segments=true)
Definition: gmshreader.hh:667
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:708
static void read(Dune::GridFactory< Grid > &factory, const std::string &fileName, bool verbose=true, bool insert_boundary_segments=true)
Definition: gmshreader.hh:699
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:680
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:257
Simple fixed size array class. This replaces std::array, if that is not available.
Definition: array.hh:40
A reference counting smart pointer.
Definition: shared_ptr.hh:64
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:244
Dune namespace.
Definition: alignment.hh:14
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 (Jul 15, 22:36, 2024)