MFEM  v3.4
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
mesh_readers.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #include "mesh_headers.hpp"
13 #include "../fem/fem.hpp"
14 #include "../general/text.hpp"
15 
16 #include <iostream>
17 #include <cstdio>
18 
19 #ifdef MFEM_USE_NETCDF
20 #include "netcdf.h"
21 #endif
22 
23 using namespace std;
24 
25 namespace mfem
26 {
27 
28 bool Mesh::remove_unused_vertices = true;
29 
30 void Mesh::ReadMFEMMesh(std::istream &input, bool mfem_v11, int &curved)
31 {
32  // Read MFEM mesh v1.0 format
33  string ident;
34 
35  // read lines beginning with '#' (comments)
36  skip_comment_lines(input, '#');
37  input >> ident; // 'dimension'
38 
39  MFEM_VERIFY(ident == "dimension", "invalid mesh file");
40  input >> Dim;
41 
42  skip_comment_lines(input, '#');
43  input >> ident; // 'elements'
44 
45  MFEM_VERIFY(ident == "elements", "invalid mesh file");
46  input >> NumOfElements;
47  elements.SetSize(NumOfElements);
48  for (int j = 0; j < NumOfElements; j++)
49  {
50  elements[j] = ReadElement(input);
51  }
52 
53  skip_comment_lines(input, '#');
54  input >> ident; // 'boundary'
55 
56  MFEM_VERIFY(ident == "boundary", "invalid mesh file");
57  input >> NumOfBdrElements;
58  boundary.SetSize(NumOfBdrElements);
59  for (int j = 0; j < NumOfBdrElements; j++)
60  {
61  boundary[j] = ReadElement(input);
62  }
63 
64  skip_comment_lines(input, '#');
65  input >> ident;
66 
67  if (mfem_v11 && ident == "vertex_parents")
68  {
69  ncmesh = new NCMesh(this, &input);
70  // NOTE: the constructor above will call LoadVertexParents
71 
72  skip_comment_lines(input, '#');
73  input >> ident;
74 
75  if (ident == "coarse_elements")
76  {
77  ncmesh->LoadCoarseElements(input);
78 
79  skip_comment_lines(input, '#');
80  input >> ident;
81  }
82  }
83 
84  MFEM_VERIFY(ident == "vertices", "invalid mesh file");
85  input >> NumOfVertices;
86  vertices.SetSize(NumOfVertices);
87 
88  input >> ws >> ident;
89  if (ident != "nodes")
90  {
91  // read the vertices
92  spaceDim = atoi(ident.c_str());
93  for (int j = 0; j < NumOfVertices; j++)
94  {
95  for (int i = 0; i < spaceDim; i++)
96  {
97  input >> vertices[j](i);
98  }
99  }
100 
101  // initialize vertex positions in NCMesh
102  if (ncmesh) { ncmesh->SetVertexPositions(vertices); }
103  }
104  else
105  {
106  // prepare to read the nodes
107  input >> ws;
108  curved = 1;
109  }
110 
111  // When visualizing solutions on non-conforming grids, PETSc
112  // may dump additional vertices
113  if (remove_unused_vertices) { RemoveUnusedVertices(); }
114 }
115 
116 void Mesh::ReadLineMesh(std::istream &input)
117 {
118  int j,p1,p2,a;
119 
120  Dim = 1;
121 
122  input >> NumOfVertices;
123  vertices.SetSize(NumOfVertices);
124  // Sets vertices and the corresponding coordinates
125  for (j = 0; j < NumOfVertices; j++)
126  {
127  input >> vertices[j](0);
128  }
129 
130  input >> NumOfElements;
131  elements.SetSize(NumOfElements);
132  // Sets elements and the corresponding indices of vertices
133  for (j = 0; j < NumOfElements; j++)
134  {
135  input >> a >> p1 >> p2;
136  elements[j] = new Segment(p1-1, p2-1, a);
137  }
138 
139  int ind[1];
140  input >> NumOfBdrElements;
141  boundary.SetSize(NumOfBdrElements);
142  for (j = 0; j < NumOfBdrElements; j++)
143  {
144  input >> a >> ind[0];
145  ind[0]--;
146  boundary[j] = new Point(ind,a);
147  }
148 }
149 
150 void Mesh::ReadNetgen2DMesh(std::istream &input, int &curved)
151 {
152  int ints[32], attr, n;
153 
154  // Read planar mesh in Netgen format.
155  Dim = 2;
156 
157  // Read the boundary elements.
158  input >> NumOfBdrElements;
159  boundary.SetSize(NumOfBdrElements);
160  for (int i = 0; i < NumOfBdrElements; i++)
161  {
162  input >> attr
163  >> ints[0] >> ints[1];
164  ints[0]--; ints[1]--;
165  boundary[i] = new Segment(ints, attr);
166  }
167 
168  // Read the elements.
169  input >> NumOfElements;
170  elements.SetSize(NumOfElements);
171  for (int i = 0; i < NumOfElements; i++)
172  {
173  input >> attr >> n;
174  for (int j = 0; j < n; j++)
175  {
176  input >> ints[j];
177  ints[j]--;
178  }
179  switch (n)
180  {
181  case 2:
182  elements[i] = new Segment(ints, attr);
183  break;
184  case 3:
185  elements[i] = new Triangle(ints, attr);
186  break;
187  case 4:
188  elements[i] = new Quadrilateral(ints, attr);
189  break;
190  }
191  }
192 
193  if (!curved)
194  {
195  // Read the vertices.
196  input >> NumOfVertices;
197  vertices.SetSize(NumOfVertices);
198  for (int i = 0; i < NumOfVertices; i++)
199  for (int j = 0; j < Dim; j++)
200  {
201  input >> vertices[i](j);
202  }
203  }
204  else
205  {
206  input >> NumOfVertices;
207  vertices.SetSize(NumOfVertices);
208  input >> ws;
209  }
210 }
211 
212 void Mesh::ReadNetgen3DMesh(std::istream &input)
213 {
214  int ints[32], attr;
215 
216  // Read a Netgen format mesh of tetrahedra.
217  Dim = 3;
218 
219  // Read the vertices
220  input >> NumOfVertices;
221 
222  vertices.SetSize(NumOfVertices);
223  for (int i = 0; i < NumOfVertices; i++)
224  for (int j = 0; j < Dim; j++)
225  {
226  input >> vertices[i](j);
227  }
228 
229  // Read the elements
230  input >> NumOfElements;
231  elements.SetSize(NumOfElements);
232  for (int i = 0; i < NumOfElements; i++)
233  {
234  input >> attr;
235  for (int j = 0; j < 4; j++)
236  {
237  input >> ints[j];
238  ints[j]--;
239  }
240 #ifdef MFEM_USE_MEMALLOC
241  Tetrahedron *tet;
242  tet = TetMemory.Alloc();
243  tet->SetVertices(ints);
244  tet->SetAttribute(attr);
245  elements[i] = tet;
246 #else
247  elements[i] = new Tetrahedron(ints, attr);
248 #endif
249  }
250 
251  // Read the boundary information.
252  input >> NumOfBdrElements;
253  boundary.SetSize(NumOfBdrElements);
254  for (int i = 0; i < NumOfBdrElements; i++)
255  {
256  input >> attr;
257  for (int j = 0; j < 3; j++)
258  {
259  input >> ints[j];
260  ints[j]--;
261  }
262  boundary[i] = new Triangle(ints, attr);
263  }
264 }
265 
266 void Mesh::ReadTrueGridMesh(std::istream &input)
267 {
268  int i, j, ints[32], attr;
269  const int buflen = 1024;
270  char buf[buflen];
271 
272  // TODO: find the actual dimension
273  Dim = 3;
274 
275  if (Dim == 2)
276  {
277  int vari;
278  double varf;
279 
280  input >> vari >> NumOfVertices >> vari >> vari >> NumOfElements;
281  input.getline(buf, buflen);
282  input.getline(buf, buflen);
283  input >> vari;
284  input.getline(buf, buflen);
285  input.getline(buf, buflen);
286  input.getline(buf, buflen);
287 
288  // Read the vertices.
289  vertices.SetSize(NumOfVertices);
290  for (i = 0; i < NumOfVertices; i++)
291  {
292  input >> vari >> varf >> vertices[i](0) >> vertices[i](1);
293  input.getline(buf, buflen);
294  }
295 
296  // Read the elements.
297  elements.SetSize(NumOfElements);
298  for (i = 0; i < NumOfElements; i++)
299  {
300  input >> vari >> attr;
301  for (j = 0; j < 4; j++)
302  {
303  input >> ints[j];
304  ints[j]--;
305  }
306  input.getline(buf, buflen);
307  input.getline(buf, buflen);
308  elements[i] = new Quadrilateral(ints, attr);
309  }
310  }
311  else if (Dim == 3)
312  {
313  int vari;
314  double varf;
315  input >> vari >> NumOfVertices >> NumOfElements;
316  input.getline(buf, buflen);
317  input.getline(buf, buflen);
318  input >> vari >> vari >> NumOfBdrElements;
319  input.getline(buf, buflen);
320  input.getline(buf, buflen);
321  input.getline(buf, buflen);
322  // Read the vertices.
323  vertices.SetSize(NumOfVertices);
324  for (i = 0; i < NumOfVertices; i++)
325  {
326  input >> vari >> varf >> vertices[i](0) >> vertices[i](1)
327  >> vertices[i](2);
328  input.getline(buf, buflen);
329  }
330  // Read the elements.
331  elements.SetSize(NumOfElements);
332  for (i = 0; i < NumOfElements; i++)
333  {
334  input >> vari >> attr;
335  for (j = 0; j < 8; j++)
336  {
337  input >> ints[j];
338  ints[j]--;
339  }
340  input.getline(buf, buflen);
341  elements[i] = new Hexahedron(ints, attr);
342  }
343  // Read the boundary elements.
344  boundary.SetSize(NumOfBdrElements);
345  for (i = 0; i < NumOfBdrElements; i++)
346  {
347  input >> attr;
348  for (j = 0; j < 4; j++)
349  {
350  input >> ints[j];
351  ints[j]--;
352  }
353  input.getline(buf, buflen);
354  boundary[i] = new Quadrilateral(ints, attr);
355  }
356  }
357 }
358 
359 // see Tetrahedron::edges
360 const int Mesh::vtk_quadratic_tet[10] =
361 { 0, 1, 2, 3, 4, 7, 5, 6, 8, 9 };
362 
363 // see Hexahedron::edges & Mesh::GenerateFaces
364 const int Mesh::vtk_quadratic_hex[27] =
365 {
366  0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19,
367  24, 22, 21, 23, 20, 25, 26
368 };
369 
370 void Mesh::ReadVTKMesh(std::istream &input, int &curved, int &read_gf)
371 {
372  int i, j, n, attr;
373 
374  string buff;
375  getline(input, buff); // comment line
376  getline(input, buff);
377  filter_dos(buff);
378  if (buff != "ASCII")
379  {
380  MFEM_ABORT("VTK mesh is not in ASCII format!");
381  return;
382  }
383  getline(input, buff);
384  filter_dos(buff);
385  if (buff != "DATASET UNSTRUCTURED_GRID")
386  {
387  MFEM_ABORT("VTK mesh is not UNSTRUCTURED_GRID!");
388  return;
389  }
390 
391  // Read the points, skipping optional sections such as the FIELD data from
392  // VisIt's VTK export (or from Mesh::PrintVTK with field_data==1).
393  do
394  {
395  input >> buff;
396  if (!input.good())
397  {
398  MFEM_ABORT("VTK mesh does not have POINTS data!");
399  }
400  }
401  while (buff != "POINTS");
402  int np = 0;
403  Vector points;
404  {
405  input >> np >> ws;
406  points.SetSize(3*np);
407  getline(input, buff); // "double"
408  for (i = 0; i < points.Size(); i++)
409  {
410  input >> points(i);
411  }
412  }
413 
414  // Read the cells
415  NumOfElements = n = 0;
416  Array<int> cells_data;
417  input >> ws >> buff;
418  if (buff == "CELLS")
419  {
420  input >> NumOfElements >> n >> ws;
421  cells_data.SetSize(n);
422  for (i = 0; i < n; i++)
423  {
424  input >> cells_data[i];
425  }
426  }
427 
428  // Read the cell types
429  Dim = 0;
430  int order = 1;
431  input >> ws >> buff;
432  if (buff == "CELL_TYPES")
433  {
434  input >> NumOfElements;
435  elements.SetSize(NumOfElements);
436  for (j = i = 0; i < NumOfElements; i++)
437  {
438  int ct;
439  input >> ct;
440  switch (ct)
441  {
442  case 5: // triangle
443  Dim = 2;
444  elements[i] = new Triangle(&cells_data[j+1]);
445  break;
446  case 9: // quadrilateral
447  Dim = 2;
448  elements[i] = new Quadrilateral(&cells_data[j+1]);
449  break;
450  case 10: // tetrahedron
451  Dim = 3;
452 #ifdef MFEM_USE_MEMALLOC
453  elements[i] = TetMemory.Alloc();
454  elements[i]->SetVertices(&cells_data[j+1]);
455 #else
456  elements[i] = new Tetrahedron(&cells_data[j+1]);
457 #endif
458  break;
459  case 12: // hexahedron
460  Dim = 3;
461  elements[i] = new Hexahedron(&cells_data[j+1]);
462  break;
463 
464  case 22: // quadratic triangle
465  Dim = 2;
466  order = 2;
467  elements[i] = new Triangle(&cells_data[j+1]);
468  break;
469  case 28: // biquadratic quadrilateral
470  Dim = 2;
471  order = 2;
472  elements[i] = new Quadrilateral(&cells_data[j+1]);
473  break;
474  case 24: // quadratic tetrahedron
475  Dim = 3;
476  order = 2;
477 #ifdef MFEM_USE_MEMALLOC
478  elements[i] = TetMemory.Alloc();
479  elements[i]->SetVertices(&cells_data[j+1]);
480 #else
481  elements[i] = new Tetrahedron(&cells_data[j+1]);
482 #endif
483  break;
484  case 29: // triquadratic hexahedron
485  Dim = 3;
486  order = 2;
487  elements[i] = new Hexahedron(&cells_data[j+1]);
488  break;
489  default:
490  MFEM_ABORT("VTK mesh : cell type " << ct << " is not supported!");
491  return;
492  }
493  j += cells_data[j] + 1;
494  }
495  }
496 
497  // Read attributes
498  streampos sp = input.tellg();
499  input >> ws >> buff;
500  if (buff == "CELL_DATA")
501  {
502  input >> n >> ws;
503  getline(input, buff);
504  filter_dos(buff);
505  // "SCALARS material dataType numComp"
506  if (!strncmp(buff.c_str(), "SCALARS material", 16))
507  {
508  getline(input, buff); // "LOOKUP_TABLE default"
509  for (i = 0; i < NumOfElements; i++)
510  {
511  input >> attr;
512  elements[i]->SetAttribute(attr);
513  }
514  }
515  else
516  {
517  input.seekg(sp);
518  }
519  }
520  else
521  {
522  input.seekg(sp);
523  }
524 
525  if (order == 1)
526  {
527  cells_data.DeleteAll();
528  NumOfVertices = np;
529  vertices.SetSize(np);
530  for (i = 0; i < np; i++)
531  {
532  vertices[i](0) = points(3*i+0);
533  vertices[i](1) = points(3*i+1);
534  vertices[i](2) = points(3*i+2);
535  }
536  points.Destroy();
537 
538  // No boundary is defined in a VTK mesh
539  NumOfBdrElements = 0;
540  }
541  else if (order == 2)
542  {
543  curved = 1;
544 
545  // generate new enumeration for the vertices
546  Array<int> pts_dof(np);
547  pts_dof = -1;
548  for (n = i = 0; i < NumOfElements; i++)
549  {
550  int *v = elements[i]->GetVertices();
551  int nv = elements[i]->GetNVertices();
552  for (j = 0; j < nv; j++)
553  if (pts_dof[v[j]] == -1)
554  {
555  pts_dof[v[j]] = n++;
556  }
557  }
558  // keep the original ordering of the vertices
559  for (n = i = 0; i < np; i++)
560  if (pts_dof[i] != -1)
561  {
562  pts_dof[i] = n++;
563  }
564  // update the element vertices
565  for (i = 0; i < NumOfElements; i++)
566  {
567  int *v = elements[i]->GetVertices();
568  int nv = elements[i]->GetNVertices();
569  for (j = 0; j < nv; j++)
570  {
571  v[j] = pts_dof[v[j]];
572  }
573  }
574  // Define the 'vertices' from the 'points' through the 'pts_dof' map
575  NumOfVertices = n;
576  vertices.SetSize(n);
577  for (i = 0; i < np; i++)
578  {
579  if ((j = pts_dof[i]) != -1)
580  {
581  vertices[j](0) = points(3*i+0);
582  vertices[j](1) = points(3*i+1);
583  vertices[j](2) = points(3*i+2);
584  }
585  }
586 
587  // No boundary is defined in a VTK mesh
588  NumOfBdrElements = 0;
589 
590  // Generate faces and edges so that we can define quadratic
591  // FE space on the mesh
592 
593  // Generate faces
594  if (Dim > 2)
595  {
596  GetElementToFaceTable();
597  GenerateFaces();
598  }
599  else
600  {
601  NumOfFaces = 0;
602  }
603 
604  // Generate edges
605  el_to_edge = new Table;
606  NumOfEdges = GetElementToEdgeTable(*el_to_edge, be_to_edge);
607  if (Dim == 2)
608  {
609  GenerateFaces(); // 'Faces' in 2D refers to the edges
610  }
611 
612  // Define quadratic FE space
614  FiniteElementSpace *fes = new FiniteElementSpace(this, fec, Dim);
615  Nodes = new GridFunction(fes);
616  Nodes->MakeOwner(fec); // Nodes will destroy 'fec' and 'fes'
617  own_nodes = 1;
618 
619  // Map vtk points to edge/face/element dofs
620  Array<int> dofs;
621  for (n = i = 0; i < NumOfElements; i++)
622  {
623  fes->GetElementDofs(i, dofs);
624  const int *vtk_mfem;
625  switch (elements[i]->GetGeometryType())
626  {
627  case Geometry::TRIANGLE:
628  case Geometry::SQUARE:
629  vtk_mfem = vtk_quadratic_hex; break; // identity map
630  case Geometry::TETRAHEDRON:
631  vtk_mfem = vtk_quadratic_tet; break;
632  case Geometry::CUBE:
633  default:
634  vtk_mfem = vtk_quadratic_hex; break;
635  }
636 
637  for (n++, j = 0; j < dofs.Size(); j++, n++)
638  {
639  if (pts_dof[cells_data[n]] == -1)
640  {
641  pts_dof[cells_data[n]] = dofs[vtk_mfem[j]];
642  }
643  else
644  {
645  if (pts_dof[cells_data[n]] != dofs[vtk_mfem[j]])
646  {
647  MFEM_ABORT("VTK mesh : inconsistent quadratic mesh!");
648  }
649  }
650  }
651  }
652 
653  // Define the 'Nodes' from the 'points' through the 'pts_dof' map
654  for (i = 0; i < np; i++)
655  {
656  dofs.SetSize(1);
657  if ((dofs[0] = pts_dof[i]) != -1)
658  {
659  fes->DofsToVDofs(dofs);
660  for (j = 0; j < dofs.Size(); j++)
661  {
662  (*Nodes)(dofs[j]) = points(3*i+j);
663  }
664  }
665  }
666 
667  read_gf = 0;
668  }
669 }
670 
671 void Mesh::ReadNURBSMesh(std::istream &input, int &curved, int &read_gf)
672 {
673  NURBSext = new NURBSExtension(input);
674 
675  Dim = NURBSext->Dimension();
676  NumOfVertices = NURBSext->GetNV();
677  NumOfElements = NURBSext->GetNE();
678  NumOfBdrElements = NURBSext->GetNBE();
679 
680  NURBSext->GetElementTopo(elements);
681  NURBSext->GetBdrElementTopo(boundary);
682 
683  vertices.SetSize(NumOfVertices);
684  curved = 1;
685  if (NURBSext->HavePatches())
686  {
687  NURBSFECollection *fec = new NURBSFECollection(NURBSext->GetOrder());
688  FiniteElementSpace *fes = new FiniteElementSpace(this, fec, Dim,
689  Ordering::byVDIM);
690  Nodes = new GridFunction(fes);
691  Nodes->MakeOwner(fec);
692  NURBSext->SetCoordsFromPatches(*Nodes);
693  own_nodes = 1;
694  read_gf = 0;
695  int vd = Nodes->VectorDim();
696  for (int i = 0; i < vd; i++)
697  {
698  Vector vert_val;
699  Nodes->GetNodalValues(vert_val, i+1);
700  for (int j = 0; j < NumOfVertices; j++)
701  {
702  vertices[j](i) = vert_val(j);
703  }
704  }
705  }
706  else
707  {
708  read_gf = 1;
709  }
710 }
711 
712 void Mesh::ReadInlineMesh(std::istream &input, int generate_edges)
713 {
714  // Initialize to negative numbers so that we know if they've been set. We're
715  // using Element::POINT as our flag, since we're not going to make a 0D mesh,
716  // ever.
717  int nx = -1;
718  int ny = -1;
719  int nz = -1;
720  double sx = -1.0;
721  double sy = -1.0;
722  double sz = -1.0;
723  Element::Type type = Element::POINT;
724 
725  while (true)
726  {
727  skip_comment_lines(input, '#');
728  // Break out if we reached the end of the file after gobbling up the
729  // whitespace and comments after the last keyword.
730  if (!input.good())
731  {
732  break;
733  }
734 
735  // Read the next keyword
736  std::string name;
737  input >> name;
738  input >> std::ws;
739  // Make sure there's an equal sign
740  MFEM_VERIFY(input.get() == '=',
741  "Inline mesh expected '=' after keyword " << name);
742  input >> std::ws;
743 
744  if (name == "nx")
745  {
746  input >> nx;
747  }
748  else if (name == "ny")
749  {
750  input >> ny;
751  }
752  else if (name == "nz")
753  {
754  input >> nz;
755  }
756  else if (name == "sx")
757  {
758  input >> sx;
759  }
760  else if (name == "sy")
761  {
762  input >> sy;
763  }
764  else if (name == "sz")
765  {
766  input >> sz;
767  }
768  else if (name == "type")
769  {
770  std::string eltype;
771  input >> eltype;
772  if (eltype == "segment")
773  {
774  type = Element::SEGMENT;
775  }
776  else if (eltype == "quad")
777  {
778  type = Element::QUADRILATERAL;
779  }
780  else if (eltype == "tri")
781  {
782  type = Element::TRIANGLE;
783  }
784  else if (eltype == "hex")
785  {
786  type = Element::HEXAHEDRON;
787  }
788  else if (eltype == "tet")
789  {
790  type = Element::TETRAHEDRON;
791  }
792  else
793  {
794  MFEM_ABORT("unrecognized element type (read '" << eltype
795  << "') in inline mesh format. "
796  "Allowed: segment, tri, tet, quad, hex");
797  }
798  }
799  else
800  {
801  MFEM_ABORT("unrecognized keyword (" << name
802  << ") in inline mesh format. "
803  "Allowed: nx, ny, nz, type, sx, sy, sz");
804  }
805 
806  input >> std::ws;
807  // Allow an optional semi-colon at the end of each line.
808  if (input.peek() == ';')
809  {
810  input.get();
811  }
812 
813  // Done reading file
814  if (!input)
815  {
816  break;
817  }
818  }
819 
820  // Now make the mesh.
821  if (type == Element::SEGMENT)
822  {
823  MFEM_VERIFY(nx > 0 && sx > 0.0,
824  "invalid 1D inline mesh format, all values must be "
825  "positive\n"
826  << " nx = " << nx << "\n"
827  << " sx = " << sx << "\n");
828  Make1D(nx, sx);
829  }
830  else if (type == Element::TRIANGLE || type == Element::QUADRILATERAL)
831  {
832  MFEM_VERIFY(nx > 0 && ny > 0 && sx > 0.0 && sy > 0.0,
833  "invalid 2D inline mesh format, all values must be "
834  "positive\n"
835  << " nx = " << nx << "\n"
836  << " ny = " << ny << "\n"
837  << " sx = " << sx << "\n"
838  << " sy = " << sy << "\n");
839  Make2D(nx, ny, type, generate_edges, sx, sy);
840  }
841  else if (type == Element::TETRAHEDRON || type == Element::HEXAHEDRON)
842  {
843  MFEM_VERIFY(nx > 0 && ny > 0 && nz > 0 &&
844  sx > 0.0 && sy > 0.0 && sz > 0.0,
845  "invalid 3D inline mesh format, all values must be "
846  "positive\n"
847  << " nx = " << nx << "\n"
848  << " ny = " << ny << "\n"
849  << " nz = " << nz << "\n"
850  << " sx = " << sx << "\n"
851  << " sy = " << sy << "\n"
852  << " sz = " << sz << "\n");
853  Make3D(nx, ny, nz, type, generate_edges, sx, sy, sz);
854  }
855  else
856  {
857  MFEM_ABORT("For inline mesh, must specify an element type ="
858  " [segment, tri, quad, tet, hex]");
859  }
860  InitBaseGeom();
861  return; // done with inline mesh construction
862 }
863 
864 void Mesh::ReadGmshMesh(std::istream &input)
865 {
866  string buff;
867  double version;
868  int binary, dsize;
869  input >> version >> binary >> dsize;
870  if (version < 2.2)
871  {
872  MFEM_ABORT("Gmsh file version < 2.2");
873  }
874  if (dsize != sizeof(double))
875  {
876  MFEM_ABORT("Gmsh file : dsize != sizeof(double)");
877  }
878  getline(input, buff);
879  // There is a number 1 in binary format
880  if (binary)
881  {
882  int one;
883  input.read(reinterpret_cast<char*>(&one), sizeof(one));
884  if (one != 1)
885  {
886  MFEM_ABORT("Gmsh file : wrong binary format");
887  }
888  }
889 
890  // A map between a serial number of the vertex and its number in the file
891  // (there may be gaps in the numbering, and also Gmsh enumerates vertices
892  // starting from 1, not 0)
893  map<int, int> vertices_map;
894  // Read the lines of the mesh file. If we face specific keyword, we'll treat
895  // the section.
896  while (input >> buff)
897  {
898  if (buff == "$Nodes") // reading mesh vertices
899  {
900  input >> NumOfVertices;
901  getline(input, buff);
902  vertices.SetSize(NumOfVertices);
903  int serial_number;
904  const int gmsh_dim = 3; // Gmsh always outputs 3 coordinates
905  double coord[gmsh_dim];
906  for (int ver = 0; ver < NumOfVertices; ++ver)
907  {
908  if (binary)
909  {
910  input.read(reinterpret_cast<char*>(&serial_number), sizeof(int));
911  input.read(reinterpret_cast<char*>(coord), gmsh_dim*sizeof(double));
912  }
913  else // ASCII
914  {
915  input >> serial_number;
916  for (int ci = 0; ci < gmsh_dim; ++ci)
917  {
918  input >> coord[ci];
919  }
920  }
921  vertices[ver] = Vertex(coord, gmsh_dim);
922  vertices_map[serial_number] = ver;
923  }
924  if (static_cast<int>(vertices_map.size()) != NumOfVertices)
925  {
926  MFEM_ABORT("Gmsh file : vertices indices are not unique");
927  }
928  } // section '$Nodes'
929  else if (buff == "$Elements") // reading mesh elements
930  {
931  int num_of_all_elements;
932  input >> num_of_all_elements;
933  // = NumOfElements + NumOfBdrElements + (maybe, PhysicalPoints)
934  getline(input, buff);
935 
936  int serial_number; // serial number of an element
937  int type_of_element; // ID describing a type of a mesh element
938  int n_tags; // number of different tags describing an element
939  int phys_domain; // element's attribute
940  int elem_domain; // another element's attribute (rarely used)
941  int n_partitions; // number of partitions where an element takes place
942 
943  // number of nodes for each type of Gmsh elements, type is the index of
944  // the array + 1
945  int nodes_of_gmsh_element[] =
946  {
947  2, // 2-node line.
948  3, // 3-node triangle.
949  4, // 4-node quadrangle.
950  4, // 4-node tetrahedron.
951  8, // 8-node hexahedron.
952  6, // 6-node prism.
953  5, // 5-node pyramid.
954  3, /* 3-node second order line (2 nodes associated with the vertices
955  and 1 with the edge). */
956  6, /* 6-node second order triangle (3 nodes associated with the
957  vertices and 3 with the edges). */
958  9, /* 9-node second order quadrangle (4 nodes associated with the
959  vertices, 4 with the edges and 1 with the face). */
960  10,/* 10-node second order tetrahedron (4 nodes associated with the
961  vertices and 6 with the edges). */
962  27,/* 27-node second order hexahedron (8 nodes associated with the
963  vertices, 12 with the edges, 6 with the faces and 1 with
964  the volume). */
965  18,/* 18-node second order prism (6 nodes associated with the
966  vertices, 9 with the edges and 3 with the quadrangular
967  faces). */
968  14,/* 14-node second order pyramid (5 nodes associated with the
969  vertices, 8 with the edges and 1 with the quadrangular
970  face). */
971  1, // 1-node point.
972  8, /* 8-node second order quadrangle (4 nodes associated with the
973  vertices and 4 with the edges). */
974  20,/* 20-node second order hexahedron (8 nodes associated with the
975  vertices and 12 with the edges). */
976  15,/* 15-node second order prism (6 nodes associated with the
977  vertices and 9 with the edges). */
978  13,/* 13-node second order pyramid (5 nodes associated with the
979  vertices and 8 with the edges). */
980  9, /* 9-node third order incomplete triangle (3 nodes associated
981  with the vertices, 6 with the edges) */
982  10,/* 10-node third order triangle (3 nodes associated with the
983  vertices, 6 with the edges, 1 with the face) */
984  12,/* 12-node fourth order incomplete triangle (3 nodes associated
985  with the vertices, 9 with the edges) */
986  15,/* 15-node fourth order triangle (3 nodes associated with the
987  vertices, 9 with the edges, 3 with the face) */
988  15,/* 15-node fifth order incomplete triangle (3 nodes associated
989  with the vertices, 12 with the edges) */
990  21,/* 21-node fifth order complete triangle (3 nodes associated with
991  the vertices, 12 with the edges, 6 with the face) */
992  4, /* 4-node third order edge (2 nodes associated with the vertices,
993  2 internal to the edge) */
994  5, /* 5-node fourth order edge (2 nodes associated with the
995  vertices, 3 internal to the edge) */
996  6, /* 6-node fifth order edge (2 nodes associated with the vertices,
997  4 internal to the edge) */
998  20 /* 20-node third order tetrahedron (4 nodes associated with the
999  vertices, 12 with the edges, 4 with the faces) */
1000  };
1001 
1002  vector<Element*> elements_0D, elements_1D, elements_2D, elements_3D;
1003  elements_0D.reserve(num_of_all_elements);
1004  elements_1D.reserve(num_of_all_elements);
1005  elements_2D.reserve(num_of_all_elements);
1006  elements_3D.reserve(num_of_all_elements);
1007 
1008  if (binary)
1009  {
1010  int n_elem_part = 0; // partial sum of elements that are read
1011  const int header_size = 3;
1012  // header consists of 3 numbers: type of the element, number of
1013  // elements of this type, and number of tags
1014  int header[header_size];
1015  int n_elem_one_type; // number of elements of a specific type
1016 
1017  while (n_elem_part < num_of_all_elements)
1018  {
1019  input.read(reinterpret_cast<char*>(header),
1020  header_size*sizeof(int));
1021  type_of_element = header[0];
1022  n_elem_one_type = header[1];
1023  n_tags = header[2];
1024 
1025  n_elem_part += n_elem_one_type;
1026 
1027  const int n_elem_nodes = nodes_of_gmsh_element[type_of_element-1];
1028  vector<int> data(1+n_tags+n_elem_nodes);
1029  for (int el = 0; el < n_elem_one_type; ++el)
1030  {
1031  input.read(reinterpret_cast<char*>(&data[0]),
1032  data.size()*sizeof(int));
1033  int dd = 0; // index for data array
1034  serial_number = data[dd++];
1035  // physical domain - the most important value (to distinguish
1036  // materials with different properties)
1037  phys_domain = (n_tags > 0) ? data[dd++] : 1;
1038  // elementary domain - to distinguish different geometrical
1039  // domains (typically, it's used rarely)
1040  elem_domain = (n_tags > 1) ? data[dd++] : 0;
1041  // the number of tags is bigger than 2 if there are some
1042  // partitions (domain decompositions)
1043  n_partitions = (n_tags > 2) ? data[dd++] : 0;
1044  // we currently just skip the partitions if they exist, and go
1045  // directly to vertices describing the mesh element
1046  vector<int> vert_indices(n_elem_nodes);
1047  for (int vi = 0; vi < n_elem_nodes; ++vi)
1048  {
1049  map<int, int>::const_iterator it =
1050  vertices_map.find(data[1+n_tags+vi]);
1051  if (it == vertices_map.end())
1052  {
1053  MFEM_ABORT("Gmsh file : vertex index doesn't exist");
1054  }
1055  vert_indices[vi] = it->second;
1056  }
1057 
1058  // non-positive attributes are not allowed in MFEM
1059  if (phys_domain <= 0)
1060  {
1061  MFEM_ABORT("Non-positive element attribute in Gmsh mesh!");
1062  }
1063 
1064  // initialize the mesh element
1065  switch (type_of_element)
1066  {
1067  case 1: // 2-node line
1068  {
1069  elements_1D.push_back(
1070  new Segment(&vert_indices[0], phys_domain));
1071  break;
1072  }
1073  case 2: // 3-node triangle
1074  {
1075  elements_2D.push_back(
1076  new Triangle(&vert_indices[0], phys_domain));
1077  break;
1078  }
1079  case 3: // 4-node quadrangle
1080  {
1081  elements_2D.push_back(
1082  new Quadrilateral(&vert_indices[0], phys_domain));
1083  break;
1084  }
1085  case 4: // 4-node tetrahedron
1086  {
1087  elements_3D.push_back(
1088  new Tetrahedron(&vert_indices[0], phys_domain));
1089  break;
1090  }
1091  case 5: // 8-node hexahedron
1092  {
1093  elements_3D.push_back(
1094  new Hexahedron(&vert_indices[0], phys_domain));
1095  break;
1096  }
1097  case 15: // 1-node point
1098  {
1099  elements_0D.push_back(
1100  new Point(&vert_indices[0], phys_domain));
1101  break;
1102  }
1103  default: // any other element
1104  MFEM_WARNING("Unsupported Gmsh element type.");
1105  break;
1106 
1107  } // switch (type_of_element)
1108  } // el (elements of one type)
1109  } // all elements
1110  } // if binary
1111  else // ASCII
1112  {
1113  for (int el = 0; el < num_of_all_elements; ++el)
1114  {
1115  input >> serial_number >> type_of_element >> n_tags;
1116  vector<int> data(n_tags);
1117  for (int i = 0; i < n_tags; ++i) { input >> data[i]; }
1118  // physical domain - the most important value (to distinguish
1119  // materials with different properties)
1120  phys_domain = (n_tags > 0) ? data[0] : 1;
1121  // elementary domain - to distinguish different geometrical
1122  // domains (typically, it's used rarely)
1123  elem_domain = (n_tags > 1) ? data[1] : 0;
1124  // the number of tags is bigger than 2 if there are some
1125  // partitions (domain decompositions)
1126  n_partitions = (n_tags > 2) ? data[2] : 0;
1127  // we currently just skip the partitions if they exist, and go
1128  // directly to vertices describing the mesh element
1129  const int n_elem_nodes = nodes_of_gmsh_element[type_of_element-1];
1130  vector<int> vert_indices(n_elem_nodes);
1131  int index;
1132  for (int vi = 0; vi < n_elem_nodes; ++vi)
1133  {
1134  input >> index;
1135  map<int, int>::const_iterator it = vertices_map.find(index);
1136  if (it == vertices_map.end())
1137  {
1138  MFEM_ABORT("Gmsh file : vertex index doesn't exist");
1139  }
1140  vert_indices[vi] = it->second;
1141  }
1142 
1143  // non-positive attributes are not allowed in MFEM
1144  if (phys_domain <= 0)
1145  {
1146  MFEM_ABORT("Non-positive element attribute in Gmsh mesh!");
1147  }
1148 
1149  // initialize the mesh element
1150  switch (type_of_element)
1151  {
1152  case 1: // 2-node line
1153  {
1154  elements_1D.push_back(
1155  new Segment(&vert_indices[0], phys_domain));
1156  break;
1157  }
1158  case 2: // 3-node triangle
1159  {
1160  elements_2D.push_back(
1161  new Triangle(&vert_indices[0], phys_domain));
1162  break;
1163  }
1164  case 3: // 4-node quadrangle
1165  {
1166  elements_2D.push_back(
1167  new Quadrilateral(&vert_indices[0], phys_domain));
1168  break;
1169  }
1170  case 4: // 4-node tetrahedron
1171  {
1172  elements_3D.push_back(
1173  new Tetrahedron(&vert_indices[0], phys_domain));
1174  break;
1175  }
1176  case 5: // 8-node hexahedron
1177  {
1178  elements_3D.push_back(
1179  new Hexahedron(&vert_indices[0], phys_domain));
1180  break;
1181  }
1182  case 15: // 1-node point
1183  {
1184  elements_0D.push_back(
1185  new Point(&vert_indices[0], phys_domain));
1186  break;
1187  }
1188  default: // any other element
1189  MFEM_WARNING("Unsupported Gmsh element type.");
1190  break;
1191 
1192  } // switch (type_of_element)
1193  } // el (all elements)
1194  } // if ASCII
1195 
1196  if (!elements_3D.empty())
1197  {
1198  Dim = 3;
1199  NumOfElements = elements_3D.size();
1200  elements.SetSize(NumOfElements);
1201  for (int el = 0; el < NumOfElements; ++el)
1202  {
1203  elements[el] = elements_3D[el];
1204  }
1205  NumOfBdrElements = elements_2D.size();
1206  boundary.SetSize(NumOfBdrElements);
1207  for (int el = 0; el < NumOfBdrElements; ++el)
1208  {
1209  boundary[el] = elements_2D[el];
1210  }
1211  // discard other elements
1212  for (size_t el = 0; el < elements_1D.size(); ++el)
1213  {
1214  delete elements_1D[el];
1215  }
1216  for (size_t el = 0; el < elements_0D.size(); ++el)
1217  {
1218  delete elements_0D[el];
1219  }
1220  }
1221  else if (!elements_2D.empty())
1222  {
1223  Dim = 2;
1224  NumOfElements = elements_2D.size();
1225  elements.SetSize(NumOfElements);
1226  for (int el = 0; el < NumOfElements; ++el)
1227  {
1228  elements[el] = elements_2D[el];
1229  }
1230  NumOfBdrElements = elements_1D.size();
1231  boundary.SetSize(NumOfBdrElements);
1232  for (int el = 0; el < NumOfBdrElements; ++el)
1233  {
1234  boundary[el] = elements_1D[el];
1235  }
1236  // discard other elements
1237  for (size_t el = 0; el < elements_0D.size(); ++el)
1238  {
1239  delete elements_0D[el];
1240  }
1241  }
1242  else if (!elements_1D.empty())
1243  {
1244  Dim = 1;
1245  NumOfElements = elements_1D.size();
1246  elements.SetSize(NumOfElements);
1247  for (int el = 0; el < NumOfElements; ++el)
1248  {
1249  elements[el] = elements_1D[el];
1250  }
1251  NumOfBdrElements = elements_0D.size();
1252  boundary.SetSize(NumOfBdrElements);
1253  for (int el = 0; el < NumOfBdrElements; ++el)
1254  {
1255  boundary[el] = elements_0D[el];
1256  }
1257  }
1258  else
1259  {
1260  MFEM_ABORT("Gmsh file : no elements found");
1261  return;
1262  }
1263 
1264  MFEM_CONTRACT_VAR(n_partitions);
1265  MFEM_CONTRACT_VAR(elem_domain);
1266 
1267  } // section '$Elements'
1268  } // we reach the end of the file
1269 }
1270 
1271 
1272 #ifdef MFEM_USE_NETCDF
1273 void Mesh::ReadCubit(const char *filename, int &curved, int &read_gf)
1274 {
1275  read_gf = 0;
1276 
1277  // curved set to zero will change if mesh is indeed curved
1278  curved = 0;
1279 
1280  const int sideMapHex8[6][4] =
1281  {
1282  {1,2,6,5},
1283  {2,3,7,6},
1284  {4,3,7,8},
1285  {1,4,8,5},
1286  {1,4,3,2},
1287  {5,8,7,6}
1288  };
1289 
1290  const int sideMapTet4[4][3] =
1291  {
1292  {1,2,4},
1293  {2,3,4},
1294  {1,4,3},
1295  {1,3,2}
1296  };
1297 
1298  // 1,2,3,4,5,6,7,8,9,10
1299  const int mfemToGenesisTet10[10] = {1,2,3,4,5,7,8,6,9,10};
1300 
1301  // 1,2,3,4,5,6,7,8,9,10,11,
1302  const int mfemToGenesisHex27[27] = {1,2,3,4,5,6,7,8,9,10,11,
1303  // 12,13,14,15,16,17,18,19
1304  12,17,18,19,20,13,14,15,
1305  // 20,21,22,23,24,25,26,27
1306  16,22,26,25,27,24,23,21
1307  };
1308 
1309  const int mfemToGenesisTri6[6] = {1,2,3,4,5,6};
1310  const int mfemToGenesisQuad9[9] = {1,2,3,4,5,6,7,8,9};
1311 
1312  const int sideMapHex27[6][9] =
1313  {
1314  {1,2,6,5,9,14,17,13,26},
1315  {2,3,7,6,10,15,18,14,25},
1316  {4,3,7,8,11,15,19,16,27},
1317  {1,4,8,5,12,16,20,13,24},
1318  {1,4,3,2,12,11,10,9,22},
1319  {5,8,7,6,20,19,18,17,23}
1320  };
1321 
1322  const int sideMapTet10[4][6] =
1323  {
1324  {1,2,4,5,9,8},
1325  {2,3,4,6,10,9},
1326  {1,4,3,8,10,7},
1327  {1,3,2,7,6,5}
1328  };
1329 
1330  // error handling.
1331  int retval;
1332 
1333  // dummy string
1334  char str_dummy[256];
1335 
1336  char temp_str[256];
1337  int temp_id;
1338 
1339  // open the file.
1340  int ncid;
1341  if ((retval = nc_open(filename, NC_NOWRITE, &ncid)))
1342  {
1343  MFEM_ABORT("Fatal NetCDF error: " << nc_strerror(retval));
1344  }
1345 
1346  // read important dimensions
1347 
1348  int id;
1349  size_t num_dim=0, num_nodes=0, num_elem=0, num_el_blk=0, num_side_sets=0;
1350 
1351  if ((retval = nc_inq_dimid(ncid, "num_dim", &id)) ||
1352  (retval = nc_inq_dim(ncid, id, str_dummy, &num_dim)) ||
1353 
1354  (retval = nc_inq_dimid(ncid, "num_nodes", &id)) ||
1355  (retval = nc_inq_dim(ncid, id, str_dummy, &num_nodes)) ||
1356 
1357  (retval = nc_inq_dimid(ncid, "num_elem", &id)) ||
1358  (retval = nc_inq_dim(ncid, id, str_dummy, &num_elem)) ||
1359 
1360  (retval = nc_inq_dimid(ncid, "num_el_blk", &id)) ||
1361  (retval = nc_inq_dim(ncid, id, str_dummy, &num_el_blk)) ||
1362 
1363  (retval = nc_inq_dimid(ncid, "num_side_sets", &id)) ||
1364  (retval = nc_inq_dim(ncid, id, str_dummy, &num_side_sets)))
1365  {
1366  MFEM_ABORT("Fatal NetCDF error: " << nc_strerror(retval));
1367  }
1368 
1369  Dim = num_dim;
1370 
1371  // create arrays for the next set of dimensions,
1372  // i.e. element blocks and side sets
1373  size_t *num_el_in_blk = new size_t[num_el_blk];
1374  size_t *num_nod_per_el = new size_t[num_el_blk];
1375  size_t *num_side_in_ss = new size_t[num_side_sets];
1376 
1377  int previous_num_nod_per_el = 0;
1378  for (int i = 0; i < (int) num_el_blk; i++)
1379  {
1380  sprintf(temp_str, "num_el_in_blk%d", i+1);
1381  if ((retval = nc_inq_dimid(ncid, temp_str, &temp_id)) ||
1382  (retval = nc_inq_dim(ncid, temp_id, str_dummy, &num_el_in_blk[i])))
1383  {
1384  MFEM_ABORT("Fatal NetCDF error: " << nc_strerror(retval));
1385  }
1386 
1387  sprintf(temp_str, "num_nod_per_el%d", i+1);
1388  if ((retval = nc_inq_dimid(ncid, temp_str, &temp_id)) ||
1389  (retval = nc_inq_dim(ncid, temp_id, str_dummy, &num_nod_per_el[i])))
1390  {
1391  MFEM_ABORT("Fatal NetCDF error: " << nc_strerror(retval));
1392  }
1393 
1394  // check for different element types in each block
1395  // which is not currently supported
1396  if (i != 0)
1397  {
1398  if ((int) num_nod_per_el[i] != previous_num_nod_per_el)
1399  {
1400  MFEM_ABORT("Element blocks of different element types not supported");
1401  }
1402  }
1403  previous_num_nod_per_el = num_nod_per_el[i];
1404  }
1405 
1406  int order = 0;
1407  if (num_nod_per_el[0] == 4 || num_nod_per_el[0] == 8)
1408  {
1409  order = 1;
1410  }
1411  else if (num_nod_per_el[0] == 10 || num_nod_per_el[0] == 27)
1412  {
1413  order = 2;
1414  }
1415  else
1416  {
1417  MFEM_ABORT("Don't know what to do with a " << num_nod_per_el[0] <<
1418  " node element\n");
1419  }
1420 
1421  for (int i = 0; i < (int) num_side_sets; i++)
1422  {
1423  sprintf(temp_str, "num_side_ss%d", i+1);
1424  if ((retval = nc_inq_dimid(ncid, temp_str, &temp_id)) ||
1425  (retval = nc_inq_dim(ncid, temp_id, str_dummy, &num_side_in_ss[i])))
1426  {
1427  MFEM_ABORT("Fatal NetCDF error: " << nc_strerror(retval));
1428  }
1429  }
1430 
1431  // read the coordinates
1432  double *coordx = new double[num_nodes];
1433  double *coordy = new double[num_nodes];
1434  double *coordz = new double[num_nodes];
1435 
1436  if ((retval = nc_inq_varid(ncid, "coordx", &id)) ||
1437  (retval = nc_get_var_double(ncid, id, coordx)) ||
1438  (retval = nc_inq_varid(ncid, "coordy", &id)) ||
1439  (retval = nc_get_var_double(ncid, id, coordy)))
1440  {
1441  MFEM_ABORT("Fatal NetCDF error: " << nc_strerror(retval));
1442  }
1443 
1444  if (num_dim == 3)
1445  {
1446  if ((retval = nc_inq_varid(ncid, "coordz", &id)) ||
1447  (retval = nc_get_var_double(ncid, id, coordz)))
1448  {
1449  MFEM_ABORT("Fatal NetCDF error: " << nc_strerror(retval));
1450  }
1451  }
1452 
1453  // read the element blocks
1454  int **elem_blk = new int*[num_el_blk];
1455  for (int i = 0; i < (int) num_el_blk; i++)
1456  {
1457  elem_blk[i] = new int[num_el_in_blk[i] * num_nod_per_el[i]];
1458  sprintf(temp_str, "connect%d", i+1);
1459  if ((retval = nc_inq_varid(ncid, temp_str, &temp_id)) ||
1460  (retval = nc_get_var_int(ncid, temp_id, elem_blk[i])))
1461  {
1462  MFEM_ABORT("Fatal NetCDF error: " << nc_strerror(retval));
1463  }
1464  }
1465  int *ebprop = new int[num_el_blk];
1466  if ((retval = nc_inq_varid(ncid, "eb_prop1", &id)) ||
1467  (retval = nc_get_var_int(ncid, id, ebprop)))
1468  {
1469  MFEM_ABORT("Fatal NetCDF error: " << nc_strerror(retval));
1470  }
1471 
1472  // read the side sets, a side is is given by (element, face) pairs
1473 
1474  int **elem_ss = new int*[num_side_sets];
1475  int **side_ss = new int*[num_side_sets];
1476 
1477  for (int i = 0; i < (int) num_side_sets; i++)
1478  {
1479  elem_ss[i] = new int[num_side_in_ss[i]];
1480  side_ss[i] = new int[num_side_in_ss[i]];
1481 
1482  sprintf(temp_str, "elem_ss%d", i+1);
1483  if ((retval = nc_inq_varid(ncid, temp_str, &temp_id)) ||
1484  (retval = nc_get_var_int(ncid, temp_id, elem_ss[i])))
1485  {
1486  MFEM_ABORT("Fatal NetCDF error: " << nc_strerror(retval));
1487  }
1488 
1489  sprintf(temp_str,"side_ss%d",i+1);
1490  if ((retval = nc_inq_varid(ncid, temp_str, &temp_id)) ||
1491  (retval = nc_get_var_int(ncid, temp_id, side_ss[i])))
1492  {
1493  MFEM_ABORT("Fatal NetCDF error: " << nc_strerror(retval));
1494  }
1495  }
1496 
1497  int *ssprop = new int[num_side_sets];
1498  if ((retval = nc_inq_varid(ncid, "ss_prop1", &id)) ||
1499  (retval = nc_get_var_int(ncid, id, ssprop)))
1500  {
1501  MFEM_ABORT("Fatal NetCDF error: " << nc_strerror(retval));
1502  }
1503 
1504  // convert (elem,side) pairs to 2D elements
1505 
1506  int num_nod_per_side = 0;
1507  if (num_nod_per_el[0] == 4)
1508  {
1509  num_nod_per_side = 3;
1510  }
1511  else if (num_nod_per_el[0] == 8)
1512  {
1513  num_nod_per_side = 4;
1514  }
1515  else if (num_nod_per_el[0] == 10)
1516  {
1517  num_nod_per_side = 6;
1518  }
1519  else if (num_nod_per_el[0] == 27)
1520  {
1521  num_nod_per_side = 9;
1522  }
1523 
1524  // given a global element number, determine the element block and local
1525  // element number
1526 
1527  int *start_of_block = new int[num_el_blk+1];
1528  start_of_block[0] = 0;
1529  for (int i = 1; i < (int) num_el_blk+1; i++)
1530  {
1531  start_of_block[i] = start_of_block[i-1] + num_el_in_blk[i-1];
1532  }
1533 
1534  int **ss_node_id = new int*[num_side_sets];
1535 
1536  for (int i = 0; i < (int) num_side_sets; i++)
1537  {
1538  ss_node_id[i] = new int[num_side_in_ss[i]*num_nod_per_side];
1539  for (int j = 0; j < (int) num_side_in_ss[i]; j++)
1540  {
1541  int glob_ind = elem_ss[i][j]-1;
1542  int iblk = 0;
1543  int loc_ind;
1544  while (iblk < (int) num_el_blk && glob_ind >= start_of_block[iblk+1])
1545  {
1546  iblk++;
1547  }
1548  if (iblk >= (int) num_el_blk)
1549  {
1550  MFEM_ABORT("Sideset element does not exist");
1551  }
1552  loc_ind = glob_ind - start_of_block[iblk];
1553  int this_side = side_ss[i][j];
1554  int ielem = loc_ind*num_nod_per_el[iblk];
1555  if (num_nod_per_side == 3)
1556  {
1557  for (int k = 0; k < num_nod_per_side; k++)
1558  {
1559  int inode = sideMapTet4[this_side-1][k];
1560  ss_node_id[i][j*num_nod_per_side+k] =
1561  elem_blk[iblk][ielem + inode - 1];
1562  }
1563  }
1564  else if (num_nod_per_side == 6)
1565  {
1566  for (int k = 0; k < num_nod_per_side; k++)
1567  {
1568  int inode = sideMapTet10[this_side-1][k];
1569  ss_node_id[i][j*num_nod_per_side+k] =
1570  elem_blk[iblk][ielem + inode - 1];
1571  }
1572  }
1573  else if (num_nod_per_side == 4)
1574  {
1575  for (int k = 0; k < num_nod_per_side; k++)
1576  {
1577  int inode = sideMapHex8[this_side-1][k];
1578  ss_node_id[i][j*num_nod_per_side+k] =
1579  elem_blk[iblk][ielem + inode - 1];
1580  }
1581  }
1582  else if (num_nod_per_side == 9)
1583  {
1584  for (int k = 0; k < num_nod_per_side; k++)
1585  {
1586  int inode = sideMapHex27[this_side-1][k];
1587  ss_node_id[i][j*num_nod_per_side+k] =
1588  elem_blk[iblk][ielem + inode - 1];
1589  }
1590  }
1591  }
1592  }
1593 
1594  // we need another node ID mapping since MFEM needs contiguous vertex IDs
1595 
1596  std::vector<int> uniqueVertexID;
1597  for (int iblk = 0; iblk < (int) num_el_blk; iblk++)
1598  {
1599  for (int i = 0; i < (int) num_el_in_blk[iblk]; i++)
1600  {
1601  int NumNodePerEl = num_nod_per_el[iblk];
1602  int NumVertPerEl = 0;
1603  if (NumNodePerEl == 8 || NumNodePerEl == 27)
1604  {
1605  NumVertPerEl = 8;
1606  }
1607  if (NumNodePerEl == 4 || NumNodePerEl == 10)
1608  {
1609  NumVertPerEl = 4;
1610  }
1611  for (int j = 0; j < NumVertPerEl; j++)
1612  {
1613  uniqueVertexID.push_back(elem_blk[iblk][i*NumNodePerEl + j]);
1614  }
1615  }
1616  }
1617  std::sort(uniqueVertexID.begin(), uniqueVertexID.end());
1618  std::vector<int>::iterator newEnd;
1619  newEnd = std::unique(uniqueVertexID.begin(), uniqueVertexID.end());
1620  uniqueVertexID.resize(std::distance(uniqueVertexID.begin(), newEnd));
1621 
1622  // OK at this point uniqueVertexID contains a list of all the nodes that are
1623  // actually used by the mesh, 1-based, and sorted. We need to invert this
1624  // list, the inverse is a map
1625 
1626  std::map<int,int> cubitToMFEMVertMap;
1627  for (int i = 0; i < (int) uniqueVertexID.size(); i++)
1628  {
1629  cubitToMFEMVertMap[uniqueVertexID[i]] = i+1;
1630  }
1631  MFEM_ASSERT(cubitToMFEMVertMap.size() == uniqueVertexID.size(),
1632  "This should never happen\n");
1633 
1634  // OK now load up the MFEM mesh structures
1635 
1636  // load up the vertices
1637 
1638  NumOfVertices = uniqueVertexID.size();
1639  vertices.SetSize(NumOfVertices);
1640  for (int i = 0; i < (int) uniqueVertexID.size(); i++)
1641  {
1642  vertices[i](0) = coordx[uniqueVertexID[i] - 1];
1643  vertices[i](1) = coordy[uniqueVertexID[i] - 1];
1644  vertices[i](2) = coordz[uniqueVertexID[i] - 1];
1645  }
1646 
1647  NumOfElements = num_elem;
1648  elements.SetSize(num_elem);
1649  int elcount = 0;
1650  int renumberedVertID[8];
1651  for (int iblk = 0; iblk < (int) num_el_blk; iblk++)
1652  {
1653  int NumNodePerEl = num_nod_per_el[iblk];
1654  if (NumNodePerEl == 4 || NumNodePerEl == 10)
1655  {
1656  for (int i = 0; i < (int) num_el_in_blk[iblk]; i++)
1657  {
1658  for (int j = 0; j < 4; j++)
1659  {
1660  renumberedVertID[j] =
1661  cubitToMFEMVertMap[elem_blk[iblk][i*NumNodePerEl+j]]-1;
1662  }
1663  elements[elcount] = new Tetrahedron(renumberedVertID,ebprop[iblk]);
1664  elcount++;
1665  }
1666  }
1667  else if (NumNodePerEl == 8 || NumNodePerEl == 27)
1668  {
1669  for (int i = 0; i < (int) num_el_in_blk[iblk]; i++)
1670  {
1671  for (int j = 0; j < 8; j++)
1672  {
1673  renumberedVertID[j] =
1674  cubitToMFEMVertMap[elem_blk[iblk][i*NumNodePerEl+j]]-1;
1675  }
1676  elements[elcount] = new Hexahedron(renumberedVertID,ebprop[iblk]);
1677  elcount++;
1678  }
1679  }
1680  }
1681 
1682  // load up the boundary elements
1683 
1684  NumOfBdrElements = 0;
1685  for (int iss = 0; iss < (int) num_side_sets; iss++)
1686  {
1687  NumOfBdrElements += num_side_in_ss[iss];
1688  }
1689  boundary.SetSize(NumOfBdrElements);
1690  int sidecount = 0;
1691  for (int iss = 0; iss < (int) num_side_sets; iss++)
1692  {
1693  if (num_nod_per_side == 3 || num_nod_per_side == 6)
1694  {
1695  for (int i = 0; i < (int) num_side_in_ss[iss]; i++)
1696  {
1697  for (int j = 0; j < 3; j++)
1698  {
1699  renumberedVertID[j] =
1700  cubitToMFEMVertMap[ss_node_id[iss][i*num_nod_per_side+j]] - 1;
1701  }
1702  boundary[sidecount] = new Triangle(renumberedVertID,ssprop[iss]);
1703  sidecount++;
1704  }
1705  }
1706  else if (num_nod_per_side == 4 || num_nod_per_side == 9)
1707  {
1708  for (int i = 0; i < (int) num_side_in_ss[iss]; i++)
1709  {
1710  for (int j = 0; j < 4; j++)
1711  {
1712  renumberedVertID[j] =
1713  cubitToMFEMVertMap[ss_node_id[iss][i*num_nod_per_side+j]] - 1;
1714  }
1715  boundary[sidecount] = new Quadrilateral(renumberedVertID,ssprop[iss]);
1716  sidecount++;
1717  }
1718  }
1719  }
1720 
1721  if (order == 2)
1722  {
1723  curved = 1;
1724  int *mymap = NULL;
1725 
1726  if (num_nod_per_el[0] == 10)
1727  {
1728  mymap = (int *) mfemToGenesisTet10;
1729  }
1730  else if (num_nod_per_el[0] == 27)
1731  {
1732  mymap = (int *) mfemToGenesisHex27;
1733  }
1734  else if (num_nod_per_el[0] == 6)
1735  {
1736  mymap = (int *) mfemToGenesisTri6;
1737  }
1738  else if (num_nod_per_el[0] == 9)
1739  {
1740  mymap = (int *) mfemToGenesisQuad9;
1741  }
1742 
1743  // Generate faces and edges so that we can define quadratic
1744  // FE space on the mesh
1745 
1746  // Generate faces
1747  if (Dim > 2)
1748  {
1749  GetElementToFaceTable();
1750  GenerateFaces();
1751  }
1752  else
1753  {
1754  NumOfFaces = 0;
1755  }
1756 
1757  // Generate edges
1758  el_to_edge = new Table;
1759  NumOfEdges = GetElementToEdgeTable(*el_to_edge, be_to_edge);
1760  if (Dim == 2)
1761  {
1762  GenerateFaces(); // 'Faces' in 2D refers to the edges
1763  }
1764 
1765  // Define quadratic FE space
1766  FiniteElementCollection *fec = new H1_FECollection(2,3);
1767  FiniteElementSpace *fes = new FiniteElementSpace(this, fec, Dim,
1768  Ordering::byVDIM);
1769  Nodes = new GridFunction(fes);
1770  Nodes->MakeOwner(fec); // Nodes will destroy 'fec' and 'fes'
1771  own_nodes = 1;
1772 
1773  // int nTotDofs = fes->GetNDofs();
1774  // int nTotVDofs = fes->GetVSize();
1775  // mfem::out << endl << "nTotDofs = " << nTotDofs << " nTotVDofs "
1776  // << nTotVDofs << endl << endl;
1777 
1778  for (int i = 0; i < NumOfElements; i++)
1779  {
1780  Array<int> dofs;
1781 
1782  fes->GetElementDofs(i, dofs);
1783  Array<int> vdofs;
1784  vdofs.SetSize(dofs.Size());
1785  for (int l = 0; l < dofs.Size(); l++) { vdofs[l] = dofs[l]; }
1786  fes->DofsToVDofs(vdofs);
1787  int iblk = 0;
1788  int loc_ind;
1789  while (iblk < (int) num_el_blk && i >= start_of_block[iblk+1]) { iblk++; }
1790  loc_ind = i - start_of_block[iblk];
1791 
1792  for (int j = 0; j < dofs.Size(); j++)
1793  {
1794  int point_id = elem_blk[iblk][loc_ind*num_nod_per_el[iblk] + mymap[j] - 1] - 1;
1795  (*Nodes)(vdofs[j]) = coordx[point_id];
1796  (*Nodes)(vdofs[j]+1) = coordy[point_id];
1797  (*Nodes)(vdofs[j]+2) = coordz[point_id];
1798  }
1799  }
1800  }
1801 
1802  // clean up all netcdf stuff
1803 
1804  nc_close(ncid);
1805 
1806  for (int i = 0; i < (int) num_side_sets; i++)
1807  {
1808  delete [] elem_ss[i];
1809  delete [] side_ss[i];
1810  }
1811 
1812  delete [] elem_ss;
1813  delete [] side_ss;
1814  delete [] num_el_in_blk;
1815  delete [] num_nod_per_el;
1816  delete [] num_side_in_ss;
1817  delete [] coordx;
1818  delete [] coordy;
1819  delete [] coordz;
1820 
1821  for (int i = 0; i < (int) num_el_blk; i++)
1822  {
1823  delete [] elem_blk[i];
1824  }
1825 
1826  delete [] elem_blk;
1827  delete [] start_of_block;
1828 
1829  for (int i = 0; i < (int) num_side_sets; i++)
1830  {
1831  delete [] ss_node_id[i];
1832  }
1833  delete [] ss_node_id;
1834  delete [] ebprop;
1835  delete [] ssprop;
1836 
1837 }
1838 #endif // #ifdef MFEM_USE_NETCDF
1839 
1840 } // namespace mfem
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Definition: fe_coll.hpp:266
int Size() const
Logical size of the array.
Definition: array.hpp:133
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:328
int Size() const
Returns the size of the vector.
Definition: vector.hpp:120
Data type for vertex.
Definition: vertex.hpp:21
virtual void SetVertices(const int *ind)
Set the vertices according to the given input.
Data type quadrilateral element.
void skip_comment_lines(std::istream &is, const char comment_char)
Definition: text.hpp:27
void DeleteAll()
Delete whole array.
Definition: array.hpp:709
Data type hexahedron element.
Definition: hexahedron.hpp:22
Type
Constants for the classes derived from Element.
Definition: element.hpp:37
virtual void GetElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:1233
Data type triangle element.
Definition: triangle.hpp:23
A class for non-conforming AMR on higher-order hexahedral, quadrilateral or triangular meshes...
Definition: ncmesh.hpp:83
Data type tetrahedron element.
Definition: tetrahedron.hpp:22
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:66
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:569
void Destroy()
Destroy a vector.
Definition: vector.hpp:347
void SetAttribute(const int attr)
Set element&#39;s attribute.
Definition: element.hpp:53
Piecewise-(bi)quadratic continuous finite elements.
Definition: fe_coll.hpp:342
void filter_dos(std::string &line)
Definition: text.hpp:41
Vector data type.
Definition: vector.hpp:48
Data type point element.
Definition: point.hpp:22
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:79
Data type line segment element.
Definition: segment.hpp:22
void DofsToVDofs(Array< int > &dofs, int ndofs=-1) const
Definition: fespace.cpp:106