13 #include "../fem/fem.hpp"
14 #include "../general/text.hpp"
19 #ifdef MFEM_USE_NETCDF
28 bool Mesh::remove_unused_vertices =
true;
30 void Mesh::ReadMFEMMesh(std::istream &input,
bool mfem_v11,
int &curved)
39 MFEM_VERIFY(ident ==
"dimension",
"invalid mesh file");
45 MFEM_VERIFY(ident ==
"elements",
"invalid mesh file");
46 input >> NumOfElements;
47 elements.SetSize(NumOfElements);
48 for (
int j = 0; j < NumOfElements; j++)
50 elements[j] = ReadElement(input);
56 MFEM_VERIFY(ident ==
"boundary",
"invalid mesh file");
57 input >> NumOfBdrElements;
58 boundary.SetSize(NumOfBdrElements);
59 for (
int j = 0; j < NumOfBdrElements; j++)
61 boundary[j] = ReadElement(input);
67 if (mfem_v11 && ident ==
"vertex_parents")
69 ncmesh =
new NCMesh(
this, &input);
75 if (ident ==
"coarse_elements")
77 ncmesh->LoadCoarseElements(input);
84 MFEM_VERIFY(ident ==
"vertices",
"invalid mesh file");
85 input >> NumOfVertices;
86 vertices.SetSize(NumOfVertices);
92 spaceDim = atoi(ident.c_str());
93 for (
int j = 0; j < NumOfVertices; j++)
95 for (
int i = 0; i < spaceDim; i++)
97 input >> vertices[j](i);
102 if (ncmesh) { ncmesh->SetVertexPositions(vertices); }
113 if (remove_unused_vertices) { RemoveUnusedVertices(); }
116 void Mesh::ReadLineMesh(std::istream &input)
122 input >> NumOfVertices;
123 vertices.SetSize(NumOfVertices);
125 for (j = 0; j < NumOfVertices; j++)
127 input >> vertices[j](0);
130 input >> NumOfElements;
131 elements.SetSize(NumOfElements);
133 for (j = 0; j < NumOfElements; j++)
135 input >> a >> p1 >> p2;
136 elements[j] =
new Segment(p1-1, p2-1, a);
140 input >> NumOfBdrElements;
141 boundary.SetSize(NumOfBdrElements);
142 for (j = 0; j < NumOfBdrElements; j++)
144 input >> a >> ind[0];
146 boundary[j] =
new Point(ind,a);
150 void Mesh::ReadNetgen2DMesh(std::istream &input,
int &curved)
152 int ints[32], attr, n;
158 input >> NumOfBdrElements;
159 boundary.SetSize(NumOfBdrElements);
160 for (
int i = 0; i < NumOfBdrElements; i++)
163 >> ints[0] >> ints[1];
164 ints[0]--; ints[1]--;
165 boundary[i] =
new Segment(ints, attr);
169 input >> NumOfElements;
170 elements.SetSize(NumOfElements);
171 for (
int i = 0; i < NumOfElements; i++)
174 for (
int j = 0; j < n; j++)
182 elements[i] =
new Segment(ints, attr);
185 elements[i] =
new Triangle(ints, attr);
196 input >> NumOfVertices;
197 vertices.SetSize(NumOfVertices);
198 for (
int i = 0; i < NumOfVertices; i++)
199 for (
int j = 0; j < Dim; j++)
201 input >> vertices[i](j);
206 input >> NumOfVertices;
207 vertices.SetSize(NumOfVertices);
212 void Mesh::ReadNetgen3DMesh(std::istream &input)
220 input >> NumOfVertices;
222 vertices.SetSize(NumOfVertices);
223 for (
int i = 0; i < NumOfVertices; i++)
224 for (
int j = 0; j < Dim; j++)
226 input >> vertices[i](j);
230 input >> NumOfElements;
231 elements.SetSize(NumOfElements);
232 for (
int i = 0; i < NumOfElements; i++)
235 for (
int j = 0; j < 4; j++)
240 #ifdef MFEM_USE_MEMALLOC
242 tet = TetMemory.Alloc();
252 input >> NumOfBdrElements;
253 boundary.SetSize(NumOfBdrElements);
254 for (
int i = 0; i < NumOfBdrElements; i++)
257 for (
int j = 0; j < 3; j++)
262 boundary[i] =
new Triangle(ints, attr);
266 void Mesh::ReadTrueGridMesh(std::istream &input)
268 int i, j, ints[32], attr;
269 const int buflen = 1024;
280 input >> vari >> NumOfVertices >> vari >> vari >> NumOfElements;
281 input.getline(buf, buflen);
282 input.getline(buf, buflen);
284 input.getline(buf, buflen);
285 input.getline(buf, buflen);
286 input.getline(buf, buflen);
289 vertices.SetSize(NumOfVertices);
290 for (i = 0; i < NumOfVertices; i++)
292 input >> vari >> varf >> vertices[i](0) >> vertices[i](1);
293 input.getline(buf, buflen);
297 elements.SetSize(NumOfElements);
298 for (i = 0; i < NumOfElements; i++)
300 input >> vari >> attr;
301 for (j = 0; j < 4; j++)
306 input.getline(buf, buflen);
307 input.getline(buf, buflen);
315 input >> vari >> NumOfVertices >> NumOfElements;
316 input.getline(buf, buflen);
317 input.getline(buf, buflen);
318 input >> vari >> vari >> NumOfBdrElements;
319 input.getline(buf, buflen);
320 input.getline(buf, buflen);
321 input.getline(buf, buflen);
323 vertices.SetSize(NumOfVertices);
324 for (i = 0; i < NumOfVertices; i++)
326 input >> vari >> varf >> vertices[i](0) >> vertices[i](1)
328 input.getline(buf, buflen);
331 elements.SetSize(NumOfElements);
332 for (i = 0; i < NumOfElements; i++)
334 input >> vari >> attr;
335 for (j = 0; j < 8; j++)
340 input.getline(buf, buflen);
344 boundary.SetSize(NumOfBdrElements);
345 for (i = 0; i < NumOfBdrElements; i++)
348 for (j = 0; j < 4; j++)
353 input.getline(buf, buflen);
360 const int Mesh::vtk_quadratic_tet[10] =
361 { 0, 1, 2, 3, 4, 7, 5, 6, 8, 9 };
364 const int Mesh::vtk_quadratic_hex[27] =
366 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19,
367 24, 22, 21, 23, 20, 25, 26
370 void Mesh::ReadVTKMesh(std::istream &input,
int &curved,
int &read_gf)
375 getline(input, buff);
376 getline(input, buff);
380 MFEM_ABORT(
"VTK mesh is not in ASCII format!");
383 getline(input, buff);
385 if (buff !=
"DATASET UNSTRUCTURED_GRID")
387 MFEM_ABORT(
"VTK mesh is not UNSTRUCTURED_GRID!");
398 MFEM_ABORT(
"VTK mesh does not have POINTS data!");
401 while (buff !=
"POINTS");
407 getline(input, buff);
408 for (i = 0; i < points.
Size(); i++)
415 NumOfElements = n = 0;
420 input >> NumOfElements >> n >> ws;
422 for (i = 0; i < n; i++)
424 input >> cells_data[i];
432 if (buff ==
"CELL_TYPES")
434 input >> NumOfElements;
435 elements.
SetSize(NumOfElements);
436 for (j = i = 0; i < NumOfElements; i++)
444 elements[i] =
new Triangle(&cells_data[j+1]);
452 #ifdef MFEM_USE_MEMALLOC
453 elements[i] = TetMemory.Alloc();
454 elements[i]->SetVertices(&cells_data[j+1]);
461 elements[i] =
new Hexahedron(&cells_data[j+1]);
467 elements[i] =
new Triangle(&cells_data[j+1]);
477 #ifdef MFEM_USE_MEMALLOC
478 elements[i] = TetMemory.Alloc();
479 elements[i]->SetVertices(&cells_data[j+1]);
487 elements[i] =
new Hexahedron(&cells_data[j+1]);
490 MFEM_ABORT(
"VTK mesh : cell type " << ct <<
" is not supported!");
493 j += cells_data[j] + 1;
498 streampos sp = input.tellg();
500 if (buff ==
"CELL_DATA")
503 getline(input, buff);
506 if (!strncmp(buff.c_str(),
"SCALARS material", 16))
508 getline(input, buff);
509 for (i = 0; i < NumOfElements; i++)
512 elements[i]->SetAttribute(attr);
529 vertices.SetSize(np);
530 for (i = 0; i < np; i++)
532 vertices[i](0) = points(3*i+0);
533 vertices[i](1) = points(3*i+1);
534 vertices[i](2) = points(3*i+2);
539 NumOfBdrElements = 0;
548 for (n = i = 0; i < NumOfElements; i++)
550 int *v = elements[i]->GetVertices();
551 int nv = elements[i]->GetNVertices();
552 for (j = 0; j < nv; j++)
553 if (pts_dof[v[j]] == -1)
559 for (n = i = 0; i < np; i++)
560 if (pts_dof[i] != -1)
565 for (i = 0; i < NumOfElements; i++)
567 int *v = elements[i]->GetVertices();
568 int nv = elements[i]->GetNVertices();
569 for (j = 0; j < nv; j++)
571 v[j] = pts_dof[v[j]];
577 for (i = 0; i < np; i++)
579 if ((j = pts_dof[i]) != -1)
581 vertices[j](0) = points(3*i+0);
582 vertices[j](1) = points(3*i+1);
583 vertices[j](2) = points(3*i+2);
588 NumOfBdrElements = 0;
596 GetElementToFaceTable();
605 el_to_edge =
new Table;
606 NumOfEdges = GetElementToEdgeTable(*el_to_edge, be_to_edge);
616 Nodes->MakeOwner(fec);
621 for (n = i = 0; i < NumOfElements; i++)
625 switch (elements[i]->GetGeometryType())
627 case Geometry::TRIANGLE:
628 case Geometry::SQUARE:
629 vtk_mfem = vtk_quadratic_hex;
break;
630 case Geometry::TETRAHEDRON:
631 vtk_mfem = vtk_quadratic_tet;
break;
634 vtk_mfem = vtk_quadratic_hex;
break;
637 for (n++, j = 0; j < dofs.
Size(); j++, n++)
639 if (pts_dof[cells_data[n]] == -1)
641 pts_dof[cells_data[n]] = dofs[vtk_mfem[j]];
645 if (pts_dof[cells_data[n]] != dofs[vtk_mfem[j]])
647 MFEM_ABORT(
"VTK mesh : inconsistent quadratic mesh!");
654 for (i = 0; i < np; i++)
657 if ((dofs[0] = pts_dof[i]) != -1)
660 for (j = 0; j < dofs.
Size(); j++)
662 (*Nodes)(dofs[j]) = points(3*i+j);
671 void Mesh::ReadNURBSMesh(std::istream &input,
int &curved,
int &read_gf)
675 Dim = NURBSext->Dimension();
676 NumOfVertices = NURBSext->GetNV();
677 NumOfElements = NURBSext->GetNE();
678 NumOfBdrElements = NURBSext->GetNBE();
680 NURBSext->GetElementTopo(elements);
681 NURBSext->GetBdrElementTopo(boundary);
683 vertices.SetSize(NumOfVertices);
685 if (NURBSext->HavePatches())
691 Nodes->MakeOwner(fec);
692 NURBSext->SetCoordsFromPatches(*Nodes);
695 int vd = Nodes->VectorDim();
696 for (
int i = 0; i < vd; i++)
699 Nodes->GetNodalValues(vert_val, i+1);
700 for (
int j = 0; j < NumOfVertices; j++)
702 vertices[j](i) = vert_val(j);
712 void Mesh::ReadInlineMesh(std::istream &input,
int generate_edges)
740 MFEM_VERIFY(input.get() ==
'=',
741 "Inline mesh expected '=' after keyword " << name);
748 else if (name ==
"ny")
752 else if (name ==
"nz")
756 else if (name ==
"sx")
760 else if (name ==
"sy")
764 else if (name ==
"sz")
768 else if (name ==
"type")
772 if (eltype ==
"segment")
774 type = Element::SEGMENT;
776 else if (eltype ==
"quad")
778 type = Element::QUADRILATERAL;
780 else if (eltype ==
"tri")
782 type = Element::TRIANGLE;
784 else if (eltype ==
"hex")
786 type = Element::HEXAHEDRON;
788 else if (eltype ==
"tet")
790 type = Element::TETRAHEDRON;
794 MFEM_ABORT(
"unrecognized element type (read '" << eltype
795 <<
"') in inline mesh format. "
796 "Allowed: segment, tri, tet, quad, hex");
801 MFEM_ABORT(
"unrecognized keyword (" << name
802 <<
") in inline mesh format. "
803 "Allowed: nx, ny, nz, type, sx, sy, sz");
808 if (input.peek() ==
';')
821 if (type == Element::SEGMENT)
823 MFEM_VERIFY(nx > 0 && sx > 0.0,
824 "invalid 1D inline mesh format, all values must be "
826 <<
" nx = " << nx <<
"\n"
827 <<
" sx = " << sx <<
"\n");
830 else if (type == Element::TRIANGLE || type == Element::QUADRILATERAL)
832 MFEM_VERIFY(nx > 0 && ny > 0 && sx > 0.0 && sy > 0.0,
833 "invalid 2D inline mesh format, all values must be "
835 <<
" nx = " << nx <<
"\n"
836 <<
" ny = " << ny <<
"\n"
837 <<
" sx = " << sx <<
"\n"
838 <<
" sy = " << sy <<
"\n");
839 Make2D(nx, ny, type, generate_edges, sx, sy);
841 else if (type == Element::TETRAHEDRON || type == Element::HEXAHEDRON)
843 MFEM_VERIFY(nx > 0 && ny > 0 && nz > 0 &&
844 sx > 0.0 && sy > 0.0 && sz > 0.0,
845 "invalid 3D inline mesh format, all values must be "
847 <<
" nx = " << nx <<
"\n"
848 <<
" ny = " << ny <<
"\n"
849 <<
" nz = " << nz <<
"\n"
850 <<
" sx = " << sx <<
"\n"
851 <<
" sy = " << sy <<
"\n"
852 <<
" sz = " << sz <<
"\n");
853 Make3D(nx, ny, nz, type, generate_edges, sx, sy, sz);
857 MFEM_ABORT(
"For inline mesh, must specify an element type ="
858 " [segment, tri, quad, tet, hex]");
864 void Mesh::ReadGmshMesh(std::istream &input)
869 input >> version >> binary >> dsize;
872 MFEM_ABORT(
"Gmsh file version < 2.2");
874 if (dsize !=
sizeof(
double))
876 MFEM_ABORT(
"Gmsh file : dsize != sizeof(double)");
878 getline(input, buff);
883 input.read(reinterpret_cast<char*>(&one),
sizeof(one));
886 MFEM_ABORT(
"Gmsh file : wrong binary format");
893 map<int, int> vertices_map;
896 while (input >> buff)
898 if (buff ==
"$Nodes")
900 input >> NumOfVertices;
901 getline(input, buff);
902 vertices.SetSize(NumOfVertices);
904 const int gmsh_dim = 3;
905 double coord[gmsh_dim];
906 for (
int ver = 0; ver < NumOfVertices; ++ver)
910 input.read(reinterpret_cast<char*>(&serial_number),
sizeof(
int));
911 input.read(reinterpret_cast<char*>(coord), gmsh_dim*
sizeof(
double));
915 input >> serial_number;
916 for (
int ci = 0; ci < gmsh_dim; ++ci)
921 vertices[ver] =
Vertex(coord, gmsh_dim);
922 vertices_map[serial_number] = ver;
924 if (static_cast<int>(vertices_map.size()) != NumOfVertices)
926 MFEM_ABORT(
"Gmsh file : vertices indices are not unique");
929 else if (buff ==
"$Elements")
931 int num_of_all_elements;
932 input >> num_of_all_elements;
934 getline(input, buff);
945 int nodes_of_gmsh_element[] =
1002 vector<Element*> elements_0D, elements_1D, elements_2D, elements_3D;
1003 elements_0D.reserve(num_of_all_elements);
1004 elements_1D.reserve(num_of_all_elements);
1005 elements_2D.reserve(num_of_all_elements);
1006 elements_3D.reserve(num_of_all_elements);
1010 int n_elem_part = 0;
1011 const int header_size = 3;
1014 int header[header_size];
1015 int n_elem_one_type;
1017 while (n_elem_part < num_of_all_elements)
1019 input.read(reinterpret_cast<char*>(header),
1020 header_size*
sizeof(
int));
1021 type_of_element = header[0];
1022 n_elem_one_type = header[1];
1025 n_elem_part += n_elem_one_type;
1027 const int n_elem_nodes = nodes_of_gmsh_element[type_of_element-1];
1028 vector<int> data(1+n_tags+n_elem_nodes);
1029 for (
int el = 0; el < n_elem_one_type; ++el)
1031 input.read(reinterpret_cast<char*>(&data[0]),
1032 data.size()*
sizeof(int));
1034 serial_number = data[dd++];
1037 phys_domain = (n_tags > 0) ? data[dd++] : 1;
1040 elem_domain = (n_tags > 1) ? data[dd++] : 0;
1043 n_partitions = (n_tags > 2) ? data[dd++] : 0;
1046 vector<int> vert_indices(n_elem_nodes);
1047 for (
int vi = 0; vi < n_elem_nodes; ++vi)
1049 map<int, int>::const_iterator it =
1050 vertices_map.find(data[1+n_tags+vi]);
1051 if (it == vertices_map.end())
1053 MFEM_ABORT(
"Gmsh file : vertex index doesn't exist");
1055 vert_indices[vi] = it->second;
1059 if (phys_domain <= 0)
1061 MFEM_ABORT(
"Non-positive element attribute in Gmsh mesh!");
1065 switch (type_of_element)
1069 elements_1D.push_back(
1070 new Segment(&vert_indices[0], phys_domain));
1075 elements_2D.push_back(
1076 new Triangle(&vert_indices[0], phys_domain));
1081 elements_2D.push_back(
1087 elements_3D.push_back(
1093 elements_3D.push_back(
1094 new Hexahedron(&vert_indices[0], phys_domain));
1099 elements_0D.push_back(
1100 new Point(&vert_indices[0], phys_domain));
1104 MFEM_WARNING(
"Unsupported Gmsh element type.");
1113 for (
int el = 0; el < num_of_all_elements; ++el)
1115 input >> serial_number >> type_of_element >> n_tags;
1116 vector<int> data(n_tags);
1117 for (
int i = 0; i < n_tags; ++i) { input >> data[i]; }
1120 phys_domain = (n_tags > 0) ? data[0] : 1;
1123 elem_domain = (n_tags > 1) ? data[1] : 0;
1126 n_partitions = (n_tags > 2) ? data[2] : 0;
1129 const int n_elem_nodes = nodes_of_gmsh_element[type_of_element-1];
1130 vector<int> vert_indices(n_elem_nodes);
1132 for (
int vi = 0; vi < n_elem_nodes; ++vi)
1135 map<int, int>::const_iterator it = vertices_map.find(index);
1136 if (it == vertices_map.end())
1138 MFEM_ABORT(
"Gmsh file : vertex index doesn't exist");
1140 vert_indices[vi] = it->second;
1144 if (phys_domain <= 0)
1146 MFEM_ABORT(
"Non-positive element attribute in Gmsh mesh!");
1150 switch (type_of_element)
1154 elements_1D.push_back(
1155 new Segment(&vert_indices[0], phys_domain));
1160 elements_2D.push_back(
1161 new Triangle(&vert_indices[0], phys_domain));
1166 elements_2D.push_back(
1172 elements_3D.push_back(
1178 elements_3D.push_back(
1179 new Hexahedron(&vert_indices[0], phys_domain));
1184 elements_0D.push_back(
1185 new Point(&vert_indices[0], phys_domain));
1189 MFEM_WARNING(
"Unsupported Gmsh element type.");
1196 if (!elements_3D.empty())
1199 NumOfElements = elements_3D.size();
1200 elements.SetSize(NumOfElements);
1201 for (
int el = 0; el < NumOfElements; ++el)
1203 elements[el] = elements_3D[el];
1205 NumOfBdrElements = elements_2D.size();
1206 boundary.SetSize(NumOfBdrElements);
1207 for (
int el = 0; el < NumOfBdrElements; ++el)
1209 boundary[el] = elements_2D[el];
1212 for (
size_t el = 0; el < elements_1D.size(); ++el)
1214 delete elements_1D[el];
1216 for (
size_t el = 0; el < elements_0D.size(); ++el)
1218 delete elements_0D[el];
1221 else if (!elements_2D.empty())
1224 NumOfElements = elements_2D.size();
1225 elements.SetSize(NumOfElements);
1226 for (
int el = 0; el < NumOfElements; ++el)
1228 elements[el] = elements_2D[el];
1230 NumOfBdrElements = elements_1D.size();
1231 boundary.SetSize(NumOfBdrElements);
1232 for (
int el = 0; el < NumOfBdrElements; ++el)
1234 boundary[el] = elements_1D[el];
1237 for (
size_t el = 0; el < elements_0D.size(); ++el)
1239 delete elements_0D[el];
1242 else if (!elements_1D.empty())
1245 NumOfElements = elements_1D.size();
1246 elements.SetSize(NumOfElements);
1247 for (
int el = 0; el < NumOfElements; ++el)
1249 elements[el] = elements_1D[el];
1251 NumOfBdrElements = elements_0D.size();
1252 boundary.SetSize(NumOfBdrElements);
1253 for (
int el = 0; el < NumOfBdrElements; ++el)
1255 boundary[el] = elements_0D[el];
1260 MFEM_ABORT(
"Gmsh file : no elements found");
1264 MFEM_CONTRACT_VAR(n_partitions);
1265 MFEM_CONTRACT_VAR(elem_domain);
1272 #ifdef MFEM_USE_NETCDF
1273 void Mesh::ReadCubit(
const char *filename,
int &curved,
int &read_gf)
1280 const int sideMapHex8[6][4] =
1290 const int sideMapTet4[4][3] =
1299 const int mfemToGenesisTet10[10] = {1,2,3,4,5,7,8,6,9,10};
1302 const int mfemToGenesisHex27[27] = {1,2,3,4,5,6,7,8,9,10,11,
1304 12,17,18,19,20,13,14,15,
1306 16,22,26,25,27,24,23,21
1309 const int mfemToGenesisTri6[6] = {1,2,3,4,5,6};
1310 const int mfemToGenesisQuad9[9] = {1,2,3,4,5,6,7,8,9};
1312 const int sideMapHex27[6][9] =
1314 {1,2,6,5,9,14,17,13,26},
1315 {2,3,7,6,10,15,18,14,25},
1316 {4,3,7,8,11,15,19,16,27},
1317 {1,4,8,5,12,16,20,13,24},
1318 {1,4,3,2,12,11,10,9,22},
1319 {5,8,7,6,20,19,18,17,23}
1322 const int sideMapTet10[4][6] =
1334 char str_dummy[256];
1341 if ((retval = nc_open(filename, NC_NOWRITE, &ncid)))
1343 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1349 size_t num_dim=0, num_nodes=0, num_elem=0, num_el_blk=0, num_side_sets=0;
1351 if ((retval = nc_inq_dimid(ncid,
"num_dim", &
id)) ||
1352 (retval = nc_inq_dim(ncid,
id, str_dummy, &num_dim)) ||
1354 (retval = nc_inq_dimid(ncid,
"num_nodes", &
id)) ||
1355 (retval = nc_inq_dim(ncid,
id, str_dummy, &num_nodes)) ||
1357 (retval = nc_inq_dimid(ncid,
"num_elem", &
id)) ||
1358 (retval = nc_inq_dim(ncid,
id, str_dummy, &num_elem)) ||
1360 (retval = nc_inq_dimid(ncid,
"num_el_blk", &
id)) ||
1361 (retval = nc_inq_dim(ncid,
id, str_dummy, &num_el_blk)) ||
1363 (retval = nc_inq_dimid(ncid,
"num_side_sets", &
id)) ||
1364 (retval = nc_inq_dim(ncid,
id, str_dummy, &num_side_sets)))
1366 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1373 size_t *num_el_in_blk =
new size_t[num_el_blk];
1374 size_t *num_nod_per_el =
new size_t[num_el_blk];
1375 size_t *num_side_in_ss =
new size_t[num_side_sets];
1377 int previous_num_nod_per_el = 0;
1378 for (
int i = 0; i < (int) num_el_blk; i++)
1380 sprintf(temp_str,
"num_el_in_blk%d", i+1);
1381 if ((retval = nc_inq_dimid(ncid, temp_str, &temp_id)) ||
1382 (retval = nc_inq_dim(ncid, temp_id, str_dummy, &num_el_in_blk[i])))
1384 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1387 sprintf(temp_str,
"num_nod_per_el%d", i+1);
1388 if ((retval = nc_inq_dimid(ncid, temp_str, &temp_id)) ||
1389 (retval = nc_inq_dim(ncid, temp_id, str_dummy, &num_nod_per_el[i])))
1391 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1398 if ((
int) num_nod_per_el[i] != previous_num_nod_per_el)
1400 MFEM_ABORT(
"Element blocks of different element types not supported");
1403 previous_num_nod_per_el = num_nod_per_el[i];
1407 if (num_nod_per_el[0] == 4 || num_nod_per_el[0] == 8)
1411 else if (num_nod_per_el[0] == 10 || num_nod_per_el[0] == 27)
1417 MFEM_ABORT(
"Don't know what to do with a " << num_nod_per_el[0] <<
1421 for (
int i = 0; i < (int) num_side_sets; i++)
1423 sprintf(temp_str,
"num_side_ss%d", i+1);
1424 if ((retval = nc_inq_dimid(ncid, temp_str, &temp_id)) ||
1425 (retval = nc_inq_dim(ncid, temp_id, str_dummy, &num_side_in_ss[i])))
1427 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1432 double *coordx =
new double[num_nodes];
1433 double *coordy =
new double[num_nodes];
1434 double *coordz =
new double[num_nodes];
1436 if ((retval = nc_inq_varid(ncid,
"coordx", &
id)) ||
1437 (retval = nc_get_var_double(ncid,
id, coordx)) ||
1438 (retval = nc_inq_varid(ncid,
"coordy", &
id)) ||
1439 (retval = nc_get_var_double(ncid,
id, coordy)))
1441 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1446 if ((retval = nc_inq_varid(ncid,
"coordz", &
id)) ||
1447 (retval = nc_get_var_double(ncid,
id, coordz)))
1449 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1454 int **elem_blk =
new int*[num_el_blk];
1455 for (
int i = 0; i < (int) num_el_blk; i++)
1457 elem_blk[i] =
new int[num_el_in_blk[i] * num_nod_per_el[i]];
1458 sprintf(temp_str,
"connect%d", i+1);
1459 if ((retval = nc_inq_varid(ncid, temp_str, &temp_id)) ||
1460 (retval = nc_get_var_int(ncid, temp_id, elem_blk[i])))
1462 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1465 int *ebprop =
new int[num_el_blk];
1466 if ((retval = nc_inq_varid(ncid,
"eb_prop1", &
id)) ||
1467 (retval = nc_get_var_int(ncid,
id, ebprop)))
1469 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1474 int **elem_ss =
new int*[num_side_sets];
1475 int **side_ss =
new int*[num_side_sets];
1477 for (
int i = 0; i < (int) num_side_sets; i++)
1479 elem_ss[i] =
new int[num_side_in_ss[i]];
1480 side_ss[i] =
new int[num_side_in_ss[i]];
1482 sprintf(temp_str,
"elem_ss%d", i+1);
1483 if ((retval = nc_inq_varid(ncid, temp_str, &temp_id)) ||
1484 (retval = nc_get_var_int(ncid, temp_id, elem_ss[i])))
1486 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1489 sprintf(temp_str,
"side_ss%d",i+1);
1490 if ((retval = nc_inq_varid(ncid, temp_str, &temp_id)) ||
1491 (retval = nc_get_var_int(ncid, temp_id, side_ss[i])))
1493 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1497 int *ssprop =
new int[num_side_sets];
1498 if ((retval = nc_inq_varid(ncid,
"ss_prop1", &
id)) ||
1499 (retval = nc_get_var_int(ncid,
id, ssprop)))
1501 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1506 int num_nod_per_side = 0;
1507 if (num_nod_per_el[0] == 4)
1509 num_nod_per_side = 3;
1511 else if (num_nod_per_el[0] == 8)
1513 num_nod_per_side = 4;
1515 else if (num_nod_per_el[0] == 10)
1517 num_nod_per_side = 6;
1519 else if (num_nod_per_el[0] == 27)
1521 num_nod_per_side = 9;
1527 int *start_of_block =
new int[num_el_blk+1];
1528 start_of_block[0] = 0;
1529 for (
int i = 1; i < (int) num_el_blk+1; i++)
1531 start_of_block[i] = start_of_block[i-1] + num_el_in_blk[i-1];
1534 int **ss_node_id =
new int*[num_side_sets];
1536 for (
int i = 0; i < (int) num_side_sets; i++)
1538 ss_node_id[i] =
new int[num_side_in_ss[i]*num_nod_per_side];
1539 for (
int j = 0; j < (int) num_side_in_ss[i]; j++)
1541 int glob_ind = elem_ss[i][j]-1;
1544 while (iblk < (
int) num_el_blk && glob_ind >= start_of_block[iblk+1])
1548 if (iblk >= (
int) num_el_blk)
1550 MFEM_ABORT(
"Sideset element does not exist");
1552 loc_ind = glob_ind - start_of_block[iblk];
1553 int this_side = side_ss[i][j];
1554 int ielem = loc_ind*num_nod_per_el[iblk];
1555 if (num_nod_per_side == 3)
1557 for (
int k = 0; k < num_nod_per_side; k++)
1559 int inode = sideMapTet4[this_side-1][k];
1560 ss_node_id[i][j*num_nod_per_side+k] =
1561 elem_blk[iblk][ielem + inode - 1];
1564 else if (num_nod_per_side == 6)
1566 for (
int k = 0; k < num_nod_per_side; k++)
1568 int inode = sideMapTet10[this_side-1][k];
1569 ss_node_id[i][j*num_nod_per_side+k] =
1570 elem_blk[iblk][ielem + inode - 1];
1573 else if (num_nod_per_side == 4)
1575 for (
int k = 0; k < num_nod_per_side; k++)
1577 int inode = sideMapHex8[this_side-1][k];
1578 ss_node_id[i][j*num_nod_per_side+k] =
1579 elem_blk[iblk][ielem + inode - 1];
1582 else if (num_nod_per_side == 9)
1584 for (
int k = 0; k < num_nod_per_side; k++)
1586 int inode = sideMapHex27[this_side-1][k];
1587 ss_node_id[i][j*num_nod_per_side+k] =
1588 elem_blk[iblk][ielem + inode - 1];
1596 std::vector<int> uniqueVertexID;
1597 for (
int iblk = 0; iblk < (int) num_el_blk; iblk++)
1599 for (
int i = 0; i < (int) num_el_in_blk[iblk]; i++)
1601 int NumNodePerEl = num_nod_per_el[iblk];
1602 int NumVertPerEl = 0;
1603 if (NumNodePerEl == 8 || NumNodePerEl == 27)
1607 if (NumNodePerEl == 4 || NumNodePerEl == 10)
1611 for (
int j = 0; j < NumVertPerEl; j++)
1613 uniqueVertexID.push_back(elem_blk[iblk][i*NumNodePerEl + j]);
1617 std::sort(uniqueVertexID.begin(), uniqueVertexID.end());
1618 std::vector<int>::iterator newEnd;
1619 newEnd = std::unique(uniqueVertexID.begin(), uniqueVertexID.end());
1620 uniqueVertexID.resize(std::distance(uniqueVertexID.begin(), newEnd));
1626 std::map<int,int> cubitToMFEMVertMap;
1627 for (
int i = 0; i < (int) uniqueVertexID.size(); i++)
1629 cubitToMFEMVertMap[uniqueVertexID[i]] = i+1;
1631 MFEM_ASSERT(cubitToMFEMVertMap.size() == uniqueVertexID.size(),
1632 "This should never happen\n");
1638 NumOfVertices = uniqueVertexID.size();
1639 vertices.SetSize(NumOfVertices);
1640 for (
int i = 0; i < (int) uniqueVertexID.size(); i++)
1642 vertices[i](0) = coordx[uniqueVertexID[i] - 1];
1643 vertices[i](1) = coordy[uniqueVertexID[i] - 1];
1644 vertices[i](2) = coordz[uniqueVertexID[i] - 1];
1647 NumOfElements = num_elem;
1648 elements.SetSize(num_elem);
1650 int renumberedVertID[8];
1651 for (
int iblk = 0; iblk < (int) num_el_blk; iblk++)
1653 int NumNodePerEl = num_nod_per_el[iblk];
1654 if (NumNodePerEl == 4 || NumNodePerEl == 10)
1656 for (
int i = 0; i < (int) num_el_in_blk[iblk]; i++)
1658 for (
int j = 0; j < 4; j++)
1660 renumberedVertID[j] =
1661 cubitToMFEMVertMap[elem_blk[iblk][i*NumNodePerEl+j]]-1;
1663 elements[elcount] =
new Tetrahedron(renumberedVertID,ebprop[iblk]);
1667 else if (NumNodePerEl == 8 || NumNodePerEl == 27)
1669 for (
int i = 0; i < (int) num_el_in_blk[iblk]; i++)
1671 for (
int j = 0; j < 8; j++)
1673 renumberedVertID[j] =
1674 cubitToMFEMVertMap[elem_blk[iblk][i*NumNodePerEl+j]]-1;
1676 elements[elcount] =
new Hexahedron(renumberedVertID,ebprop[iblk]);
1684 NumOfBdrElements = 0;
1685 for (
int iss = 0; iss < (int) num_side_sets; iss++)
1687 NumOfBdrElements += num_side_in_ss[iss];
1689 boundary.SetSize(NumOfBdrElements);
1691 for (
int iss = 0; iss < (int) num_side_sets; iss++)
1693 if (num_nod_per_side == 3 || num_nod_per_side == 6)
1695 for (
int i = 0; i < (int) num_side_in_ss[iss]; i++)
1697 for (
int j = 0; j < 3; j++)
1699 renumberedVertID[j] =
1700 cubitToMFEMVertMap[ss_node_id[iss][i*num_nod_per_side+j]] - 1;
1702 boundary[sidecount] =
new Triangle(renumberedVertID,ssprop[iss]);
1706 else if (num_nod_per_side == 4 || num_nod_per_side == 9)
1708 for (
int i = 0; i < (int) num_side_in_ss[iss]; i++)
1710 for (
int j = 0; j < 4; j++)
1712 renumberedVertID[j] =
1713 cubitToMFEMVertMap[ss_node_id[iss][i*num_nod_per_side+j]] - 1;
1715 boundary[sidecount] =
new Quadrilateral(renumberedVertID,ssprop[iss]);
1726 if (num_nod_per_el[0] == 10)
1728 mymap = (
int *) mfemToGenesisTet10;
1730 else if (num_nod_per_el[0] == 27)
1732 mymap = (
int *) mfemToGenesisHex27;
1734 else if (num_nod_per_el[0] == 6)
1736 mymap = (
int *) mfemToGenesisTri6;
1738 else if (num_nod_per_el[0] == 9)
1740 mymap = (
int *) mfemToGenesisQuad9;
1749 GetElementToFaceTable();
1758 el_to_edge =
new Table;
1759 NumOfEdges = GetElementToEdgeTable(*el_to_edge, be_to_edge);
1770 Nodes->MakeOwner(fec);
1778 for (
int i = 0; i < NumOfElements; i++)
1785 for (
int l = 0; l < dofs.
Size(); l++) { vdofs[l] = dofs[l]; }
1789 while (iblk < (
int) num_el_blk && i >= start_of_block[iblk+1]) { iblk++; }
1790 loc_ind = i - start_of_block[iblk];
1792 for (
int j = 0; j < dofs.
Size(); j++)
1794 int point_id = elem_blk[iblk][loc_ind*num_nod_per_el[iblk] + mymap[j] - 1] - 1;
1795 (*Nodes)(vdofs[j]) = coordx[point_id];
1796 (*Nodes)(vdofs[j]+1) = coordy[point_id];
1797 (*Nodes)(vdofs[j]+2) = coordz[point_id];
1806 for (
int i = 0; i < (int) num_side_sets; i++)
1808 delete [] elem_ss[i];
1809 delete [] side_ss[i];
1814 delete [] num_el_in_blk;
1815 delete [] num_nod_per_el;
1816 delete [] num_side_in_ss;
1821 for (
int i = 0; i < (int) num_el_blk; i++)
1823 delete [] elem_blk[i];
1827 delete [] start_of_block;
1829 for (
int i = 0; i < (int) num_side_sets; i++)
1831 delete [] ss_node_id[i];
1833 delete [] ss_node_id;
1838 #endif // #ifdef MFEM_USE_NETCDF
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
int Size() const
Logical size of the array.
Class for grid function - Vector with associated FE space.
void SetSize(int s)
Resize the vector to size s.
int Size() const
Returns the size of the vector.
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)
void DeleteAll()
Delete whole array.
Data type hexahedron element.
Type
Constants for the classes derived from Element.
virtual void GetElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom in array dofs for i'th element.
Data type triangle element.
A class for non-conforming AMR on higher-order hexahedral, quadrilateral or triangular meshes...
Data type tetrahedron element.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
void Destroy()
Destroy a vector.
void SetAttribute(const int attr)
Set element's attribute.
Piecewise-(bi)quadratic continuous finite elements.
void filter_dos(std::string &line)
Arbitrary order H1-conforming (continuous) finite elements.
Data type line segment element.
void DofsToVDofs(Array< int > &dofs, int ndofs=-1) const