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