13 #include "../fem/fem.hpp"
14 #include "../general/text.hpp"
19 #ifdef MFEM_USE_NETCDF
28 void Mesh::ReadMFEMMesh(std::istream &input,
bool mfem_v11,
int &curved)
37 MFEM_VERIFY(ident ==
"dimension",
"invalid mesh file");
43 MFEM_VERIFY(ident ==
"elements",
"invalid mesh file");
44 input >> NumOfElements;
45 elements.SetSize(NumOfElements);
46 for (
int j = 0; j < NumOfElements; j++)
48 elements[j] = ReadElement(input);
54 MFEM_VERIFY(ident ==
"boundary",
"invalid mesh file");
55 input >> NumOfBdrElements;
56 boundary.SetSize(NumOfBdrElements);
57 for (
int j = 0; j < NumOfBdrElements; j++)
59 boundary[j] = ReadElement(input);
65 if (mfem_v11 && ident ==
"vertex_parents")
67 ncmesh =
new NCMesh(
this, &input);
73 if (ident ==
"coarse_elements")
75 ncmesh->LoadCoarseElements(input);
82 MFEM_VERIFY(ident ==
"vertices",
"invalid mesh file");
83 input >> NumOfVertices;
84 vertices.SetSize(NumOfVertices);
90 spaceDim = atoi(ident.c_str());
91 for (
int j = 0; j < NumOfVertices; j++)
93 for (
int i = 0; i < spaceDim; i++)
95 input >> vertices[j](i);
100 if (ncmesh) { ncmesh->SetVertexPositions(vertices); }
110 void Mesh::ReadLineMesh(std::istream &input)
116 input >> NumOfVertices;
117 vertices.SetSize(NumOfVertices);
119 for (j = 0; j < NumOfVertices; j++)
121 input >> vertices[j](0);
124 input >> NumOfElements;
125 elements.SetSize(NumOfElements);
127 for (j = 0; j < NumOfElements; j++)
129 input >> a >> p1 >> p2;
130 elements[j] =
new Segment(p1-1, p2-1, a);
134 input >> NumOfBdrElements;
135 boundary.SetSize(NumOfBdrElements);
136 for (j = 0; j < NumOfBdrElements; j++)
138 input >> a >> ind[0];
140 boundary[j] =
new Point(ind,a);
144 void Mesh::ReadNetgen2DMesh(std::istream &input,
int &curved)
146 int ints[32], attr, n;
152 input >> NumOfBdrElements;
153 boundary.SetSize(NumOfBdrElements);
154 for (
int i = 0; i < NumOfBdrElements; i++)
157 >> ints[0] >> ints[1];
158 ints[0]--; ints[1]--;
159 boundary[i] =
new Segment(ints, attr);
163 input >> NumOfElements;
164 elements.SetSize(NumOfElements);
165 for (
int i = 0; i < NumOfElements; i++)
168 for (
int j = 0; j < n; j++)
176 elements[i] =
new Segment(ints, attr);
179 elements[i] =
new Triangle(ints, attr);
190 input >> NumOfVertices;
191 vertices.SetSize(NumOfVertices);
192 for (
int i = 0; i < NumOfVertices; i++)
193 for (
int j = 0; j < Dim; j++)
195 input >> vertices[i](j);
200 input >> NumOfVertices;
201 vertices.SetSize(NumOfVertices);
206 void Mesh::ReadNetgen3DMesh(std::istream &input)
214 input >> NumOfVertices;
216 vertices.SetSize(NumOfVertices);
217 for (
int i = 0; i < NumOfVertices; i++)
218 for (
int j = 0; j < Dim; j++)
220 input >> vertices[i](j);
224 input >> NumOfElements;
225 elements.SetSize(NumOfElements);
226 for (
int i = 0; i < NumOfElements; i++)
229 for (
int j = 0; j < 4; j++)
234 #ifdef MFEM_USE_MEMALLOC
236 tet = TetMemory.Alloc();
246 input >> NumOfBdrElements;
247 boundary.SetSize(NumOfBdrElements);
248 for (
int i = 0; i < NumOfBdrElements; i++)
251 for (
int j = 0; j < 3; j++)
256 boundary[i] =
new Triangle(ints, attr);
260 void Mesh::ReadTrueGridMesh(std::istream &input)
262 int i, j, ints[32], attr;
263 const int buflen = 1024;
274 input >> vari >> NumOfVertices >> vari >> vari >> NumOfElements;
275 input.getline(buf, buflen);
276 input.getline(buf, buflen);
278 input.getline(buf, buflen);
279 input.getline(buf, buflen);
280 input.getline(buf, buflen);
283 vertices.SetSize(NumOfVertices);
284 for (i = 0; i < NumOfVertices; i++)
286 input >> vari >> varf >> vertices[i](0) >> vertices[i](1);
287 input.getline(buf, buflen);
291 elements.SetSize(NumOfElements);
292 for (i = 0; i < NumOfElements; i++)
294 input >> vari >> attr;
295 for (j = 0; j < 4; j++)
300 input.getline(buf, buflen);
301 input.getline(buf, buflen);
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);
317 vertices.SetSize(NumOfVertices);
318 for (i = 0; i < NumOfVertices; i++)
320 input >> vari >> varf >> vertices[i](0) >> vertices[i](1)
322 input.getline(buf, buflen);
325 elements.SetSize(NumOfElements);
326 for (i = 0; i < NumOfElements; i++)
328 input >> vari >> attr;
329 for (j = 0; j < 8; j++)
334 input.getline(buf, buflen);
338 boundary.SetSize(NumOfBdrElements);
339 for (i = 0; i < NumOfBdrElements; i++)
342 for (j = 0; j < 4; j++)
347 input.getline(buf, buflen);
354 const int Mesh::vtk_quadratic_tet[10] =
355 { 0, 1, 2, 3, 4, 7, 5, 6, 8, 9 };
358 const int Mesh::vtk_quadratic_hex[27] =
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
364 void Mesh::ReadVTKMesh(std::istream &input,
int &curved,
int &read_gf)
369 getline(input, buff);
370 getline(input, buff);
374 MFEM_ABORT(
"VTK mesh is not in ASCII format!");
377 getline(input, buff);
379 if (buff !=
"DATASET UNSTRUCTURED_GRID")
381 MFEM_ABORT(
"VTK mesh is not UNSTRUCTURED_GRID!");
392 MFEM_ABORT(
"VTK mesh does not have POINTS data!");
395 while (buff !=
"POINTS");
401 getline(input, buff);
402 for (i = 0; i < points.
Size(); i++)
409 NumOfElements = n = 0;
414 input >> NumOfElements >> n >> ws;
416 for (i = 0; i < n; i++)
418 input >> cells_data[i];
426 if (buff ==
"CELL_TYPES")
428 input >> NumOfElements;
429 elements.
SetSize(NumOfElements);
430 for (j = i = 0; i < NumOfElements; i++)
438 elements[i] =
new Triangle(&cells_data[j+1]);
446 #ifdef MFEM_USE_MEMALLOC
447 elements[i] = TetMemory.Alloc();
448 elements[i]->SetVertices(&cells_data[j+1]);
455 elements[i] =
new Hexahedron(&cells_data[j+1]);
461 elements[i] =
new Triangle(&cells_data[j+1]);
471 #ifdef MFEM_USE_MEMALLOC
472 elements[i] = TetMemory.Alloc();
473 elements[i]->SetVertices(&cells_data[j+1]);
481 elements[i] =
new Hexahedron(&cells_data[j+1]);
484 MFEM_ABORT(
"VTK mesh : cell type " << ct <<
" is not supported!");
487 j += cells_data[j] + 1;
492 streampos sp = input.tellg();
494 if (buff ==
"CELL_DATA")
497 getline(input, buff);
500 if (!strncmp(buff.c_str(),
"SCALARS material", 16))
502 getline(input, buff);
503 for (i = 0; i < NumOfElements; i++)
506 elements[i]->SetAttribute(attr);
523 vertices.SetSize(np);
524 for (i = 0; i < np; i++)
526 vertices[i](0) = points(3*i+0);
527 vertices[i](1) = points(3*i+1);
528 vertices[i](2) = points(3*i+2);
533 NumOfBdrElements = 0;
542 for (n = i = 0; i < NumOfElements; i++)
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)
553 for (n = i = 0; i < np; i++)
554 if (pts_dof[i] != -1)
559 for (i = 0; i < NumOfElements; i++)
561 int *v = elements[i]->GetVertices();
562 int nv = elements[i]->GetNVertices();
563 for (j = 0; j < nv; j++)
565 v[j] = pts_dof[v[j]];
571 for (i = 0; i < np; i++)
573 if ((j = pts_dof[i]) != -1)
575 vertices[j](0) = points(3*i+0);
576 vertices[j](1) = points(3*i+1);
577 vertices[j](2) = points(3*i+2);
582 NumOfBdrElements = 0;
590 GetElementToFaceTable();
599 el_to_edge =
new Table;
600 NumOfEdges = GetElementToEdgeTable(*el_to_edge, be_to_edge);
610 Nodes->MakeOwner(fec);
615 for (n = i = 0; i < NumOfElements; i++)
619 switch (elements[i]->GetGeometryType())
621 case Geometry::TRIANGLE:
622 case Geometry::SQUARE:
623 vtk_mfem = vtk_quadratic_hex;
break;
624 case Geometry::TETRAHEDRON:
625 vtk_mfem = vtk_quadratic_tet;
break;
628 vtk_mfem = vtk_quadratic_hex;
break;
631 for (n++, j = 0; j < dofs.
Size(); j++, n++)
633 if (pts_dof[cells_data[n]] == -1)
635 pts_dof[cells_data[n]] = dofs[vtk_mfem[j]];
639 if (pts_dof[cells_data[n]] != dofs[vtk_mfem[j]])
641 MFEM_ABORT(
"VTK mesh : inconsistent quadratic mesh!");
648 for (i = 0; i < np; i++)
651 if ((dofs[0] = pts_dof[i]) != -1)
654 for (j = 0; j < dofs.
Size(); j++)
656 (*Nodes)(dofs[j]) = points(3*i+j);
665 void Mesh::ReadNURBSMesh(std::istream &input,
int &curved,
int &read_gf)
669 Dim = NURBSext->Dimension();
670 NumOfVertices = NURBSext->GetNV();
671 NumOfElements = NURBSext->GetNE();
672 NumOfBdrElements = NURBSext->GetNBE();
674 NURBSext->GetElementTopo(elements);
675 NURBSext->GetBdrElementTopo(boundary);
677 vertices.SetSize(NumOfVertices);
679 if (NURBSext->HavePatches())
685 Nodes->MakeOwner(fec);
686 NURBSext->SetCoordsFromPatches(*Nodes);
689 int vd = Nodes->VectorDim();
690 for (
int i = 0; i < vd; i++)
693 Nodes->GetNodalValues(vert_val, i+1);
694 for (
int j = 0; j < NumOfVertices; j++)
696 vertices[j](i) = vert_val(j);
706 void Mesh::ReadInlineMesh(std::istream &input,
int generate_edges)
734 MFEM_VERIFY(input.get() ==
'=',
735 "Inline mesh expected '=' after keyword " << name);
742 else if (name ==
"ny")
746 else if (name ==
"nz")
750 else if (name ==
"sx")
754 else if (name ==
"sy")
758 else if (name ==
"sz")
762 else if (name ==
"type")
766 if (eltype ==
"segment")
768 type = Element::SEGMENT;
770 else if (eltype ==
"quad")
772 type = Element::QUADRILATERAL;
774 else if (eltype ==
"tri")
776 type = Element::TRIANGLE;
778 else if (eltype ==
"hex")
780 type = Element::HEXAHEDRON;
782 else if (eltype ==
"tet")
784 type = Element::TETRAHEDRON;
788 MFEM_ABORT(
"unrecognized element type (read '" << eltype
789 <<
"') in inline mesh format. "
790 "Allowed: segment, tri, tet, quad, hex");
795 MFEM_ABORT(
"unrecognized keyword (" << name
796 <<
") in inline mesh format. "
797 "Allowed: nx, ny, nz, type, sx, sy, sz");
802 if (input.peek() ==
';')
815 if (type == Element::SEGMENT)
817 MFEM_VERIFY(nx > 0 && sx > 0.0,
818 "invalid 1D inline mesh format, all values must be "
820 <<
" nx = " << nx <<
"\n"
821 <<
" sx = " << sx <<
"\n");
824 else if (type == Element::TRIANGLE || type == Element::QUADRILATERAL)
826 MFEM_VERIFY(nx > 0 && ny > 0 && sx > 0.0 && sy > 0.0,
827 "invalid 2D inline mesh format, all values must be "
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);
835 else if (type == Element::TETRAHEDRON || type == Element::HEXAHEDRON)
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 "
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);
851 MFEM_ABORT(
"For inline mesh, must specify an element type ="
852 " [segment, tri, quad, tet, hex]");
858 void Mesh::ReadGmshMesh(std::istream &input)
863 input >> version >> binary >> dsize;
866 MFEM_ABORT(
"Gmsh file version < 2.2");
868 if (dsize !=
sizeof(
double))
870 MFEM_ABORT(
"Gmsh file : dsize != sizeof(double)");
872 getline(input, buff);
877 input.read(reinterpret_cast<char*>(&one),
sizeof(one));
880 MFEM_ABORT(
"Gmsh file : wrong binary format");
887 map<int, int> vertices_map;
890 while (input >> buff)
892 if (buff ==
"$Nodes")
894 input >> NumOfVertices;
895 getline(input, buff);
896 vertices.SetSize(NumOfVertices);
898 const int gmsh_dim = 3;
899 double coord[gmsh_dim];
900 for (
int ver = 0; ver < NumOfVertices; ++ver)
904 input.read(reinterpret_cast<char*>(&serial_number),
sizeof(
int));
905 input.read(reinterpret_cast<char*>(coord), gmsh_dim*
sizeof(
double));
909 input >> serial_number;
910 for (
int ci = 0; ci < gmsh_dim; ++ci)
915 vertices[ver] =
Vertex(coord, gmsh_dim);
916 vertices_map[serial_number] = ver;
918 if (static_cast<int>(vertices_map.size()) != NumOfVertices)
920 MFEM_ABORT(
"Gmsh file : vertices indices are not unique");
923 else if (buff ==
"$Elements")
925 int num_of_all_elements;
926 input >> num_of_all_elements;
928 getline(input, buff);
939 int nodes_of_gmsh_element[] =
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);
1004 int n_elem_part = 0;
1005 const int header_size = 3;
1008 int header[header_size];
1009 int n_elem_one_type;
1011 while (n_elem_part < num_of_all_elements)
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];
1019 n_elem_part += n_elem_one_type;
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)
1025 input.read(reinterpret_cast<char*>(&data[0]),
1026 data.size()*
sizeof(int));
1028 serial_number = data[dd++];
1031 phys_domain = (n_tags > 0) ? data[dd++] : 0;
1034 elem_domain = (n_tags > 1) ? data[dd++] : 0;
1037 n_partitions = (n_tags > 2) ? data[dd++] : 0;
1040 vector<int> vert_indices(n_elem_nodes);
1041 for (
int vi = 0; vi < n_elem_nodes; ++vi)
1043 map<int, int>::const_iterator it =
1044 vertices_map.find(data[1+n_tags+vi]);
1045 if (it == vertices_map.end())
1047 MFEM_ABORT(
"Gmsh file : vertex index doesn't exist");
1049 vert_indices[vi] = it->second;
1053 if (phys_domain <= 0)
1055 MFEM_ABORT(
"Non-positive element attribute in Gmsh mesh!");
1059 switch (type_of_element)
1063 elements_1D.push_back(
1064 new Segment(&vert_indices[0], phys_domain));
1069 elements_2D.push_back(
1070 new Triangle(&vert_indices[0], phys_domain));
1075 elements_2D.push_back(
1081 elements_3D.push_back(
1087 elements_3D.push_back(
1088 new Hexahedron(&vert_indices[0], phys_domain));
1093 elements_0D.push_back(
1094 new Point(&vert_indices[0], phys_domain));
1098 MFEM_WARNING(
"Unsupported Gmsh element type.");
1107 for (
int el = 0; el < num_of_all_elements; ++el)
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]; }
1114 phys_domain = (n_tags > 0) ? data[0] : 0;
1117 elem_domain = (n_tags > 1) ? data[1] : 0;
1120 n_partitions = (n_tags > 2) ? data[2] : 0;
1123 const int n_elem_nodes = nodes_of_gmsh_element[type_of_element-1];
1124 vector<int> vert_indices(n_elem_nodes);
1126 for (
int vi = 0; vi < n_elem_nodes; ++vi)
1129 map<int, int>::const_iterator it = vertices_map.find(index);
1130 if (it == vertices_map.end())
1132 MFEM_ABORT(
"Gmsh file : vertex index doesn't exist");
1134 vert_indices[vi] = it->second;
1138 if (phys_domain <= 0)
1140 MFEM_ABORT(
"Non-positive element attribute in Gmsh mesh!");
1144 switch (type_of_element)
1148 elements_1D.push_back(
1149 new Segment(&vert_indices[0], phys_domain));
1154 elements_2D.push_back(
1155 new Triangle(&vert_indices[0], phys_domain));
1160 elements_2D.push_back(
1166 elements_3D.push_back(
1172 elements_3D.push_back(
1173 new Hexahedron(&vert_indices[0], phys_domain));
1178 elements_0D.push_back(
1179 new Point(&vert_indices[0], phys_domain));
1183 MFEM_WARNING(
"Unsupported Gmsh element type.");
1190 if (!elements_3D.empty())
1193 NumOfElements = elements_3D.size();
1194 elements.SetSize(NumOfElements);
1195 for (
int el = 0; el < NumOfElements; ++el)
1197 elements[el] = elements_3D[el];
1199 NumOfBdrElements = elements_2D.size();
1200 boundary.SetSize(NumOfBdrElements);
1201 for (
int el = 0; el < NumOfBdrElements; ++el)
1203 boundary[el] = elements_2D[el];
1206 for (
size_t el = 0; el < elements_1D.size(); ++el)
1208 delete elements_1D[el];
1210 for (
size_t el = 0; el < elements_0D.size(); ++el)
1212 delete elements_0D[el];
1215 else if (!elements_2D.empty())
1218 NumOfElements = elements_2D.size();
1219 elements.SetSize(NumOfElements);
1220 for (
int el = 0; el < NumOfElements; ++el)
1222 elements[el] = elements_2D[el];
1224 NumOfBdrElements = elements_1D.size();
1225 boundary.SetSize(NumOfBdrElements);
1226 for (
int el = 0; el < NumOfBdrElements; ++el)
1228 boundary[el] = elements_1D[el];
1231 for (
size_t el = 0; el < elements_0D.size(); ++el)
1233 delete elements_0D[el];
1236 else if (!elements_1D.empty())
1239 NumOfElements = elements_1D.size();
1240 elements.SetSize(NumOfElements);
1241 for (
int el = 0; el < NumOfElements; ++el)
1243 elements[el] = elements_1D[el];
1245 NumOfBdrElements = elements_0D.size();
1246 boundary.SetSize(NumOfBdrElements);
1247 for (
int el = 0; el < NumOfBdrElements; ++el)
1249 boundary[el] = elements_0D[el];
1254 MFEM_ABORT(
"Gmsh file : no elements found");
1258 MFEM_CONTRACT_VAR(n_partitions);
1259 MFEM_CONTRACT_VAR(elem_domain);
1266 #ifdef MFEM_USE_NETCDF
1267 void Mesh::ReadCubit(
const char *filename,
int &curved,
int &read_gf)
1274 const int sideMapHex8[6][4] =
1284 const int sideMapTet4[4][3] =
1293 const int mfemToGenesisTet10[10] = {1,2,3,4,5,7,8,6,9,10};
1296 const int mfemToGenesisHex27[27] = {1,2,3,4,5,6,7,8,9,10,11,
1298 12,17,18,19,20,13,14,15,
1300 16,22,26,25,27,24,23,21
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};
1306 const int sideMapHex27[6][9] =
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}
1316 const int sideMapTet10[4][6] =
1328 char str_dummy[256];
1335 if ((retval = nc_open(filename, NC_NOWRITE, &ncid)))
1337 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1343 size_t num_dim=0, num_nodes=0, num_elem=0, num_el_blk=0, num_side_sets=0;
1345 if ((retval = nc_inq_dimid(ncid,
"num_dim", &
id)) ||
1346 (retval = nc_inq_dim(ncid,
id, str_dummy, &num_dim)) ||
1348 (retval = nc_inq_dimid(ncid,
"num_nodes", &
id)) ||
1349 (retval = nc_inq_dim(ncid,
id, str_dummy, &num_nodes)) ||
1351 (retval = nc_inq_dimid(ncid,
"num_elem" ,&
id)) ||
1352 (retval = nc_inq_dim(ncid,
id, str_dummy, &num_elem)) ||
1354 (retval = nc_inq_dimid(ncid,
"num_el_blk", &
id)) ||
1355 (retval = nc_inq_dim(ncid,
id, str_dummy, &num_el_blk)) ||
1357 (retval = nc_inq_dimid(ncid,
"num_side_sets", &
id)) ||
1358 (retval = nc_inq_dim(ncid,
id, str_dummy, &num_side_sets)))
1360 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
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];
1371 int previous_num_nod_per_el = 0;
1372 for (
int i = 0; i < (int) num_el_blk; i++)
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])))
1378 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
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])))
1385 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1392 if ((
int) num_nod_per_el[i] != previous_num_nod_per_el)
1394 MFEM_ABORT(
"Element blocks of different element types not supported");
1397 previous_num_nod_per_el = num_nod_per_el[i];
1401 if (num_nod_per_el[0] == 4 || num_nod_per_el[0] == 8)
1405 else if (num_nod_per_el[0] == 10 || num_nod_per_el[0] == 27)
1411 MFEM_ABORT(
"Don't know what to do with a " << num_nod_per_el[0] <<
1415 for (
int i = 0; i < (int) num_side_sets; i++)
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])))
1421 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1426 double *coordx =
new double[num_nodes];
1427 double *coordy =
new double[num_nodes];
1428 double *coordz =
new double[num_nodes];
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)))
1435 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1440 if ((retval = nc_inq_varid(ncid,
"coordz", &
id)) ||
1441 (retval = nc_get_var_double(ncid,
id, coordz)))
1443 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1448 int **elem_blk =
new int*[num_el_blk];
1449 for (
int i = 0; i < (int) num_el_blk; i++)
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])))
1456 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
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)))
1463 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1468 int **elem_ss =
new int*[num_side_sets];
1469 int **side_ss =
new int*[num_side_sets];
1471 for (
int i = 0; i < (int) num_side_sets; i++)
1473 elem_ss[i] =
new int[num_side_in_ss[i]];
1474 side_ss[i] =
new int[num_side_in_ss[i]];
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])))
1480 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
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])))
1487 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
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)))
1495 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1500 int num_nod_per_side = 0;
1501 if (num_nod_per_el[0] == 4)
1503 num_nod_per_side = 3;
1505 else if (num_nod_per_el[0] == 8)
1507 num_nod_per_side = 4;
1509 else if (num_nod_per_el[0] == 10)
1511 num_nod_per_side = 6;
1513 else if (num_nod_per_el[0] == 27)
1515 num_nod_per_side = 9;
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++)
1525 start_of_block[i] = start_of_block[i-1] + num_el_in_blk[i-1];
1528 int **ss_node_id =
new int*[num_side_sets];
1530 for (
int i = 0; i < (int) num_side_sets; i++)
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++)
1535 int glob_ind = elem_ss[i][j]-1;
1538 while (iblk < (
int) num_el_blk && glob_ind >= start_of_block[iblk+1])
1542 if (iblk >= (
int) num_el_blk)
1544 MFEM_ABORT(
"Sideset element does not exist");
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)
1551 for (
int k = 0; k < num_nod_per_side; k++)
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];
1558 else if (num_nod_per_side == 6)
1560 for (
int k = 0; k < num_nod_per_side; k++)
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];
1567 else if (num_nod_per_side == 4)
1569 for (
int k = 0; k < num_nod_per_side; k++)
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];
1576 else if (num_nod_per_side == 9)
1578 for (
int k = 0; k < num_nod_per_side; k++)
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];
1590 std::vector<int> uniqueVertexID;
1591 for (
int iblk = 0; iblk < (int) num_el_blk; iblk++)
1593 for (
int i = 0; i < (int) num_el_in_blk[iblk]; i++)
1595 int NumNodePerEl = num_nod_per_el[iblk];
1596 int NumVertPerEl = 0;
1597 if (NumNodePerEl == 8 || NumNodePerEl == 27)
1601 if (NumNodePerEl == 4 || NumNodePerEl == 10)
1605 for (
int j = 0; j < NumVertPerEl; j++)
1607 uniqueVertexID.push_back(elem_blk[iblk][i*NumNodePerEl + j]);
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));
1620 std::map<int,int> cubitToMFEMVertMap;
1621 for (
int i = 0; i < (int) uniqueVertexID.size(); i++)
1623 cubitToMFEMVertMap[uniqueVertexID[i]] = i+1;
1625 MFEM_ASSERT(cubitToMFEMVertMap.size() == uniqueVertexID.size(),
1626 "This should never happen\n");
1632 NumOfVertices = uniqueVertexID.size();
1633 vertices.SetSize(NumOfVertices);
1634 for (
int i = 0; i < (int) uniqueVertexID.size(); i++)
1636 vertices[i](0) = coordx[uniqueVertexID[i] - 1];
1637 vertices[i](1) = coordy[uniqueVertexID[i] - 1];
1638 vertices[i](2) = coordz[uniqueVertexID[i] - 1];
1641 NumOfElements = num_elem;
1642 elements.SetSize(num_elem);
1644 int renumberedVertID[8];
1645 for (
int iblk = 0; iblk < (int) num_el_blk; iblk++)
1647 int NumNodePerEl = num_nod_per_el[iblk];
1648 if (NumNodePerEl == 4 || NumNodePerEl == 10)
1650 for (
int i = 0; i < (int) num_el_in_blk[iblk]; i++)
1652 for (
int j = 0; j < 4; j++)
1654 renumberedVertID[j] =
1655 cubitToMFEMVertMap[elem_blk[iblk][i*NumNodePerEl+j]]-1;
1657 elements[elcount] =
new Tetrahedron(renumberedVertID,ebprop[iblk]);
1661 else if (NumNodePerEl == 8 || NumNodePerEl == 27)
1663 for (
int i = 0; i < (int) num_el_in_blk[iblk]; i++)
1665 for (
int j = 0; j < 8; j++)
1667 renumberedVertID[j] =
1668 cubitToMFEMVertMap[elem_blk[iblk][i*NumNodePerEl+j]]-1;
1670 elements[elcount] =
new Hexahedron(renumberedVertID,ebprop[iblk]);
1678 NumOfBdrElements = 0;
1679 for (
int iss = 0; iss < (int) num_side_sets; iss++)
1681 NumOfBdrElements += num_side_in_ss[iss];
1683 boundary.SetSize(NumOfBdrElements);
1685 for (
int iss = 0; iss < (int) num_side_sets; iss++)
1687 if (num_nod_per_side == 3 || num_nod_per_side == 6)
1689 for (
int i = 0; i < (int) num_side_in_ss[iss]; i++)
1691 for (
int j = 0; j < 3; j++)
1693 renumberedVertID[j] =
1694 cubitToMFEMVertMap[ss_node_id[iss][i*num_nod_per_side+j]] - 1;
1696 boundary[sidecount] =
new Triangle(renumberedVertID,ssprop[iss]);
1700 else if (num_nod_per_side == 4 || num_nod_per_side == 9)
1702 for (
int i = 0; i < (int) num_side_in_ss[iss]; i++)
1704 for (
int j = 0; j < 4; j++)
1706 renumberedVertID[j] =
1707 cubitToMFEMVertMap[ss_node_id[iss][i*num_nod_per_side+j]] - 1;
1709 boundary[sidecount] =
new Quadrilateral(renumberedVertID,ssprop[iss]);
1720 if (num_nod_per_el[0] == 10)
1722 mymap = (
int *) mfemToGenesisTet10;
1724 else if (num_nod_per_el[0] == 27)
1726 mymap = (
int *) mfemToGenesisHex27;
1728 else if (num_nod_per_el[0] == 6)
1730 mymap = (
int *) mfemToGenesisTri6;
1732 else if (num_nod_per_el[0] == 9)
1734 mymap = (
int *) mfemToGenesisQuad9;
1743 GetElementToFaceTable();
1752 el_to_edge =
new Table;
1753 NumOfEdges = GetElementToEdgeTable(*el_to_edge, be_to_edge);
1764 Nodes->MakeOwner(fec);
1772 for (
int i = 0; i < NumOfElements; i++)
1779 for (
int l = 0; l < dofs.
Size(); l++) { vdofs[l] = dofs[l]; }
1783 while (iblk < (
int) num_el_blk && i >= start_of_block[iblk+1]) { iblk++; }
1784 loc_ind = i - start_of_block[iblk];
1786 for (
int j = 0; j < dofs.
Size(); j++)
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];
1797 delete [] num_el_in_blk;
1798 delete [] num_nod_per_el;
1799 delete [] num_side_in_ss;
1803 for (
int i = 0; i < (int) num_el_blk; i++)
1808 delete [] start_of_block;
1809 for (
int i = 0; i < (int) num_side_sets; i++)
1811 delete [] ss_node_id[i];
1813 delete [] ss_node_id;
1818 #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.
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