13 #include "../fem/fem.hpp"
14 #include "../general/text.hpp"
18 #ifdef MFEM_USE_NETCDF
27 void Mesh::ReadMFEMMesh(std::istream &input,
bool mfem_v11,
int &curved)
36 MFEM_VERIFY(ident ==
"dimension",
"invalid mesh file");
42 MFEM_VERIFY(ident ==
"elements",
"invalid mesh file");
43 input >> NumOfElements;
44 elements.SetSize(NumOfElements);
45 for (
int j = 0; j < NumOfElements; j++)
47 elements[j] = ReadElement(input);
53 MFEM_VERIFY(ident ==
"boundary",
"invalid mesh file");
54 input >> NumOfBdrElements;
55 boundary.SetSize(NumOfBdrElements);
56 for (
int j = 0; j < NumOfBdrElements; j++)
58 boundary[j] = ReadElement(input);
64 if (mfem_v11 && ident ==
"vertex_parents")
66 ncmesh =
new NCMesh(
this, &input);
72 if (ident ==
"coarse_elements")
74 ncmesh->LoadCoarseElements(input);
81 MFEM_VERIFY(ident ==
"vertices",
"invalid mesh file");
82 input >> NumOfVertices;
83 vertices.SetSize(NumOfVertices);
89 spaceDim = atoi(ident.c_str());
90 for (
int j = 0; j < NumOfVertices; j++)
92 for (
int i = 0; i < spaceDim; i++)
94 input >> vertices[j](i);
99 if (ncmesh) { ncmesh->SetVertexPositions(vertices); }
109 void Mesh::ReadLineMesh(std::istream &input)
115 input >> NumOfVertices;
116 vertices.SetSize(NumOfVertices);
118 for (j = 0; j < NumOfVertices; j++)
120 input >> vertices[j](0);
123 input >> NumOfElements;
124 elements.SetSize(NumOfElements);
126 for (j = 0; j < NumOfElements; j++)
128 input >> a >> p1 >> p2;
129 elements[j] =
new Segment(p1-1, p2-1, a);
133 input >> NumOfBdrElements;
134 boundary.SetSize(NumOfBdrElements);
135 for (j = 0; j < NumOfBdrElements; j++)
137 input >> a >> ind[0];
139 boundary[j] =
new Point(ind,a);
143 void Mesh::ReadNetgen2DMesh(std::istream &input,
int &curved)
145 int ints[32], attr, n;
151 input >> NumOfBdrElements;
152 boundary.SetSize(NumOfBdrElements);
153 for (
int i = 0; i < NumOfBdrElements; i++)
156 >> ints[0] >> ints[1];
157 ints[0]--; ints[1]--;
158 boundary[i] =
new Segment(ints, attr);
162 input >> NumOfElements;
163 elements.SetSize(NumOfElements);
164 for (
int i = 0; i < NumOfElements; i++)
167 for (
int j = 0; j < n; j++)
175 elements[i] =
new Segment(ints, attr);
178 elements[i] =
new Triangle(ints, attr);
189 input >> NumOfVertices;
190 vertices.SetSize(NumOfVertices);
191 for (
int i = 0; i < NumOfVertices; i++)
192 for (
int j = 0; j < Dim; j++)
194 input >> vertices[i](j);
199 input >> NumOfVertices;
200 vertices.SetSize(NumOfVertices);
205 void Mesh::ReadNetgen3DMesh(std::istream &input)
213 input >> NumOfVertices;
215 vertices.SetSize(NumOfVertices);
216 for (
int i = 0; i < NumOfVertices; i++)
217 for (
int j = 0; j < Dim; j++)
219 input >> vertices[i](j);
223 input >> NumOfElements;
224 elements.SetSize(NumOfElements);
225 for (
int i = 0; i < NumOfElements; i++)
228 for (
int j = 0; j < 4; j++)
233 #ifdef MFEM_USE_MEMALLOC
235 tet = TetMemory.Alloc();
245 input >> NumOfBdrElements;
246 boundary.SetSize(NumOfBdrElements);
247 for (
int i = 0; i < NumOfBdrElements; i++)
250 for (
int j = 0; j < 3; j++)
255 boundary[i] =
new Triangle(ints, attr);
259 void Mesh::ReadTrueGridMesh(std::istream &input)
261 int i, j, ints[32], attr;
262 const int buflen = 1024;
273 input >> vari >> NumOfVertices >> vari >> vari >> NumOfElements;
274 input.getline(buf, buflen);
275 input.getline(buf, buflen);
277 input.getline(buf, buflen);
278 input.getline(buf, buflen);
279 input.getline(buf, buflen);
282 vertices.SetSize(NumOfVertices);
283 for (i = 0; i < NumOfVertices; i++)
285 input >> vari >> varf >> vertices[i](0) >> vertices[i](1);
286 input.getline(buf, buflen);
290 elements.SetSize(NumOfElements);
291 for (i = 0; i < NumOfElements; i++)
293 input >> vari >> attr;
294 for (j = 0; j < 4; j++)
299 input.getline(buf, buflen);
300 input.getline(buf, buflen);
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);
316 vertices.SetSize(NumOfVertices);
317 for (i = 0; i < NumOfVertices; i++)
319 input >> vari >> varf >> vertices[i](0) >> vertices[i](1)
321 input.getline(buf, buflen);
324 elements.SetSize(NumOfElements);
325 for (i = 0; i < NumOfElements; i++)
327 input >> vari >> attr;
328 for (j = 0; j < 8; j++)
333 input.getline(buf, buflen);
337 boundary.SetSize(NumOfBdrElements);
338 for (i = 0; i < NumOfBdrElements; i++)
341 for (j = 0; j < 4; j++)
346 input.getline(buf, buflen);
353 const int Mesh::vtk_quadratic_tet[10] =
354 { 0, 1, 2, 3, 4, 7, 5, 6, 8, 9 };
357 const int Mesh::vtk_quadratic_hex[27] =
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
363 void Mesh::ReadVTKMesh(std::istream &input,
int &curved,
int &read_gf)
368 getline(input, buff);
369 getline(input, buff);
373 MFEM_ABORT(
"VTK mesh is not in ASCII format!");
376 getline(input, buff);
378 if (buff !=
"DATASET UNSTRUCTURED_GRID")
380 MFEM_ABORT(
"VTK mesh is not UNSTRUCTURED_GRID!");
391 MFEM_ABORT(
"VTK mesh does not have POINTS data!");
394 while (buff !=
"POINTS");
400 getline(input, buff);
401 for (i = 0; i < points.
Size(); i++)
408 NumOfElements = n = 0;
413 input >> NumOfElements >> n >> ws;
415 for (i = 0; i < n; i++)
417 input >> cells_data[i];
425 if (buff ==
"CELL_TYPES")
427 input >> NumOfElements;
428 elements.
SetSize(NumOfElements);
429 for (j = i = 0; i < NumOfElements; i++)
437 elements[i] =
new Triangle(&cells_data[j+1]);
445 #ifdef MFEM_USE_MEMALLOC
446 elements[i] = TetMemory.Alloc();
447 elements[i]->SetVertices(&cells_data[j+1]);
454 elements[i] =
new Hexahedron(&cells_data[j+1]);
460 elements[i] =
new Triangle(&cells_data[j+1]);
470 #ifdef MFEM_USE_MEMALLOC
471 elements[i] = TetMemory.Alloc();
472 elements[i]->SetVertices(&cells_data[j+1]);
480 elements[i] =
new Hexahedron(&cells_data[j+1]);
483 MFEM_ABORT(
"VTK mesh : cell type " << ct <<
" is not supported!");
486 j += cells_data[j] + 1;
491 streampos sp = input.tellg();
493 if (buff ==
"CELL_DATA")
496 getline(input, buff);
499 if (!strncmp(buff.c_str(),
"SCALARS material", 16))
501 getline(input, buff);
502 for (i = 0; i < NumOfElements; i++)
505 elements[i]->SetAttribute(attr);
522 vertices.SetSize(np);
523 for (i = 0; i < np; i++)
525 vertices[i](0) = points(3*i+0);
526 vertices[i](1) = points(3*i+1);
527 vertices[i](2) = points(3*i+2);
532 NumOfBdrElements = 0;
541 for (n = i = 0; i < NumOfElements; i++)
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)
552 for (n = i = 0; i < np; i++)
553 if (pts_dof[i] != -1)
558 for (i = 0; i < NumOfElements; i++)
560 int *v = elements[i]->GetVertices();
561 int nv = elements[i]->GetNVertices();
562 for (j = 0; j < nv; j++)
564 v[j] = pts_dof[v[j]];
570 for (i = 0; i < np; i++)
572 if ((j = pts_dof[i]) != -1)
574 vertices[j](0) = points(3*i+0);
575 vertices[j](1) = points(3*i+1);
576 vertices[j](2) = points(3*i+2);
581 NumOfBdrElements = 0;
589 GetElementToFaceTable();
598 el_to_edge =
new Table;
599 NumOfEdges = GetElementToEdgeTable(*el_to_edge, be_to_edge);
609 Nodes->MakeOwner(fec);
614 for (n = i = 0; i < NumOfElements; i++)
618 switch (elements[i]->GetGeometryType())
620 case Geometry::TRIANGLE:
621 case Geometry::SQUARE:
622 vtk_mfem = vtk_quadratic_hex;
break;
623 case Geometry::TETRAHEDRON:
624 vtk_mfem = vtk_quadratic_tet;
break;
627 vtk_mfem = vtk_quadratic_hex;
break;
630 for (n++, j = 0; j < dofs.
Size(); j++, n++)
632 if (pts_dof[cells_data[n]] == -1)
634 pts_dof[cells_data[n]] = dofs[vtk_mfem[j]];
638 if (pts_dof[cells_data[n]] != dofs[vtk_mfem[j]])
640 MFEM_ABORT(
"VTK mesh : inconsistent quadratic mesh!");
647 for (i = 0; i < np; i++)
650 if ((dofs[0] = pts_dof[i]) != -1)
653 for (j = 0; j < dofs.
Size(); j++)
655 (*Nodes)(dofs[j]) = points(3*i+j);
664 void Mesh::ReadNURBSMesh(std::istream &input,
int &curved,
int &read_gf)
668 Dim = NURBSext->Dimension();
669 NumOfVertices = NURBSext->GetNV();
670 NumOfElements = NURBSext->GetNE();
671 NumOfBdrElements = NURBSext->GetNBE();
673 NURBSext->GetElementTopo(elements);
674 NURBSext->GetBdrElementTopo(boundary);
676 vertices.SetSize(NumOfVertices);
678 if (NURBSext->HavePatches())
684 Nodes->MakeOwner(fec);
685 NURBSext->SetCoordsFromPatches(*Nodes);
688 int vd = Nodes->VectorDim();
689 for (
int i = 0; i < vd; i++)
692 Nodes->GetNodalValues(vert_val, i+1);
693 for (
int j = 0; j < NumOfVertices; j++)
695 vertices[j](i) = vert_val(j);
705 void Mesh::ReadInlineMesh(std::istream &input,
int generate_edges)
733 MFEM_VERIFY(input.get() ==
'=',
734 "Inline mesh expected '=' after keyword " << name);
741 else if (name ==
"ny")
745 else if (name ==
"nz")
749 else if (name ==
"sx")
753 else if (name ==
"sy")
757 else if (name ==
"sz")
761 else if (name ==
"type")
765 if (eltype ==
"segment")
767 type = Element::SEGMENT;
769 else if (eltype ==
"quad")
771 type = Element::QUADRILATERAL;
773 else if (eltype ==
"tri")
775 type = Element::TRIANGLE;
777 else if (eltype ==
"hex")
779 type = Element::HEXAHEDRON;
781 else if (eltype ==
"tet")
783 type = Element::TETRAHEDRON;
787 MFEM_ABORT(
"unrecognized element type (read '" << eltype
788 <<
"') in inline mesh format. "
789 "Allowed: segment, tri, tet, quad, hex");
794 MFEM_ABORT(
"unrecognized keyword (" << name
795 <<
") in inline mesh format. "
796 "Allowed: nx, ny, nz, type, sx, sy, sz");
801 if (input.peek() ==
';')
814 if (type == Element::SEGMENT)
816 MFEM_VERIFY(nx > 0 && sx > 0.0,
817 "invalid 1D inline mesh format, all values must be "
819 <<
" nx = " << nx <<
"\n"
820 <<
" sx = " << sx <<
"\n");
823 else if (type == Element::TRIANGLE || type == Element::QUADRILATERAL)
825 MFEM_VERIFY(nx > 0 && ny > 0 && sx > 0.0 && sy > 0.0,
826 "invalid 2D inline mesh format, all values must be "
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);
834 else if (type == Element::TETRAHEDRON || type == Element::HEXAHEDRON)
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 "
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);
850 MFEM_ABORT(
"For inline mesh, must specify an element type ="
851 " [segment, tri, quad, tet, hex]");
857 void Mesh::ReadGmshMesh(std::istream &input)
862 input >> version >> binary >> dsize;
865 MFEM_ABORT(
"Gmsh file version < 2.2");
867 if (dsize !=
sizeof(
double))
869 MFEM_ABORT(
"Gmsh file : dsize != sizeof(double)");
871 getline(input, buff);
876 input.read(reinterpret_cast<char*>(&one),
sizeof(one));
879 MFEM_ABORT(
"Gmsh file : wrong binary format");
886 map<int, int> vertices_map;
889 while (input >> buff)
891 if (buff ==
"$Nodes")
893 input >> NumOfVertices;
894 getline(input, buff);
895 vertices.SetSize(NumOfVertices);
897 const int gmsh_dim = 3;
898 double coord[gmsh_dim];
899 for (
int ver = 0; ver < NumOfVertices; ++ver)
903 input.read(reinterpret_cast<char*>(&serial_number),
sizeof(
int));
904 input.read(reinterpret_cast<char*>(coord), gmsh_dim*
sizeof(
double));
908 input >> serial_number;
909 for (
int ci = 0; ci < gmsh_dim; ++ci)
914 vertices[ver] =
Vertex(coord, gmsh_dim);
915 vertices_map[serial_number] = ver;
917 if (static_cast<int>(vertices_map.size()) != NumOfVertices)
919 MFEM_ABORT(
"Gmsh file : vertices indices are not unique");
922 else if (buff ==
"$Elements")
924 int num_of_all_elements;
925 input >> num_of_all_elements;
927 getline(input, buff);
938 int nodes_of_gmsh_element[] =
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);
1003 int n_elem_part = 0;
1004 const int header_size = 3;
1007 int header[header_size];
1008 int n_elem_one_type;
1010 while (n_elem_part < num_of_all_elements)
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];
1018 n_elem_part += n_elem_one_type;
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)
1024 input.read(reinterpret_cast<char*>(&data[0]),
1025 data.size()*
sizeof(int));
1027 serial_number = data[dd++];
1030 phys_domain = (n_tags > 0) ? data[dd++] : 0;
1033 elem_domain = (n_tags > 1) ? data[dd++] : 0;
1036 n_partitions = (n_tags > 2) ? data[dd++] : 0;
1039 vector<int> vert_indices(n_elem_nodes);
1040 for (
int vi = 0; vi < n_elem_nodes; ++vi)
1042 map<int, int>::const_iterator it =
1043 vertices_map.find(data[1+n_tags+vi]);
1044 if (it == vertices_map.end())
1046 MFEM_ABORT(
"Gmsh file : vertex index doesn't exist");
1048 vert_indices[vi] = it->second;
1052 if (phys_domain <= 0)
1054 MFEM_ABORT(
"Non-positive element attribute in Gmsh mesh!");
1058 switch (type_of_element)
1062 elements_1D.push_back(
1063 new Segment(&vert_indices[0], phys_domain));
1068 elements_2D.push_back(
1069 new Triangle(&vert_indices[0], phys_domain));
1074 elements_2D.push_back(
1080 elements_3D.push_back(
1086 elements_3D.push_back(
1087 new Hexahedron(&vert_indices[0], phys_domain));
1092 elements_0D.push_back(
1093 new Point(&vert_indices[0], phys_domain));
1097 MFEM_WARNING(
"Unsupported Gmsh element type.");
1106 for (
int el = 0; el < num_of_all_elements; ++el)
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]; }
1113 phys_domain = (n_tags > 0) ? data[0] : 0;
1116 elem_domain = (n_tags > 1) ? data[1] : 0;
1119 n_partitions = (n_tags > 2) ? data[2] : 0;
1122 const int n_elem_nodes = nodes_of_gmsh_element[type_of_element-1];
1123 vector<int> vert_indices(n_elem_nodes);
1125 for (
int vi = 0; vi < n_elem_nodes; ++vi)
1128 map<int, int>::const_iterator it = vertices_map.find(index);
1129 if (it == vertices_map.end())
1131 MFEM_ABORT(
"Gmsh file : vertex index doesn't exist");
1133 vert_indices[vi] = it->second;
1137 if (phys_domain <= 0)
1139 MFEM_ABORT(
"Non-positive element attribute in Gmsh mesh!");
1143 switch (type_of_element)
1147 elements_1D.push_back(
1148 new Segment(&vert_indices[0], phys_domain));
1153 elements_2D.push_back(
1154 new Triangle(&vert_indices[0], phys_domain));
1159 elements_2D.push_back(
1165 elements_3D.push_back(
1171 elements_3D.push_back(
1172 new Hexahedron(&vert_indices[0], phys_domain));
1177 elements_0D.push_back(
1178 new Point(&vert_indices[0], phys_domain));
1182 MFEM_WARNING(
"Unsupported Gmsh element type.");
1189 if (!elements_3D.empty())
1192 NumOfElements = elements_3D.size();
1193 elements.SetSize(NumOfElements);
1194 for (
int el = 0; el < NumOfElements; ++el)
1196 elements[el] = elements_3D[el];
1198 NumOfBdrElements = elements_2D.size();
1199 boundary.SetSize(NumOfBdrElements);
1200 for (
int el = 0; el < NumOfBdrElements; ++el)
1202 boundary[el] = elements_2D[el];
1205 for (
size_t el = 0; el < elements_1D.size(); ++el)
1207 delete elements_1D[el];
1209 for (
size_t el = 0; el < elements_0D.size(); ++el)
1211 delete elements_0D[el];
1214 else if (!elements_2D.empty())
1217 NumOfElements = elements_2D.size();
1218 elements.SetSize(NumOfElements);
1219 for (
int el = 0; el < NumOfElements; ++el)
1221 elements[el] = elements_2D[el];
1223 NumOfBdrElements = elements_1D.size();
1224 boundary.SetSize(NumOfBdrElements);
1225 for (
int el = 0; el < NumOfBdrElements; ++el)
1227 boundary[el] = elements_1D[el];
1230 for (
size_t el = 0; el < elements_0D.size(); ++el)
1232 delete elements_0D[el];
1235 else if (!elements_1D.empty())
1238 NumOfElements = elements_1D.size();
1239 elements.SetSize(NumOfElements);
1240 for (
int el = 0; el < NumOfElements; ++el)
1242 elements[el] = elements_1D[el];
1244 NumOfBdrElements = elements_0D.size();
1245 boundary.SetSize(NumOfBdrElements);
1246 for (
int el = 0; el < NumOfBdrElements; ++el)
1248 boundary[el] = elements_0D[el];
1253 MFEM_ABORT(
"Gmsh file : no elements found");
1257 MFEM_CONTRACT_VAR(n_partitions);
1258 MFEM_CONTRACT_VAR(elem_domain);
1265 #ifdef MFEM_USE_NETCDF
1266 void Mesh::ReadCubit(
const char *filename,
int &curved,
int &read_gf)
1273 const int sideMapHex8[6][4] =
1283 const int sideMapTet4[4][3] =
1292 const int mfemToGenesisTet10[10] = {1,2,3,4,5,7,8,6,9,10};
1295 const int mfemToGenesisHex27[27] = {1,2,3,4,5,6,7,8,9,10,11,
1297 12,17,18,19,20,13,14,15,
1299 16,22,26,25,27,24,23,21
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};
1305 const int sideMapHex27[6][9] =
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}
1315 const int sideMapTet10[4][6] =
1327 char str_dummy[256];
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 = 0;
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 = 0;
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];
1595 int NumVertPerEl = 0;
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 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