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