Dune Core Modules (2.7.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 <memory>
14#include <string>
15#include <vector>
16
19#include <dune/common/to_unique_ptr.hh>
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 public:
52 // empty function since this class does not implement anything
53 static void registerFactory() {}
54 };
55
56 // quadratic boundary segments in 1d
57 /*
58 Note the points
59
60 (0) (alpha) (1)
61
62 are mapped to the points in global coordinates
63
64 p0 p2 p1
65
66 alpha is determined automatically from the given points.
67 */
68 template< int dimWorld >
69 struct GmshReaderQuadraticBoundarySegment< 2, dimWorld >
70 : public Dune::BoundarySegment< 2, dimWorld >
71 {
72 typedef GmshReaderQuadraticBoundarySegment< 2, dimWorld > ThisType;
73 typedef typename Dune::BoundarySegment< 2, dimWorld > :: ObjectStreamType ObjectStreamType;
74 typedef Dune::FieldVector< double, dimWorld > GlobalVector;
75
76 GmshReaderQuadraticBoundarySegment ( const GlobalVector &p0_, const GlobalVector &p1_, const GlobalVector &p2_)
77 : p0(p0_), p1(p1_), p2(p2_)
78 {
79 init();
80 }
81
82 GmshReaderQuadraticBoundarySegment( ObjectStreamType& in )
83 {
84 // key is read before by the factory
85 const int bytes = sizeof(double)*dimWorld;
86 in.read( (char *) &p0[ 0 ], bytes );
87 in.read( (char *) &p1[ 0 ], bytes );
88 in.read( (char *) &p2[ 0 ], bytes );
89 init();
90 }
91
92 static void registerFactory()
93 {
94 if( key() < 0 )
95 {
96 key() = Dune::BoundarySegment< 2, dimWorld >::template registerFactory< ThisType >();
97 }
98 }
99
100 virtual GlobalVector operator() ( const Dune::FieldVector<double,1> &local ) const
101 {
102 GlobalVector y;
103 y = 0.0;
104 y.axpy((local[0]-alpha)*(local[0]-1.0)/alpha,p0);
105 y.axpy(local[0]*(local[0]-1.0)/(alpha*(alpha-1.0)),p1);
106 y.axpy(local[0]*(local[0]-alpha)/(1.0-alpha),p2);
107 return y;
108 }
109
110 void backup( ObjectStreamType& out ) const
111 {
112 // backup key to identify object
113 out.write( (const char *) &key(), sizeof( int ) );
114 // backup data
115 const int bytes = sizeof(double)*dimWorld;
116 out.write( (const char*) &p0[ 0 ], bytes );
117 out.write( (const char*) &p1[ 0 ], bytes );
118 out.write( (const char*) &p2[ 0 ], bytes );
119 }
120
121 protected:
122 void init()
123 {
124 GlobalVector d1 = p1;
125 d1 -= p0;
126 GlobalVector d2 = p2;
127 d2 -= p1;
128
129 alpha=d1.two_norm()/(d1.two_norm()+d2.two_norm());
130 if (alpha<1E-6 || alpha>1-1E-6)
131 DUNE_THROW(Dune::IOError, "ration in quadratic boundary segment bad");
132 }
133
134 static int& key() {
135 static int k = -1;
136 return k;
137 }
138
139 private:
140 GlobalVector p0,p1,p2;
141 double alpha;
142 };
143
144
145 // quadratic boundary segments in 2d
146 /* numbering of points corresponding to gmsh:
147
148 2
149
150 5 4
151
152 0 3 1
153
154 Note: The vertices 3, 4, 5 are not necessarily at the edge midpoints but can
155 be placed with parameters alpha, beta , gamma at the following positions
156 in local coordinates:
157
158
159 2 = (0,1)
160
161 5 = (0,beta) 4 = (1-gamma/sqrt(2),gamma/sqrt(2))
162
163 0 = (0,0) 3 = (alpha,0) 1 = (1,0)
164
165 The parameters alpha, beta, gamma are determined from the given vertices in
166 global coordinates.
167 */
168 template<>
169 class GmshReaderQuadraticBoundarySegment< 3, 3 >
170 : public Dune::BoundarySegment< 3 >
171 {
172 typedef GmshReaderQuadraticBoundarySegment< 3, 3 > ThisType;
173 typedef typename Dune::BoundarySegment< 3 > :: ObjectStreamType ObjectStreamType;
174 public:
175 GmshReaderQuadraticBoundarySegment (Dune::FieldVector<double,3> p0_, Dune::FieldVector<double,3> p1_,
178 : p0(p0_), p1(p1_), p2(p2_), p3(p3_), p4(p4_), p5(p5_)
179 {
180 init();
181 }
182
183 GmshReaderQuadraticBoundarySegment( ObjectStreamType& in )
184 {
185 const int bytes = sizeof(double)*3;
186 in.read( (char *) &p0[ 0 ], bytes );
187 in.read( (char *) &p1[ 0 ], bytes );
188 in.read( (char *) &p2[ 0 ], bytes );
189 in.read( (char *) &p3[ 0 ], bytes );
190 in.read( (char *) &p4[ 0 ], bytes );
191 in.read( (char *) &p5[ 0 ], bytes );
192 init();
193 }
194
195 static void registerFactory()
196 {
197 if( key() < 0 )
198 {
199 key() = Dune::BoundarySegment< 3 >::template registerFactory< ThisType >();
200 }
201 }
202
203 virtual Dune::FieldVector<double,3> operator() (const Dune::FieldVector<double,2>& local) const
204 {
206 y = 0.0;
207 y.axpy(phi0(local),p0);
208 y.axpy(phi1(local),p1);
209 y.axpy(phi2(local),p2);
210 y.axpy(phi3(local),p3);
211 y.axpy(phi4(local),p4);
212 y.axpy(phi5(local),p5);
213 return y;
214 }
215
216 void backup( ObjectStreamType& out ) const
217 {
218 // backup key to identify object in factory
219 out.write( (const char*) &key(), sizeof( int ) );
220 // backup data
221 const int bytes = sizeof(double)*3;
222 out.write( (const char*) &p0[ 0 ], bytes );
223 out.write( (const char*) &p1[ 0 ], bytes );
224 out.write( (const char*) &p2[ 0 ], bytes );
225 out.write( (const char*) &p3[ 0 ], bytes );
226 out.write( (const char*) &p4[ 0 ], bytes );
227 out.write( (const char*) &p5[ 0 ], bytes );
228 }
229
230 protected:
231 void init()
232 {
233 sqrt2 = sqrt(2.0);
235
236 d1 = p3; d1 -= p0;
237 d2 = p1; d2 -= p3;
238 alpha=d1.two_norm()/(d1.two_norm()+d2.two_norm());
239 if (alpha<1E-6 || alpha>1-1E-6)
240 DUNE_THROW(Dune::IOError, "alpha in quadratic boundary segment bad");
241
242 d1 = p5; d1 -= p0;
243 d2 = p2; d2 -= p5;
244 beta=d1.two_norm()/(d1.two_norm()+d2.two_norm());
245 if (beta<1E-6 || beta>1-1E-6)
246 DUNE_THROW(Dune::IOError, "beta in quadratic boundary segment bad");
247
248 d1 = p4; d1 -= p1;
249 d2 = p2; d2 -= p4;
250 gamma=sqrt2*(d1.two_norm()/(d1.two_norm()+d2.two_norm()));
251 if (gamma<1E-6 || gamma>1-1E-6)
252 DUNE_THROW(Dune::IOError, "gamma in quadratic boundary segment bad");
253 }
254
255 static int& key() {
256 static int k = -1;
257 return k;
258 }
259
260 private:
261 // The six Lagrange basis function on the reference element
262 // for the points given above
263
264 double phi0 (const Dune::FieldVector<double,2>& local) const
265 {
266 return (alpha*beta-beta*local[0]-alpha*local[1])*(1-local[0]-local[1])/(alpha*beta);
267 }
268 double phi3 (const Dune::FieldVector<double,2>& local) const
269 {
270 return local[0]*(1-local[0]-local[1])/(alpha*(1-alpha));
271 }
272 double phi1 (const Dune::FieldVector<double,2>& local) const
273 {
274 return local[0]*(gamma*local[0]-(sqrt2-gamma-sqrt2*alpha)*local[1]-alpha*gamma)/(gamma*(1-alpha));
275 }
276 double phi5 (const Dune::FieldVector<double,2>& local) const
277 {
278 return local[1]*(1-local[0]-local[1])/(beta*(1-beta));
279 }
280 double phi4 (const Dune::FieldVector<double,2>& local) const
281 {
282 return local[0]*local[1]/((1-gamma/sqrt2)*gamma/sqrt2);
283 }
284 double phi2 (const Dune::FieldVector<double,2>& local) const
285 {
286 return local[1]*(beta*(1-gamma/sqrt2)-local[0]*(beta-gamma/sqrt2)-local[1]*(1-gamma/sqrt2))/((1-gamma/sqrt2)*(beta-1));
287 }
288
289 Dune::FieldVector<double,3> p0,p1,p2,p3,p4,p5;
290 double alpha,beta,gamma,sqrt2;
291 };
292
293 } // end empty namespace
294
296 template<typename GridType>
298 {
299 protected:
300 // private data
302 bool verbose;
303 bool insert_boundary_segments;
304 unsigned int number_of_real_vertices;
305 int boundary_element_count;
306 int element_count;
307 // read buffer
308 char buf[512];
309 std::string fileName;
310 // exported data
311 std::vector<int> boundary_id_to_physical_entity;
312 std::vector<int> element_index_to_physical_entity;
313
314 // static data
315 static const int dim = GridType::dimension;
316 static const int dimWorld = GridType::dimensionworld;
317 static_assert( (dimWorld <= 3), "GmshReader requires dimWorld <= 3." );
318
319 // typedefs
321
322 // don't use something like
323 // readfile(file, 1, "%s\n", buf);
324 // to skip the rest of of the line -- that will only skip the next
325 // whitespace-separated word! Use skipline() instead.
326 void readfile(FILE * file, int cnt, const char * format,
327 void* t1, void* t2 = 0, void* t3 = 0, void* t4 = 0,
328 void* t5 = 0, void* t6 = 0, void* t7 = 0, void* t8 = 0,
329 void* t9 = 0, void* t10 = 0)
330 {
331 off_t pos = ftello(file);
332 int c = fscanf(file, format, t1, t2, t3, t4, t5, t6, t7, t8, t9, t10);
333 if (c != cnt)
334 DUNE_THROW(Dune::IOError, "Error parsing " << fileName << " "
335 "file pos " << pos
336 << ": Expected '" << format << "', only read " << c << " entries instead of " << cnt << ".");
337 }
338
339 // skip over the rest of the line, including the terminating newline
340 void skipline(FILE * file)
341 {
342 int c;
343 do {
344 c = std::fgetc(file);
345 } while(c != '\n' && c != EOF);
346 }
347
348 public:
349
350 GmshReaderParser(Dune::GridFactory<GridType>& _factory, bool v, bool i) :
351 factory(_factory), verbose(v), insert_boundary_segments(i) {}
352
353 std::vector<int> & boundaryIdMap()
354 {
355 return boundary_id_to_physical_entity;
356 }
357
358 std::vector<int> & elementIndexMap()
359 {
360 return element_index_to_physical_entity;
361 }
362
363 void read (const std::string& f)
364 {
365 if (verbose) std::cout << "Reading " << dim << "d Gmsh grid..." << std::endl;
366
367 // open file name, we use C I/O
368 fileName = f;
369 FILE* file = fopen(fileName.c_str(),"rb");
370 if (file==0)
371 DUNE_THROW(Dune::IOError, "Could not open " << fileName);
372
373 //=========================================
374 // Header: Read vertices into vector
375 // Check vertices that are needed
376 //=========================================
377
378 number_of_real_vertices = 0;
379 boundary_element_count = 0;
380 element_count = 0;
381
382 // process header
383 double version_number;
384 int file_type, data_size;
385
386 readfile(file,1,"%s\n",buf);
387 if (strcmp(buf,"$MeshFormat")!=0)
388 DUNE_THROW(Dune::IOError, "expected $MeshFormat in first line");
389 readfile(file,3,"%lg %d %d\n",&version_number,&file_type,&data_size);
390 if( (version_number < 2.0) || (version_number > 2.2) )
391 DUNE_THROW(Dune::IOError, "can only read Gmsh version 2 files");
392 if (verbose) std::cout << "version " << version_number << " Gmsh file detected" << std::endl;
393 readfile(file,1,"%s\n",buf);
394 if (strcmp(buf,"$EndMeshFormat")!=0)
395 DUNE_THROW(Dune::IOError, "expected $EndMeshFormat");
396
397 // node section
398 int number_of_nodes;
399
400 readfile(file,1,"%s\n",buf);
401 if (strcmp(buf,"$Nodes")!=0)
402 DUNE_THROW(Dune::IOError, "expected $Nodes");
403 readfile(file,1,"%d\n",&number_of_nodes);
404 if (verbose) std::cout << "file contains " << number_of_nodes << " nodes" << std::endl;
405
406 // read nodes
407 // The '+1' is due to the fact that gmsh numbers node starting from 1 rather than from 0
408 std::vector< GlobalVector > nodes( number_of_nodes+1 );
409 {
410 int id;
411 double x[ 3 ];
412 for( int i = 1; i <= number_of_nodes; ++i )
413 {
414 readfile(file,4, "%d %lg %lg %lg\n", &id, &x[ 0 ], &x[ 1 ], &x[ 2 ] );
415
416 if (id > number_of_nodes) {
418 "Only dense sequences of node indices are currently supported (node index "
419 << id << " is invalid).");
420 }
421
422 // just store node position
423 for( int j = 0; j < dimWorld; ++j )
424 nodes[ id ][ j ] = x[ j ];
425 }
426 readfile(file,1,"%s\n",buf);
427 if (strcmp(buf,"$EndNodes")!=0)
428 DUNE_THROW(Dune::IOError, "expected $EndNodes");
429 }
430
431 // element section
432 readfile(file,1,"%s\n",buf);
433 if (strcmp(buf,"$Elements")!=0)
434 DUNE_THROW(Dune::IOError, "expected $Elements");
435 int number_of_elements;
436 readfile(file,1,"%d\n",&number_of_elements);
437 if (verbose) std::cout << "file contains " << number_of_elements << " elements" << std::endl;
438
439 //=========================================
440 // Pass 1: Select and insert those vertices in the file that
441 // actually occur as corners of an element.
442 //=========================================
443
444 off_t section_element_offset = ftello(file);
445 std::map<int,unsigned int> renumber;
446 for (int i=1; i<=number_of_elements; i++)
447 {
448 int id, elm_type, number_of_tags;
449 readfile(file,3,"%d %d %d ",&id,&elm_type,&number_of_tags);
450 for (int k=1; k<=number_of_tags; k++)
451 {
452 int blub;
453 readfile(file,1,"%d ",&blub);
454 // k == 1: physical entity (not used here)
455 // k == 2: elementary entity (not used here either)
456 // if version_number < 2.2:
457 // k == 3: mesh partition 0
458 // else
459 // k == 3: number of mesh partitions
460 // k => 4: mesh partition k-4
461 }
462 pass1HandleElement(file, elm_type, renumber, nodes);
463 }
464 if (verbose) std::cout << "number of real vertices = " << number_of_real_vertices << std::endl;
465 if (verbose) std::cout << "number of boundary elements = " << boundary_element_count << std::endl;
466 if (verbose) std::cout << "number of elements = " << element_count << std::endl;
467 readfile(file,1,"%s\n",buf);
468 if (strcmp(buf,"$EndElements")!=0)
469 DUNE_THROW(Dune::IOError, "expected $EndElements");
470 boundary_id_to_physical_entity.resize(boundary_element_count);
471 element_index_to_physical_entity.resize(element_count);
472
473 //==============================================
474 // Pass 2: Insert boundary segments and elements
475 //==============================================
476
477 fseeko(file, section_element_offset, SEEK_SET);
478 boundary_element_count = 0;
479 element_count = 0;
480 for (int i=1; i<=number_of_elements; i++)
481 {
482 int id, elm_type, number_of_tags;
483 readfile(file,3,"%d %d %d ",&id,&elm_type,&number_of_tags);
484 int physical_entity = -1;
485
486 for (int k=1; k<=number_of_tags; k++)
487 {
488 int blub;
489 readfile(file,1,"%d ",&blub);
490 if (k==1) physical_entity = blub;
491 }
492 pass2HandleElement(file, elm_type, renumber, nodes, physical_entity);
493 }
494 readfile(file,1,"%s\n",buf);
495 if (strcmp(buf,"$EndElements")!=0)
496 DUNE_THROW(Dune::IOError, "expected $EndElements");
497
498 fclose(file);
499 }
500
506 void pass1HandleElement(FILE* file, const int elm_type,
507 std::map<int,unsigned int> & renumber,
508 const std::vector< GlobalVector > & nodes)
509 {
510 // some data about gmsh elements
511 const int nDofs[16] = {-1, 2, 3, 4, 4, 8, 6, 5, 3, 6, -1, 10, -1, -1, -1, 1};
512 const int nVertices[16] = {-1, 2, 3, 4, 4, 8, 6, 5, 2, 3, -1, 4, -1, -1, -1, 1};
513 const int elementDim[16] = {-1, 1, 2, 2, 3, 3, 3, 3, 1, 2, -1, 3, -1, -1, -1, 0};
514
515 // test whether we support the element type
516 if ( not (elm_type > 0 && elm_type <= 15 // 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 // insert each vertex if it hasn't been inserted already
539 for (int i=0; i<nVertices[elm_type]; i++)
540 if (renumber.find(elementDofs[i])==renumber.end())
541 {
542 renumber[elementDofs[i]] = number_of_real_vertices++;
543 factory.insertVertex(nodes[elementDofs[i]]);
544 }
545
546 // count elements and boundary elements
547 if (elementDim[elm_type] == dim)
548 element_count++;
549 else
550 boundary_element_count++;
551
552 }
553
554
555
556 // generic-case: This is not supposed to be used at runtime.
557 template <class E, class V, class V2>
558 void boundarysegment_insert(
559 const V& nodes,
560 const E& elementDofs,
561 const V2& vertices
562 )
563 {
564 DUNE_THROW(Dune::IOError, "tried to create a 3D boundary segment in a non-3D Grid");
565 }
566
567 // 3d-case:
568 template <class E, class V>
569 void boundarysegment_insert(
570 const std::vector<FieldVector<double, 3> >& nodes,
571 const E& elementDofs,
572 const V& vertices
573 )
574 {
575 std::array<FieldVector<double,dimWorld>, 6> v;
576 for (int i=0; i<6; i++)
577 for (int j=0; j<dimWorld; j++)
578 v[i][j] = nodes[elementDofs[i]][j];
579
580 BoundarySegment<dim,dimWorld>* newBoundarySegment
581 = (BoundarySegment<dim,dimWorld>*) new GmshReaderQuadraticBoundarySegment< 3, 3 >( v[0], v[1], v[2],
582 v[3], v[4], v[5] );
583
584 factory.insertBoundarySegment( vertices,
585 std::shared_ptr<BoundarySegment<dim,dimWorld> >(newBoundarySegment) );
586 }
587
588
589
594 virtual void pass2HandleElement(FILE* file, const int elm_type,
595 std::map<int,unsigned int> & renumber,
596 const std::vector< GlobalVector > & nodes,
597 const int physical_entity)
598 {
599 // some data about gmsh elements
600 const int nDofs[16] = {-1, 2, 3, 4, 4, 8, 6, 5, 3, 6, -1, 10, -1, -1, -1, 1};
601 const int nVertices[16] = {-1, 2, 3, 4, 4, 8, 6, 5, 2, 3, -1, 4, -1, -1, -1, 1};
602 const int elementDim[16] = {-1, 1, 2, 2, 3, 3, 3, 3, 1, 2, -1, 3, -1, -1, -1, 0};
603
604 // test whether we support the element type
605 if ( not (elm_type > 0 && elm_type <= 15 // index in suitable range?
606 && (elementDim[elm_type] == dim || elementDim[elm_type] == (dim-1) ) ) ) // real element or boundary element?
607 {
608 skipline(file); // skip rest of line if element is unknown
609 return;
610 }
611
612 // The format string for parsing is n times '%d' in a row
613 std::string formatString = "%d";
614 for (int i=1; i<nDofs[elm_type]; i++)
615 formatString += " %d";
616 formatString += "\n";
617
618 // '10' is the largest number of dofs we may encounter in a .msh file
619 std::vector<int> elementDofs(10);
620
621 readfile(file,nDofs[elm_type], formatString.c_str(),
622 &(elementDofs[0]),&(elementDofs[1]),&(elementDofs[2]),
623 &(elementDofs[3]),&(elementDofs[4]),&(elementDofs[5]),
624 &(elementDofs[6]),&(elementDofs[7]),&(elementDofs[8]),
625 &(elementDofs[9]));
626
627 // correct differences between gmsh and Dune in the local vertex numbering
628 switch (elm_type)
629 {
630 case 3 : // 4-node quadrilateral
631 std::swap(elementDofs[2],elementDofs[3]);
632 break;
633 case 5 : // 8-node hexahedron
634 std::swap(elementDofs[2],elementDofs[3]);
635 std::swap(elementDofs[6],elementDofs[7]);
636 break;
637 case 7 : // 5-node pyramid
638 std::swap(elementDofs[2],elementDofs[3]);
639 break;
640 }
641
642 // renumber corners to account for the explicitly given vertex
643 // numbering in the file
644 std::vector<unsigned int> vertices(nVertices[elm_type]);
645
646 for (int i=0; i<nVertices[elm_type]; i++)
647 vertices[i] = renumber[elementDofs[i]];
648
649 // If it is an element, insert it as such
650 if (elementDim[elm_type] == dim) {
651
652 switch (elm_type)
653 {
654 case 1 : // 2-node line
656 break;
657 case 2 : // 3-node triangle
659 break;
660 case 3 : // 4-node quadrilateral
662 break;
663 case 4 : // 4-node tetrahedron
665 break;
666 case 5 : // 8-node hexahedron
668 break;
669 case 6 : // 6-node prism
671 break;
672 case 7 : // 5-node pyramid
674 break;
675 case 9 : // 6-node triangle
677 break;
678 case 11 : // 10-node tetrahedron
680 break;
681 }
682
683 } else {
684 // it must be a boundary segment then
685 if (insert_boundary_segments) {
686
687 switch (elm_type)
688 {
689 case 1 : // 2-node line
690 factory.insertBoundarySegment(vertices);
691 break;
692
693 case 2 : // 3-node triangle
694 factory.insertBoundarySegment(vertices);
695 break;
696
697 case 3 : // 4-node quadrilateral
698 factory.insertBoundarySegment(vertices);
699 break;
700
701 case 15 : // 1-node point
702 factory.insertBoundarySegment(vertices);
703 break;
704
705 case 8 : { // 3-node line
706 std::array<FieldVector<double,dimWorld>, 3> v;
707 for (int i=0; i<dimWorld; i++) {
708 v[0][i] = nodes[elementDofs[0]][i];
709 v[1][i] = nodes[elementDofs[2]][i]; // yes, the renumbering is intended!
710 v[2][i] = nodes[elementDofs[1]][i];
711 }
712 BoundarySegment<dim,dimWorld>* newBoundarySegment
713 = (BoundarySegment<dim,dimWorld>*) new GmshReaderQuadraticBoundarySegment< 2, dimWorld >(v[0], v[1], v[2]);
714 factory.insertBoundarySegment(vertices,
715 std::shared_ptr<BoundarySegment<dim,dimWorld> >(newBoundarySegment));
716 break;
717 }
718 case 9 : { // 6-node triangle
719 boundarysegment_insert(nodes, elementDofs, vertices);
720 break;
721 }
722
723 }
724
725 }
726 }
727
728 // count elements and boundary elements
729 if (elementDim[elm_type] == dim) {
730 element_index_to_physical_entity[element_count] = physical_entity;
731 element_count++;
732 } else {
733 boundary_id_to_physical_entity[boundary_element_count] = physical_entity;
734 boundary_element_count++;
735 }
736
737 }
738
739 };
740
765 template<typename GridType>
767 {
768 static void registerFactory(Dune::GridFactory<GridType>& factory)
769 {
770 // register boundary segment to boundary segment factory for possible load balancing
771 // this needs to be done on all cores since the type might not be known otherwise
772 GmshReaderQuadraticBoundarySegment< Grid::dimension, Grid::dimensionworld >::registerFactory();
773
774#ifndef NDEBUG
775 // check that this method is called on all cores
776 factory.comm().barrier();
777#endif
778 }
779 public:
780 typedef GridType Grid;
781
788 static ToUniquePtr<Grid> read (const std::string& fileName, bool verbose = true, bool insertBoundarySegments=true)
789 {
790 // make a grid factory
792
793 // register create methods for boundary segments
794 registerFactory( factory );
795
796 // create parse object and read grid on process 0
797 if (factory.comm().rank() == 0)
798 {
799 GmshReaderParser<Grid> parser(factory,verbose,insertBoundarySegments);
800 parser.read(fileName);
801 }
802
803 return factory.createGrid();
804 }
805
812 static ToUniquePtr<Grid> read (const std::string& fileName,
813 std::vector<int>& boundarySegmentToPhysicalEntity,
814 std::vector<int>& elementToPhysicalEntity,
815 bool verbose = true, bool insertBoundarySegments=true)
816 {
817 // make a grid factory
819
820 // register create methods for boundary segments
821 registerFactory( factory );
822
823 // create parse object and read grid on process 0
824 if (factory.comm().rank() == 0)
825 {
826 GmshReaderParser<Grid> parser(factory,verbose,insertBoundarySegments);
827 parser.read(fileName);
828
829 boundarySegmentToPhysicalEntity.swap(parser.boundaryIdMap());
830 elementToPhysicalEntity.swap(parser.elementIndexMap());
831 }
832
833 return factory.createGrid();
834 }
835
837 static void read (Dune::GridFactory<Grid>& factory, const std::string& fileName,
838 bool verbose = true, bool insertBoundarySegments=true)
839 {
840 // register create methods for boundary segments
841 registerFactory( factory );
842
843 // create parse object and read grid on process 0
844 if (factory.comm().rank() == 0)
845 {
846 GmshReaderParser<Grid> parser(factory,verbose,insertBoundarySegments);
847 parser.read(fileName);
848 }
849 }
850
852 static void read (Dune::GridFactory<Grid>& factory,
853 const std::string& fileName,
854 std::vector<int>& boundarySegmentToPhysicalEntity,
855 std::vector<int>& elementToPhysicalEntity,
856 bool verbose = true, bool insertBoundarySegments=true)
857 {
858 // register create methods for boundary segments
859 registerFactory( factory );
860
861 // create parse object and read grid on process 0
862 if (factory.comm().rank() == 0)
863 {
864 GmshReaderParser<Grid> parser(factory,verbose,insertBoundarySegments);
865 parser.read(fileName);
866
867 boundarySegmentToPhysicalEntity.swap(parser.boundaryIdMap());
868 elementToPhysicalEntity.swap(parser.elementIndexMap());
869 }
870 }
871 };
872
875} // namespace Dune
876
877#endif
Base class for grid boundary segments of arbitrary geometry.
int rank() const
Return rank, is between 0 and size()-1.
Definition: communication.hh:95
int barrier() const
Wait until all processes have arrived at this point in the program.
Definition: communication.hh:251
FieldTraits< value_type >::real_type two_norm() const
two norm sqrt(sum over squared values of entries)
Definition: densevector.hh:642
derived_type & axpy(const field_type &a, const DenseVector< Other > &x)
vector space axpy operation ( *this += a x )
Definition: densevector.hh:576
vector space out of a tensor product of fields.
Definition: fvector.hh:96
dimension independent parts for GmshReaderParser
Definition: gmshreader.hh:298
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:506
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:594
Read Gmsh mesh file.
Definition: gmshreader.hh:767
static ToUniquePtr< Grid > read(const std::string &fileName, bool verbose=true, bool insertBoundarySegments=true)
Definition: gmshreader.hh:788
static void read(Dune::GridFactory< Grid > &factory, const std::string &fileName, bool verbose=true, bool insertBoundarySegments=true)
Definition: gmshreader.hh:837
static ToUniquePtr< Grid > read(const std::string &fileName, std::vector< int > &boundarySegmentToPhysicalEntity, std::vector< int > &elementToPhysicalEntity, bool verbose=true, bool insertBoundarySegments=true)
Definition: gmshreader.hh:812
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:852
Communication comm() const
Return the Communication used by the grid factory.
Definition: gridfactory.hh:252
Provide a generic factory class for unstructured grids.
Definition: gridfactory.hh:269
virtual void insertElement(const GeometryType &type, const std::vector< unsigned int > &vertices)
Insert an element into the coarse grid.
Definition: gridfactory.hh:301
virtual ToUniquePtr< GridType > createGrid()
Finalize grid creation and hand over the grid.
Definition: gridfactory.hh:327
virtual void insertVertex(const FieldVector< ctype, dimworld > &pos)
Insert a vertex into the coarse grid.
Definition: gridfactory.hh:290
virtual void insertBoundarySegment(const std::vector< unsigned int > &vertices)
insert a boundary segment
Definition: gridfactory.hh:319
Default exception class for I/O errors.
Definition: exceptions.hh:229
An owning pointer wrapper that can be assigned to (smart) pointers. Cannot be copied....
Definition: to_unique_ptr.hh:37
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:803
constexpr GeometryType prism
GeometryType representing a 3D prism.
Definition: type.hh:833
constexpr GeometryType triangle
GeometryType representing a triangle.
Definition: type.hh:809
constexpr GeometryType quadrilateral
GeometryType representing a quadrilateral (a square).
Definition: type.hh:815
constexpr GeometryType hexahedron
GeometryType representing a hexahedron.
Definition: type.hh:839
constexpr GeometryType pyramid
GeometryType representing a 3D pyramid.
Definition: type.hh:827
constexpr GeometryType tetrahedron
GeometryType representing a tetrahedron.
Definition: type.hh:821
static std::string formatString(const std::string &s, const T &... args)
Format values according to printf format string.
Definition: stringutility.hh:71
Dune namespace.
Definition: alignedallocator.hh:14
Base class for classes implementing geometries of boundary segments.
Definition: boundarysegment.hh:92
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 (Nov 12, 23:30, 2024)