MFEM  v3.3
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 
18 #ifdef MFEM_USE_NETCDF
19 #include "netcdf.h"
20 #endif
21 
22 using namespace std;
23 
24 namespace mfem
25 {
26 
27 void Mesh::ReadMFEMMesh(std::istream &input, bool mfem_v11, int &curved)
28 {
29  // Read MFEM mesh v1.0 format
30  string ident;
31 
32  // read lines beginning with '#' (comments)
33  skip_comment_lines(input, '#');
34  input >> ident; // 'dimension'
35 
36  MFEM_VERIFY(ident == "dimension", "invalid mesh file");
37  input >> Dim;
38 
39  skip_comment_lines(input, '#');
40  input >> ident; // 'elements'
41 
42  MFEM_VERIFY(ident == "elements", "invalid mesh file");
43  input >> NumOfElements;
44  elements.SetSize(NumOfElements);
45  for (int j = 0; j < NumOfElements; j++)
46  {
47  elements[j] = ReadElement(input);
48  }
49 
50  skip_comment_lines(input, '#');
51  input >> ident; // 'boundary'
52 
53  MFEM_VERIFY(ident == "boundary", "invalid mesh file");
54  input >> NumOfBdrElements;
55  boundary.SetSize(NumOfBdrElements);
56  for (int j = 0; j < NumOfBdrElements; j++)
57  {
58  boundary[j] = ReadElement(input);
59  }
60 
61  skip_comment_lines(input, '#');
62  input >> ident;
63 
64  if (mfem_v11 && ident == "vertex_parents")
65  {
66  ncmesh = new NCMesh(this, &input);
67  // NOTE: the constructor above will call LoadVertexParents
68 
69  skip_comment_lines(input, '#');
70  input >> ident;
71 
72  if (ident == "coarse_elements")
73  {
74  ncmesh->LoadCoarseElements(input);
75 
76  skip_comment_lines(input, '#');
77  input >> ident;
78  }
79  }
80 
81  MFEM_VERIFY(ident == "vertices", "invalid mesh file");
82  input >> NumOfVertices;
83  vertices.SetSize(NumOfVertices);
84 
85  input >> ws >> ident;
86  if (ident != "nodes")
87  {
88  // read the vertices
89  spaceDim = atoi(ident.c_str());
90  for (int j = 0; j < NumOfVertices; j++)
91  {
92  for (int i = 0; i < spaceDim; i++)
93  {
94  input >> vertices[j](i);
95  }
96  }
97 
98  // initialize vertex positions in NCMesh
99  if (ncmesh) { ncmesh->SetVertexPositions(vertices); }
100  }
101  else
102  {
103  // prepare to read the nodes
104  input >> ws;
105  curved = 1;
106  }
107 }
108 
109 void Mesh::ReadLineMesh(std::istream &input)
110 {
111  int j,p1,p2,a;
112 
113  Dim = 1;
114 
115  input >> NumOfVertices;
116  vertices.SetSize(NumOfVertices);
117  // Sets vertices and the corresponding coordinates
118  for (j = 0; j < NumOfVertices; j++)
119  {
120  input >> vertices[j](0);
121  }
122 
123  input >> NumOfElements;
124  elements.SetSize(NumOfElements);
125  // Sets elements and the corresponding indices of vertices
126  for (j = 0; j < NumOfElements; j++)
127  {
128  input >> a >> p1 >> p2;
129  elements[j] = new Segment(p1-1, p2-1, a);
130  }
131 
132  int ind[1];
133  input >> NumOfBdrElements;
134  boundary.SetSize(NumOfBdrElements);
135  for (j = 0; j < NumOfBdrElements; j++)
136  {
137  input >> a >> ind[0];
138  ind[0]--;
139  boundary[j] = new Point(ind,a);
140  }
141 }
142 
143 void Mesh::ReadNetgen2DMesh(std::istream &input, int &curved)
144 {
145  int ints[32], attr, n;
146 
147  // Read planar mesh in Netgen format.
148  Dim = 2;
149 
150  // Read the boundary elements.
151  input >> NumOfBdrElements;
152  boundary.SetSize(NumOfBdrElements);
153  for (int i = 0; i < NumOfBdrElements; i++)
154  {
155  input >> attr
156  >> ints[0] >> ints[1];
157  ints[0]--; ints[1]--;
158  boundary[i] = new Segment(ints, attr);
159  }
160 
161  // Read the elements.
162  input >> NumOfElements;
163  elements.SetSize(NumOfElements);
164  for (int i = 0; i < NumOfElements; i++)
165  {
166  input >> attr >> n;
167  for (int j = 0; j < n; j++)
168  {
169  input >> ints[j];
170  ints[j]--;
171  }
172  switch (n)
173  {
174  case 2:
175  elements[i] = new Segment(ints, attr);
176  break;
177  case 3:
178  elements[i] = new Triangle(ints, attr);
179  break;
180  case 4:
181  elements[i] = new Quadrilateral(ints, attr);
182  break;
183  }
184  }
185 
186  if (!curved)
187  {
188  // Read the vertices.
189  input >> NumOfVertices;
190  vertices.SetSize(NumOfVertices);
191  for (int i = 0; i < NumOfVertices; i++)
192  for (int j = 0; j < Dim; j++)
193  {
194  input >> vertices[i](j);
195  }
196  }
197  else
198  {
199  input >> NumOfVertices;
200  vertices.SetSize(NumOfVertices);
201  input >> ws;
202  }
203 }
204 
205 void Mesh::ReadNetgen3DMesh(std::istream &input)
206 {
207  int ints[32], attr;
208 
209  // Read a Netgen format mesh of tetrahedra.
210  Dim = 3;
211 
212  // Read the vertices
213  input >> NumOfVertices;
214 
215  vertices.SetSize(NumOfVertices);
216  for (int i = 0; i < NumOfVertices; i++)
217  for (int j = 0; j < Dim; j++)
218  {
219  input >> vertices[i](j);
220  }
221 
222  // Read the elements
223  input >> NumOfElements;
224  elements.SetSize(NumOfElements);
225  for (int i = 0; i < NumOfElements; i++)
226  {
227  input >> attr;
228  for (int j = 0; j < 4; j++)
229  {
230  input >> ints[j];
231  ints[j]--;
232  }
233 #ifdef MFEM_USE_MEMALLOC
234  Tetrahedron *tet;
235  tet = TetMemory.Alloc();
236  tet->SetVertices(ints);
237  tet->SetAttribute(attr);
238  elements[i] = tet;
239 #else
240  elements[i] = new Tetrahedron(ints, attr);
241 #endif
242  }
243 
244  // Read the boundary information.
245  input >> NumOfBdrElements;
246  boundary.SetSize(NumOfBdrElements);
247  for (int i = 0; i < NumOfBdrElements; i++)
248  {
249  input >> attr;
250  for (int j = 0; j < 3; j++)
251  {
252  input >> ints[j];
253  ints[j]--;
254  }
255  boundary[i] = new Triangle(ints, attr);
256  }
257 }
258 
259 void Mesh::ReadTrueGridMesh(std::istream &input)
260 {
261  int i, j, ints[32], attr;
262  const int buflen = 1024;
263  char buf[buflen];
264 
265  // TODO: find the actual dimension
266  Dim = 3;
267 
268  if (Dim == 2)
269  {
270  int vari;
271  double varf;
272 
273  input >> vari >> NumOfVertices >> vari >> vari >> NumOfElements;
274  input.getline(buf, buflen);
275  input.getline(buf, buflen);
276  input >> vari;
277  input.getline(buf, buflen);
278  input.getline(buf, buflen);
279  input.getline(buf, buflen);
280 
281  // Read the vertices.
282  vertices.SetSize(NumOfVertices);
283  for (i = 0; i < NumOfVertices; i++)
284  {
285  input >> vari >> varf >> vertices[i](0) >> vertices[i](1);
286  input.getline(buf, buflen);
287  }
288 
289  // Read the elements.
290  elements.SetSize(NumOfElements);
291  for (i = 0; i < NumOfElements; i++)
292  {
293  input >> vari >> attr;
294  for (j = 0; j < 4; j++)
295  {
296  input >> ints[j];
297  ints[j]--;
298  }
299  input.getline(buf, buflen);
300  input.getline(buf, buflen);
301  elements[i] = new Quadrilateral(ints, attr);
302  }
303  }
304  else if (Dim == 3)
305  {
306  int vari;
307  double varf;
308  input >> vari >> NumOfVertices >> NumOfElements;
309  input.getline(buf, buflen);
310  input.getline(buf, buflen);
311  input >> vari >> vari >> NumOfBdrElements;
312  input.getline(buf, buflen);
313  input.getline(buf, buflen);
314  input.getline(buf, buflen);
315  // Read the vertices.
316  vertices.SetSize(NumOfVertices);
317  for (i = 0; i < NumOfVertices; i++)
318  {
319  input >> vari >> varf >> vertices[i](0) >> vertices[i](1)
320  >> vertices[i](2);
321  input.getline(buf, buflen);
322  }
323  // Read the elements.
324  elements.SetSize(NumOfElements);
325  for (i = 0; i < NumOfElements; i++)
326  {
327  input >> vari >> attr;
328  for (j = 0; j < 8; j++)
329  {
330  input >> ints[j];
331  ints[j]--;
332  }
333  input.getline(buf, buflen);
334  elements[i] = new Hexahedron(ints, attr);
335  }
336  // Read the boundary elements.
337  boundary.SetSize(NumOfBdrElements);
338  for (i = 0; i < NumOfBdrElements; i++)
339  {
340  input >> attr;
341  for (j = 0; j < 4; j++)
342  {
343  input >> ints[j];
344  ints[j]--;
345  }
346  input.getline(buf, buflen);
347  boundary[i] = new Quadrilateral(ints, attr);
348  }
349  }
350 }
351 
352 // see Tetrahedron::edges
353 const int Mesh::vtk_quadratic_tet[10] =
354 { 0, 1, 2, 3, 4, 7, 5, 6, 8, 9 };
355 
356 // see Hexahedron::edges & Mesh::GenerateFaces
357 const int Mesh::vtk_quadratic_hex[27] =
358 {
359  0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19,
360  24, 22, 21, 23, 20, 25, 26
361 };
362 
363 void Mesh::ReadVTKMesh(std::istream &input, int &curved, int &read_gf)
364 {
365  int i, j, n, attr;
366 
367  string buff;
368  getline(input, buff); // comment line
369  getline(input, buff);
370  filter_dos(buff);
371  if (buff != "ASCII")
372  {
373  MFEM_ABORT("VTK mesh is not in ASCII format!");
374  return;
375  }
376  getline(input, buff);
377  filter_dos(buff);
378  if (buff != "DATASET UNSTRUCTURED_GRID")
379  {
380  MFEM_ABORT("VTK mesh is not UNSTRUCTURED_GRID!");
381  return;
382  }
383 
384  // Read the points, skipping optional sections such as the FIELD data from
385  // VisIt's VTK export (or from Mesh::PrintVTK with field_data==1).
386  do
387  {
388  input >> buff;
389  if (!input.good())
390  {
391  MFEM_ABORT("VTK mesh does not have POINTS data!");
392  }
393  }
394  while (buff != "POINTS");
395  int np = 0;
396  Vector points;
397  {
398  input >> np >> ws;
399  points.SetSize(3*np);
400  getline(input, buff); // "double"
401  for (i = 0; i < points.Size(); i++)
402  {
403  input >> points(i);
404  }
405  }
406 
407  // Read the cells
408  NumOfElements = n = 0;
409  Array<int> cells_data;
410  input >> ws >> buff;
411  if (buff == "CELLS")
412  {
413  input >> NumOfElements >> n >> ws;
414  cells_data.SetSize(n);
415  for (i = 0; i < n; i++)
416  {
417  input >> cells_data[i];
418  }
419  }
420 
421  // Read the cell types
422  Dim = 0;
423  int order = 1;
424  input >> ws >> buff;
425  if (buff == "CELL_TYPES")
426  {
427  input >> NumOfElements;
428  elements.SetSize(NumOfElements);
429  for (j = i = 0; i < NumOfElements; i++)
430  {
431  int ct;
432  input >> ct;
433  switch (ct)
434  {
435  case 5: // triangle
436  Dim = 2;
437  elements[i] = new Triangle(&cells_data[j+1]);
438  break;
439  case 9: // quadrilateral
440  Dim = 2;
441  elements[i] = new Quadrilateral(&cells_data[j+1]);
442  break;
443  case 10: // tetrahedron
444  Dim = 3;
445 #ifdef MFEM_USE_MEMALLOC
446  elements[i] = TetMemory.Alloc();
447  elements[i]->SetVertices(&cells_data[j+1]);
448 #else
449  elements[i] = new Tetrahedron(&cells_data[j+1]);
450 #endif
451  break;
452  case 12: // hexahedron
453  Dim = 3;
454  elements[i] = new Hexahedron(&cells_data[j+1]);
455  break;
456 
457  case 22: // quadratic triangle
458  Dim = 2;
459  order = 2;
460  elements[i] = new Triangle(&cells_data[j+1]);
461  break;
462  case 28: // biquadratic quadrilateral
463  Dim = 2;
464  order = 2;
465  elements[i] = new Quadrilateral(&cells_data[j+1]);
466  break;
467  case 24: // quadratic tetrahedron
468  Dim = 3;
469  order = 2;
470 #ifdef MFEM_USE_MEMALLOC
471  elements[i] = TetMemory.Alloc();
472  elements[i]->SetVertices(&cells_data[j+1]);
473 #else
474  elements[i] = new Tetrahedron(&cells_data[j+1]);
475 #endif
476  break;
477  case 29: // triquadratic hexahedron
478  Dim = 3;
479  order = 2;
480  elements[i] = new Hexahedron(&cells_data[j+1]);
481  break;
482  default:
483  MFEM_ABORT("VTK mesh : cell type " << ct << " is not supported!");
484  return;
485  }
486  j += cells_data[j] + 1;
487  }
488  }
489 
490  // Read attributes
491  streampos sp = input.tellg();
492  input >> ws >> buff;
493  if (buff == "CELL_DATA")
494  {
495  input >> n >> ws;
496  getline(input, buff);
497  filter_dos(buff);
498  // "SCALARS material dataType numComp"
499  if (!strncmp(buff.c_str(), "SCALARS material", 16))
500  {
501  getline(input, buff); // "LOOKUP_TABLE default"
502  for (i = 0; i < NumOfElements; i++)
503  {
504  input >> attr;
505  elements[i]->SetAttribute(attr);
506  }
507  }
508  else
509  {
510  input.seekg(sp);
511  }
512  }
513  else
514  {
515  input.seekg(sp);
516  }
517 
518  if (order == 1)
519  {
520  cells_data.DeleteAll();
521  NumOfVertices = np;
522  vertices.SetSize(np);
523  for (i = 0; i < np; i++)
524  {
525  vertices[i](0) = points(3*i+0);
526  vertices[i](1) = points(3*i+1);
527  vertices[i](2) = points(3*i+2);
528  }
529  points.Destroy();
530 
531  // No boundary is defined in a VTK mesh
532  NumOfBdrElements = 0;
533  }
534  else if (order == 2)
535  {
536  curved = 1;
537 
538  // generate new enumeration for the vertices
539  Array<int> pts_dof(np);
540  pts_dof = -1;
541  for (n = i = 0; i < NumOfElements; i++)
542  {
543  int *v = elements[i]->GetVertices();
544  int nv = elements[i]->GetNVertices();
545  for (j = 0; j < nv; j++)
546  if (pts_dof[v[j]] == -1)
547  {
548  pts_dof[v[j]] = n++;
549  }
550  }
551  // keep the original ordering of the vertices
552  for (n = i = 0; i < np; i++)
553  if (pts_dof[i] != -1)
554  {
555  pts_dof[i] = n++;
556  }
557  // update the element vertices
558  for (i = 0; i < NumOfElements; i++)
559  {
560  int *v = elements[i]->GetVertices();
561  int nv = elements[i]->GetNVertices();
562  for (j = 0; j < nv; j++)
563  {
564  v[j] = pts_dof[v[j]];
565  }
566  }
567  // Define the 'vertices' from the 'points' through the 'pts_dof' map
568  NumOfVertices = n;
569  vertices.SetSize(n);
570  for (i = 0; i < np; i++)
571  {
572  if ((j = pts_dof[i]) != -1)
573  {
574  vertices[j](0) = points(3*i+0);
575  vertices[j](1) = points(3*i+1);
576  vertices[j](2) = points(3*i+2);
577  }
578  }
579 
580  // No boundary is defined in a VTK mesh
581  NumOfBdrElements = 0;
582 
583  // Generate faces and edges so that we can define quadratic
584  // FE space on the mesh
585 
586  // Generate faces
587  if (Dim > 2)
588  {
589  GetElementToFaceTable();
590  GenerateFaces();
591  }
592  else
593  {
594  NumOfFaces = 0;
595  }
596 
597  // Generate edges
598  el_to_edge = new Table;
599  NumOfEdges = GetElementToEdgeTable(*el_to_edge, be_to_edge);
600  if (Dim == 2)
601  {
602  GenerateFaces(); // 'Faces' in 2D refers to the edges
603  }
604 
605  // Define quadratic FE space
607  FiniteElementSpace *fes = new FiniteElementSpace(this, fec, Dim);
608  Nodes = new GridFunction(fes);
609  Nodes->MakeOwner(fec); // Nodes will destroy 'fec' and 'fes'
610  own_nodes = 1;
611 
612  // Map vtk points to edge/face/element dofs
613  Array<int> dofs;
614  for (n = i = 0; i < NumOfElements; i++)
615  {
616  fes->GetElementDofs(i, dofs);
617  const int *vtk_mfem;
618  switch (elements[i]->GetGeometryType())
619  {
620  case Geometry::TRIANGLE:
621  case Geometry::SQUARE:
622  vtk_mfem = vtk_quadratic_hex; break; // identity map
623  case Geometry::TETRAHEDRON:
624  vtk_mfem = vtk_quadratic_tet; break;
625  case Geometry::CUBE:
626  default:
627  vtk_mfem = vtk_quadratic_hex; break;
628  }
629 
630  for (n++, j = 0; j < dofs.Size(); j++, n++)
631  {
632  if (pts_dof[cells_data[n]] == -1)
633  {
634  pts_dof[cells_data[n]] = dofs[vtk_mfem[j]];
635  }
636  else
637  {
638  if (pts_dof[cells_data[n]] != dofs[vtk_mfem[j]])
639  {
640  MFEM_ABORT("VTK mesh : inconsistent quadratic mesh!");
641  }
642  }
643  }
644  }
645 
646  // Define the 'Nodes' from the 'points' through the 'pts_dof' map
647  for (i = 0; i < np; i++)
648  {
649  dofs.SetSize(1);
650  if ((dofs[0] = pts_dof[i]) != -1)
651  {
652  fes->DofsToVDofs(dofs);
653  for (j = 0; j < dofs.Size(); j++)
654  {
655  (*Nodes)(dofs[j]) = points(3*i+j);
656  }
657  }
658  }
659 
660  read_gf = 0;
661  }
662 }
663 
664 void Mesh::ReadNURBSMesh(std::istream &input, int &curved, int &read_gf)
665 {
666  NURBSext = new NURBSExtension(input);
667 
668  Dim = NURBSext->Dimension();
669  NumOfVertices = NURBSext->GetNV();
670  NumOfElements = NURBSext->GetNE();
671  NumOfBdrElements = NURBSext->GetNBE();
672 
673  NURBSext->GetElementTopo(elements);
674  NURBSext->GetBdrElementTopo(boundary);
675 
676  vertices.SetSize(NumOfVertices);
677  curved = 1;
678  if (NURBSext->HavePatches())
679  {
680  NURBSFECollection *fec = new NURBSFECollection(NURBSext->GetOrder());
681  FiniteElementSpace *fes = new FiniteElementSpace(this, fec, Dim,
682  Ordering::byVDIM);
683  Nodes = new GridFunction(fes);
684  Nodes->MakeOwner(fec);
685  NURBSext->SetCoordsFromPatches(*Nodes);
686  own_nodes = 1;
687  read_gf = 0;
688  int vd = Nodes->VectorDim();
689  for (int i = 0; i < vd; i++)
690  {
691  Vector vert_val;
692  Nodes->GetNodalValues(vert_val, i+1);
693  for (int j = 0; j < NumOfVertices; j++)
694  {
695  vertices[j](i) = vert_val(j);
696  }
697  }
698  }
699  else
700  {
701  read_gf = 1;
702  }
703 }
704 
705 void Mesh::ReadInlineMesh(std::istream &input, int generate_edges)
706 {
707  // Initialize to negative numbers so that we know if they've been set. We're
708  // using Element::POINT as our flag, since we're not going to make a 0D mesh,
709  // ever.
710  int nx = -1;
711  int ny = -1;
712  int nz = -1;
713  double sx = -1.0;
714  double sy = -1.0;
715  double sz = -1.0;
716  Element::Type type = Element::POINT;
717 
718  while (true)
719  {
720  skip_comment_lines(input, '#');
721  // Break out if we reached the end of the file after gobbling up the
722  // whitespace and comments after the last keyword.
723  if (!input.good())
724  {
725  break;
726  }
727 
728  // Read the next keyword
729  std::string name;
730  input >> name;
731  input >> std::ws;
732  // Make sure there's an equal sign
733  MFEM_VERIFY(input.get() == '=',
734  "Inline mesh expected '=' after keyword " << name);
735  input >> std::ws;
736 
737  if (name == "nx")
738  {
739  input >> nx;
740  }
741  else if (name == "ny")
742  {
743  input >> ny;
744  }
745  else if (name == "nz")
746  {
747  input >> nz;
748  }
749  else if (name == "sx")
750  {
751  input >> sx;
752  }
753  else if (name == "sy")
754  {
755  input >> sy;
756  }
757  else if (name == "sz")
758  {
759  input >> sz;
760  }
761  else if (name == "type")
762  {
763  std::string eltype;
764  input >> eltype;
765  if (eltype == "segment")
766  {
767  type = Element::SEGMENT;
768  }
769  else if (eltype == "quad")
770  {
771  type = Element::QUADRILATERAL;
772  }
773  else if (eltype == "tri")
774  {
775  type = Element::TRIANGLE;
776  }
777  else if (eltype == "hex")
778  {
779  type = Element::HEXAHEDRON;
780  }
781  else if (eltype == "tet")
782  {
783  type = Element::TETRAHEDRON;
784  }
785  else
786  {
787  MFEM_ABORT("unrecognized element type (read '" << eltype
788  << "') in inline mesh format. "
789  "Allowed: segment, tri, tet, quad, hex");
790  }
791  }
792  else
793  {
794  MFEM_ABORT("unrecognized keyword (" << name
795  << ") in inline mesh format. "
796  "Allowed: nx, ny, nz, type, sx, sy, sz");
797  }
798 
799  input >> std::ws;
800  // Allow an optional semi-colon at the end of each line.
801  if (input.peek() == ';')
802  {
803  input.get();
804  }
805 
806  // Done reading file
807  if (!input)
808  {
809  break;
810  }
811  }
812 
813  // Now make the mesh.
814  if (type == Element::SEGMENT)
815  {
816  MFEM_VERIFY(nx > 0 && sx > 0.0,
817  "invalid 1D inline mesh format, all values must be "
818  "positive\n"
819  << " nx = " << nx << "\n"
820  << " sx = " << sx << "\n");
821  Make1D(nx, sx);
822  }
823  else if (type == Element::TRIANGLE || type == Element::QUADRILATERAL)
824  {
825  MFEM_VERIFY(nx > 0 && ny > 0 && sx > 0.0 && sy > 0.0,
826  "invalid 2D inline mesh format, all values must be "
827  "positive\n"
828  << " nx = " << nx << "\n"
829  << " ny = " << ny << "\n"
830  << " sx = " << sx << "\n"
831  << " sy = " << sy << "\n");
832  Make2D(nx, ny, type, generate_edges, sx, sy);
833  }
834  else if (type == Element::TETRAHEDRON || type == Element::HEXAHEDRON)
835  {
836  MFEM_VERIFY(nx > 0 && ny > 0 && nz > 0 &&
837  sx > 0.0 && sy > 0.0 && sz > 0.0,
838  "invalid 3D inline mesh format, all values must be "
839  "positive\n"
840  << " nx = " << nx << "\n"
841  << " ny = " << ny << "\n"
842  << " nz = " << nz << "\n"
843  << " sx = " << sx << "\n"
844  << " sy = " << sy << "\n"
845  << " sz = " << sz << "\n");
846  Make3D(nx, ny, nz, type, generate_edges, sx, sy, sz);
847  }
848  else
849  {
850  MFEM_ABORT("For inline mesh, must specify an element type ="
851  " [segment, tri, quad, tet, hex]");
852  }
853  InitBaseGeom();
854  return; // done with inline mesh construction
855 }
856 
857 void Mesh::ReadGmshMesh(std::istream &input)
858 {
859  string buff;
860  double version;
861  int binary, dsize;
862  input >> version >> binary >> dsize;
863  if (version < 2.2)
864  {
865  MFEM_ABORT("Gmsh file version < 2.2");
866  }
867  if (dsize != sizeof(double))
868  {
869  MFEM_ABORT("Gmsh file : dsize != sizeof(double)");
870  }
871  getline(input, buff);
872  // There is a number 1 in binary format
873  if (binary)
874  {
875  int one;
876  input.read(reinterpret_cast<char*>(&one), sizeof(one));
877  if (one != 1)
878  {
879  MFEM_ABORT("Gmsh file : wrong binary format");
880  }
881  }
882 
883  // A map between a serial number of the vertex and its number in the file
884  // (there may be gaps in the numbering, and also Gmsh enumerates vertices
885  // starting from 1, not 0)
886  map<int, int> vertices_map;
887  // Read the lines of the mesh file. If we face specific keyword, we'll treat
888  // the section.
889  while (input >> buff)
890  {
891  if (buff == "$Nodes") // reading mesh vertices
892  {
893  input >> NumOfVertices;
894  getline(input, buff);
895  vertices.SetSize(NumOfVertices);
896  int serial_number;
897  const int gmsh_dim = 3; // Gmsh always outputs 3 coordinates
898  double coord[gmsh_dim];
899  for (int ver = 0; ver < NumOfVertices; ++ver)
900  {
901  if (binary)
902  {
903  input.read(reinterpret_cast<char*>(&serial_number), sizeof(int));
904  input.read(reinterpret_cast<char*>(coord), gmsh_dim*sizeof(double));
905  }
906  else // ASCII
907  {
908  input >> serial_number;
909  for (int ci = 0; ci < gmsh_dim; ++ci)
910  {
911  input >> coord[ci];
912  }
913  }
914  vertices[ver] = Vertex(coord, gmsh_dim);
915  vertices_map[serial_number] = ver;
916  }
917  if (static_cast<int>(vertices_map.size()) != NumOfVertices)
918  {
919  MFEM_ABORT("Gmsh file : vertices indices are not unique");
920  }
921  } // section '$Nodes'
922  else if (buff == "$Elements") // reading mesh elements
923  {
924  int num_of_all_elements;
925  input >> num_of_all_elements;
926  // = NumOfElements + NumOfBdrElements + (maybe, PhysicalPoints)
927  getline(input, buff);
928 
929  int serial_number; // serial number of an element
930  int type_of_element; // ID describing a type of a mesh element
931  int n_tags; // number of different tags describing an element
932  int phys_domain; // element's attribute
933  int elem_domain; // another element's attribute (rarely used)
934  int n_partitions; // number of partitions where an element takes place
935 
936  // number of nodes for each type of Gmsh elements, type is the index of
937  // the array + 1
938  int nodes_of_gmsh_element[] =
939  {
940  2, // 2-node line.
941  3, // 3-node triangle.
942  4, // 4-node quadrangle.
943  4, // 4-node tetrahedron.
944  8, // 8-node hexahedron.
945  6, // 6-node prism.
946  5, // 5-node pyramid.
947  3, /* 3-node second order line (2 nodes associated with the vertices
948  and 1 with the edge). */
949  6, /* 6-node second order triangle (3 nodes associated with the
950  vertices and 3 with the edges). */
951  9, /* 9-node second order quadrangle (4 nodes associated with the
952  vertices, 4 with the edges and 1 with the face). */
953  10,/* 10-node second order tetrahedron (4 nodes associated with the
954  vertices and 6 with the edges). */
955  27,/* 27-node second order hexahedron (8 nodes associated with the
956  vertices, 12 with the edges, 6 with the faces and 1 with
957  the volume). */
958  18,/* 18-node second order prism (6 nodes associated with the
959  vertices, 9 with the edges and 3 with the quadrangular
960  faces). */
961  14,/* 14-node second order pyramid (5 nodes associated with the
962  vertices, 8 with the edges and 1 with the quadrangular
963  face). */
964  1, // 1-node point.
965  8, /* 8-node second order quadrangle (4 nodes associated with the
966  vertices and 4 with the edges). */
967  20,/* 20-node second order hexahedron (8 nodes associated with the
968  vertices and 12 with the edges). */
969  15,/* 15-node second order prism (6 nodes associated with the
970  vertices and 9 with the edges). */
971  13,/* 13-node second order pyramid (5 nodes associated with the
972  vertices and 8 with the edges). */
973  9, /* 9-node third order incomplete triangle (3 nodes associated
974  with the vertices, 6 with the edges) */
975  10,/* 10-node third order triangle (3 nodes associated with the
976  vertices, 6 with the edges, 1 with the face) */
977  12,/* 12-node fourth order incomplete triangle (3 nodes associated
978  with the vertices, 9 with the edges) */
979  15,/* 15-node fourth order triangle (3 nodes associated with the
980  vertices, 9 with the edges, 3 with the face) */
981  15,/* 15-node fifth order incomplete triangle (3 nodes associated
982  with the vertices, 12 with the edges) */
983  21,/* 21-node fifth order complete triangle (3 nodes associated with
984  the vertices, 12 with the edges, 6 with the face) */
985  4, /* 4-node third order edge (2 nodes associated with the vertices,
986  2 internal to the edge) */
987  5, /* 5-node fourth order edge (2 nodes associated with the
988  vertices, 3 internal to the edge) */
989  6, /* 6-node fifth order edge (2 nodes associated with the vertices,
990  4 internal to the edge) */
991  20 /* 20-node third order tetrahedron (4 nodes associated with the
992  vertices, 12 with the edges, 4 with the faces) */
993  };
994 
995  vector<Element*> elements_0D, elements_1D, elements_2D, elements_3D;
996  elements_0D.reserve(num_of_all_elements);
997  elements_1D.reserve(num_of_all_elements);
998  elements_2D.reserve(num_of_all_elements);
999  elements_3D.reserve(num_of_all_elements);
1000 
1001  if (binary)
1002  {
1003  int n_elem_part = 0; // partial sum of elements that are read
1004  const int header_size = 3;
1005  // header consists of 3 numbers: type of the element, number of
1006  // elements of this type, and number of tags
1007  int header[header_size];
1008  int n_elem_one_type; // number of elements of a specific type
1009 
1010  while (n_elem_part < num_of_all_elements)
1011  {
1012  input.read(reinterpret_cast<char*>(header),
1013  header_size*sizeof(int));
1014  type_of_element = header[0];
1015  n_elem_one_type = header[1];
1016  n_tags = header[2];
1017 
1018  n_elem_part += n_elem_one_type;
1019 
1020  const int n_elem_nodes = nodes_of_gmsh_element[type_of_element-1];
1021  vector<int> data(1+n_tags+n_elem_nodes);
1022  for (int el = 0; el < n_elem_one_type; ++el)
1023  {
1024  input.read(reinterpret_cast<char*>(&data[0]),
1025  data.size()*sizeof(int));
1026  int dd = 0; // index for data array
1027  serial_number = data[dd++];
1028  // physical domain - the most important value (to distinguish
1029  // materials with different properties)
1030  phys_domain = (n_tags > 0) ? data[dd++] : 0;
1031  // elementary domain - to distinguish different geometrical
1032  // domains (typically, it's used rarely)
1033  elem_domain = (n_tags > 1) ? data[dd++] : 0;
1034  // the number of tags is bigger than 2 if there are some
1035  // partitions (domain decompositions)
1036  n_partitions = (n_tags > 2) ? data[dd++] : 0;
1037  // we currently just skip the partitions if they exist, and go
1038  // directly to vertices describing the mesh element
1039  vector<int> vert_indices(n_elem_nodes);
1040  for (int vi = 0; vi < n_elem_nodes; ++vi)
1041  {
1042  map<int, int>::const_iterator it =
1043  vertices_map.find(data[1+n_tags+vi]);
1044  if (it == vertices_map.end())
1045  {
1046  MFEM_ABORT("Gmsh file : vertex index doesn't exist");
1047  }
1048  vert_indices[vi] = it->second;
1049  }
1050 
1051  // non-positive attributes are not allowed in MFEM
1052  if (phys_domain <= 0)
1053  {
1054  MFEM_ABORT("Non-positive element attribute in Gmsh mesh!");
1055  }
1056 
1057  // initialize the mesh element
1058  switch (type_of_element)
1059  {
1060  case 1: // 2-node line
1061  {
1062  elements_1D.push_back(
1063  new Segment(&vert_indices[0], phys_domain));
1064  break;
1065  }
1066  case 2: // 3-node triangle
1067  {
1068  elements_2D.push_back(
1069  new Triangle(&vert_indices[0], phys_domain));
1070  break;
1071  }
1072  case 3: // 4-node quadrangle
1073  {
1074  elements_2D.push_back(
1075  new Quadrilateral(&vert_indices[0], phys_domain));
1076  break;
1077  }
1078  case 4: // 4-node tetrahedron
1079  {
1080  elements_3D.push_back(
1081  new Tetrahedron(&vert_indices[0], phys_domain));
1082  break;
1083  }
1084  case 5: // 8-node hexahedron
1085  {
1086  elements_3D.push_back(
1087  new Hexahedron(&vert_indices[0], phys_domain));
1088  break;
1089  }
1090  case 15: // 1-node point
1091  {
1092  elements_0D.push_back(
1093  new Point(&vert_indices[0], phys_domain));
1094  break;
1095  }
1096  default: // any other element
1097  MFEM_WARNING("Unsupported Gmsh element type.");
1098  break;
1099 
1100  } // switch (type_of_element)
1101  } // el (elements of one type)
1102  } // all elements
1103  } // if binary
1104  else // ASCII
1105  {
1106  for (int el = 0; el < num_of_all_elements; ++el)
1107  {
1108  input >> serial_number >> type_of_element >> n_tags;
1109  vector<int> data(n_tags);
1110  for (int i = 0; i < n_tags; ++i) { input >> data[i]; }
1111  // physical domain - the most important value (to distinguish
1112  // materials with different properties)
1113  phys_domain = (n_tags > 0) ? data[0] : 0;
1114  // elementary domain - to distinguish different geometrical
1115  // domains (typically, it's used rarely)
1116  elem_domain = (n_tags > 1) ? data[1] : 0;
1117  // the number of tags is bigger than 2 if there are some
1118  // partitions (domain decompositions)
1119  n_partitions = (n_tags > 2) ? data[2] : 0;
1120  // we currently just skip the partitions if they exist, and go
1121  // directly to vertices describing the mesh element
1122  const int n_elem_nodes = nodes_of_gmsh_element[type_of_element-1];
1123  vector<int> vert_indices(n_elem_nodes);
1124  int index;
1125  for (int vi = 0; vi < n_elem_nodes; ++vi)
1126  {
1127  input >> index;
1128  map<int, int>::const_iterator it = vertices_map.find(index);
1129  if (it == vertices_map.end())
1130  {
1131  MFEM_ABORT("Gmsh file : vertex index doesn't exist");
1132  }
1133  vert_indices[vi] = it->second;
1134  }
1135 
1136  // non-positive attributes are not allowed in MFEM
1137  if (phys_domain <= 0)
1138  {
1139  MFEM_ABORT("Non-positive element attribute in Gmsh mesh!");
1140  }
1141 
1142  // initialize the mesh element
1143  switch (type_of_element)
1144  {
1145  case 1: // 2-node line
1146  {
1147  elements_1D.push_back(
1148  new Segment(&vert_indices[0], phys_domain));
1149  break;
1150  }
1151  case 2: // 3-node triangle
1152  {
1153  elements_2D.push_back(
1154  new Triangle(&vert_indices[0], phys_domain));
1155  break;
1156  }
1157  case 3: // 4-node quadrangle
1158  {
1159  elements_2D.push_back(
1160  new Quadrilateral(&vert_indices[0], phys_domain));
1161  break;
1162  }
1163  case 4: // 4-node tetrahedron
1164  {
1165  elements_3D.push_back(
1166  new Tetrahedron(&vert_indices[0], phys_domain));
1167  break;
1168  }
1169  case 5: // 8-node hexahedron
1170  {
1171  elements_3D.push_back(
1172  new Hexahedron(&vert_indices[0], phys_domain));
1173  break;
1174  }
1175  case 15: // 1-node point
1176  {
1177  elements_0D.push_back(
1178  new Point(&vert_indices[0], phys_domain));
1179  break;
1180  }
1181  default: // any other element
1182  MFEM_WARNING("Unsupported Gmsh element type.");
1183  break;
1184 
1185  } // switch (type_of_element)
1186  } // el (all elements)
1187  } // if ASCII
1188 
1189  if (!elements_3D.empty())
1190  {
1191  Dim = 3;
1192  NumOfElements = elements_3D.size();
1193  elements.SetSize(NumOfElements);
1194  for (int el = 0; el < NumOfElements; ++el)
1195  {
1196  elements[el] = elements_3D[el];
1197  }
1198  NumOfBdrElements = elements_2D.size();
1199  boundary.SetSize(NumOfBdrElements);
1200  for (int el = 0; el < NumOfBdrElements; ++el)
1201  {
1202  boundary[el] = elements_2D[el];
1203  }
1204  // discard other elements
1205  for (size_t el = 0; el < elements_1D.size(); ++el)
1206  {
1207  delete elements_1D[el];
1208  }
1209  for (size_t el = 0; el < elements_0D.size(); ++el)
1210  {
1211  delete elements_0D[el];
1212  }
1213  }
1214  else if (!elements_2D.empty())
1215  {
1216  Dim = 2;
1217  NumOfElements = elements_2D.size();
1218  elements.SetSize(NumOfElements);
1219  for (int el = 0; el < NumOfElements; ++el)
1220  {
1221  elements[el] = elements_2D[el];
1222  }
1223  NumOfBdrElements = elements_1D.size();
1224  boundary.SetSize(NumOfBdrElements);
1225  for (int el = 0; el < NumOfBdrElements; ++el)
1226  {
1227  boundary[el] = elements_1D[el];
1228  }
1229  // discard other elements
1230  for (size_t el = 0; el < elements_0D.size(); ++el)
1231  {
1232  delete elements_0D[el];
1233  }
1234  }
1235  else if (!elements_1D.empty())
1236  {
1237  Dim = 1;
1238  NumOfElements = elements_1D.size();
1239  elements.SetSize(NumOfElements);
1240  for (int el = 0; el < NumOfElements; ++el)
1241  {
1242  elements[el] = elements_1D[el];
1243  }
1244  NumOfBdrElements = elements_0D.size();
1245  boundary.SetSize(NumOfBdrElements);
1246  for (int el = 0; el < NumOfBdrElements; ++el)
1247  {
1248  boundary[el] = elements_0D[el];
1249  }
1250  }
1251  else
1252  {
1253  MFEM_ABORT("Gmsh file : no elements found");
1254  return;
1255  }
1256 
1257  MFEM_CONTRACT_VAR(n_partitions);
1258  MFEM_CONTRACT_VAR(elem_domain);
1259 
1260  } // section '$Elements'
1261  } // we reach the end of the file
1262 }
1263 
1264 
1265 #ifdef MFEM_USE_NETCDF
1266 void Mesh::ReadCubit(const char *filename, int &curved, int &read_gf)
1267 {
1268  read_gf = 0;
1269 
1270  // curved set to zero will change if mesh is indeed curved
1271  curved = 0;
1272 
1273  const int sideMapHex8[6][4] =
1274  {
1275  {1,2,6,5},
1276  {2,3,7,6},
1277  {4,3,7,8},
1278  {1,4,8,5},
1279  {1,4,3,2},
1280  {5,8,7,6}
1281  };
1282 
1283  const int sideMapTet4[4][3] =
1284  {
1285  {1,2,4},
1286  {2,3,4},
1287  {1,4,3},
1288  {1,3,2}
1289  };
1290 
1291  // 1,2,3,4,5,6,7,8,9,10
1292  const int mfemToGenesisTet10[10] = {1,2,3,4,5,7,8,6,9,10};
1293 
1294  // 1,2,3,4,5,6,7,8,9,10,11,
1295  const int mfemToGenesisHex27[27] = {1,2,3,4,5,6,7,8,9,10,11,
1296  // 12,13,14,15,16,17,18,19
1297  12,17,18,19,20,13,14,15,
1298  // 20,21,22,23,24,25,26,27
1299  16,22,26,25,27,24,23,21
1300  };
1301 
1302  const int mfemToGenesisTri6[6] = {1,2,3,4,5,6};
1303  const int mfemToGenesisQuad9[9] = {1,2,3,4,5,6,7,8,9};
1304 
1305  const int sideMapHex27[6][9] =
1306  {
1307  {1,2,6,5,9,14,17,13,26},
1308  {2,3,7,6,10,15,18,14,25},
1309  {4,3,7,8,11,15,19,16,27},
1310  {1,4,8,5,12,16,20,13,24},
1311  {1,4,3,2,12,11,10,9,22},
1312  {5,8,7,6,20,19,18,17,23}
1313  };
1314 
1315  const int sideMapTet10[4][6] =
1316  {
1317  {1,2,4,5,9,8},
1318  {2,3,4,6,10,9},
1319  {1,4,3,8,10,7},
1320  {1,3,2,7,6,5}
1321  };
1322 
1323  // error handling.
1324  int retval;
1325 
1326  // dummy string
1327  char str_dummy[256];
1328 
1329  char temp_str[256];
1330  int temp_id;
1331 
1332  // open the file.
1333  int ncid;
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 = 0;
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 = 0;
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 = 0;
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 = 0;
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 = NULL;
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:333
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 to size s.
Definition: vector.hpp:310
int Size() const
Returns the size of the vector.
Definition: vector.hpp:106
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:481
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:1035
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:80
Data type tetrahedron element.
Definition: tetrahedron.hpp:22
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:349
void Destroy()
Destroy a vector.
Definition: vector.hpp:329
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:36
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