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