13 #include "../fem/fem.hpp"
17 #ifdef MFEM_USE_NETCDF
26 void Mesh::ReadMFEMMesh(std::istream &input,
bool mfem_v11,
int &curved)
32 skip_comment_lines(input,
'#');
35 MFEM_VERIFY(ident ==
"dimension",
"invalid mesh file");
38 skip_comment_lines(input,
'#');
41 MFEM_VERIFY(ident ==
"elements",
"invalid mesh file");
42 input >> NumOfElements;
43 elements.SetSize(NumOfElements);
44 for (
int j = 0; j < NumOfElements; j++)
46 elements[j] = ReadElement(input);
49 skip_comment_lines(input,
'#');
52 MFEM_VERIFY(ident ==
"boundary",
"invalid mesh file");
53 input >> NumOfBdrElements;
54 boundary.SetSize(NumOfBdrElements);
55 for (
int j = 0; j < NumOfBdrElements; j++)
57 boundary[j] = ReadElement(input);
60 skip_comment_lines(input,
'#');
63 if (ident ==
"vertex_parents" && mfem_v11)
65 ncmesh =
new NCMesh(
this, &input);
68 skip_comment_lines(input,
'#');
71 if (ident ==
"coarse_elements")
73 ncmesh->LoadCoarseElements(input);
75 skip_comment_lines(input,
'#');
80 MFEM_VERIFY(ident ==
"vertices",
"invalid mesh file");
81 input >> NumOfVertices;
82 vertices.SetSize(NumOfVertices);
88 spaceDim = atoi(ident.c_str());
89 for (
int j = 0; j < NumOfVertices; j++)
91 for (
int i = 0; i < spaceDim; i++)
93 input >> vertices[j](i);
98 if (ncmesh) { ncmesh->SetVertexPositions(vertices); }
108 void Mesh::ReadLineMesh(std::istream &input)
114 input >> NumOfVertices;
115 vertices.SetSize(NumOfVertices);
117 for (j = 0; j < NumOfVertices; j++)
119 input >> vertices[j](0);
122 input >> NumOfElements;
123 elements.SetSize(NumOfElements);
125 for (j = 0; j < NumOfElements; j++)
127 input >> a >> p1 >> p2;
128 elements[j] =
new Segment(p1-1, p2-1, a);
132 input >> NumOfBdrElements;
133 boundary.SetSize(NumOfBdrElements);
134 for (j = 0; j < NumOfBdrElements; j++)
136 input >> a >> ind[0];
138 boundary[j] =
new Point(ind,a);
142 void Mesh::ReadNetgen2DMesh(std::istream &input,
int &curved)
144 int ints[32], attr, n;
150 input >> NumOfBdrElements;
151 boundary.SetSize(NumOfBdrElements);
152 for (
int i = 0; i < NumOfBdrElements; i++)
155 >> ints[0] >> ints[1];
156 ints[0]--; ints[1]--;
157 boundary[i] =
new Segment(ints, attr);
161 input >> NumOfElements;
162 elements.SetSize(NumOfElements);
163 for (
int i = 0; i < NumOfElements; i++)
166 for (
int j = 0; j < n; j++)
174 elements[i] =
new Segment(ints, attr);
177 elements[i] =
new Triangle(ints, attr);
188 input >> NumOfVertices;
189 vertices.SetSize(NumOfVertices);
190 for (
int i = 0; i < NumOfVertices; i++)
191 for (
int j = 0; j < Dim; j++)
193 input >> vertices[i](j);
198 input >> NumOfVertices;
199 vertices.SetSize(NumOfVertices);
204 void Mesh::ReadNetgen3DMesh(std::istream &input)
212 input >> NumOfVertices;
214 vertices.SetSize(NumOfVertices);
215 for (
int i = 0; i < NumOfVertices; i++)
216 for (
int j = 0; j < Dim; j++)
218 input >> vertices[i](j);
222 input >> NumOfElements;
223 elements.SetSize(NumOfElements);
224 for (
int i = 0; i < NumOfElements; i++)
227 for (
int j = 0; j < 4; j++)
232 #ifdef MFEM_USE_MEMALLOC
234 tet = TetMemory.Alloc();
244 input >> NumOfBdrElements;
245 boundary.SetSize(NumOfBdrElements);
246 for (
int i = 0; i < NumOfBdrElements; i++)
249 for (
int j = 0; j < 3; j++)
254 boundary[i] =
new Triangle(ints, attr);
258 void Mesh::ReadTrueGridMesh(std::istream &input)
260 int i, j, ints[32], attr;
261 const int buflen = 1024;
272 input >> vari >> NumOfVertices >> vari >> vari >> NumOfElements;
273 input.getline(buf, buflen);
274 input.getline(buf, buflen);
276 input.getline(buf, buflen);
277 input.getline(buf, buflen);
278 input.getline(buf, buflen);
281 vertices.SetSize(NumOfVertices);
282 for (i = 0; i < NumOfVertices; i++)
284 input >> vari >> varf >> vertices[i](0) >> vertices[i](1);
285 input.getline(buf, buflen);
289 elements.SetSize(NumOfElements);
290 for (i = 0; i < NumOfElements; i++)
292 input >> vari >> attr;
293 for (j = 0; j < 4; j++)
298 input.getline(buf, buflen);
299 input.getline(buf, buflen);
307 input >> vari >> NumOfVertices >> NumOfElements;
308 input.getline(buf, buflen);
309 input.getline(buf, buflen);
310 input >> vari >> vari >> NumOfBdrElements;
311 input.getline(buf, buflen);
312 input.getline(buf, buflen);
313 input.getline(buf, buflen);
315 vertices.SetSize(NumOfVertices);
316 for (i = 0; i < NumOfVertices; i++)
318 input >> vari >> varf >> vertices[i](0) >> vertices[i](1)
320 input.getline(buf, buflen);
323 elements.SetSize(NumOfElements);
324 for (i = 0; i < NumOfElements; i++)
326 input >> vari >> attr;
327 for (j = 0; j < 8; j++)
332 input.getline(buf, buflen);
336 boundary.SetSize(NumOfBdrElements);
337 for (i = 0; i < NumOfBdrElements; i++)
340 for (j = 0; j < 4; j++)
345 input.getline(buf, buflen);
352 const int Mesh::vtk_quadratic_tet[10] =
353 { 0, 1, 2, 3, 4, 7, 5, 6, 8, 9 };
356 const int Mesh::vtk_quadratic_hex[27] =
358 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19,
359 24, 22, 21, 23, 20, 25, 26
362 void Mesh::ReadVTKMesh(std::istream &input,
int &curved,
int &read_gf)
367 getline(input, buff);
368 getline(input, buff);
372 MFEM_ABORT(
"VTK mesh is not in ASCII format!");
375 getline(input, buff);
377 if (buff !=
"DATASET UNSTRUCTURED_GRID")
379 MFEM_ABORT(
"VTK mesh is not UNSTRUCTURED_GRID!");
390 MFEM_ABORT(
"VTK mesh does not have POINTS data!");
393 while (buff !=
"POINTS");
399 getline(input, buff);
400 for (i = 0; i < points.
Size(); i++)
407 NumOfElements = n = 0;
412 input >> NumOfElements >> n >> ws;
414 for (i = 0; i < n; i++)
416 input >> cells_data[i];
424 if (buff ==
"CELL_TYPES")
426 input >> NumOfElements;
427 elements.
SetSize(NumOfElements);
428 for (j = i = 0; i < NumOfElements; i++)
436 elements[i] =
new Triangle(&cells_data[j+1]);
444 #ifdef MFEM_USE_MEMALLOC
445 elements[i] = TetMemory.Alloc();
446 elements[i]->SetVertices(&cells_data[j+1]);
453 elements[i] =
new Hexahedron(&cells_data[j+1]);
459 elements[i] =
new Triangle(&cells_data[j+1]);
469 #ifdef MFEM_USE_MEMALLOC
470 elements[i] = TetMemory.Alloc();
471 elements[i]->SetVertices(&cells_data[j+1]);
479 elements[i] =
new Hexahedron(&cells_data[j+1]);
482 MFEM_ABORT(
"VTK mesh : cell type " << ct <<
" is not supported!");
485 j += cells_data[j] + 1;
490 streampos sp = input.tellg();
492 if (buff ==
"CELL_DATA")
495 getline(input, buff);
498 if (!strncmp(buff.c_str(),
"SCALARS material", 16))
500 getline(input, buff);
501 for (i = 0; i < NumOfElements; i++)
504 elements[i]->SetAttribute(attr);
521 vertices.SetSize(np);
522 for (i = 0; i < np; i++)
524 vertices[i](0) = points(3*i+0);
525 vertices[i](1) = points(3*i+1);
526 vertices[i](2) = points(3*i+2);
531 NumOfBdrElements = 0;
540 for (n = i = 0; i < NumOfElements; i++)
542 int *v = elements[i]->GetVertices();
543 int nv = elements[i]->GetNVertices();
544 for (j = 0; j < nv; j++)
545 if (pts_dof[v[j]] == -1)
551 for (n = i = 0; i < np; i++)
552 if (pts_dof[i] != -1)
557 for (i = 0; i < NumOfElements; i++)
559 int *v = elements[i]->GetVertices();
560 int nv = elements[i]->GetNVertices();
561 for (j = 0; j < nv; j++)
563 v[j] = pts_dof[v[j]];
569 for (i = 0; i < np; i++)
571 if ((j = pts_dof[i]) != -1)
573 vertices[j](0) = points(3*i+0);
574 vertices[j](1) = points(3*i+1);
575 vertices[j](2) = points(3*i+2);
580 NumOfBdrElements = 0;
588 GetElementToFaceTable();
597 el_to_edge =
new Table;
598 NumOfEdges = GetElementToEdgeTable(*el_to_edge, be_to_edge);
608 Nodes->MakeOwner(fec);
613 for (n = i = 0; i < NumOfElements; i++)
617 switch (elements[i]->GetGeometryType())
619 case Geometry::TRIANGLE:
620 case Geometry::SQUARE:
621 vtk_mfem = vtk_quadratic_hex;
break;
622 case Geometry::TETRAHEDRON:
623 vtk_mfem = vtk_quadratic_tet;
break;
626 vtk_mfem = vtk_quadratic_hex;
break;
629 for (n++, j = 0; j < dofs.
Size(); j++, n++)
631 if (pts_dof[cells_data[n]] == -1)
633 pts_dof[cells_data[n]] = dofs[vtk_mfem[j]];
637 if (pts_dof[cells_data[n]] != dofs[vtk_mfem[j]])
639 MFEM_ABORT(
"VTK mesh : inconsistent quadratic mesh!");
646 for (i = 0; i < np; i++)
649 if ((dofs[0] = pts_dof[i]) != -1)
652 for (j = 0; j < dofs.
Size(); j++)
654 (*Nodes)(dofs[j]) = points(3*i+j);
663 void Mesh::ReadNURBSMesh(std::istream &input,
int &curved,
int &read_gf)
667 Dim = NURBSext->Dimension();
668 NumOfVertices = NURBSext->GetNV();
669 NumOfElements = NURBSext->GetNE();
670 NumOfBdrElements = NURBSext->GetNBE();
672 NURBSext->GetElementTopo(elements);
673 NURBSext->GetBdrElementTopo(boundary);
675 vertices.SetSize(NumOfVertices);
677 if (NURBSext->HavePatches())
683 Nodes->MakeOwner(fec);
684 NURBSext->SetCoordsFromPatches(*Nodes);
687 int vd = Nodes->VectorDim();
688 for (
int i = 0; i < vd; i++)
691 Nodes->GetNodalValues(vert_val, i+1);
692 for (
int j = 0; j < NumOfVertices; j++)
694 vertices[j](i) = vert_val(j);
704 void Mesh::ReadInlineMesh(std::istream &input,
int generate_edges)
719 skip_comment_lines(input,
'#');
732 MFEM_VERIFY(input.get() ==
'=',
733 "Inline mesh expected '=' after keyword " << name);
740 else if (name ==
"ny")
744 else if (name ==
"nz")
748 else if (name ==
"sx")
752 else if (name ==
"sy")
756 else if (name ==
"sz")
760 else if (name ==
"type")
764 if (eltype ==
"segment")
766 type = Element::SEGMENT;
768 else if (eltype ==
"quad")
770 type = Element::QUADRILATERAL;
772 else if (eltype ==
"tri")
774 type = Element::TRIANGLE;
776 else if (eltype ==
"hex")
778 type = Element::HEXAHEDRON;
780 else if (eltype ==
"tet")
782 type = Element::TETRAHEDRON;
786 MFEM_ABORT(
"unrecognized element type (read '" << eltype
787 <<
"') in inline mesh format. "
788 "Allowed: segment, tri, tet, quad, hex");
793 MFEM_ABORT(
"unrecognized keyword (" << name
794 <<
") in inline mesh format. "
795 "Allowed: nx, ny, nz, type, sx, sy, sz");
800 if (input.peek() ==
';')
813 if (type == Element::SEGMENT)
815 MFEM_VERIFY(nx > 0 && sx > 0.0,
816 "invalid 1D inline mesh format, all values must be "
818 <<
" nx = " << nx <<
"\n"
819 <<
" sx = " << sx <<
"\n");
822 else if (type == Element::TRIANGLE || type == Element::QUADRILATERAL)
824 MFEM_VERIFY(nx > 0 && ny > 0 && sx > 0.0 && sy > 0.0,
825 "invalid 2D inline mesh format, all values must be "
827 <<
" nx = " << nx <<
"\n"
828 <<
" ny = " << ny <<
"\n"
829 <<
" sx = " << sx <<
"\n"
830 <<
" sy = " << sy <<
"\n");
831 Make2D(nx, ny, type, generate_edges, sx, sy);
833 else if (type == Element::TETRAHEDRON || type == Element::HEXAHEDRON)
835 MFEM_VERIFY(nx > 0 && ny > 0 && nz > 0 &&
836 sx > 0.0 && sy > 0.0 && sz > 0.0,
837 "invalid 3D inline mesh format, all values must be "
839 <<
" nx = " << nx <<
"\n"
840 <<
" ny = " << ny <<
"\n"
841 <<
" nz = " << nz <<
"\n"
842 <<
" sx = " << sx <<
"\n"
843 <<
" sy = " << sy <<
"\n"
844 <<
" sz = " << sz <<
"\n");
845 Make3D(nx, ny, nz, type, generate_edges, sx, sy, sz);
849 MFEM_ABORT(
"For inline mesh, must specify an element type ="
850 " [segment, tri, quad, tet, hex]");
856 void Mesh::ReadGmshMesh(std::istream &input)
861 input >> version >> binary >> dsize;
864 MFEM_ABORT(
"Gmsh file version < 2.2");
866 if (dsize !=
sizeof(
double))
868 MFEM_ABORT(
"Gmsh file : dsize != sizeof(double)");
870 getline(input, buff);
875 input.read(reinterpret_cast<char*>(&one),
sizeof(one));
878 MFEM_ABORT(
"Gmsh file : wrong binary format");
885 map<int, int> vertices_map;
888 while (input >> buff)
890 if (buff ==
"$Nodes")
892 input >> NumOfVertices;
893 getline(input, buff);
894 vertices.SetSize(NumOfVertices);
896 const int gmsh_dim = 3;
897 double coord[gmsh_dim];
898 for (
int ver = 0; ver < NumOfVertices; ++ver)
902 input.read(reinterpret_cast<char*>(&serial_number),
sizeof(
int));
903 input.read(reinterpret_cast<char*>(coord), gmsh_dim*
sizeof(
double));
907 input >> serial_number;
908 for (
int ci = 0; ci < gmsh_dim; ++ci)
913 vertices[ver] =
Vertex(coord, gmsh_dim);
914 vertices_map[serial_number] = ver;
916 if (static_cast<int>(vertices_map.size()) != NumOfVertices)
918 MFEM_ABORT(
"Gmsh file : vertices indices are not unique");
921 else if (buff ==
"$Elements")
923 int num_of_all_elements;
924 input >> num_of_all_elements;
926 getline(input, buff);
937 int nodes_of_gmsh_element[] =
994 vector<Element*> elements_0D, elements_1D, elements_2D, elements_3D;
995 elements_0D.reserve(num_of_all_elements);
996 elements_1D.reserve(num_of_all_elements);
997 elements_2D.reserve(num_of_all_elements);
998 elements_3D.reserve(num_of_all_elements);
1002 int n_elem_part = 0;
1003 const int header_size = 3;
1006 int header[header_size];
1007 int n_elem_one_type;
1009 while (n_elem_part < num_of_all_elements)
1011 input.read(reinterpret_cast<char*>(header),
1012 header_size*
sizeof(
int));
1013 type_of_element = header[0];
1014 n_elem_one_type = header[1];
1017 n_elem_part += n_elem_one_type;
1019 const int n_elem_nodes = nodes_of_gmsh_element[type_of_element-1];
1020 vector<int> data(1+n_tags+n_elem_nodes);
1021 for (
int el = 0; el < n_elem_one_type; ++el)
1023 input.read(reinterpret_cast<char*>(&data[0]),
1024 data.size()*
sizeof(int));
1026 serial_number = data[dd++];
1029 phys_domain = (n_tags > 0) ? data[dd++] : 0;
1032 elem_domain = (n_tags > 1) ? data[dd++] : 0;
1035 n_partitions = (n_tags > 2) ? data[dd++] : 0;
1038 vector<int> vert_indices(n_elem_nodes);
1039 for (
int vi = 0; vi < n_elem_nodes; ++vi)
1041 map<int, int>::const_iterator it =
1042 vertices_map.find(data[1+n_tags+vi]);
1043 if (it == vertices_map.end())
1045 MFEM_ABORT(
"Gmsh file : vertex index doesn't exist");
1047 vert_indices[vi] = it->second;
1051 if (phys_domain <= 0)
1053 MFEM_ABORT(
"Non-positive element attribute in Gmsh mesh!");
1057 switch (type_of_element)
1061 elements_1D.push_back(
1062 new Segment(&vert_indices[0], phys_domain));
1067 elements_2D.push_back(
1068 new Triangle(&vert_indices[0], phys_domain));
1073 elements_2D.push_back(
1079 elements_3D.push_back(
1085 elements_3D.push_back(
1086 new Hexahedron(&vert_indices[0], phys_domain));
1091 elements_0D.push_back(
1092 new Point(&vert_indices[0], phys_domain));
1096 MFEM_WARNING(
"Unsupported Gmsh element type.");
1105 for (
int el = 0; el < num_of_all_elements; ++el)
1107 input >> serial_number >> type_of_element >> n_tags;
1108 vector<int> data(n_tags);
1109 for (
int i = 0; i < n_tags; ++i) { input >> data[i]; }
1112 phys_domain = (n_tags > 0) ? data[0] : 0;
1115 elem_domain = (n_tags > 1) ? data[1] : 0;
1118 n_partitions = (n_tags > 2) ? data[2] : 0;
1121 const int n_elem_nodes = nodes_of_gmsh_element[type_of_element-1];
1122 vector<int> vert_indices(n_elem_nodes);
1124 for (
int vi = 0; vi < n_elem_nodes; ++vi)
1127 map<int, int>::const_iterator it = vertices_map.find(index);
1128 if (it == vertices_map.end())
1130 MFEM_ABORT(
"Gmsh file : vertex index doesn't exist");
1132 vert_indices[vi] = it->second;
1136 if (phys_domain <= 0)
1138 MFEM_ABORT(
"Non-positive element attribute in Gmsh mesh!");
1142 switch (type_of_element)
1146 elements_1D.push_back(
1147 new Segment(&vert_indices[0], phys_domain));
1152 elements_2D.push_back(
1153 new Triangle(&vert_indices[0], phys_domain));
1158 elements_2D.push_back(
1164 elements_3D.push_back(
1170 elements_3D.push_back(
1171 new Hexahedron(&vert_indices[0], phys_domain));
1176 elements_0D.push_back(
1177 new Point(&vert_indices[0], phys_domain));
1181 MFEM_WARNING(
"Unsupported Gmsh element type.");
1188 if (!elements_3D.empty())
1191 NumOfElements = elements_3D.size();
1192 elements.SetSize(NumOfElements);
1193 for (
int el = 0; el < NumOfElements; ++el)
1195 elements[el] = elements_3D[el];
1197 NumOfBdrElements = elements_2D.size();
1198 boundary.SetSize(NumOfBdrElements);
1199 for (
int el = 0; el < NumOfBdrElements; ++el)
1201 boundary[el] = elements_2D[el];
1204 for (
size_t el = 0; el < elements_1D.size(); ++el)
1206 delete elements_1D[el];
1208 for (
size_t el = 0; el < elements_0D.size(); ++el)
1210 delete elements_0D[el];
1213 else if (!elements_2D.empty())
1216 NumOfElements = elements_2D.size();
1217 elements.SetSize(NumOfElements);
1218 for (
int el = 0; el < NumOfElements; ++el)
1220 elements[el] = elements_2D[el];
1222 NumOfBdrElements = elements_1D.size();
1223 boundary.SetSize(NumOfBdrElements);
1224 for (
int el = 0; el < NumOfBdrElements; ++el)
1226 boundary[el] = elements_1D[el];
1229 for (
size_t el = 0; el < elements_0D.size(); ++el)
1231 delete elements_0D[el];
1234 else if (!elements_1D.empty())
1237 NumOfElements = elements_1D.size();
1238 elements.SetSize(NumOfElements);
1239 for (
int el = 0; el < NumOfElements; ++el)
1241 elements[el] = elements_1D[el];
1243 NumOfBdrElements = elements_0D.size();
1244 boundary.SetSize(NumOfBdrElements);
1245 for (
int el = 0; el < NumOfBdrElements; ++el)
1247 boundary[el] = elements_0D[el];
1252 MFEM_ABORT(
"Gmsh file : no elements found");
1256 MFEM_CONTRACT_VAR(n_partitions);
1257 MFEM_CONTRACT_VAR(elem_domain);
1264 #ifdef MFEM_USE_NETCDF
1272 const int sideMapHex8[6][4] =
1282 const int sideMapTet4[4][3] =
1291 const int mfemToGenesisTet10[10] = {1,2,3,4,5,7,8,6,9,10};
1294 const int mfemToGenesisHex27[27] = {1,2,3,4,5,6,7,8,9,10,11,
1296 12,17,18,19,20,13,14,15,
1298 16,22,26,25,27,24,23,21
1301 const int mfemToGenesisTri6[6] = {1,2,3,4,5,6};
1302 const int mfemToGenesisQuad9[9] = {1,2,3,4,5,6,7,8,9};
1304 const int sideMapHex27[6][9] =
1306 {1,2,6,5,9,14,17,13,26},
1307 {2,3,7,6,10,15,18,14,25},
1308 {4,3,7,8,11,15,19,16,27},
1309 {1,4,8,5,12,16,20,13,24},
1310 {1,4,3,2,12,11,10,9,22},
1311 {5,8,7,6,20,19,18,17,23}
1314 const int sideMapTet10[4][6] =
1326 char str_dummy[256];
1333 const char* filename = input.
filename;
1334 if ((retval = nc_open(filename, NC_NOWRITE, &ncid)))
1336 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1342 size_t num_dim=0, num_nodes=0, num_elem=0, num_el_blk=0, num_side_sets=0;
1344 if ((retval = nc_inq_dimid(ncid,
"num_dim", &
id)) ||
1345 (retval = nc_inq_dim(ncid,
id, str_dummy, &num_dim)) ||
1347 (retval = nc_inq_dimid(ncid,
"num_nodes", &
id)) ||
1348 (retval = nc_inq_dim(ncid,
id, str_dummy, &num_nodes)) ||
1350 (retval = nc_inq_dimid(ncid,
"num_elem" ,&
id)) ||
1351 (retval = nc_inq_dim(ncid,
id, str_dummy, &num_elem)) ||
1353 (retval = nc_inq_dimid(ncid,
"num_el_blk", &
id)) ||
1354 (retval = nc_inq_dim(ncid,
id, str_dummy, &num_el_blk)) ||
1356 (retval = nc_inq_dimid(ncid,
"num_side_sets", &
id)) ||
1357 (retval = nc_inq_dim(ncid,
id, str_dummy, &num_side_sets)))
1359 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
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];
1370 int previous_num_nod_per_el;
1371 for (
int i = 0; i < (int) num_el_blk; i++)
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])))
1377 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
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])))
1384 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1391 if ((
int) num_nod_per_el[i] != previous_num_nod_per_el)
1393 MFEM_ABORT(
"Element blocks of different element types not supported");
1396 previous_num_nod_per_el = num_nod_per_el[i];
1400 if (num_nod_per_el[0] == 4 || num_nod_per_el[0] == 8)
1404 else if (num_nod_per_el[0] == 10 || num_nod_per_el[0] == 27)
1410 MFEM_ABORT(
"Don't know what to do with a " << num_nod_per_el[0] <<
1414 for (
int i = 0; i < (int) num_side_sets; i++)
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])))
1420 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1425 double *coordx =
new double[num_nodes];
1426 double *coordy =
new double[num_nodes];
1427 double *coordz =
new double[num_nodes];
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)))
1434 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1439 if ((retval = nc_inq_varid(ncid,
"coordz", &
id)) ||
1440 (retval = nc_get_var_double(ncid,
id, coordz)))
1442 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1447 int **elem_blk =
new int*[num_el_blk];
1448 for (
int i = 0; i < (int) num_el_blk; i++)
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])))
1455 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
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)))
1462 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1467 int **elem_ss =
new int*[num_side_sets];
1468 int **side_ss =
new int*[num_side_sets];
1470 for (
int i = 0; i < (int) num_side_sets; i++)
1472 elem_ss[i] =
new int[num_side_in_ss[i]];
1473 side_ss[i] =
new int[num_side_in_ss[i]];
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])))
1479 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
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])))
1486 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
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)))
1494 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
1499 int num_nod_per_side;
1500 if (num_nod_per_el[0] == 4)
1502 num_nod_per_side = 3;
1504 else if (num_nod_per_el[0] == 8)
1506 num_nod_per_side = 4;
1508 else if (num_nod_per_el[0] == 10)
1510 num_nod_per_side = 6;
1512 else if (num_nod_per_el[0] == 27)
1514 num_nod_per_side = 9;
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++)
1524 start_of_block[i] = start_of_block[i-1] + num_el_in_blk[i-1];
1527 int **ss_node_id =
new int*[num_side_sets];
1529 for (
int i = 0; i < (int) num_side_sets; i++)
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++)
1534 int glob_ind = elem_ss[i][j]-1;
1537 while (iblk < (
int) num_el_blk && glob_ind >= start_of_block[iblk+1])
1541 if (iblk >= (
int) num_el_blk)
1543 MFEM_ABORT(
"Sideset element does not exist");
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)
1550 for (
int k = 0; k < num_nod_per_side; k++)
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];
1557 else if (num_nod_per_side == 6)
1559 for (
int k = 0; k < num_nod_per_side; k++)
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];
1566 else if (num_nod_per_side == 4)
1568 for (
int k = 0; k < num_nod_per_side; k++)
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];
1575 else if (num_nod_per_side == 9)
1577 for (
int k = 0; k < num_nod_per_side; k++)
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];
1589 std::vector<int> uniqueVertexID;
1590 for (
int iblk = 0; iblk < (int) num_el_blk; iblk++)
1592 for (
int i = 0; i < (int) num_el_in_blk[iblk]; i++)
1594 int NumNodePerEl = num_nod_per_el[iblk];
1596 if (NumNodePerEl == 8 || NumNodePerEl == 27)
1600 if (NumNodePerEl == 4 || NumNodePerEl == 10)
1604 for (
int j = 0; j < NumVertPerEl; j++)
1606 uniqueVertexID.push_back(elem_blk[iblk][i*NumNodePerEl + j]);
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));
1619 std::map<int,int> cubitToMFEMVertMap;
1620 for (
int i = 0; i < (int) uniqueVertexID.size(); i++)
1622 cubitToMFEMVertMap[uniqueVertexID[i]] = i+1;
1624 MFEM_ASSERT(cubitToMFEMVertMap.size() == uniqueVertexID.size(),
1625 "This should never happen\n");
1631 NumOfVertices = uniqueVertexID.size();
1632 vertices.SetSize(NumOfVertices);
1633 for (
int i = 0; i < (int) uniqueVertexID.size(); i++)
1635 vertices[i](0) = coordx[uniqueVertexID[i] - 1];
1636 vertices[i](1) = coordy[uniqueVertexID[i] - 1];
1637 vertices[i](2) = coordz[uniqueVertexID[i] - 1];
1640 NumOfElements = num_elem;
1641 elements.SetSize(num_elem);
1643 int renumberedVertID[8];
1644 for (
int iblk = 0; iblk < (int) num_el_blk; iblk++)
1646 int NumNodePerEl = num_nod_per_el[iblk];
1647 if (NumNodePerEl == 4 || NumNodePerEl == 10)
1649 for (
int i = 0; i < (int) num_el_in_blk[iblk]; i++)
1651 for (
int j = 0; j < 4; j++)
1653 renumberedVertID[j] =
1654 cubitToMFEMVertMap[elem_blk[iblk][i*NumNodePerEl+j]]-1;
1656 elements[elcount] =
new Tetrahedron(renumberedVertID,ebprop[iblk]);
1660 else if (NumNodePerEl == 8 || NumNodePerEl == 27)
1662 for (
int i = 0; i < (int) num_el_in_blk[iblk]; i++)
1664 for (
int j = 0; j < 8; j++)
1666 renumberedVertID[j] =
1667 cubitToMFEMVertMap[elem_blk[iblk][i*NumNodePerEl+j]]-1;
1669 elements[elcount] =
new Hexahedron(renumberedVertID,ebprop[iblk]);
1677 NumOfBdrElements = 0;
1678 for (
int iss = 0; iss < (int) num_side_sets; iss++)
1680 NumOfBdrElements += num_side_in_ss[iss];
1682 boundary.SetSize(NumOfBdrElements);
1684 for (
int iss = 0; iss < (int) num_side_sets; iss++)
1686 if (num_nod_per_side == 3 || num_nod_per_side == 6)
1688 for (
int i = 0; i < (int) num_side_in_ss[iss]; i++)
1690 for (
int j = 0; j < 3; j++)
1692 renumberedVertID[j] =
1693 cubitToMFEMVertMap[ss_node_id[iss][i*num_nod_per_side+j]] - 1;
1695 boundary[sidecount] =
new Triangle(renumberedVertID,ssprop[iss]);
1699 else if (num_nod_per_side == 4 || num_nod_per_side == 9)
1701 for (
int i = 0; i < (int) num_side_in_ss[iss]; i++)
1703 for (
int j = 0; j < 4; j++)
1705 renumberedVertID[j] =
1706 cubitToMFEMVertMap[ss_node_id[iss][i*num_nod_per_side+j]] - 1;
1708 boundary[sidecount] =
new Quadrilateral(renumberedVertID,ssprop[iss]);
1719 if (num_nod_per_el[0] == 10)
1721 mymap = (
int *) mfemToGenesisTet10;
1723 else if (num_nod_per_el[0] == 27)
1725 mymap = (
int *) mfemToGenesisHex27;
1727 else if (num_nod_per_el[0] == 6)
1729 mymap = (
int *) mfemToGenesisTri6;
1731 else if (num_nod_per_el[0] == 9)
1733 mymap = (
int *) mfemToGenesisQuad9;
1742 GetElementToFaceTable();
1751 el_to_edge =
new Table;
1752 NumOfEdges = GetElementToEdgeTable(*el_to_edge, be_to_edge);
1763 Nodes->MakeOwner(fec);
1771 for (
int i = 0; i < NumOfElements; i++)
1778 for (
int l = 0; l < dofs.
Size(); l++) { vdofs[l] = dofs[l]; }
1782 while (iblk < (
int) num_el_blk && i >= start_of_block[iblk+1]) { iblk++; }
1783 loc_ind = i - start_of_block[iblk];
1785 for (
int j = 0; j < dofs.
Size(); j++)
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];
1796 delete [] num_el_in_blk;
1797 delete [] num_nod_per_el;
1798 delete [] num_side_in_ss;
1802 for (
int i = 0; i < (int) num_el_blk; i++)
1807 delete [] start_of_block;
1808 for (
int i = 0; i < (int) num_side_sets; i++)
1810 delete [] ss_node_id[i];
1812 delete [] ss_node_id;
1817 #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 if the new size is different.
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 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.
Arbitrary order H1-conforming (continuous) finite elements.
Data type line segment element.
void DofsToVDofs(Array< int > &dofs, int ndofs=-1) const