17#include "../general/tinyxml2.h"
44 MFEM_VERIFY(version == 10 || version == 12 || version == 13,
45 "unknown MFEM mesh version");
53 MFEM_VERIFY(ident ==
"dimension",
"invalid mesh file");
59 MFEM_VERIFY(ident ==
"elements",
"invalid mesh file");
72 MFEM_VERIFY(ident ==
"attribute_sets",
"invalid mesh file");
82 MFEM_VERIFY(ident ==
"boundary",
"invalid mesh file");
95 MFEM_VERIFY(ident ==
"bdr_attribute_sets",
"invalid mesh file");
105 MFEM_VERIFY(ident ==
"vertices",
"invalid mesh file");
109 input >> ws >> ident;
110 if (ident !=
"nodes")
153 input >>
a >> p1 >> p2;
162 input >>
a >> ind[0];
170 int ints[32], attr, n;
181 >> ints[0] >> ints[1];
182 ints[0]--; ints[1]--;
192 for (
int j = 0; j < n; j++)
217 for (
int j = 0; j <
Dim; j++)
242 for (
int j = 0; j <
Dim; j++)
253 for (
int j = 0; j < 4; j++)
258#ifdef MFEM_USE_MEMALLOC
275 for (
int j = 0; j < 3; j++)
286 int i, j, ints[32], attr;
287 const int buflen = 1024;
299 input.getline(buf, buflen);
300 input.getline(buf, buflen);
302 input.getline(buf, buflen);
303 input.getline(buf, buflen);
304 input.getline(buf, buflen);
311 input.getline(buf, buflen);
318 input >> vari >> attr;
319 for (j = 0; j < 4; j++)
324 input.getline(buf, buflen);
325 input.getline(buf, buflen);
334 input.getline(buf, buflen);
335 input.getline(buf, buflen);
337 input.getline(buf, buflen);
338 input.getline(buf, buflen);
339 input.getline(buf, buflen);
346 input.getline(buf, buflen);
352 input >> vari >> attr;
353 for (j = 0; j < 8; j++)
358 input.getline(buf, buflen);
366 for (j = 0; j < 4; j++)
371 input.getline(buf, buflen);
379{ 0, 1, 2, 3, 4, 7, 5, 6, 8, 9 };
384{ 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12};
389{ 0, 2, 1, 3, 5, 4, 8, 7, 6, 11, 10, 9, 12, 14, 13, 17, 16, 15};
394 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19,
395 24, 22, 21, 23, 20, 25, 26
402 int &curved,
int &read_gf,
bool &finalize_topo)
404 int np = points.
Size()/3;
410 bool legacy_elem =
false, lagrange_elem =
false;
414 int j = (i > 0) ? cell_offsets[i-1] : 0;
415 int ct = cell_types[i];
418 if (cell_attributes.
Size() > 0)
420 elements[i]->SetAttribute(cell_attributes[i]);
426 int prism_vertices[6];
427 for (
int k=0; k<6; ++k)
431 elements[i]->SetVertices(prism_vertices);
435 elements[i]->SetVertices(&cell_data[j]);
442 else { legacy_elem =
true; }
444 MFEM_VERIFY(
Dim == -1 ||
Dim == elem_dim,
445 "Elements with different dimensions are not supported");
446 MFEM_VERIFY(order == -1 || order == elem_order,
447 "Elements with different orders are not supported");
448 MFEM_VERIFY(legacy_elem != lagrange_elem,
449 "Mixing of legacy and Lagrange cell types is not supported");
458 real_t min_value, max_value;
459 for (
int d = 3; d > 0; --d)
461 min_value = max_value = points(3*0 + d-1);
462 for (
int i = 1; i < np; i++)
464 min_value = std::min(min_value, points(3*i + d-1));
465 max_value = std::max(max_value, points(3*i + d-1));
466 if (min_value != max_value)
476 if (order == 1 && !lagrange_elem)
480 for (
int i = 0; i < np; i++)
503 int *v =
elements[i]->GetVertices();
504 int nv =
elements[i]->GetNVertices();
505 for (
int j = 0; j < nv; j++)
507 if (pts_dof[v[j]] == -1) { pts_dof[v[j]] = 0; }
516 for (
int i = 0; i < np; i++)
518 if (pts_dof[i] != -1)
526 int *v =
elements[i]->GetVertices();
527 int nv =
elements[i]->GetNVertices();
528 for (
int j = 0; j < nv; j++)
530 v[j] = pts_dof[v[j]];
535 for (
int i = 0; i < np; i++)
569 switch (
elements[i]->GetGeometryType())
587 int offset = (i == 0) ? 0 : cell_offsets[i-1];
588 for (
int j = 0; j < dofs.
Size(); j++)
590 if (pts_dof[cell_data[offset+j]] == -1)
592 pts_dof[cell_data[offset+j]] = dofs[vtk_mfem[j]];
596 if (pts_dof[cell_data[offset+j]] != dofs[vtk_mfem[j]])
598 MFEM_ABORT(
"VTK mesh: inconsistent quadratic mesh!");
614 std::map<Geometry::Type,Array<int>> vtk_inv_maps;
615 std::map<Geometry::Type,const Array<int>*> lex_orderings;
624 if (vtk_inv_map.
Size() == 0)
629 for (
int j=0; j<vtk_map.
Size(); ++j)
631 vtk_inv_map[vtk_map[j]] = j;
634 const Array<int> *&lex_ordering = lex_orderings[geom];
640 MFEM_ASSERT(nodal_fe != NULL,
"Unsupported element type");
644 for (
int lex_idx = 0; lex_idx < dofs.
Size(); lex_idx++)
646 int mfem_idx = (*lex_ordering)[lex_idx];
647 int vtk_idx = vtk_inv_map[lex_idx];
648 int pt_idx = cell_data[n + vtk_idx];
649 if (pts_dof[pt_idx] == -1)
651 pts_dof[pt_idx] = dofs[mfem_idx];
655 if (pts_dof[pt_idx] != dofs[mfem_idx])
657 MFEM_ABORT(
"VTK mesh: inconsistent Lagrange mesh!");
666 for (
int i = 0; i < np; i++)
669 if (pts_dof[i] != -1)
671 dofs[0] = pts_dof[i];
673 for (
int d = 0; d < dofs.
Size(); d++)
675 (*Nodes)(dofs[d]) = points(3*i+d);
686using namespace tinyxml2;
692 if (s1 == NULL || s2 == NULL) {
return false; }
693 return strcmp(s1, s2) == 0;
701struct BufferReaderBase
703 enum HeaderType { UINT32_HEADER, UINT64_HEADER };
704 virtual void ReadBinary(
const char *buf,
void *dest,
int n)
const = 0;
705 virtual void ReadBase64(
const char *txt,
void *dest,
int n)
const = 0;
706 virtual ~BufferReaderBase() { }
719template <
typename T,
typename F>
720struct BufferReader : BufferReaderBase
723 HeaderType header_type;
724 BufferReader(
bool compressed_, HeaderType header_type_)
725 : compressed(compressed_), header_type(header_type_) { }
728 size_t HeaderEntrySize()
const
730 return header_type == UINT64_HEADER ?
sizeof(uint64_t) : sizeof(uint32_t);
736 uint64_t ReadHeaderEntry(
const char *header_buf)
const
739 : bin_io::
read<uint32_t>(header_buf);
748 int NumHeaderBytes(
const char *header_buf)
const
750 if (!compressed) {
return static_cast<int>(HeaderEntrySize()); }
751 return (3 + ReadHeaderEntry(header_buf))*HeaderEntrySize();
759 void ReadBinaryWithHeader(
const char *header_buf,
const char *buf,
760 void *dest_void,
int n)
const
762 std::vector<char> uncompressed_data;
763 T *dest =
static_cast<T*
>(dest_void);
773 int header_entry_size = HeaderEntrySize();
774 int nblocks = ReadHeaderEntry(header_buf);
775 header_buf += header_entry_size;
776 std::vector<size_t> header(nblocks + 2);
777 for (
int i=0; i<nblocks+2; ++i)
779 header[i] = ReadHeaderEntry(header_buf);
780 header_buf += header_entry_size;
782 uncompressed_data.resize((nblocks-1)*header[0] + header[1]);
783 Bytef *dest_ptr = (Bytef *)uncompressed_data.data();
784 Bytef *dest_start = dest_ptr;
785 const Bytef *source_ptr = (
const Bytef *)buf;
786 for (
int i=0; i<nblocks; ++i)
788 uLongf source_len = header[i+2];
789 uLong dest_len = (i == nblocks-1) ? header[1] : header[0];
790 int res = uncompress(dest_ptr, &dest_len, source_ptr, source_len);
791 MFEM_VERIFY(res == Z_OK,
"Error uncompressing");
792 dest_ptr += dest_len;
793 source_ptr += source_len;
795 MFEM_VERIFY(
size_t(
sizeof(F)*n) ==
size_t(dest_ptr - dest_start),
796 "AppendedData: wrong data size");
797 buf = uncompressed_data.data();
799 MFEM_ABORT(
"MFEM must be compiled with zlib enabled to uncompress.")
806 MFEM_VERIFY(
sizeof(F)*n == ReadHeaderEntry(header_buf),
807 "AppendedData: wrong data size");
810 if (std::is_same_v<T, F>)
813 memcpy(dest, buf,
sizeof(T)*n);
817 for (
int i=0; i<n; ++i)
828 void ReadBinary(
const char *buf,
void *dest,
int n)
const override
830 ReadBinaryWithHeader(buf, buf + NumHeaderBytes(buf), dest, n);
837 void ReadBase64(
const char *txt,
void *dest,
int n)
const override
842 if (*txt !=
' ' && *txt !=
'\n') {
break; }
849 std::vector<char> nblocks_buf;
852 std::vector<char> data, header;
856 nblocks_buf.data())));
860 ReadBinaryWithHeader(header.data(), data.data(), dest, n);
864 std::vector<char> data;
866 ReadBinary(data.data(), dest, n);
877 const char *appended_data, *byte_order, *compressor;
878 enum AppendedDataEncoding { RAW, BASE64 };
879 map<string,BufferReaderBase*> type_map;
880 AppendedDataEncoding encoding;
886 XMLDataReader(
const XMLElement *vtk,
const XMLElement *vtu)
889 BufferReaderBase::HeaderType htype;
892 htype = BufferReaderBase::UINT64_HEADER;
896 htype = BufferReaderBase::UINT32_HEADER;
900 byte_order = vtk->Attribute(
"byte_order");
904 compressor = vtk->Attribute(
"compressor");
905 bool compressed = (compressor != NULL);
908 appended_data = NULL;
909 for (
const XMLElement *xml_elem = vtu->NextSiblingElement();
911 xml_elem = xml_elem->NextSiblingElement())
915 const char *encoding_str = xml_elem->Attribute(
"encoding");
918 appended_data = xml_elem->GetAppendedData();
923 appended_data = xml_elem->GetText();
926 MFEM_VERIFY(appended_data != NULL,
"Invalid AppendedData");
928 bool found_leading_underscore =
false;
929 while (*appended_data)
932 if (*appended_data ==
'_')
934 found_leading_underscore =
true;
939 MFEM_VERIFY(found_leading_underscore,
"Invalid AppendedData");
944 type_map[
"Int8"] =
new BufferReader<int, int8_t>(compressed, htype);
945 type_map[
"Int16"] =
new BufferReader<int, int16_t>(compressed, htype);
946 type_map[
"Int32"] =
new BufferReader<int, int32_t>(compressed, htype);
947 type_map[
"Int64"] =
new BufferReader<int, int64_t>(compressed, htype);
948 type_map[
"UInt8"] =
new BufferReader<int, uint8_t>(compressed, htype);
949 type_map[
"UInt16"] =
new BufferReader<int, uint16_t>(compressed, htype);
950 type_map[
"UInt32"] =
new BufferReader<int, uint32_t>(compressed, htype);
951 type_map[
"UInt64"] =
new BufferReader<int, uint64_t>(compressed, htype);
952 type_map[
"Float32"] =
new BufferReader<double, float>(compressed, htype);
953 type_map[
"Float64"] =
new BufferReader<double, double>(compressed, htype);
959 template <
typename T>
960 void Read(
const XMLElement *xml_elem, T *dest,
int n)
962 static const char *erstr =
"Error reading XML DataArray";
963 MFEM_VERIFY(
StringCompare(xml_elem->Name(),
"DataArray"), erstr);
964 const char *format = xml_elem->Attribute(
"format");
967 const char *txt = xml_elem->GetText();
968 MFEM_VERIFY(txt != NULL, erstr);
969 std::istringstream data_stream(txt);
970 for (
int i=0; i<n; ++i) { data_stream >> dest[i]; }
974 VerifyBinaryOptions();
975 int offset = xml_elem->IntAttribute(
"offset");
976 const char *
type = xml_elem->Attribute(
"type");
977 MFEM_VERIFY(type != NULL, erstr);
978 BufferReaderBase *reader = type_map[
type];
979 MFEM_VERIFY(reader != NULL, erstr);
980 MFEM_VERIFY(appended_data != NULL,
"No AppendedData found");
983 reader->ReadBinary(appended_data + offset, dest, n);
987 reader->ReadBase64(appended_data + offset, dest, n);
992 VerifyBinaryOptions();
993 const char *txt = xml_elem->GetText();
994 MFEM_VERIFY(txt != NULL, erstr);
995 const char *
type = xml_elem->Attribute(
"type");
996 if (type == NULL) { MFEM_ABORT(erstr); }
997 BufferReaderBase *reader = type_map[
type];
998 if (reader == NULL) { MFEM_ABORT(erstr); }
999 reader->ReadBase64(txt, dest, n);
1003 MFEM_ABORT(
"Invalid XML VTK DataArray format");
1011 void VerifyByteOrder()
const
1016 MFEM_ABORT(
"Converting between different byte orders is unsupported.");
1023 void VerifyCompressor()
const
1025 if (compressor && !
StringCompare(compressor,
"vtkZLibDataCompressor"))
1027 MFEM_ABORT(
"Unsupported compressor. Only zlib is supported.")
1029#ifndef MFEM_USE_ZLIB
1030 MFEM_VERIFY(compressor == NULL,
"MFEM must be compiled with zlib enabled "
1031 "to support reading compressed data.");
1037 void VerifyBinaryOptions()
const
1045 for (
auto &x : type_map) {
delete x.second; }
1052 bool &finalize_topo,
const std::string &xml_prefix)
1054 using namespace vtk_xml;
1056 static const char *erstr =
"XML parsing error";
1059 std::vector<char> buf(xml_prefix.begin(), xml_prefix.end());
1060 std::istreambuf_iterator<char> eos;
1061 buf.insert(buf.end(), std::istreambuf_iterator<char>(input), eos);
1062 buf.push_back(
'\0');
1065 xml.Parse(buf.data(), buf.size());
1066 if (xml.ErrorID() != XML_SUCCESS)
1068 MFEM_ABORT(
"Error parsing XML VTK file.\n" << xml.ErrorStr());
1071 const XMLElement *vtkfile = xml.FirstChildElement();
1072 MFEM_VERIFY(vtkfile, erstr);
1073 MFEM_VERIFY(StringCompare(vtkfile->Name(),
"VTKFile"), erstr);
1074 const XMLElement *vtu = vtkfile->FirstChildElement();
1075 MFEM_VERIFY(vtu, erstr);
1076 MFEM_VERIFY(StringCompare(vtu->Name(),
"UnstructuredGrid"), erstr);
1078 XMLDataReader data_reader(vtkfile, vtu);
1081 const XMLElement *piece = vtu->FirstChildElement();
1082 MFEM_VERIFY(StringCompare(piece->Name(),
"Piece"), erstr);
1083 MFEM_VERIFY(piece->NextSiblingElement() == NULL,
1084 "XML VTK meshes with more than one Piece are not supported");
1085 int npts = piece->IntAttribute(
"NumberOfPoints");
1086 int ncells = piece->IntAttribute(
"NumberOfCells");
1090 const XMLElement *pts_xml;
1091 for (pts_xml = piece->FirstChildElement();
1093 pts_xml = pts_xml->NextSiblingElement())
1095 if (StringCompare(pts_xml->Name(),
"Points"))
1097 const XMLElement *pts_data = pts_xml->FirstChildElement();
1098 MFEM_VERIFY(pts_data->IntAttribute(
"NumberOfComponents") == 3,
1099 "XML VTK Points DataArray must have 3 components");
1100 data_reader.Read(pts_data, points.
GetData(), points.
Size());
1104 if (pts_xml == NULL) { MFEM_ABORT(erstr); }
1107 Array<int> cell_data, cell_offsets(ncells), cell_types(ncells);
1108 const XMLElement *cells_xml;
1109 for (cells_xml = piece->FirstChildElement();
1111 cells_xml = cells_xml->NextSiblingElement())
1113 if (StringCompare(cells_xml->Name(),
"Cells"))
1115 const XMLElement *cell_data_xml = NULL;
1116 for (
const XMLElement *data_xml = cells_xml->FirstChildElement();
1118 data_xml = data_xml->NextSiblingElement())
1120 const char *data_name = data_xml->Attribute(
"Name");
1121 if (StringCompare(data_name,
"offsets"))
1123 data_reader.Read(data_xml, cell_offsets.GetData(), ncells);
1125 else if (StringCompare(data_name,
"types"))
1127 data_reader.Read(data_xml, cell_types.
GetData(), ncells);
1129 else if (StringCompare(data_name,
"connectivity"))
1135 cell_data_xml = data_xml;
1138 MFEM_VERIFY(cell_data_xml != NULL, erstr);
1139 int cell_data_size = cell_offsets.Last();
1140 cell_data.
SetSize(cell_data_size);
1141 data_reader.Read(cell_data_xml, cell_data.
GetData(), cell_data_size);
1145 if (cells_xml == NULL) { MFEM_ABORT(erstr); }
1151 bool found_attributes =
false;
1152 for (
const XMLElement *cell_data_xml = piece->FirstChildElement();
1153 cell_data_xml != NULL;
1154 cell_data_xml = cell_data_xml->NextSiblingElement())
1156 const bool is_cell_data =
1157 StringCompare(cell_data_xml->Name(),
"CellData");
1158 const bool is_material =
1159 StringCompare(cell_data_xml->Attribute(
"Scalars"),
"material");
1160 const bool is_attribute =
1161 StringCompare(cell_data_xml->Attribute(
"Scalars"),
"attribute");
1162 if (is_cell_data && (is_material || (is_attribute && !found_attributes)))
1164 found_attributes =
true;
1165 const XMLElement *data_xml = cell_data_xml->FirstChildElement();
1166 if (data_xml != NULL && StringCompare(data_xml->Name(),
"DataArray"))
1168 cell_attributes.
SetSize(ncells);
1169 data_reader.Read(data_xml, cell_attributes.
GetData(), ncells);
1174 CreateVTKMesh(points, cell_data, cell_offsets, cell_types, cell_attributes,
1175 curved, read_gf, finalize_topo);
1179 bool &finalize_topo)
1188 getline(input, buff);
1189 getline(input, buff);
1191 if (buff !=
"ASCII")
1193 MFEM_ABORT(
"VTK mesh is not in ASCII format!");
1198 getline(input, buff);
1200 if (!input.good()) { MFEM_ABORT(
"VTK mesh is not UNSTRUCTURED_GRID!"); }
1202 while (buff !=
"DATASET UNSTRUCTURED_GRID");
1211 MFEM_ABORT(
"VTK mesh does not have POINTS data!");
1214 while (buff !=
"POINTS");
1219 getline(input, buff);
1220 points.
Load(input, 3*np);
1235 MFEM_ABORT(
"VTK mesh does not have CELLS data!");
1238 while (buff !=
"CELLS");
1242 if (buff ==
"CELLS")
1245 input >> ncells >> n >> ws;
1247 cell_data.
SetSize(n - ncells);
1249 for (
int i=0; i<ncells; ++i)
1253 cell_offsets[i] = offset + nv;
1254 for (
int j=0; j<nv; ++j)
1256 input >> cell_data[offset + j];
1263 input >> ws >> buff;
1266 MFEM_VERIFY(buff ==
"CELL_TYPES",
"CELL_TYPES not provided in VTK mesh.")
1268 cell_types.
Load(ncells, input);
1270 while ((input.good()) && (buff !=
"CELL_DATA"))
1274 getline(input, buff);
1279 bool found_attributes =
false;
1280 while ((input.good()))
1282 getline(input, buff);
1283 if (buff.rfind(
"POINT_DATA") == 0)
1287 else if (buff.rfind(
"SCALARS material") == 0 ||
1288 (buff.rfind(
"SCALARS attribute") == 0 && !found_attributes))
1290 found_attributes =
true;
1291 getline(input, buff);
1292 if (buff.rfind(
"LOOKUP_TABLE default") != 0)
1294 MFEM_ABORT(
"Invalid LOOKUP_TABLE for material array in VTK file.");
1296 cell_attributes.
Load(ncells, input);
1308 CreateVTKMesh(points, cell_data, cell_offsets, cell_types, cell_attributes,
1309 curved, read_gf, finalize_topo);
1313 bool spacing,
bool nc)
1383 MFEM_VERIFY(input.get() ==
'=',
1384 "Inline mesh expected '=' after keyword " << name);
1391 else if (name ==
"ny")
1395 else if (name ==
"nz")
1399 else if (name ==
"sx")
1403 else if (name ==
"sy")
1407 else if (name ==
"sz")
1411 else if (name ==
"type")
1415 if (eltype ==
"segment")
1419 else if (eltype ==
"quad")
1423 else if (eltype ==
"tri")
1427 else if (eltype ==
"hex")
1431 else if (eltype ==
"wedge")
1435 else if (eltype ==
"pyramid")
1439 else if (eltype ==
"tet")
1445 MFEM_ABORT(
"unrecognized element type (read '" << eltype
1446 <<
"') in inline mesh format. "
1447 "Allowed: segment, tri, quad, tet, hex, wedge");
1452 MFEM_ABORT(
"unrecognized keyword (" << name
1453 <<
") in inline mesh format. "
1454 "Allowed: nx, ny, nz, type, sx, sy, sz");
1459 if (input.peek() ==
';')
1474 MFEM_VERIFY(nx > 0 && sx > 0.0,
1475 "invalid 1D inline mesh format, all values must be "
1477 <<
" nx = " << nx <<
"\n"
1478 <<
" sx = " << sx <<
"\n");
1483 MFEM_VERIFY(nx > 0 && ny > 0 && sx > 0.0 && sy > 0.0,
1484 "invalid 2D inline mesh format, all values must be "
1486 <<
" nx = " << nx <<
"\n"
1487 <<
" ny = " << ny <<
"\n"
1488 <<
" sx = " << sx <<
"\n"
1489 <<
" sy = " << sy <<
"\n");
1490 Make2D(nx, ny, type, sx, sy, generate_edges,
true);
1495 MFEM_VERIFY(nx > 0 && ny > 0 && nz > 0 &&
1496 sx > 0.0 && sy > 0.0 && sz > 0.0,
1497 "invalid 3D inline mesh format, all values must be "
1499 <<
" nx = " << nx <<
"\n"
1500 <<
" ny = " << ny <<
"\n"
1501 <<
" nz = " << nz <<
"\n"
1502 <<
" sx = " << sx <<
"\n"
1503 <<
" sy = " << sy <<
"\n"
1504 <<
" sz = " << sz <<
"\n");
1505 Make3D(nx, ny, nz, type, sx, sy, sz,
true);
1510 MFEM_ABORT(
"For inline mesh, must specify an element type ="
1511 " [segment, tri, quad, tet, hex, wedge]");
1520 input >> version >> binary >> dsize;
1523 MFEM_ABORT(
"Gmsh file version < 2.2");
1525 if (dsize !=
sizeof(
double))
1527 MFEM_ABORT(
"Gmsh file : dsize != sizeof(double)");
1529 getline(input, buff);
1534 input.read(
reinterpret_cast<char*
>(&one),
sizeof(one));
1537 MFEM_ABORT(
"Gmsh file : wrong binary format");
1544 map<int, int> vertices_map;
1550 map<int,map<int,std::string> > phys_names_by_dim;
1570 bool periodic =
false;
1577 while (input >> buff)
1579 if (buff ==
"$Nodes")
1582 getline(input, buff);
1585 const int gmsh_dim = 3;
1591 input.read(
reinterpret_cast<char*
>(&serial_number),
sizeof(
int));
1592 input.read(
reinterpret_cast<char*
>(coord), gmsh_dim*
sizeof(
double));
1596 input >> serial_number;
1597 for (
int ci = 0; ci < gmsh_dim; ++ci)
1603 vertices_map[serial_number] = ver;
1605 for (
int ci = 0; ci < gmsh_dim; ++ci)
1607 bb_min[ci] = (ver == 0) ? coord[ci] :
1608 std::min(
bb_min[ci], coord[ci]);
1609 bb_max[ci] = (ver == 0) ? coord[ci] :
1610 std::max(
bb_max[ci], coord[ci]);
1628 MFEM_ABORT(
"Gmsh file : vertices indices are not unique");
1631 else if (buff ==
"$Elements")
1633 int num_of_all_elements;
1634 input >> num_of_all_elements;
1636 getline(input, buff);
1639 int type_of_element;
1647 int nodes_of_gmsh_element[] =
1796 -1,-1,-1,-1,-1,-1,-1,
1812 -1,-1,-1,-1,-1,-1,-1,
1857 int lin3[] = {0,2,1};
1858 int lin4[] = {0,2,3,1};
1859 int tri6[] = {0,3,1,5,4,2};
1860 int tri10[] = {0,3,4,1,8,9,5,7,6,2};
1861 int quad9[] = {0,4,1,7,8,5,3,6,2};
1862 int quad16[] = {0,4,5,1,11,12,13,6,
1865 int tet10[] {0,4,1,6,5,2,7,9,8,3};
1866 int tet20[] = {0,4,5,1,9,16,6,8,7,2,
1867 11,17,15,18,19,13,10,14,12,3
1869 int hex27[] {0,8,1,9,20,11,3,13,2,
1870 10,21,12,22,26,23,15,24,14,
1871 4,16,5,17,25,18,7,19,6
1873 int hex64[] {0,8,9,1,10,32,35,14,
1874 11,33,34,15,3,19,18,2,
1875 12,36,37,16,40,56,57,44,
1876 43,59,58,45,22,49,48,20,
1877 13,39,38,17,41,60,61,47,
1878 42,63,62,46,23,50,51,21,
1879 4,24,25,5,26,52,53,28,
1880 27,55,54,29,7,31,30,6
1883 int wdg18[] = {0,6,1,7,9,2,8,15,10,
1884 16,17,11,3,12,4,13,14,5
1886 int wdg40[] = {0,6,7,1,8,24,12,9,13,2,
1887 10,26,27,14,30,38,34,33,35,16,
1888 11,29,28,15,31,39,37,32,36,17,
1889 3,18,19,4,20,25,22,21,23,5
1892 int pyr14[] = {0,5,1,6,13,8,3,
1895 int pyr30[] = {0,5,6,1,7,25,28,11,8,26,
1896 27,12,3,16,15,2,9,21,13,22,
1897 29,23,19,24,17,10,14,20,18,4
1900 vector<Element*> elements_0D, elements_1D, elements_2D, elements_3D;
1901 elements_0D.reserve(num_of_all_elements);
1902 elements_1D.reserve(num_of_all_elements);
1903 elements_2D.reserve(num_of_all_elements);
1904 elements_3D.reserve(num_of_all_elements);
1907 vector<Array<int>*> ho_verts_1D, ho_verts_2D, ho_verts_3D;
1908 ho_verts_1D.reserve(num_of_all_elements);
1909 ho_verts_2D.reserve(num_of_all_elements);
1910 ho_verts_3D.reserve(num_of_all_elements);
1913 vector<int> ho_el_order_1D, ho_el_order_2D, ho_el_order_3D;
1914 ho_el_order_1D.reserve(num_of_all_elements);
1915 ho_el_order_2D.reserve(num_of_all_elements);
1916 ho_el_order_3D.reserve(num_of_all_elements);
1928 ho_lin[2] = lin3; ho_lin[3] = lin4;
1929 ho_tri[2] = tri6; ho_tri[3] = tri10;
1930 ho_sqr[2] = quad9; ho_sqr[3] = quad16;
1931 ho_tet[2] = tet10; ho_tet[3] = tet20;
1932 ho_hex[2] = hex27; ho_hex[3] = hex64;
1933 ho_wdg[2] = wdg18; ho_wdg[3] = wdg40;
1934 ho_pyr[2] = pyr14; ho_pyr[3] = pyr30;
1936 bool has_nonpositive_phys_domain =
false;
1937 bool has_positive_phys_domain =
false;
1941 int n_elem_part = 0;
1942 const int header_size = 3;
1945 int header[header_size];
1946 int n_elem_one_type;
1948 while (n_elem_part < num_of_all_elements)
1950 input.read(
reinterpret_cast<char*
>(header),
1951 header_size*
sizeof(
int));
1952 type_of_element = header[0];
1953 n_elem_one_type = header[1];
1956 n_elem_part += n_elem_one_type;
1958 const int n_elem_nodes = nodes_of_gmsh_element[type_of_element-1];
1959 vector<int> data(1+n_tags+n_elem_nodes);
1960 for (
int el = 0; el < n_elem_one_type; ++el)
1962 input.read(
reinterpret_cast<char*
>(&data[0]),
1963 data.size()*
sizeof(
int));
1965 serial_number = data[dd++];
1968 phys_domain = (n_tags > 0) ? data[dd++] : 1;
1971 elem_domain = (n_tags > 1) ? data[dd++] : 0;
1974 n_partitions = (n_tags > 2) ? data[dd++] : 0;
1977 vector<int> vert_indices(n_elem_nodes);
1978 for (
int vi = 0; vi < n_elem_nodes; ++vi)
1980 map<int, int>::const_iterator it =
1981 vertices_map.find(data[1+n_tags+vi]);
1982 if (it == vertices_map.end())
1984 MFEM_ABORT(
"Gmsh file : vertex index doesn't exist");
1986 vert_indices[vi] = it->second;
1994 if (phys_domain <= 0)
1996 has_nonpositive_phys_domain =
true;
2001 has_positive_phys_domain =
true;
2006 switch (type_of_element)
2019 elements_1D.push_back(
2020 new Segment(&vert_indices[0], phys_domain));
2021 if (type_of_element != 1)
2023 el_order = n_elem_nodes - 1;
2025 hov->
Append(&vert_indices[0], n_elem_nodes);
2026 ho_verts_1D.push_back(hov);
2027 ho_el_order_1D.push_back(el_order);
2033 case 21: el_order--;
2034 case 23: el_order--;
2035 case 25: el_order--;
2036 case 42: el_order--;
2037 case 43: el_order--;
2038 case 44: el_order--;
2039 case 45: el_order--;
2040 case 46: el_order--;
2042 elements_2D.push_back(
2043 new Triangle(&vert_indices[0], phys_domain));
2047 hov->
Append(&vert_indices[0], n_elem_nodes);
2048 ho_verts_2D.push_back(hov);
2049 ho_el_order_2D.push_back(el_order);
2054 case 10: el_order--;
2055 case 36: el_order--;
2056 case 37: el_order--;
2057 case 38: el_order--;
2058 case 47: el_order--;
2059 case 48: el_order--;
2060 case 49: el_order--;
2061 case 50: el_order--;
2062 case 51: el_order--;
2064 elements_2D.push_back(
2069 hov->
Append(&vert_indices[0], n_elem_nodes);
2070 ho_verts_2D.push_back(hov);
2071 ho_el_order_2D.push_back(el_order);
2076 case 11: el_order--;
2077 case 29: el_order--;
2078 case 30: el_order--;
2079 case 31: el_order--;
2080 case 71: el_order--;
2081 case 72: el_order--;
2082 case 73: el_order--;
2083 case 74: el_order--;
2084 case 75: el_order--;
2086#ifdef MFEM_USE_MEMALLOC
2087 elements_3D.push_back(
TetMemory.Alloc());
2088 elements_3D.back()->SetVertices(&vert_indices[0]);
2089 elements_3D.back()->SetAttribute(phys_domain);
2091 elements_3D.push_back(
2097 hov->
Append(&vert_indices[0], n_elem_nodes);
2098 ho_verts_3D.push_back(hov);
2099 ho_el_order_3D.push_back(el_order);
2104 case 12: el_order--;
2105 case 92: el_order--;
2106 case 93: el_order--;
2107 case 94: el_order--;
2108 case 95: el_order--;
2109 case 96: el_order--;
2110 case 97: el_order--;
2111 case 98: el_order--;
2114 elements_3D.push_back(
2115 new Hexahedron(&vert_indices[0], phys_domain));
2119 hov->
Append(&vert_indices[0], n_elem_nodes);
2120 ho_verts_3D.push_back(hov);
2121 ho_el_order_3D.push_back(el_order);
2126 case 13: el_order--;
2127 case 90: el_order--;
2128 case 91: el_order--;
2129 case 106: el_order--;
2130 case 107: el_order--;
2131 case 108: el_order--;
2132 case 109: el_order--;
2133 case 110: el_order--;
2136 elements_3D.push_back(
2137 new Wedge(&vert_indices[0], phys_domain));
2141 hov->
Append(&vert_indices[0], n_elem_nodes);
2142 ho_verts_3D.push_back(hov);
2143 ho_el_order_3D.push_back(el_order);
2148 case 14: el_order--;
2149 case 118: el_order--;
2150 case 119: el_order--;
2151 case 120: el_order--;
2152 case 121: el_order--;
2153 case 122: el_order--;
2154 case 123: el_order--;
2155 case 124: el_order--;
2158 elements_3D.push_back(
2159 new Pyramid(&vert_indices[0], phys_domain));
2163 hov->
Append(&vert_indices[0], n_elem_nodes);
2164 ho_verts_3D.push_back(hov);
2165 ho_el_order_3D.push_back(el_order);
2171 elements_0D.push_back(
2172 new Point(&vert_indices[0], phys_domain));
2176 MFEM_WARNING(
"Unsupported Gmsh element type.");
2185 for (
int el = 0; el < num_of_all_elements; ++el)
2187 input >> serial_number >> type_of_element >> n_tags;
2188 vector<int> data(n_tags);
2189 for (
int i = 0; i < n_tags; ++i) { input >> data[i]; }
2192 phys_domain = (n_tags > 0) ? data[0] : 1;
2195 elem_domain = (n_tags > 1) ? data[1] : 0;
2198 n_partitions = (n_tags > 2) ? data[2] : 0;
2201 const int n_elem_nodes = nodes_of_gmsh_element[type_of_element-1];
2202 vector<int> vert_indices(n_elem_nodes);
2204 for (
int vi = 0; vi < n_elem_nodes; ++vi)
2207 map<int, int>::const_iterator it = vertices_map.find(
index);
2208 if (it == vertices_map.end())
2210 MFEM_ABORT(
"Gmsh file : vertex index doesn't exist");
2212 vert_indices[vi] = it->second;
2220 if (phys_domain <= 0)
2222 has_nonpositive_phys_domain =
true;
2227 has_positive_phys_domain =
true;
2232 switch (type_of_element)
2245 elements_1D.push_back(
2246 new Segment(&vert_indices[0], phys_domain));
2247 if (type_of_element != 1)
2250 hov->
Append(&vert_indices[0], n_elem_nodes);
2251 ho_verts_1D.push_back(hov);
2252 el_order = n_elem_nodes - 1;
2253 ho_el_order_1D.push_back(el_order);
2259 case 21: el_order--;
2260 case 23: el_order--;
2261 case 25: el_order--;
2262 case 42: el_order--;
2263 case 43: el_order--;
2264 case 44: el_order--;
2265 case 45: el_order--;
2266 case 46: el_order--;
2268 elements_2D.push_back(
2269 new Triangle(&vert_indices[0], phys_domain));
2273 hov->
Append(&vert_indices[0], n_elem_nodes);
2274 ho_verts_2D.push_back(hov);
2275 ho_el_order_2D.push_back(el_order);
2280 case 10: el_order--;
2281 case 36: el_order--;
2282 case 37: el_order--;
2283 case 38: el_order--;
2284 case 47: el_order--;
2285 case 48: el_order--;
2286 case 49: el_order--;
2287 case 50: el_order--;
2288 case 51: el_order--;
2290 elements_2D.push_back(
2295 hov->
Append(&vert_indices[0], n_elem_nodes);
2296 ho_verts_2D.push_back(hov);
2297 ho_el_order_2D.push_back(el_order);
2302 case 11: el_order--;
2303 case 29: el_order--;
2304 case 30: el_order--;
2305 case 31: el_order--;
2306 case 71: el_order--;
2307 case 72: el_order--;
2308 case 73: el_order--;
2309 case 74: el_order--;
2310 case 75: el_order--;
2312#ifdef MFEM_USE_MEMALLOC
2313 elements_3D.push_back(
TetMemory.Alloc());
2314 elements_3D.back()->SetVertices(&vert_indices[0]);
2315 elements_3D.back()->SetAttribute(phys_domain);
2317 elements_3D.push_back(
2323 hov->
Append(&vert_indices[0], n_elem_nodes);
2324 ho_verts_3D.push_back(hov);
2325 ho_el_order_3D.push_back(el_order);
2330 case 12: el_order--;
2331 case 92: el_order--;
2332 case 93: el_order--;
2333 case 94: el_order--;
2334 case 95: el_order--;
2335 case 96: el_order--;
2336 case 97: el_order--;
2337 case 98: el_order--;
2340 elements_3D.push_back(
2341 new Hexahedron(&vert_indices[0], phys_domain));
2345 hov->
Append(&vert_indices[0], n_elem_nodes);
2346 ho_verts_3D.push_back(hov);
2347 ho_el_order_3D.push_back(el_order);
2352 case 13: el_order--;
2353 case 90: el_order--;
2354 case 91: el_order--;
2355 case 106: el_order--;
2356 case 107: el_order--;
2357 case 108: el_order--;
2358 case 109: el_order--;
2359 case 110: el_order--;
2362 elements_3D.push_back(
2363 new Wedge(&vert_indices[0], phys_domain));
2367 hov->
Append(&vert_indices[0], n_elem_nodes);
2368 ho_verts_3D.push_back(hov);
2369 ho_el_order_3D.push_back(el_order);
2374 case 14: el_order--;
2375 case 118: el_order--;
2376 case 119: el_order--;
2377 case 120: el_order--;
2378 case 121: el_order--;
2379 case 122: el_order--;
2380 case 123: el_order--;
2381 case 124: el_order--;
2384 elements_3D.push_back(
2385 new Pyramid(&vert_indices[0], phys_domain));
2389 hov->
Append(&vert_indices[0], n_elem_nodes);
2390 ho_verts_3D.push_back(hov);
2391 ho_el_order_3D.push_back(el_order);
2397 elements_0D.push_back(
2398 new Point(&vert_indices[0], phys_domain));
2402 MFEM_WARNING(
"Unsupported Gmsh element type.");
2409 if (has_positive_phys_domain && has_nonpositive_phys_domain)
2411 MFEM_ABORT(
"Non-positive element attribute in Gmsh mesh!\n"
2412 "By default Gmsh sets element tags (attributes)"
2413 " to '0' but MFEM requires that they be"
2414 " positive integers.\n"
2415 "Use \"Physical Curve\", \"Physical Surface\","
2416 " or \"Physical Volume\" to set tags/attributes"
2417 " for all curves, surfaces, or volumes in your"
2418 " Gmsh geometry to values which are >= 1.");
2420 else if (has_nonpositive_phys_domain)
2422 mfem::out <<
"\nGmsh reader: all element attributes were zero.\n"
2423 <<
"MFEM only supports positive element attributes.\n"
2424 <<
"Setting element attributes to 1.\n\n";
2427 if (!elements_3D.empty())
2442 for (
size_t el = 0; el < ho_el_order_3D.size(); el++)
2444 mesh_order = max(mesh_order, ho_el_order_3D[el]);
2447 for (
size_t el = 0; el < elements_1D.size(); ++el)
2449 delete elements_1D[el];
2451 for (
size_t el = 0; el < elements_0D.size(); ++el)
2453 delete elements_0D[el];
2456 else if (!elements_2D.empty())
2471 for (
size_t el = 0; el < ho_el_order_2D.size(); el++)
2473 mesh_order = max(mesh_order, ho_el_order_2D[el]);
2476 for (
size_t el = 0; el < elements_0D.size(); ++el)
2478 delete elements_0D[el];
2481 else if (!elements_1D.empty())
2496 for (
size_t el = 0; el < ho_el_order_1D.size(); el++)
2498 mesh_order = max(mesh_order, ho_el_order_1D[el]);
2503 MFEM_ABORT(
"Gmsh file : no elements found");
2533 const int * vm = NULL;
2538 ho_verts = ho_verts_1D[el];
2539 el_order = ho_el_order_1D[el];
2540 if (!ho_lin[el_order])
2542 ho_lin[el_order] =
new int[ho_verts->
Size()];
2545 vm = ho_lin[el_order];
2548 ho_verts = ho_verts_2D[el];
2549 el_order = ho_el_order_2D[el];
2550 if (!ho_tri[el_order])
2552 ho_tri[el_order] =
new int[ho_verts->
Size()];
2555 vm = ho_tri[el_order];
2558 ho_verts = ho_verts_2D[el];
2559 el_order = ho_el_order_2D[el];
2560 if (!ho_sqr[el_order])
2562 ho_sqr[el_order] =
new int[ho_verts->
Size()];
2565 vm = ho_sqr[el_order];
2568 ho_verts = ho_verts_3D[el];
2569 el_order = ho_el_order_3D[el];
2570 if (!ho_tet[el_order])
2572 ho_tet[el_order] =
new int[ho_verts->
Size()];
2575 vm = ho_tet[el_order];
2578 ho_verts = ho_verts_3D[el];
2579 el_order = ho_el_order_3D[el];
2580 if (!ho_hex[el_order])
2582 ho_hex[el_order] =
new int[ho_verts->
Size()];
2585 vm = ho_hex[el_order];
2588 ho_verts = ho_verts_3D[el];
2589 el_order = ho_el_order_3D[el];
2590 if (!ho_wdg[el_order])
2592 ho_wdg[el_order] =
new int[ho_verts->
Size()];
2595 vm = ho_wdg[el_order];
2598 ho_verts = ho_verts_3D[el];
2599 el_order = ho_el_order_3D[el];
2600 if (!ho_pyr[el_order])
2602 ho_pyr[el_order] =
new int[ho_verts->
Size()];
2605 vm = ho_pyr[el_order];
2608 MFEM_WARNING(
"Unsupported Gmsh element type.");
2611 int nv = (ho_verts) ? ho_verts->
Size() : 0;
2613 for (
int v = 0; v<nv; v++)
2618 Nodes_gf(
spaceDim * (o + v) + d) = c[d];
2626 for (
size_t el=0; el<ho_verts_1D.size(); el++)
2628 delete ho_verts_1D[el];
2630 for (
size_t el=0; el<ho_verts_2D.size(); el++)
2632 delete ho_verts_2D[el];
2634 for (
size_t el=0; el<ho_verts_3D.size(); el++)
2636 delete ho_verts_3D[el];
2640 for (
int ord=4; ord<ho_lin.
Size(); ord++)
2642 if (ho_lin[ord] != NULL) {
delete [] ho_lin[ord]; }
2644 for (
int ord=4; ord<ho_tri.
Size(); ord++)
2646 if (ho_tri[ord] != NULL) {
delete [] ho_tri[ord]; }
2648 for (
int ord=4; ord<ho_sqr.
Size(); ord++)
2650 if (ho_sqr[ord] != NULL) {
delete [] ho_sqr[ord]; }
2652 for (
int ord=4; ord<ho_tet.
Size(); ord++)
2654 if (ho_tet[ord] != NULL) {
delete [] ho_tet[ord]; }
2656 for (
int ord=4; ord<ho_hex.
Size(); ord++)
2658 if (ho_hex[ord] != NULL) {
delete [] ho_hex[ord]; }
2660 for (
int ord=4; ord<ho_wdg.
Size(); ord++)
2662 if (ho_wdg[ord] != NULL) {
delete [] ho_wdg[ord]; }
2664 for (
int ord=4; ord<ho_pyr.
Size(); ord++)
2666 if (ho_pyr[ord] != NULL) {
delete [] ho_pyr[ord]; }
2672 MFEM_CONTRACT_VAR(n_partitions);
2673 MFEM_CONTRACT_VAR(elem_domain);
2676 else if (buff ==
"$PhysicalNames")
2682 for (
int i=0; i < num_names; i++)
2684 input >> mdim >> num;
2685 getline(input, name);
2688 while (!name.empty() &&
2689 (*name.begin() ==
' ' || *name.begin() ==
'\t'))
2693 while (!name.empty() &&
2694 (*name.rbegin() ==
' ' || *name.rbegin() ==
'\t' ||
2695 *name.rbegin() ==
'\n' || *name.rbegin() ==
'\r'))
2696 { name.resize(name.length()-1);}
2699 if ( (*name.begin() ==
'"' || *name.begin() ==
'\'') &&
2700 (*name.rbegin() ==
'"' || *name.rbegin() ==
'\''))
2702 name = name.substr(1,name.length()-2);
2705 phys_names_by_dim[mdim][num] = name;
2708 else if (buff ==
"$Periodic")
2715 for (
int i = 0; i < v2v.
Size(); i++)
2721 input >> num_per_ent;
2722 getline(input, buff);
2723 for (
int i = 0; i < num_per_ent; i++)
2725 getline(input, buff);
2726 getline(input, buff);
2727 if (!strncmp(buff.c_str(),
"Affine", 6))
2733 num_nodes = atoi(buff.c_str());
2735 for (
int j=0; j<num_nodes; j++)
2738 input >> slave >> master;
2739 v2v[slave - 1] = master - 1;
2741 getline(input, buff);
2748 for (
int slave = 0; slave < v2v.
Size(); slave++)
2750 int master = v2v[slave];
2751 if (master != slave)
2754 while (v2v[master] != master && master != slave)
2756 master = v2v[master];
2758 if (master == slave)
2767 v2v[slave] = master;
2773 if (mesh_order == 1)
2782 for (
int i = 0; i < this->
GetNE(); i++)
2787 for (
int j = 0; j < nv; j++)
2794 for (
int i = 0; i < this->
GetNBE(); i++)
2799 for (
int j = 0; j < nv; j++)
2808 if (phys_names_by_dim.size() > 0)
2811 for (
auto const &bdr_attr : phys_names_by_dim[
Dim-1])
2821 for (
auto const &attr : phys_names_by_dim[
Dim])
2845#ifdef MFEM_USE_NETCDF
2854 1,2,3,4,5,7,8,6,9,10
2860 1,2,3,4,5,6,7,8,9,10,11,
2863 12,17,18,19,20,13,14,15,
2866 16,22,26,25,27,24,23,21
2871 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14
2876 1, 2, 3, 4, 5, 6, 7, 8, 9, 13, 14, 15, 10, 11, 12, 16, 17, 18
2977 CubitElement() =
delete;
2980 ~CubitElement() =
default;
2989 inline size_t GetNumFaces()
const {
return _num_faces; }
2992 inline size_t GetNumVertices()
const {
return _num_vertices; }
2995 inline size_t GetNumNodes()
const {
return _num_nodes; }
2998 size_t GetNumFaceVertices(
size_t iface = 1)
const;
3001 inline uint8_t GetOrder()
const {
return _order; }
3004 Element * BuildElement(Mesh & mesh,
const int * vertex_ids,
3005 const int block_id)
const;
3008 Element * BuildBoundaryElement(Mesh & mesh,
const int iface,
3009 const int * vertex_ids,
const int sideset_id)
const;
3022 Element * NewElement(Mesh & mesh,
Geometry::Type geom,
const int *vertices,
3023 const int attribute)
const;
3030 size_t _num_vertices;
3037 _element_type = element_type;
3039 switch (element_type)
3114 MFEM_ABORT(
"Unsupported Cubit element type " << element_type <<
".");
3140 MFEM_ABORT(
"Unsupported 3D element with " << num_nodes <<
" nodes.");
3157 MFEM_ABORT(
"Unsupported 2D element with " << num_nodes <<
" nodes.");
3166 return Get2DElementType(num_nodes);
3170 return Get3DElementType(num_nodes);
3174 MFEM_ABORT(
"Unsupported Cubit dimension " <<
dimension <<
".");
3178CubitFaceType CubitElement::GetFaceType(
size_t side_id)
const
3181 bool valid_id = (side_id >= 1 &&
3182 side_id <= GetNumFaces());
3185 MFEM_ABORT(
"Encountered invalid side ID: " << side_id <<
".");
3188 switch (_element_type)
3215 MFEM_ABORT(
"Unknown element type: " << _element_type <<
".");
3220size_t CubitElement::GetNumFaceVertices(
size_t side_id)
const
3222 switch (GetFaceType(side_id))
3234 MFEM_ABORT(
"Unrecognized Cubit face type " << GetFaceType(side_id) <<
".");
3240 const int *vertices,
3241 const int attribute)
const
3243 Element *new_element = mesh.NewElement(geom);
3245 new_element->SetAttribute(attribute);
3251 const int *vertex_ids,
3252 const int block_id)
const
3254 switch (GetElementType())
3275 MFEM_ABORT(
"Unsupported Cubit element type encountered.");
3280mfem::Element * CubitElement::BuildBoundaryElement(Mesh &mesh,
3282 const int *vertex_ids,
3283 const int sideset_id)
const
3285 switch (GetFaceType(face_id))
3297 MFEM_ABORT(
"Unsupported Cubit face type encountered.");
3310 CubitBlock() =
delete;
3311 ~CubitBlock() =
default;
3321 const CubitElement & GetBlockElement(
int block_id)
const;
3331 uint8_t GetOrder()
const;
3332 inline uint8_t GetDimension()
const {
return _dimension; }
3334 inline size_t GetNumBlocks()
const {
return BlockIDs().size(); }
3335 inline bool HasBlocks()
const {
return !BlockIDs().empty(); }
3342 void CheckElementBlockIsCompatible(
const CubitElement & new_block_element)
3348 void ClearBlockElements();
3353 inline const std::set<int> & BlockIDs()
const {
return _block_ids; }
3355 bool HasBlockID(
int block_id)
const;
3356 bool ValidBlockID(
int block_id)
const;
3357 bool ValidDimension(
int dimension)
const;
3363 std::set<int> _block_ids;
3368 std::map<int, CubitElement> _block_element_for_block_id;
3381 MFEM_ABORT(
"Invalid dimension '" <<
dimension <<
"' specified.");
3386 ClearBlockElements();
3392 if (HasBlockID(block_id))
3394 MFEM_ABORT(
"Block with ID '" << block_id <<
"' has already been added.");
3396 else if (!ValidBlockID(block_id))
3398 MFEM_ABORT(
"Illegal block ID '" << block_id <<
"'.");
3401 CubitElement block_element = CubitElement(element_type);
3406 CheckElementBlockIsCompatible(block_element);
3410 _order = block_element.GetOrder();
3413 _block_ids.insert(block_id);
3414 _block_element_for_block_id.emplace(block_id,
3419CubitBlock::GetOrder()
const
3423 MFEM_ABORT(
"No elements have been added.");
3430CubitBlock::ClearBlockElements()
3434 _block_element_for_block_id.clear();
3438CubitBlock::HasBlockID(
int block_id)
const
3440 return (_block_ids.count(block_id) > 0);
3444CubitBlock::ValidBlockID(
int block_id)
const
3446 return (block_id > 0);
3450CubitBlock::ValidDimension(
int dimension)
const
3456CubitBlock::GetBlockElement(
int block_id)
const
3458 if (!HasBlockID(block_id))
3460 MFEM_ABORT(
"No element info for block ID '" << block_id <<
"'.");
3463 return _block_element_for_block_id.at(block_id);
3467CubitBlock::CheckElementBlockIsCompatible(
const CubitElement &
3468 new_block_element)
const
3476 if (GetOrder() != new_block_element.GetOrder())
3478 MFEM_ABORT(
"All block elements must be of the same order.");
3488 NetCDFReader() =
delete;
3489 NetCDFReader(
const std::string fname);
3494 bool HasVariable(
const char * name);
3497 void ReadVariable(
const char * name,
int * data);
3500 void ReadVariable(
const char * name,
double * data);
3503 bool HasDimension(
const char * name);
3506 void ReadDimension(
const char * name,
size_t *
dimension);
3509 void BuildIDToNameMap(
const vector<int> & ids,
3510 unordered_map<int, string> & ids_to_names,
3511 const string & quantity_name);
3515 void CheckForNetCDFError();
3518 int ReadVariableID(
const char * name);
3521 int ReadDimensionID(
const char * name);
3525 void HandleNetCDFError(
const int error_code);
3527 int _netcdf_status{NC_NOERR};
3528 int _netcdf_descriptor;
3531 char *_name_buffer{NULL};
3535NetCDFReader::NetCDFReader(
const std::string fname)
3537 _netcdf_status = nc_open(fname.c_str(), NC_NOWRITE, &_netcdf_descriptor);
3538 CheckForNetCDFError();
3541 _name_buffer =
new char[NC_MAX_NAME + 1];
3544NetCDFReader::~NetCDFReader()
3546 _netcdf_status = nc_close(_netcdf_descriptor);
3547 CheckForNetCDFError();
3551 delete[] _name_buffer;
3555void NetCDFReader::CheckForNetCDFError()
3557 if (_netcdf_status != NC_NOERR)
3559 HandleNetCDFError(_netcdf_status);
3563void NetCDFReader::HandleNetCDFError(
const int error_code)
3565 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(error_code));
3568int NetCDFReader::ReadVariableID(
const char * var_name)
3572 _netcdf_status = nc_inq_varid(_netcdf_descriptor, var_name,
3574 CheckForNetCDFError();
3579int NetCDFReader::ReadDimensionID(
const char * name)
3583 _netcdf_status = nc_inq_dimid(_netcdf_descriptor, name, &dim_id);
3584 CheckForNetCDFError();
3589void NetCDFReader::ReadDimension(
const char * name,
size_t *
dimension)
3591 const int dimension_id = ReadDimensionID(name);
3594 _netcdf_status = nc_inq_dim(_netcdf_descriptor, dimension_id, _name_buffer,
3596 CheckForNetCDFError();
3599bool NetCDFReader::HasVariable(
const char * name)
3602 const int status = nc_inq_varid(_netcdf_descriptor, name, &var_id);
3611 HandleNetCDFError(status);
3616bool NetCDFReader::HasDimension(
const char * name)
3619 const int status = nc_inq_dimid(_netcdf_descriptor, name, &dim_id);
3628 HandleNetCDFError(status);
3633void NetCDFReader::ReadVariable(
const char * name,
int * data)
3635 const int variable_id = ReadVariableID(name);
3637 _netcdf_status = nc_get_var_int(_netcdf_descriptor, variable_id, data);
3638 CheckForNetCDFError();
3642void NetCDFReader::ReadVariable(
const char * name,
double * data)
3644 const int variable_id = ReadVariableID(name);
3646 _netcdf_status = nc_get_var_double(_netcdf_descriptor, variable_id, data);
3647 CheckForNetCDFError();
3651void NetCDFReader::BuildIDToNameMap(
const vector<int> & ids,
3652 unordered_map<int, string> & ids_to_names,
3653 const string & quantity_name)
3658 _netcdf_status = nc_inq_varid(_netcdf_descriptor, quantity_name.c_str(),
3662 if (_netcdf_status == NC_ENOTVAR)
3668 CheckForNetCDFError();
3673 _netcdf_status = nc_inq_vartype(_netcdf_descriptor, varid_names,
3675 CheckForNetCDFError();
3677 if (var_type == NC_CHAR)
3679 int dimids_names[2], names_ndim;
3680 size_t num_names, name_len;
3682 _netcdf_status = nc_inq_varndims(_netcdf_descriptor, varid_names,
3684 CheckForNetCDFError();
3685 MFEM_ASSERT(names_ndim == 2,
"This variable should have two dimensions");
3687 _netcdf_status = nc_inq_vardimid(_netcdf_descriptor, varid_names,
3689 CheckForNetCDFError();
3691 _netcdf_status = nc_inq_dimlen(_netcdf_descriptor, dimids_names[0], &num_names);
3692 CheckForNetCDFError();
3693 MFEM_ASSERT(num_names == ids.size(),
3694 "The block id and block name lengths don't match");
3696 _netcdf_status = nc_inq_dimlen(_netcdf_descriptor, dimids_names[1], &name_len);
3697 CheckForNetCDFError();
3700 vector<char> names(ids.size() * name_len);
3701 _netcdf_status = nc_get_var_text(_netcdf_descriptor, varid_names,
3703 CheckForNetCDFError();
3705 for (
size_t i = 0; i < ids.size(); ++i)
3707 string name(&names[i * name_len], name_len);
3709 name.resize(name.find(
'\0'));
3710 ids_to_names[ids[i]] = name;
3715 mfem_error(
"Unexpected netcdf variable type");
3721static void ReadCubitNodeCoordinates(NetCDFReader & cubit_reader,
3726 cubit_reader.ReadVariable(
"coordx", coordx);
3727 cubit_reader.ReadVariable(
"coordy", coordy);
3731 cubit_reader.ReadVariable(
"coordz", coordz);
3737static void ReadCubitNumElementsInBlock(NetCDFReader & cubit_reader,
3738 const vector<int> & block_ids,
3739 map<int, size_t> &num_elements_for_block_id)
3741 num_elements_for_block_id.clear();
3744 const int buffer_size = NC_MAX_NAME + 1;
3745 char string_buffer[buffer_size];
3748 for (
const auto block_id : block_ids)
3751 snprintf(string_buffer, buffer_size,
"num_el_in_blk%d", iblock++);
3753 size_t num_elements_for_block = 0;
3754 cubit_reader.ReadDimension(string_buffer, &num_elements_for_block);
3756 num_elements_for_block_id[block_id] = num_elements_for_block;
3762static void BuildElementIDsForBlockID(
3763 const vector<int> & block_ids,
3764 const map<int, size_t> & num_elements_for_block_id,
3765 map<
int, vector<int>> & element_ids_for_block_id,
3766 map<int, int> & block_id_for_element_id)
3768 element_ids_for_block_id.clear();
3769 block_id_for_element_id.clear();
3774 for (
int block_id : block_ids)
3776 const int num_elements_for_block = num_elements_for_block_id.at(block_id);
3778 vector<int> element_ids(num_elements_for_block);
3780 for (
int i = 0; i < num_elements_for_block; i++, element_id++)
3782 element_ids[i] = element_id;
3783 block_id_for_element_id[element_id] = block_id;
3786 element_ids_for_block_id[block_id] = std::move(element_ids);
3791static void ReadCubitBlocks(NetCDFReader & cubit_reader,
3792 const vector<int> block_ids,
3793 CubitBlock & cubit_blocks)
3795 const int buffer_size = NC_MAX_NAME + 1;
3796 char string_buffer[buffer_size];
3798 size_t num_nodes_per_element;
3801 for (
int block_id : block_ids)
3804 snprintf(string_buffer, buffer_size,
"num_nod_per_el%d", iblock++);
3806 cubit_reader.ReadDimension(string_buffer, &num_nodes_per_element);
3810 num_nodes_per_element, cubit_blocks.GetDimension());
3811 cubit_blocks.AddBlockElement(block_id, element_type);
3817static void ReadCubitDimensions(NetCDFReader & cubit_reader,
3822 size_t &num_side_sets)
3824 cubit_reader.ReadDimension(
"num_dim", &num_dim);
3825 cubit_reader.ReadDimension(
"num_nodes", &num_nodes);
3826 cubit_reader.ReadDimension(
"num_elem", &num_elem);
3827 cubit_reader.ReadDimension(
"num_el_blk", &num_el_blk);
3830 if (cubit_reader.HasDimension(
"num_side_sets"))
3832 cubit_reader.ReadDimension(
"num_side_sets", &num_side_sets);
3842static void ReadCubitBoundaries(NetCDFReader & cubit_reader,
3843 const vector<int> & boundary_ids,
3844 map<
int, vector<int>> & element_ids_for_boundary_id,
3845 map<
int, vector<int>> & side_ids_for_boundary_id)
3847 const int buffer_size = NC_MAX_NAME + 1;
3848 char string_buffer[buffer_size];
3851 for (
int boundary_id : boundary_ids)
3854 size_t num_sides = 0;
3856 snprintf(string_buffer, buffer_size,
"num_side_ss%d", ibdr);
3857 cubit_reader.ReadDimension(string_buffer, &num_sides);
3860 vector<int> boundary_element_ids(num_sides);
3861 vector<int> boundary_side_ids(num_sides);
3864 snprintf(string_buffer, buffer_size,
"elem_ss%d", ibdr);
3865 cubit_reader.ReadVariable(string_buffer, boundary_element_ids.data());
3868 snprintf(string_buffer, buffer_size,
"side_ss%d", ibdr++);
3869 cubit_reader.ReadVariable(string_buffer, boundary_side_ids.data());
3872 element_ids_for_boundary_id[boundary_id] = std::move(boundary_element_ids);
3873 side_ids_for_boundary_id[boundary_id] = std::move(boundary_side_ids);
3878static void BuildCubitBlockIDs(NetCDFReader & cubit_reader,
3879 const int num_element_blocks,
3880 vector<int> & block_ids)
3882 block_ids.resize(num_element_blocks);
3883 cubit_reader.ReadVariable(
"eb_prop1", block_ids.data());
3887static void ReadCubitBoundaryIDs(NetCDFReader & cubit_reader,
3888 const int num_boundaries, vector<int> & boundary_ids)
3890 boundary_ids.clear();
3892 if (num_boundaries < 1) {
return; }
3894 boundary_ids.resize(num_boundaries);
3895 cubit_reader.ReadVariable(
"ss_prop1", boundary_ids.data());
3899static void ReadCubitElementBlocks(NetCDFReader & cubit_reader,
3900 const CubitBlock & cubit_blocks,
3901 const vector<int> & block_ids,
3902 const map<
int, vector<int>> & element_ids_for_block_id,
3903 map<
int, vector<int>> &node_ids_for_element_id)
3905 const int buffer_size = NC_MAX_NAME + 1;
3906 char string_buffer[buffer_size];
3909 for (
const int block_id : block_ids)
3911 const CubitElement & block_element = cubit_blocks.GetBlockElement(block_id);
3913 const vector<int> & block_element_ids = element_ids_for_block_id.at(block_id);
3915 const size_t num_nodes_for_block = block_element_ids.size() *
3916 block_element.GetNumNodes();
3918 vector<int> node_ids_for_block(num_nodes_for_block);
3921 snprintf(string_buffer, buffer_size,
"connect%d", iblock++);
3923 cubit_reader.ReadVariable(string_buffer, node_ids_for_block.data());
3927 for (
int element_id : block_element_ids)
3929 vector<int> element_node_ids(block_element.GetNumNodes());
3931 for (
int i = 0; i < (int)block_element.GetNumNodes(); i++)
3933 element_node_ids[i] = node_ids_for_block[ielement * block_element.GetNumNodes()
3939 node_ids_for_element_id[element_id] = std::move(element_node_ids);
3945static void BuildBoundaryNodeIDs(
const vector<int> & boundary_ids,
3946 const CubitBlock & blocks,
3947 const map<
int, vector<int>> & node_ids_for_element_id,
3948 const map<
int, vector<int>> & element_ids_for_boundary_id,
3949 const map<
int, vector<int>> & side_ids_for_boundary_id,
3950 const map<int, int> & block_id_for_element_id,
3951 map<
int, vector<vector<int>>> & node_ids_for_boundary_id)
3953 for (
int boundary_id : boundary_ids)
3956 auto & boundary_element_ids = element_ids_for_boundary_id.at(
3958 auto & boundary_element_sides = side_ids_for_boundary_id.at(
3962 vector<vector<int>> boundary_node_ids(
3963 boundary_element_ids.size());
3966 for (
int jelement = 0; jelement < (int)boundary_element_ids.size(); jelement++)
3969 const int boundary_element_global_id = boundary_element_ids[jelement];
3970 const int boundary_side = boundary_element_sides[jelement];
3973 const int block_id = block_id_for_element_id.at(boundary_element_global_id);
3974 const CubitElement & block_element = blocks.GetBlockElement(block_id);
3976 const int num_face_vertices = block_element.GetNumFaceVertices(boundary_side);
3977 vector<int> nodes_of_element_on_side(num_face_vertices);
3980 const vector<int> & element_node_ids =
3981 node_ids_for_element_id.at(boundary_element_global_id);
3985 for (
int knode = 0; knode < num_face_vertices; knode++)
3989 switch (block_element.GetElementType())
4016 MFEM_ABORT(
"Unsupported element type encountered.\n");
4020 nodes_of_element_on_side[knode] = element_node_ids[inode - 1];
4023 boundary_node_ids[jelement] = std::move(nodes_of_element_on_side);
4027 node_ids_for_boundary_id[boundary_id] = std::move(boundary_node_ids);
4032static void BuildUniqueVertexIDs(
const vector<int> & unique_block_ids,
4033 const CubitBlock & blocks,
4034 const map<
int, vector<int>> & element_ids_for_block_id,
4035 const map<
int, vector<int>> & node_ids_for_element_id,
4036 vector<int> & unique_vertex_ids)
4039 for (
int block_id : unique_block_ids)
4041 auto & element_ids = element_ids_for_block_id.at(block_id);
4043 auto & block_element = blocks.GetBlockElement(block_id);
4045 for (
int element_id : element_ids)
4047 auto & node_ids = node_ids_for_element_id.at(element_id);
4049 for (
size_t knode = 0; knode < block_element.GetNumVertices(); knode++)
4051 unique_vertex_ids.push_back(node_ids[knode]);
4057 std::sort(unique_vertex_ids.begin(), unique_vertex_ids.end());
4059 auto new_end = std::unique(unique_vertex_ids.begin(), unique_vertex_ids.end());
4061 unique_vertex_ids.resize(std::distance(unique_vertex_ids.begin(), new_end));
4067static void BuildCubitToMFEMVertexMap(
const vector<int> & unique_vertex_ids,
4068 map<int, int> & cubit_to_mfem_vertex_map)
4070 cubit_to_mfem_vertex_map.clear();
4073 for (
int vertex_id : unique_vertex_ids)
4075 cubit_to_mfem_vertex_map[vertex_id] = ivertex++;
4083static void FinalizeCubitSecondOrderMesh(Mesh &mesh,
4084 const vector<int> & unique_block_ids,
4085 const CubitBlock & blocks,
4086 const map<
int, vector<int>> & element_ids_for_block_id,
4087 const map<
int, vector<int>> & node_ids_for_element_id,
4088 const double *coordx,
4089 const double *coordy,
4090 const double *coordz)
4092 mesh.FinalizeTopology();
4095 const int Dim = mesh.Dimension();
4096 FiniteElementCollection *fec =
new H1_FECollection(2, Dim);
4097 FiniteElementSpace *fes =
new FiniteElementSpace(&mesh, fec, Dim,
4099 GridFunction *Nodes =
new GridFunction(fes);
4100 Nodes->MakeOwner(fec);
4101 mesh.SetNodalGridFunction(Nodes,
true);
4103 for (
int block_id : unique_block_ids)
4105 const CubitElement & block_element = blocks.GetBlockElement(block_id);
4107 int *mfem_to_genesis_map = NULL;
4109 switch (block_element.GetElementType())
4130 MFEM_ABORT(
"Something went wrong. Linear elements detected when order is 2.");
4133 auto & element_ids = element_ids_for_block_id.at(block_id);
4135 for (
int element_id : element_ids)
4139 fes->GetElementDofs(element_id - 1, dofs);
4141 Array<int> vdofs = dofs;
4142 fes->DofsToVDofs(vdofs);
4144 const vector<int> & element_node_ids = node_ids_for_element_id.at(element_id);
4146 for (
int jnode = 0; jnode < dofs.Size(); jnode++)
4148 const int node_index = element_node_ids[mfem_to_genesis_map[jnode] - 1] - 1;
4150 (*Nodes)(vdofs[jnode]) = coordx[node_index];
4151 (*Nodes)(vdofs[jnode] + 1) = coordy[node_index];
4155 (*Nodes)(vdofs[jnode] + 2) = coordz[node_index];
4166 const vector<double> & coordx,
4167 const vector<double> & coordy,
4168 const vector<double> & coordz)
4175 const int original_1based_id = unique_vertex_ids[ivertex];
4177 vertices[ivertex](0) = coordx[original_1based_id - 1];
4178 vertices[ivertex](1) = coordy[original_1based_id - 1];
4182 vertices[ivertex](2) = coordz[original_1based_id - 1];
4189 const cubit::CubitBlock * blocks,
4190 const vector<int> & block_ids,
4191 const map<
int, vector<int>> & element_ids_for_block_id,
4192 const map<
int, vector<int>> & node_ids_for_element_id,
4193 const map<int, int> & cubit_to_mfem_vertex_map)
4195 using namespace cubit;
4200 int element_counter = 0;
4203 for (
int block_id : block_ids)
4205 const CubitElement & block_element = blocks->GetBlockElement(block_id);
4207 vector<int> renumbered_vertex_ids(block_element.GetNumVertices());
4209 const vector<int> &block_element_ids = element_ids_for_block_id.at(block_id);
4212 for (
int element_id : block_element_ids)
4214 const vector<int> & element_node_ids = node_ids_for_element_id.at(element_id);
4217 for (
size_t knode = 0; knode < block_element.GetNumVertices(); knode++)
4219 const int node_id = element_node_ids[knode];
4222 renumbered_vertex_ids[knode] = cubit_to_mfem_vertex_map.at(node_id) - 1;
4226 elements[element_counter++] = block_element.BuildElement(*
this,
4227 renumbered_vertex_ids.data(),
4235 const cubit::CubitBlock * blocks,
4236 const vector<int> & boundary_ids,
4237 const map<
int, vector<int>> & element_ids_for_boundary_id,
4238 const map<
int, vector<vector<int>>> & node_ids_for_boundary_id,
4239 const map<
int, vector<int>> & side_ids_for_boundary_id,
4240 const map<int, int> & block_id_for_element_id,
4241 const map<int, int> & cubit_to_mfem_vertex_map)
4243 using namespace cubit;
4246 for (
int boundary_id : boundary_ids)
4253 array<int, 8> renumbered_vertex_ids;
4256 int boundary_counter = 0;
4257 for (
int boundary_id : boundary_ids)
4259 const vector<int> &elements_on_boundary = element_ids_for_boundary_id.at(
4262 const vector<vector<int>> &nodes_on_boundary = node_ids_for_boundary_id.at(
4266 for (
int side_id : side_ids_for_boundary_id.at(boundary_id))
4269 const int element_id = elements_on_boundary.at(jelement);
4270 const int element_block = block_id_for_element_id.at(element_id);
4271 const CubitElement & block_element = blocks->GetBlockElement(element_block);
4273 const vector<int> & element_nodes_on_side = nodes_on_boundary.at(jelement);
4276 for (
size_t knode = 0; knode < element_nodes_on_side.size(); knode++)
4278 const int node_id = element_nodes_on_side[knode];
4281 renumbered_vertex_ids[knode] = cubit_to_mfem_vertex_map.at(node_id) - 1;
4285 boundary[boundary_counter++] = block_element.BuildBoundaryElement(*
this,
4287 renumbered_vertex_ids.data(),
4297 using namespace cubit;
4305 NetCDFReader cubit_reader(filename);
4310 size_t num_dimensions, num_nodes, num_elements, num_element_blocks,
4313 ReadCubitDimensions(cubit_reader, num_dimensions, num_nodes, num_elements,
4314 num_element_blocks, num_boundaries);
4316 Dim = num_dimensions;
4321 vector<int> block_ids;
4322 BuildCubitBlockIDs(cubit_reader, num_element_blocks, block_ids);
4323 unordered_map<int, string> blk_ids_to_names;
4324 cubit_reader.BuildIDToNameMap(block_ids, blk_ids_to_names,
"eb_names");
4325 for (
const auto & pr : blk_ids_to_names)
4327 const auto blk_id = pr.first;
4328 const auto & blk_name = pr.second;
4329 if (!blk_name.empty())
4339 map<int, size_t> num_elements_for_block_id;
4340 ReadCubitNumElementsInBlock(cubit_reader, block_ids,
4341 num_elements_for_block_id);
4343 map<int, vector<int>> element_ids_for_block_id;
4344 map<int, int> block_id_for_element_id;
4345 BuildElementIDsForBlockID(
4346 block_ids, num_elements_for_block_id, element_ids_for_block_id,
4347 block_id_for_element_id);
4351 CubitBlock blocks(num_dimensions);
4352 ReadCubitBlocks(cubit_reader, block_ids, blocks);
4355 map<int, vector<int>> node_ids_for_element_id;
4356 ReadCubitElementBlocks(cubit_reader,
4359 element_ids_for_block_id,
4360 node_ids_for_element_id);
4365 vector<int> boundary_ids;
4366 ReadCubitBoundaryIDs(cubit_reader, num_boundaries, boundary_ids);
4367 unordered_map<int, string> bnd_ids_to_names;
4368 cubit_reader.BuildIDToNameMap(boundary_ids, bnd_ids_to_names,
"ss_names");
4369 for (
const auto & pr : bnd_ids_to_names)
4371 const auto bnd_id = pr.first;
4372 const auto & bnd_name = pr.second;
4373 if (!bnd_name.empty())
4387 map<int, vector<int>> element_ids_for_boundary_id;
4388 map<int, vector<int>> side_ids_for_boundary_id;
4390 ReadCubitBoundaries(cubit_reader, boundary_ids,
4391 element_ids_for_boundary_id, side_ids_for_boundary_id);
4393 map<int, vector<vector<int>>> node_ids_for_boundary_id;
4395 BuildBoundaryNodeIDs(boundary_ids, blocks, node_ids_for_element_id,
4396 element_ids_for_boundary_id, side_ids_for_boundary_id,
4397 block_id_for_element_id,
4398 node_ids_for_boundary_id);
4403 vector<double> coordx(num_nodes);
4404 vector<double> coordy(num_nodes);
4405 vector<double> coordz(num_dimensions == 3 ? num_nodes : 0);
4407 ReadCubitNodeCoordinates(cubit_reader, coordx.data(), coordy.data(),
4413 vector<int> unique_vertex_ids;
4414 BuildUniqueVertexIDs(block_ids, blocks, element_ids_for_block_id,
4415 node_ids_for_element_id, unique_vertex_ids);
4423 map<int, int> cubit_to_mfem_vertex_map;
4424 BuildCubitToMFEMVertexMap(unique_vertex_ids, cubit_to_mfem_vertex_map);
4435 element_ids_for_block_id,
4436 node_ids_for_element_id, cubit_to_mfem_vertex_map);
4442 element_ids_for_boundary_id, node_ids_for_boundary_id, side_ids_for_boundary_id,
4443 block_id_for_element_id,
4444 cubit_to_mfem_vertex_map);
4449 if (blocks.GetOrder() == 2)
4453 FinalizeCubitSecondOrderMesh(*
this,
4456 element_ids_for_block_id,
4457 node_ids_for_element_id,
void Load(std::istream &in, int fmt=0)
Read an Array from the stream in using format fmt. The format fmt can be:
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
int Size() const
Return the logical size of the array.
int Append(const T &el)
Append element 'el' to array, resize if necessary.
T * GetData()
Returns the data.
void SortAll()
Sort each named array in the container.
void UniqueAll()
Remove duplicates from each, previously sorted, named array.
void Load(std::istream &in)
Load the contents of the container from an input stream.
bool AttributeSetExists(const std::string &name) const
Return true is the named attribute set is present.
void AddToAttributeSet(const std::string &set_name, int attr)
Add a single entry to an existing attribute set.
ArraysByName< int > attr_sets
Named sets of attributes.
Array< int > & CreateAttributeSet(const std::string &set_name)
Create an empty named attribute set.
@ ClosedUniform
Nodes: x_i = i/(n-1), i=0,...,n-1.
Abstract data type element.
virtual void GetVertices(Array< int > &v) const =0
Get the indices defining the vertices.
void SetAttribute(const int attr)
Set element's attribute.
Type
Constants for the classes derived from Element.
virtual int GetNVertices() const =0
virtual void SetVertices(const Array< int > &v)=0
Set the indices defining the vertices.
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
void DofsToVDofs(Array< int > &dofs, int ndofs=-1) const
Compute the full set of vdofs corresponding to each entry in dofs.
DofTransformation * GetElementDofs(int elem, Array< int > &dofs) const
Returns indices of degrees of freedom of element 'elem'. The returned indices are offsets into an ldo...
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
Abstract class for all finite elements.
static const int Dimension[NumGeom]
Class for grid function - Vector with associated FE space.
void MakeOwner(FiniteElementCollection *fec_)
Make the GridFunction the owner of fec_owned and fes.
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
int VectorDim() const
Shortcut for calling FiniteElementSpace::GetVectorDim() on the underlying fes.
void GetNodalValues(int i, Array< real_t > &nval, int vdim=1) const
Returns the values at the vertices of element i for the 1-based dimension vdim.
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the GridFunction.
Arbitrary order H1-conforming (continuous) finite elements.
Data type hexahedron element.
Arbitrary order "L2-conforming" discontinuous finite elements.
int CheckElementOrientation(bool fix_it=true)
Check (and optionally attempt to fix) the orientation of the elements.
void ReadGmshMesh(std::istream &input, int &curved, int &read_gf)
Element * NewElement(int geom)
friend class NCNURBSExtension
void BuildCubitElements(const int num_elements, const cubit::CubitBlock *blocks, const std::vector< int > &block_ids, const std::map< int, std::vector< int > > &element_ids_for_block_id, const std::map< int, std::vector< int > > &node_ids_for_element_id, const std::map< int, int > &cubit_to_mfem_vertex_map)
Called internally in ReadCubit. This method builds the mesh elements.
MemAlloc< Tetrahedron, 1024 > TetMemory
NURBSExtension * NURBSext
Optional NURBS mesh extension.
void ReadTrueGridMesh(std::istream &input)
static const int vtk_quadratic_tet[10]
Element::Type GetElementType(int i) const
Returns the type of element i.
void ReadNetgen3DMesh(std::istream &input)
void ReadInlineMesh(std::istream &input, bool generate_edges=false)
void ReadLineMesh(std::istream &input)
static const int vtk_quadratic_wedge[18]
AttributeSets bdr_attribute_sets
Named sets of boundary element attributes.
void Make1D(int n, real_t sx=1.0)
const Element * GetElement(int i) const
Return pointer to the i'th element object.
void FinalizeTopology(bool generate_bdr=true)
Finalize the construction of the secondary topology (connectivity) data of a Mesh.
void ReadXML_VTKMesh(std::istream &input, int &curved, int &read_gf, bool &finalize_topo, const std::string &xml_prefix="")
int GetNE() const
Returns number of elements.
void Make3D(int nx, int ny, int nz, Element::Type type, real_t sx, real_t sy, real_t sz, bool sfc_ordering)
Creates a mesh for the parallelepiped [0,sx]x[0,sy]x[0,sz], divided into nx*ny*nz hexahedra if type =...
const Element * GetBdrElement(int i) const
Return pointer to the i'th boundary element object.
friend class NURBSExtension
static bool remove_unused_vertices
void ReadCubit(const std::string &filename, int &curved, int &read_gf)
Load a mesh from a Genesis file.
static const int vtk_quadratic_pyramid[13]
AttributeSets attribute_sets
Named sets of element attributes.
void ReadVTKMesh(std::istream &input, int &curved, int &read_gf, bool &finalize_topo)
static const int vtk_quadratic_hex[27]
int GetNBE() const
Returns number of boundary elements.
void ReadMFEMMesh(std::istream &input, int version, int &curved)
void CreateVTKMesh(const Vector &points, const Array< int > &cell_data, const Array< int > &cell_offsets, const Array< int > &cell_types, const Array< int > &cell_attributes, int &curved, int &read_gf, bool &finalize_topo)
void BuildCubitBoundaries(const cubit::CubitBlock *blocks, const std::vector< int > &boundary_ids, const std::map< int, std::vector< int > > &element_ids_for_boundary_id, const std::map< int, std::vector< std::vector< int > > > &node_ids_for_boundary_id, const std::map< int, std::vector< int > > &side_ids_for_boundary_id, const std::map< int, int > &block_id_for_element_id, const std::map< int, int > &cubit_to_mfem_vertex_map)
Called internally in ReadCubit. This method adds the mesh boundary elements.
Array< Element * > boundary
void ReadNURBSMesh(std::istream &input, int &curved, int &read_gf, bool spacing=false, bool nc=false)
void SetMeshGen()
Determine the mesh generator bitmask meshgen, see MeshGenerator().
void Make2D(int nx, int ny, Element::Type type, real_t sx, real_t sy, bool generate_edges, bool sfc_ordering)
Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny quadrilaterals if type = QUADRILATER...
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Set the curvature of the mesh nodes using the given polynomial degree.
Element * ReadElement(std::istream &input)
Geometry::Type GetElementBaseGeometry(int i) const
void BuildCubitVertices(const std::vector< int > &unique_vertex_ids, const std::vector< double > &coordx, const std::vector< double > &coordy, const std::vector< double > &coordz)
Called internally in ReadCubit. This method creates the vertices.
void ReadNetgen2DMesh(std::istream &input, int &curved)
Array< Element * > elements
const real_t * GetVertex(int i) const
Return pointer to vertex i's coordinates.
void RemoveUnusedVertices()
Remove unused vertices and rebuild mesh connectivity.
int GetNBE() const
Return the number of active boundary elements.
void SetCoordsFromPatches(Vector &Nodes)
Set FE coordinates in Nodes, using data from patches, and erase patches.
void GetElementTopo(Array< Element * > &elements) const
Generate the active mesh elements and return them in elements.
bool HavePatches() const
Return true if at least 1 patch is defined, false otherwise.
void GetBdrElementTopo(Array< Element * > &boundary) const
Generate the active mesh boundary elements and return them in boundary.
int GetOrder() const
If all KnotVector orders are identical, return that number. Otherwise, return NURBSFECollection::Vari...
int GetNV() const
Return the local number of active vertices.
int Dimension() const
Return the dimension of the reference space (not physical space).
int GetNE() const
Return the number of active elements.
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Class for standard nodal finite elements.
const Array< int > & GetLexicographicOrdering() const
Get an Array<int> that maps lexicographically ordered indices to the indices of the respective nodes/...
Data type Pyramid element.
Piecewise-(bi)quadratic continuous finite elements.
Data type quadrilateral element.
Data type line segment element.
Data type tetrahedron element.
void SetVertices(const Array< int > &v) override
Set the indices defining the vertices.
Data type triangle element.
Vector coefficient defined by a vector GridFunction.
void Load(std::istream **in, int np, int *dim)
Reads a vector from multiple files.
int Size() const
Returns the size of the vector.
real_t * GetData() const
Return a pointer to the beginning of the Vector data.
constexpr int dimension
This example only works in 3D. Kernels for 2D are not implemented.
int index(int i, int j, int nx, int ny)
void DecodeBase64(const char *src, size_t len, std::vector< char > &buf)
Decode len base-64 encoded characters in the buffer src, and store the resulting decoded data in buf....
size_t NumBase64Chars(size_t nbytes)
Return the number of characters needed to encode nbytes in base-64.
T read(std::istream &is)
Read a value from the stream and return it.
const int mfem_to_genesis_tri6[6]
const int mfem_to_genesis_pyramid14[14]
const int cubit_side_map_tri3[3][2]
const int mfem_to_genesis_wedge18[18]
const int cubit_side_map_hex8[6][4]
const int cubit_side_map_quad4[4][2]
const int mfem_to_genesis_quad9[9]
const int cubit_side_map_tet4[4][3]
const int cubit_side_map_wedge6[5][4]
const int mfem_to_genesis_hex27[27]
const int cubit_side_map_pyramid5[5][4]
const int mfem_to_genesis_tet10[10]
MFEM_HOST_DEVICE constexpr auto type(const tuple< T... > &values)
a function intended to be used for extracting the ith type from a tuple.
bool StringCompare(const char *s1, const char *s2)
void mfem_error(const char *msg)
void GmshHOWedgeMapping(int order, int *map)
Generate Gmsh vertex mapping for a Wedge.
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
void GmshHOTriangleMapping(int order, int *map)
Generate Gmsh vertex mapping for a Triangle.
void GmshHOPyramidMapping(int order, int *map)
Generate Gmsh vertex mapping for a Pyramid.
void filter_dos(std::string &line)
Check for, and remove, a trailing '\r' from and std::string.
const char * VTKByteOrder()
Determine the byte order and return either "BigEndian" or "LittleEndian".
void GmshHOTetrahedronMapping(int order, int *map)
Generate Gmsh vertex mapping for a Tetrahedron.
void CreateVTKElementConnectivity(Array< int > &con, Geometry::Type geom, int ref)
Create the VTK element connectivity array for a given element geometry and refinement level.
void GmshHOSegmentMapping(int order, int *map)
Generate Gmsh vertex mapping for a Segment.
void GmshHOQuadrilateralMapping(int order, int *map)
Generate Gmsh vertex mapping for a Quadrilateral.
void GmshHOHexahedronMapping(int order, int *map)
Generate Gmsh vertex mapping for a Hexahedron.
void skip_comment_lines(std::istream &is, const char comment_char)
Check if the stream starts with comment_char. If so skip it.
static const int LAGRANGE_PRISM
static const int PrismMap[6]
Permutation from MFEM's prism ordering to VTK's prism ordering.
static bool IsLagrange(int vtk_geom)
Does the given VTK geometry type describe an arbitrary-order Lagrange element?
static Geometry::Type GetMFEMGeometry(int vtk_geom)
Given a VTK geometry type, return the corresponding MFEM Geometry::Type.
static int GetOrder(int vtk_geom, int npoints)
For the given VTK geometry type and number of points, return the order of the element.