13 #include "../fem/fem.hpp"
14 #include "../general/binaryio.hpp"
15 #include "../general/text.hpp"
16 #include "../general/tinyxml2.h"
24 #ifdef MFEM_USE_NETCDF
37 bool Mesh::remove_unused_vertices =
true;
39 void Mesh::ReadMFEMMesh(std::istream &input,
int version,
int &curved)
42 MFEM_VERIFY(version == 10 || version == 12,
43 "unknown MFEM mesh version");
51 MFEM_VERIFY(ident ==
"dimension",
"invalid mesh file");
57 MFEM_VERIFY(ident ==
"elements",
"invalid mesh file");
58 input >> NumOfElements;
59 elements.SetSize(NumOfElements);
60 for (
int j = 0; j < NumOfElements; j++)
62 elements[j] = ReadElement(input);
68 MFEM_VERIFY(ident ==
"boundary",
"invalid mesh file");
69 input >> NumOfBdrElements;
70 boundary.SetSize(NumOfBdrElements);
71 for (
int j = 0; j < NumOfBdrElements; j++)
73 boundary[j] = ReadElement(input);
79 MFEM_VERIFY(ident ==
"vertices",
"invalid mesh file");
80 input >> NumOfVertices;
81 vertices.SetSize(NumOfVertices);
87 spaceDim = atoi(ident.c_str());
88 for (
int j = 0; j < NumOfVertices; j++)
90 for (
int i = 0; i < spaceDim; i++)
92 input >> vertices[j](i);
105 if (remove_unused_vertices) { RemoveUnusedVertices(); }
108 void Mesh::ReadLineMesh(std::istream &input)
114 input >> NumOfVertices;
115 vertices.SetSize(NumOfVertices);
117 for (j = 0; j < NumOfVertices; j++)
119 input >> vertices[j](0);
122 input >> NumOfElements;
123 elements.SetSize(NumOfElements);
125 for (j = 0; j < NumOfElements; j++)
127 input >> a >> p1 >> p2;
128 elements[j] =
new Segment(p1-1, p2-1, a);
132 input >> NumOfBdrElements;
133 boundary.SetSize(NumOfBdrElements);
134 for (j = 0; j < NumOfBdrElements; j++)
136 input >> a >> ind[0];
138 boundary[j] =
new Point(ind,a);
142 void Mesh::ReadNetgen2DMesh(std::istream &input,
int &curved)
144 int ints[32], attr, n;
150 input >> NumOfBdrElements;
151 boundary.SetSize(NumOfBdrElements);
152 for (
int i = 0; i < NumOfBdrElements; i++)
155 >> ints[0] >> ints[1];
156 ints[0]--; ints[1]--;
157 boundary[i] =
new Segment(ints, attr);
161 input >> NumOfElements;
162 elements.SetSize(NumOfElements);
163 for (
int i = 0; i < NumOfElements; i++)
166 for (
int j = 0; j < n; j++)
174 elements[i] =
new Segment(ints, attr);
177 elements[i] =
new Triangle(ints, attr);
188 input >> NumOfVertices;
189 vertices.SetSize(NumOfVertices);
190 for (
int i = 0; i < NumOfVertices; i++)
191 for (
int j = 0; j < Dim; j++)
193 input >> vertices[i](j);
198 input >> NumOfVertices;
199 vertices.SetSize(NumOfVertices);
204 void Mesh::ReadNetgen3DMesh(std::istream &input)
212 input >> NumOfVertices;
214 vertices.SetSize(NumOfVertices);
215 for (
int i = 0; i < NumOfVertices; i++)
216 for (
int j = 0; j < Dim; j++)
218 input >> vertices[i](j);
222 input >> NumOfElements;
223 elements.SetSize(NumOfElements);
224 for (
int i = 0; i < NumOfElements; i++)
227 for (
int j = 0; j < 4; j++)
232 #ifdef MFEM_USE_MEMALLOC
234 tet = TetMemory.Alloc();
244 input >> NumOfBdrElements;
245 boundary.SetSize(NumOfBdrElements);
246 for (
int i = 0; i < NumOfBdrElements; i++)
249 for (
int j = 0; j < 3; j++)
254 boundary[i] =
new Triangle(ints, attr);
258 void Mesh::ReadTrueGridMesh(std::istream &input)
260 int i, j, ints[32], attr;
261 const int buflen = 1024;
272 input >> vari >> NumOfVertices >> vari >> vari >> NumOfElements;
273 input.getline(buf, buflen);
274 input.getline(buf, buflen);
276 input.getline(buf, buflen);
277 input.getline(buf, buflen);
278 input.getline(buf, buflen);
281 vertices.SetSize(NumOfVertices);
282 for (i = 0; i < NumOfVertices; i++)
284 input >> vari >> varf >> vertices[i](0) >> vertices[i](1);
285 input.getline(buf, buflen);
289 elements.SetSize(NumOfElements);
290 for (i = 0; i < NumOfElements; i++)
292 input >> vari >> attr;
293 for (j = 0; j < 4; j++)
298 input.getline(buf, buflen);
299 input.getline(buf, buflen);
307 input >> vari >> NumOfVertices >> NumOfElements;
308 input.getline(buf, buflen);
309 input.getline(buf, buflen);
310 input >> vari >> vari >> NumOfBdrElements;
311 input.getline(buf, buflen);
312 input.getline(buf, buflen);
313 input.getline(buf, buflen);
315 vertices.SetSize(NumOfVertices);
316 for (i = 0; i < NumOfVertices; i++)
318 input >> vari >> varf >> vertices[i](0) >> vertices[i](1)
320 input.getline(buf, buflen);
323 elements.SetSize(NumOfElements);
324 for (i = 0; i < NumOfElements; i++)
326 input >> vari >> attr;
327 for (j = 0; j < 8; j++)
332 input.getline(buf, buflen);
336 boundary.SetSize(NumOfBdrElements);
337 for (i = 0; i < NumOfBdrElements; i++)
340 for (j = 0; j < 4; j++)
345 input.getline(buf, buflen);
352 const int Mesh::vtk_quadratic_tet[10] =
353 { 0, 1, 2, 3, 4, 7, 5, 6, 8, 9 };
357 const int Mesh::vtk_quadratic_pyramid[13] =
358 { 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12};
362 const int Mesh::vtk_quadratic_wedge[18] =
363 { 0, 2, 1, 3, 5, 4, 8, 7, 6, 11, 10, 9, 12, 14, 13, 17, 16, 15};
366 const int Mesh::vtk_quadratic_hex[27] =
368 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19,
369 24, 22, 21, 23, 20, 25, 26
376 int &curved,
int &read_gf,
bool &finalize_topo)
378 int np = points.
Size()/3;
380 NumOfElements = cell_types.
Size();
381 elements.SetSize(NumOfElements);
384 bool legacy_elem =
false, lagrange_elem =
false;
386 for (
int i = 0; i < NumOfElements; i++)
388 int j = (i > 0) ? cell_offsets[i-1] : 0;
389 int ct = cell_types[i];
391 elements[i] = NewElement(geom);
392 if (cell_attributes.
Size() > 0)
394 elements[i]->SetAttribute(cell_attributes[i]);
398 if (geom == Geometry::PRISM && ct != VTKGeometry::LAGRANGE_PRISM)
400 int prism_vertices[6];
401 for (
int k=0; k<6; ++k)
403 prism_vertices[k] = cell_data[j+VTKGeometry::PrismMap[k]];
405 elements[i]->SetVertices(prism_vertices);
409 elements[i]->SetVertices(&cell_data[j]);
412 int elem_dim = Geometry::Dimension[geom];
413 int elem_order = VTKGeometry::GetOrder(ct, cell_offsets[i] - j);
415 if (VTKGeometry::IsLagrange(ct)) { lagrange_elem =
true; }
416 else { legacy_elem =
true; }
418 MFEM_VERIFY(Dim == -1 || Dim == elem_dim,
419 "Elements with different dimensions are not supported");
420 MFEM_VERIFY(order == -1 || order == elem_order,
421 "Elements with different orders are not supported");
422 MFEM_VERIFY(legacy_elem != lagrange_elem,
423 "Mixing of legacy and Lagrange cell types is not supported");
432 double min_value, max_value;
433 for (
int d = 3; d > 0; --d)
435 min_value = max_value = points(3*0 + d-1);
436 for (
int i = 1; i < np; i++)
438 min_value = std::min(min_value, points(3*i + d-1));
439 max_value = std::max(max_value, points(3*i + d-1));
440 if (min_value != max_value)
446 if (spaceDim > 0) {
break; }
450 if (order == 1 && !lagrange_elem)
453 vertices.SetSize(np);
454 for (
int i = 0; i < np; i++)
456 vertices[i](0) = points(3*i+0);
457 vertices[i](1) = points(3*i+1);
458 vertices[i](2) = points(3*i+2);
461 NumOfBdrElements = 0;
463 CheckElementOrientation(
true);
475 for (
int i = 0; i < NumOfElements; i++)
477 int *v = elements[i]->GetVertices();
478 int nv = elements[i]->GetNVertices();
479 for (
int j = 0; j < nv; j++)
481 if (pts_dof[v[j]] == -1) { pts_dof[v[j]] = 0; }
490 for (
int i = 0; i < np; i++)
492 if (pts_dof[i] != -1)
494 pts_dof[i] = NumOfVertices++;
498 for (
int i = 0; i < NumOfElements; i++)
500 int *v = elements[i]->GetVertices();
501 int nv = elements[i]->GetNVertices();
502 for (
int j = 0; j < nv; j++)
504 v[j] = pts_dof[v[j]];
508 vertices.
SetSize(NumOfVertices);
509 for (
int i = 0; i < np; i++)
514 vertices[j](0) = points(3*i+0);
515 vertices[j](1) = points(3*i+1);
516 vertices[j](2) = points(3*i+2);
521 NumOfBdrElements = 0;
534 Nodes->MakeOwner(fec);
539 for (
int i = 0; i < NumOfElements; i++)
543 switch (elements[i]->GetGeometryType())
545 case Geometry::TRIANGLE:
546 case Geometry::SQUARE:
547 vtk_mfem = vtk_quadratic_hex;
break;
548 case Geometry::TETRAHEDRON:
549 vtk_mfem = vtk_quadratic_tet;
break;
551 vtk_mfem = vtk_quadratic_hex;
break;
552 case Geometry::PRISM:
553 vtk_mfem = vtk_quadratic_wedge;
break;
554 case Geometry::PYRAMID:
555 vtk_mfem = vtk_quadratic_pyramid;
break;
561 int offset = (i == 0) ? 0 : cell_offsets[i-1];
562 for (
int j = 0; j < dofs.
Size(); j++)
564 if (pts_dof[cell_data[offset+j]] == -1)
566 pts_dof[cell_data[offset+j]] = dofs[vtk_mfem[j]];
570 if (pts_dof[cell_data[offset+j]] != dofs[vtk_mfem[j]])
572 MFEM_ABORT(
"VTK mesh: inconsistent quadratic mesh!");
584 Nodes->MakeOwner(fec);
588 std::map<Geometry::Type,Array<int>> vtk_inv_maps;
589 std::map<Geometry::Type,const Array<int>*> lex_orderings;
592 for (n = i = 0; i < NumOfElements; i++)
598 if (vtk_inv_map.
Size() == 0)
603 for (
int j=0; j<vtk_map.
Size(); ++j)
605 vtk_inv_map[vtk_map[j]] = j;
608 const Array<int> *&lex_ordering = lex_orderings[geom];
614 MFEM_ASSERT(nodal_fe != NULL,
"Unsupported element type");
618 for (
int lex_idx = 0; lex_idx < dofs.
Size(); lex_idx++)
620 int mfem_idx = (*lex_ordering)[lex_idx];
621 int vtk_idx = vtk_inv_map[lex_idx];
622 int pt_idx = cell_data[n + vtk_idx];
623 if (pts_dof[pt_idx] == -1)
625 pts_dof[pt_idx] = dofs[mfem_idx];
629 if (pts_dof[pt_idx] != dofs[mfem_idx])
631 MFEM_ABORT(
"VTK mesh: inconsistent Lagrange mesh!");
640 for (
int i = 0; i < np; i++)
643 if (pts_dof[i] != -1)
645 dofs[0] = pts_dof[i];
647 for (
int d = 0; d < dofs.
Size(); d++)
649 (*Nodes)(dofs[d]) = points(3*i+d);
660 using namespace tinyxml2;
666 if (s1 == NULL || s2 == NULL) {
return false; }
667 return strcmp(s1, s2) == 0;
675 struct BufferReaderBase
677 enum HeaderType { UINT32_HEADER, UINT64_HEADER };
678 virtual void ReadBinary(
const char *buf,
void *dest,
int n)
const = 0;
679 virtual void ReadBase64(
const char *txt,
void *dest,
int n)
const = 0;
680 virtual ~BufferReaderBase() { }
693 template <
typename T,
typename F>
694 struct BufferReader : BufferReaderBase
697 HeaderType header_type;
698 BufferReader(
bool compressed_, HeaderType header_type_)
699 : compressed(compressed_), header_type(header_type_) { }
702 size_t HeaderEntrySize()
const
704 return header_type == UINT64_HEADER ?
sizeof(uint64_t) :
sizeof(uint32_t);
710 uint64_t ReadHeaderEntry(
const char *header_buf)
const
712 return (header_type == UINT64_HEADER) ? bin_io::read<uint64_t>(header_buf)
713 : bin_io::read<uint32_t>(header_buf);
722 int NumHeaderBytes(
const char *header_buf)
const
724 if (!compressed) {
return HeaderEntrySize(); }
725 return (3 + ReadHeaderEntry(header_buf))*HeaderEntrySize();
733 void ReadBinaryWithHeader(
const char *header_buf,
const char *buf,
734 void *dest_void,
int n)
const
736 std::vector<char> uncompressed_data;
737 T *dest =
static_cast<T*
>(dest_void);
747 int header_entry_size = HeaderEntrySize();
748 int nblocks = ReadHeaderEntry(header_buf);
749 header_buf += header_entry_size;
750 std::vector<int> header(nblocks + 2);
751 for (
int i=0; i<nblocks+2; ++i)
753 header[i] = ReadHeaderEntry(header_buf);
754 header_buf += header_entry_size;
756 uncompressed_data.resize((nblocks-1)*header[0] + header[1]);
757 Bytef *dest_ptr = (Bytef *)uncompressed_data.data();
758 Bytef *dest_start = dest_ptr;
759 const Bytef *source_ptr = (
const Bytef *)buf;
760 for (
int i=0; i<nblocks; ++i)
762 uLongf source_len = header[i+2];
763 uLong dest_len = (i == nblocks-1) ? header[1] : header[0];
764 int res = uncompress(dest_ptr, &dest_len, source_ptr, source_len);
765 MFEM_VERIFY(res == Z_OK,
"Error uncompressing");
766 dest_ptr += dest_len;
767 source_ptr += source_len;
769 MFEM_VERIFY(
int(
sizeof(F)*n) == (dest_ptr - dest_start),
770 "AppendedData: wrong data size");
771 buf = uncompressed_data.data();
773 MFEM_ABORT(
"MFEM must be compiled with zlib enabled to uncompress.")
781 if (header_type == UINT32_HEADER)
783 uint32_t *data_size_32 = (uint32_t *)header_buf;
784 data_size = *data_size_32;
788 uint64_t *data_size_64 = (uint64_t *)header_buf;
789 data_size = *data_size_64;
791 MFEM_VERIFY(
sizeof(F)*n == data_size,
"AppendedData: wrong data size");
794 if (std::is_same<T, F>::value)
797 memcpy(dest, buf,
sizeof(T)*n);
801 for (
int i=0; i<n; ++i)
804 dest[i] = bin_io::read<F>(buf + i*
sizeof(F));
812 void ReadBinary(
const char *buf,
void *dest,
int n)
const override
814 ReadBinaryWithHeader(buf, buf + NumHeaderBytes(buf), dest, n);
821 void ReadBase64(
const char *txt,
void *dest,
int n)
const override
826 if (*txt !=
' ' && *txt !=
'\n') {
break; }
833 std::vector<char> nblocks_buf;
836 std::vector<char> data, header;
843 ReadBinaryWithHeader(header.data(), data.data(), dest, n);
847 std::vector<char> data;
849 ReadBinary(data.data(), dest, n);
860 const char *appended_data, *byte_order, *compressor;
861 enum AppendedDataEncoding { RAW, BASE64 };
862 map<string,BufferReaderBase*> type_map;
863 AppendedDataEncoding encoding;
869 XMLDataReader(
const XMLElement *vtk,
const XMLElement *vtu)
872 BufferReaderBase::HeaderType htype;
875 htype = BufferReaderBase::UINT64_HEADER;
879 htype = BufferReaderBase::UINT32_HEADER;
883 byte_order = vtk->Attribute(
"byte_order");
887 compressor = vtk->Attribute(
"compressor");
888 bool compressed = (compressor != NULL);
891 appended_data = NULL;
892 for (
const XMLElement *xml_elem = vtu->NextSiblingElement();
894 xml_elem = xml_elem->NextSiblingElement())
898 const char *encoding_str = xml_elem->Attribute(
"encoding");
901 appended_data = xml_elem->GetAppendedData();
906 appended_data = xml_elem->GetText();
909 MFEM_VERIFY(appended_data != NULL,
"Invalid AppendedData");
911 bool found_leading_underscore =
false;
912 while (*appended_data)
915 if (*appended_data ==
'_')
917 found_leading_underscore =
true;
922 MFEM_VERIFY(found_leading_underscore,
"Invalid AppendedData");
927 type_map[
"Int8"] =
new BufferReader<int, int8_t>(compressed, htype);
928 type_map[
"Int16"] =
new BufferReader<int, int16_t>(compressed, htype);
929 type_map[
"Int32"] =
new BufferReader<int, int32_t>(compressed, htype);
930 type_map[
"Int64"] =
new BufferReader<int, int64_t>(compressed, htype);
931 type_map[
"UInt8"] =
new BufferReader<int, uint8_t>(compressed, htype);
932 type_map[
"UInt16"] =
new BufferReader<int, uint16_t>(compressed, htype);
933 type_map[
"UInt32"] =
new BufferReader<int, uint32_t>(compressed, htype);
934 type_map[
"UInt64"] =
new BufferReader<int, uint64_t>(compressed, htype);
935 type_map[
"Float32"] =
new BufferReader<double, float>(compressed, htype);
936 type_map[
"Float64"] =
new BufferReader<double, double>(compressed, htype);
942 template <
typename T>
943 void Read(
const XMLElement *xml_elem, T *dest,
int n)
945 static const char *erstr =
"Error reading XML DataArray";
946 MFEM_VERIFY(
StringCompare(xml_elem->Name(),
"DataArray"), erstr);
947 const char *format = xml_elem->Attribute(
"format");
950 const char *txt = xml_elem->GetText();
951 MFEM_VERIFY(txt != NULL, erstr);
952 std::istringstream data_stream(txt);
953 for (
int i=0; i<n; ++i) { data_stream >> dest[i]; }
957 VerifyBinaryOptions();
958 int offset = xml_elem->IntAttribute(
"offset");
959 const char *type = xml_elem->Attribute(
"type");
960 MFEM_VERIFY(type != NULL, erstr);
961 BufferReaderBase *reader = type_map[type];
962 MFEM_VERIFY(reader != NULL, erstr);
963 MFEM_VERIFY(appended_data != NULL,
"No AppendedData found");
966 reader->ReadBinary(appended_data + offset, dest, n);
970 reader->ReadBase64(appended_data + offset, dest, n);
975 VerifyBinaryOptions();
976 const char *txt = xml_elem->GetText();
977 MFEM_VERIFY(txt != NULL, erstr);
978 const char *type = xml_elem->Attribute(
"type");
979 if (type == NULL) { MFEM_ABORT(erstr); }
980 BufferReaderBase *reader = type_map[type];
981 if (reader == NULL) { MFEM_ABORT(erstr); }
982 reader->ReadBase64(txt, dest, n);
986 MFEM_ABORT(
"Invalid XML VTK DataArray format");
994 void VerifyByteOrder()
const
999 MFEM_ABORT(
"Converting between different byte orders is unsupported.");
1006 void VerifyCompressor()
const
1008 if (compressor && !
StringCompare(compressor,
"vtkZLibDataCompressor"))
1010 MFEM_ABORT(
"Unsupported compressor. Only zlib is supported.")
1012 #ifndef MFEM_USE_ZLIB
1013 MFEM_VERIFY(compressor == NULL,
"MFEM must be compiled with zlib enabled "
1014 "to support reading compressed data.");
1020 void VerifyBinaryOptions()
const
1028 for (
auto &x : type_map) {
delete x.second; }
1034 void Mesh::ReadXML_VTKMesh(std::istream &input,
int &curved,
int &read_gf,
1035 bool &finalize_topo,
const std::string &xml_prefix)
1037 using namespace vtk_xml;
1039 static const char *erstr =
"XML parsing error";
1042 std::vector<char> buf(xml_prefix.begin(), xml_prefix.end());
1043 std::istreambuf_iterator<char> eos;
1044 buf.insert(buf.end(), std::istreambuf_iterator<char>(input), eos);
1045 buf.push_back(
'\0');
1048 xml.Parse(buf.data(), buf.size());
1049 if (xml.ErrorID() != XML_SUCCESS)
1051 MFEM_ABORT(
"Error parsing XML VTK file.\n" << xml.ErrorStr());
1054 const XMLElement *vtkfile = xml.FirstChildElement();
1055 MFEM_VERIFY(vtkfile, erstr);
1056 MFEM_VERIFY(
StringCompare(vtkfile->Name(),
"VTKFile"), erstr);
1057 const XMLElement *vtu = vtkfile->FirstChildElement();
1058 MFEM_VERIFY(vtu, erstr);
1059 MFEM_VERIFY(
StringCompare(vtu->Name(),
"UnstructuredGrid"), erstr);
1061 XMLDataReader data_reader(vtkfile, vtu);
1064 const XMLElement *piece = vtu->FirstChildElement();
1066 MFEM_VERIFY(piece->NextSiblingElement() == NULL,
1067 "XML VTK meshes with more than one Piece are not supported");
1068 int npts = piece->IntAttribute(
"NumberOfPoints");
1069 int ncells = piece->IntAttribute(
"NumberOfCells");
1073 const XMLElement *pts_xml;
1074 for (pts_xml = piece->FirstChildElement();
1076 pts_xml = pts_xml->NextSiblingElement())
1080 const XMLElement *pts_data = pts_xml->FirstChildElement();
1081 MFEM_VERIFY(pts_data->IntAttribute(
"NumberOfComponents") == 3,
1082 "XML VTK Points DataArray must have 3 components");
1083 data_reader.Read(pts_data, points.
GetData(), points.
Size());
1087 if (pts_xml == NULL) { MFEM_ABORT(erstr); }
1090 Array<int> cell_data, cell_offsets(ncells), cell_types(ncells);
1091 const XMLElement *cells_xml;
1092 for (cells_xml = piece->FirstChildElement();
1094 cells_xml = cells_xml->NextSiblingElement())
1098 const XMLElement *cell_data_xml = NULL;
1099 for (
const XMLElement *data_xml = cells_xml->FirstChildElement();
1101 data_xml = data_xml->NextSiblingElement())
1103 const char *data_name = data_xml->Attribute(
"Name");
1106 data_reader.Read(data_xml, cell_offsets.GetData(), ncells);
1110 data_reader.Read(data_xml, cell_types.
GetData(), ncells);
1118 cell_data_xml = data_xml;
1121 MFEM_VERIFY(cell_data_xml != NULL, erstr);
1122 int cell_data_size = cell_offsets.Last();
1123 cell_data.SetSize(cell_data_size);
1124 data_reader.Read(cell_data_xml, cell_data.GetData(), cell_data_size);
1128 if (cells_xml == NULL) { MFEM_ABORT(erstr); }
1132 for (
const XMLElement *cell_data_xml = piece->FirstChildElement();
1133 cell_data_xml != NULL;
1134 cell_data_xml = cell_data_xml->NextSiblingElement())
1137 &&
StringCompare(cell_data_xml->Attribute(
"Scalars"),
"material"))
1139 const XMLElement *data_xml = cell_data_xml->FirstChildElement();
1140 if (data_xml != NULL &&
StringCompare(data_xml->Name(),
"DataArray"))
1142 cell_attributes.
SetSize(ncells);
1143 data_reader.Read(data_xml, cell_attributes.GetData(), ncells);
1148 CreateVTKMesh(points, cell_data, cell_offsets, cell_types, cell_attributes,
1149 curved, read_gf, finalize_topo);
1152 void Mesh::ReadVTKMesh(std::istream &input,
int &curved,
int &read_gf,
1153 bool &finalize_topo)
1162 getline(input, buff);
1163 getline(input, buff);
1165 if (buff !=
"ASCII")
1167 MFEM_ABORT(
"VTK mesh is not in ASCII format!");
1172 getline(input, buff);
1174 if (!input.good()) { MFEM_ABORT(
"VTK mesh is not UNSTRUCTURED_GRID!"); }
1176 while (buff !=
"DATASET UNSTRUCTURED_GRID");
1185 MFEM_ABORT(
"VTK mesh does not have POINTS data!");
1188 while (buff !=
"POINTS");
1193 getline(input, buff);
1194 points.
Load(input, 3*np);
1209 MFEM_ABORT(
"VTK mesh does not have CELLS data!");
1212 while (buff !=
"CELLS");
1216 if (buff ==
"CELLS")
1219 input >> ncells >> n >> ws;
1221 cell_data.
SetSize(n - ncells);
1223 for (
int i=0; i<ncells; ++i)
1227 cell_offsets[i] = offset + nv;
1228 for (
int j=0; j<nv; ++j)
1230 input >> cell_data[offset + j];
1237 input >> ws >> buff;
1240 MFEM_VERIFY(buff ==
"CELL_TYPES",
"CELL_TYPES not provided in VTK mesh.")
1242 cell_types.
Load(ncells, input);
1244 while ((input.good()) && (buff !=
"CELL_DATA"))
1248 getline(input, buff);
1253 while ((input.good()))
1255 getline(input, buff);
1256 if (buff.rfind(
"POINT_DATA") == 0)
1260 else if (buff.rfind(
"SCALARS material") == 0)
1262 getline(input, buff);
1263 if (buff.rfind(
"LOOKUP_TABLE default") != 0)
1265 MFEM_ABORT(
"Invalid LOOKUP_TABLE for material array in VTK file.");
1267 cell_attributes.
Load(ncells, input);
1279 CreateVTKMesh(points, cell_data, cell_offsets, cell_types, cell_attributes,
1280 curved, read_gf, finalize_topo);
1283 void Mesh::ReadNURBSMesh(std::istream &input,
int &curved,
int &read_gf)
1287 Dim = NURBSext->Dimension();
1288 NumOfVertices = NURBSext->GetNV();
1289 NumOfElements = NURBSext->GetNE();
1290 NumOfBdrElements = NURBSext->GetNBE();
1292 NURBSext->GetElementTopo(elements);
1293 NURBSext->GetBdrElementTopo(boundary);
1295 vertices.SetSize(NumOfVertices);
1297 if (NURBSext->HavePatches())
1303 Nodes->MakeOwner(fec);
1304 NURBSext->SetCoordsFromPatches(*Nodes);
1307 spaceDim = Nodes->VectorDim();
1308 for (
int i = 0; i < spaceDim; i++)
1311 Nodes->GetNodalValues(vert_val, i+1);
1312 for (
int j = 0; j < NumOfVertices; j++)
1314 vertices[j](i) = vert_val(j);
1324 void Mesh::ReadInlineMesh(std::istream &input,
bool generate_edges)
1352 MFEM_VERIFY(input.get() ==
'=',
1353 "Inline mesh expected '=' after keyword " << name);
1360 else if (name ==
"ny")
1364 else if (name ==
"nz")
1368 else if (name ==
"sx")
1372 else if (name ==
"sy")
1376 else if (name ==
"sz")
1380 else if (name ==
"type")
1384 if (eltype ==
"segment")
1386 type = Element::SEGMENT;
1388 else if (eltype ==
"quad")
1390 type = Element::QUADRILATERAL;
1392 else if (eltype ==
"tri")
1394 type = Element::TRIANGLE;
1396 else if (eltype ==
"hex")
1398 type = Element::HEXAHEDRON;
1400 else if (eltype ==
"wedge")
1402 type = Element::WEDGE;
1404 else if (eltype ==
"pyramid")
1406 type = Element::PYRAMID;
1408 else if (eltype ==
"tet")
1410 type = Element::TETRAHEDRON;
1414 MFEM_ABORT(
"unrecognized element type (read '" << eltype
1415 <<
"') in inline mesh format. "
1416 "Allowed: segment, tri, quad, tet, hex, wedge");
1421 MFEM_ABORT(
"unrecognized keyword (" << name
1422 <<
") in inline mesh format. "
1423 "Allowed: nx, ny, nz, type, sx, sy, sz");
1428 if (input.peek() ==
';')
1441 if (type == Element::SEGMENT)
1443 MFEM_VERIFY(nx > 0 && sx > 0.0,
1444 "invalid 1D inline mesh format, all values must be "
1446 <<
" nx = " << nx <<
"\n"
1447 <<
" sx = " << sx <<
"\n");
1450 else if (type == Element::TRIANGLE || type == Element::QUADRILATERAL)
1452 MFEM_VERIFY(nx > 0 && ny > 0 && sx > 0.0 && sy > 0.0,
1453 "invalid 2D inline mesh format, all values must be "
1455 <<
" nx = " << nx <<
"\n"
1456 <<
" ny = " << ny <<
"\n"
1457 <<
" sx = " << sx <<
"\n"
1458 <<
" sy = " << sy <<
"\n");
1459 Make2D(nx, ny, type, sx, sy, generate_edges,
true);
1461 else if (type == Element::TETRAHEDRON || type == Element::WEDGE ||
1462 type == Element::HEXAHEDRON || type == Element::PYRAMID)
1464 MFEM_VERIFY(nx > 0 && ny > 0 && nz > 0 &&
1465 sx > 0.0 && sy > 0.0 && sz > 0.0,
1466 "invalid 3D inline mesh format, all values must be "
1468 <<
" nx = " << nx <<
"\n"
1469 <<
" ny = " << ny <<
"\n"
1470 <<
" nz = " << nz <<
"\n"
1471 <<
" sx = " << sx <<
"\n"
1472 <<
" sy = " << sy <<
"\n"
1473 <<
" sz = " << sz <<
"\n");
1474 Make3D(nx, ny, nz, type, sx, sy, sz,
true);
1479 MFEM_ABORT(
"For inline mesh, must specify an element type ="
1480 " [segment, tri, quad, tet, hex, wedge]");
1484 void Mesh::ReadGmshMesh(std::istream &input,
int &curved,
int &read_gf)
1489 input >> version >> binary >> dsize;
1492 MFEM_ABORT(
"Gmsh file version < 2.2");
1494 if (dsize !=
sizeof(
double))
1496 MFEM_ABORT(
"Gmsh file : dsize != sizeof(double)");
1498 getline(input, buff);
1503 input.read(reinterpret_cast<char*>(&one),
sizeof(one));
1506 MFEM_ABORT(
"Gmsh file : wrong binary format");
1513 map<int, int> vertices_map;
1525 double bb_tol = 1e-14;
1533 bool periodic =
false;
1540 while (input >> buff)
1542 if (buff ==
"$Nodes")
1544 input >> NumOfVertices;
1545 getline(input, buff);
1546 vertices.SetSize(NumOfVertices);
1548 const int gmsh_dim = 3;
1549 double coord[gmsh_dim];
1550 for (
int ver = 0; ver < NumOfVertices; ++ver)
1554 input.read(reinterpret_cast<char*>(&serial_number),
sizeof(
int));
1555 input.read(reinterpret_cast<char*>(coord), gmsh_dim*
sizeof(
double));
1559 input >> serial_number;
1560 for (
int ci = 0; ci < gmsh_dim; ++ci)
1565 vertices[ver] =
Vertex(coord, gmsh_dim);
1566 vertices_map[serial_number] = ver;
1568 for (
int ci = 0; ci < gmsh_dim; ++ci)
1570 bb_min[ci] = (ver == 0) ? coord[ci] :
1571 std::min(bb_min[ci], coord[ci]);
1572 bb_max[ci] = (ver == 0) ? coord[ci] :
1573 std::max(bb_max[ci], coord[ci]);
1576 double bb_size = std::max(bb_max[0] - bb_min[0],
1577 std::max(bb_max[1] - bb_min[1],
1578 bb_max[2] - bb_min[2]));
1580 if (bb_max[1] - bb_min[1] > bb_size * bb_tol)
1584 if (bb_max[2] - bb_min[2] > bb_size * bb_tol)
1589 if (static_cast<int>(vertices_map.size()) != NumOfVertices)
1591 MFEM_ABORT(
"Gmsh file : vertices indices are not unique");
1594 else if (buff ==
"$Elements")
1596 int num_of_all_elements;
1597 input >> num_of_all_elements;
1599 getline(input, buff);
1602 int type_of_element;
1610 int nodes_of_gmsh_element[] =
1759 -1,-1,-1,-1,-1,-1,-1,
1775 -1,-1,-1,-1,-1,-1,-1,
1820 int lin3[] = {0,2,1};
1821 int lin4[] = {0,2,3,1};
1822 int tri6[] = {0,3,1,5,4,2};
1823 int tri10[] = {0,3,4,1,8,9,5,7,6,2};
1824 int quad9[] = {0,4,1,7,8,5,3,6,2};
1825 int quad16[] = {0,4,5,1,11,12,13,6,
1828 int tet10[] {0,4,1,6,5,2,7,9,8,3};
1829 int tet20[] = {0,4,5,1,9,16,6,8,7,2,
1830 11,17,15,18,19,13,10,14,12,3
1832 int hex27[] {0,8,1,9,20,11,3,13,2,
1833 10,21,12,22,26,23,15,24,14,
1834 4,16,5,17,25,18,7,19,6
1836 int hex64[] {0,8,9,1,10,32,35,14,
1837 11,33,34,15,3,19,18,2,
1838 12,36,37,16,40,56,57,44,
1839 43,59,58,45,22,49,48,20,
1840 13,39,38,17,41,60,61,47,
1841 42,63,62,46,23,50,51,21,
1842 4,24,25,5,26,52,53,28,
1843 27,55,54,29,7,31,30,6
1846 int wdg18[] = {0,6,1,7,9,2,8,15,10,
1847 16,17,11,3,12,4,13,14,5
1849 int wdg40[] = {0,6,7,1,8,24,12,9,13,2,
1850 10,26,27,14,30,38,34,33,35,16,
1851 11,29,28,15,31,39,37,32,36,17,
1852 3,18,19,4,20,25,22,21,23,5
1855 int pyr14[] = {0,5,1,6,13,8,3,
1858 int pyr30[] = {0,5,6,1,7,25,28,11,8,26,
1859 27,12,3,16,15,2,9,21,13,22,
1860 29,23,19,24,17,10,14,20,18,4
1863 vector<Element*> elements_0D, elements_1D, elements_2D, elements_3D;
1864 elements_0D.reserve(num_of_all_elements);
1865 elements_1D.reserve(num_of_all_elements);
1866 elements_2D.reserve(num_of_all_elements);
1867 elements_3D.reserve(num_of_all_elements);
1870 vector<Array<int>*> ho_verts_1D, ho_verts_2D, ho_verts_3D;
1871 ho_verts_1D.reserve(num_of_all_elements);
1872 ho_verts_2D.reserve(num_of_all_elements);
1873 ho_verts_3D.reserve(num_of_all_elements);
1876 vector<int> ho_el_order_1D, ho_el_order_2D, ho_el_order_3D;
1877 ho_el_order_1D.reserve(num_of_all_elements);
1878 ho_el_order_2D.reserve(num_of_all_elements);
1879 ho_el_order_3D.reserve(num_of_all_elements);
1891 ho_lin[2] = lin3; ho_lin[3] = lin4;
1892 ho_tri[2] = tri6; ho_tri[3] = tri10;
1893 ho_sqr[2] = quad9; ho_sqr[3] = quad16;
1894 ho_tet[2] = tet10; ho_tet[3] = tet20;
1895 ho_hex[2] = hex27; ho_hex[3] = hex64;
1896 ho_wdg[2] = wdg18; ho_wdg[3] = wdg40;
1897 ho_pyr[2] = pyr14; ho_pyr[3] = pyr30;
1899 bool has_nonpositive_phys_domain =
false;
1900 bool has_positive_phys_domain =
false;
1904 int n_elem_part = 0;
1905 const int header_size = 3;
1908 int header[header_size];
1909 int n_elem_one_type;
1911 while (n_elem_part < num_of_all_elements)
1913 input.read(reinterpret_cast<char*>(header),
1914 header_size*
sizeof(
int));
1915 type_of_element = header[0];
1916 n_elem_one_type = header[1];
1919 n_elem_part += n_elem_one_type;
1921 const int n_elem_nodes = nodes_of_gmsh_element[type_of_element-1];
1922 vector<int> data(1+n_tags+n_elem_nodes);
1923 for (
int el = 0; el < n_elem_one_type; ++el)
1925 input.read(reinterpret_cast<char*>(&data[0]),
1926 data.size()*
sizeof(int));
1928 serial_number = data[dd++];
1931 phys_domain = (n_tags > 0) ? data[dd++] : 1;
1934 elem_domain = (n_tags > 1) ? data[dd++] : 0;
1937 n_partitions = (n_tags > 2) ? data[dd++] : 0;
1940 vector<int> vert_indices(n_elem_nodes);
1941 for (
int vi = 0; vi < n_elem_nodes; ++vi)
1943 map<int, int>::const_iterator it =
1944 vertices_map.find(data[1+n_tags+vi]);
1945 if (it == vertices_map.end())
1947 MFEM_ABORT(
"Gmsh file : vertex index doesn't exist");
1949 vert_indices[vi] = it->second;
1957 if (phys_domain <= 0)
1959 has_nonpositive_phys_domain =
true;
1964 has_positive_phys_domain =
true;
1969 switch (type_of_element)
1982 elements_1D.push_back(
1983 new Segment(&vert_indices[0], phys_domain));
1984 if (type_of_element != 1)
1986 el_order = n_elem_nodes - 1;
1988 hov->
Append(&vert_indices[0], n_elem_nodes);
1989 ho_verts_1D.push_back(hov);
1990 ho_el_order_1D.push_back(el_order);
1996 case 21: el_order--;
1997 case 23: el_order--;
1998 case 25: el_order--;
1999 case 42: el_order--;
2000 case 43: el_order--;
2001 case 44: el_order--;
2002 case 45: el_order--;
2003 case 46: el_order--;
2005 elements_2D.push_back(
2006 new Triangle(&vert_indices[0], phys_domain));
2010 hov->
Append(&vert_indices[0], n_elem_nodes);
2011 ho_verts_2D.push_back(hov);
2012 ho_el_order_2D.push_back(el_order);
2017 case 10: el_order--;
2018 case 36: el_order--;
2019 case 37: el_order--;
2020 case 38: el_order--;
2021 case 47: el_order--;
2022 case 48: el_order--;
2023 case 49: el_order--;
2024 case 50: el_order--;
2025 case 51: el_order--;
2027 elements_2D.push_back(
2032 hov->
Append(&vert_indices[0], n_elem_nodes);
2033 ho_verts_2D.push_back(hov);
2034 ho_el_order_2D.push_back(el_order);
2039 case 11: el_order--;
2040 case 29: el_order--;
2041 case 30: el_order--;
2042 case 31: el_order--;
2043 case 71: el_order--;
2044 case 72: el_order--;
2045 case 73: el_order--;
2046 case 74: el_order--;
2047 case 75: el_order--;
2049 #ifdef MFEM_USE_MEMALLOC
2050 elements_3D.push_back(TetMemory.Alloc());
2051 elements_3D.back()->SetVertices(&vert_indices[0]);
2052 elements_3D.back()->SetAttribute(phys_domain);
2054 elements_3D.push_back(
2060 hov->
Append(&vert_indices[0], n_elem_nodes);
2061 ho_verts_3D.push_back(hov);
2062 ho_el_order_3D.push_back(el_order);
2067 case 12: el_order--;
2068 case 92: el_order--;
2069 case 93: el_order--;
2070 case 94: el_order--;
2071 case 95: el_order--;
2072 case 96: el_order--;
2073 case 97: el_order--;
2074 case 98: el_order--;
2077 elements_3D.push_back(
2078 new Hexahedron(&vert_indices[0], phys_domain));
2082 hov->
Append(&vert_indices[0], n_elem_nodes);
2083 ho_verts_3D.push_back(hov);
2084 ho_el_order_3D.push_back(el_order);
2089 case 13: el_order--;
2090 case 90: el_order--;
2091 case 91: el_order--;
2092 case 106: el_order--;
2093 case 107: el_order--;
2094 case 108: el_order--;
2095 case 109: el_order--;
2096 case 110: el_order--;
2099 elements_3D.push_back(
2100 new Wedge(&vert_indices[0], phys_domain));
2104 hov->
Append(&vert_indices[0], n_elem_nodes);
2105 ho_verts_3D.push_back(hov);
2106 ho_el_order_3D.push_back(el_order);
2137 elements_0D.push_back(
2138 new Point(&vert_indices[0], phys_domain));
2142 MFEM_WARNING(
"Unsupported Gmsh element type.");
2151 for (
int el = 0; el < num_of_all_elements; ++el)
2153 input >> serial_number >> type_of_element >> n_tags;
2154 vector<int> data(n_tags);
2155 for (
int i = 0; i < n_tags; ++i) { input >> data[i]; }
2158 phys_domain = (n_tags > 0) ? data[0] : 1;
2161 elem_domain = (n_tags > 1) ? data[1] : 0;
2164 n_partitions = (n_tags > 2) ? data[2] : 0;
2167 const int n_elem_nodes = nodes_of_gmsh_element[type_of_element-1];
2168 vector<int> vert_indices(n_elem_nodes);
2170 for (
int vi = 0; vi < n_elem_nodes; ++vi)
2173 map<int, int>::const_iterator it = vertices_map.find(index);
2174 if (it == vertices_map.end())
2176 MFEM_ABORT(
"Gmsh file : vertex index doesn't exist");
2178 vert_indices[vi] = it->second;
2186 if (phys_domain <= 0)
2188 has_nonpositive_phys_domain =
true;
2193 has_positive_phys_domain =
true;
2198 switch (type_of_element)
2211 elements_1D.push_back(
2212 new Segment(&vert_indices[0], phys_domain));
2213 if (type_of_element != 1)
2216 hov->
Append(&vert_indices[0], n_elem_nodes);
2217 ho_verts_1D.push_back(hov);
2218 el_order = n_elem_nodes - 1;
2219 ho_el_order_1D.push_back(el_order);
2225 case 21: el_order--;
2226 case 23: el_order--;
2227 case 25: el_order--;
2228 case 42: el_order--;
2229 case 43: el_order--;
2230 case 44: el_order--;
2231 case 45: el_order--;
2232 case 46: el_order--;
2234 elements_2D.push_back(
2235 new Triangle(&vert_indices[0], phys_domain));
2239 hov->
Append(&vert_indices[0], n_elem_nodes);
2240 ho_verts_2D.push_back(hov);
2241 ho_el_order_2D.push_back(el_order);
2246 case 10: el_order--;
2247 case 36: el_order--;
2248 case 37: el_order--;
2249 case 38: el_order--;
2250 case 47: el_order--;
2251 case 48: el_order--;
2252 case 49: el_order--;
2253 case 50: el_order--;
2254 case 51: el_order--;
2256 elements_2D.push_back(
2261 hov->
Append(&vert_indices[0], n_elem_nodes);
2262 ho_verts_2D.push_back(hov);
2263 ho_el_order_2D.push_back(el_order);
2268 case 11: el_order--;
2269 case 29: el_order--;
2270 case 30: el_order--;
2271 case 31: el_order--;
2272 case 71: el_order--;
2273 case 72: el_order--;
2274 case 73: el_order--;
2275 case 74: el_order--;
2276 case 75: el_order--;
2278 #ifdef MFEM_USE_MEMALLOC
2279 elements_3D.push_back(TetMemory.Alloc());
2280 elements_3D.back()->SetVertices(&vert_indices[0]);
2281 elements_3D.back()->SetAttribute(phys_domain);
2283 elements_3D.push_back(
2289 hov->
Append(&vert_indices[0], n_elem_nodes);
2290 ho_verts_3D.push_back(hov);
2291 ho_el_order_3D.push_back(el_order);
2296 case 12: el_order--;
2297 case 92: el_order--;
2298 case 93: el_order--;
2299 case 94: el_order--;
2300 case 95: el_order--;
2301 case 96: el_order--;
2302 case 97: el_order--;
2303 case 98: el_order--;
2306 elements_3D.push_back(
2307 new Hexahedron(&vert_indices[0], phys_domain));
2311 hov->
Append(&vert_indices[0], n_elem_nodes);
2312 ho_verts_3D.push_back(hov);
2313 ho_el_order_3D.push_back(el_order);
2318 case 13: el_order--;
2319 case 90: el_order--;
2320 case 91: el_order--;
2321 case 106: el_order--;
2322 case 107: el_order--;
2323 case 108: el_order--;
2324 case 109: el_order--;
2325 case 110: el_order--;
2328 elements_3D.push_back(
2329 new Wedge(&vert_indices[0], phys_domain));
2333 hov->
Append(&vert_indices[0], n_elem_nodes);
2334 ho_verts_3D.push_back(hov);
2335 ho_el_order_3D.push_back(el_order);
2366 elements_0D.push_back(
2367 new Point(&vert_indices[0], phys_domain));
2371 MFEM_WARNING(
"Unsupported Gmsh element type.");
2378 if (has_positive_phys_domain && has_nonpositive_phys_domain)
2380 MFEM_ABORT(
"Non-positive element attribute in Gmsh mesh!\n"
2381 "By default Gmsh sets element tags (attributes)"
2382 " to '0' but MFEM requires that they be"
2383 " positive integers.\n"
2384 "Use \"Physical Curve\", \"Physical Surface\","
2385 " or \"Physical Volume\" to set tags/attributes"
2386 " for all curves, surfaces, or volumes in your"
2387 " Gmsh geometry to values which are >= 1.");
2389 else if (has_nonpositive_phys_domain)
2391 mfem::out <<
"\nGmsh reader: all element attributes were zero.\n"
2392 <<
"MFEM only supports positive element attributes.\n"
2393 <<
"Setting element attributes to 1.\n\n";
2396 if (!elements_3D.empty())
2399 NumOfElements = elements_3D.size();
2400 elements.SetSize(NumOfElements);
2401 for (
int el = 0; el < NumOfElements; ++el)
2403 elements[el] = elements_3D[el];
2405 NumOfBdrElements = elements_2D.size();
2406 boundary.SetSize(NumOfBdrElements);
2407 for (
int el = 0; el < NumOfBdrElements; ++el)
2409 boundary[el] = elements_2D[el];
2411 for (
size_t el = 0; el < ho_el_order_3D.size(); el++)
2413 mesh_order = max(mesh_order, ho_el_order_3D[el]);
2416 for (
size_t el = 0; el < elements_1D.size(); ++el)
2418 delete elements_1D[el];
2420 for (
size_t el = 0; el < elements_0D.size(); ++el)
2422 delete elements_0D[el];
2425 else if (!elements_2D.empty())
2428 NumOfElements = elements_2D.size();
2429 elements.SetSize(NumOfElements);
2430 for (
int el = 0; el < NumOfElements; ++el)
2432 elements[el] = elements_2D[el];
2434 NumOfBdrElements = elements_1D.size();
2435 boundary.SetSize(NumOfBdrElements);
2436 for (
int el = 0; el < NumOfBdrElements; ++el)
2438 boundary[el] = elements_1D[el];
2440 for (
size_t el = 0; el < ho_el_order_2D.size(); el++)
2442 mesh_order = max(mesh_order, ho_el_order_2D[el]);
2445 for (
size_t el = 0; el < elements_0D.size(); ++el)
2447 delete elements_0D[el];
2450 else if (!elements_1D.empty())
2453 NumOfElements = elements_1D.size();
2454 elements.SetSize(NumOfElements);
2455 for (
int el = 0; el < NumOfElements; ++el)
2457 elements[el] = elements_1D[el];
2459 NumOfBdrElements = elements_0D.size();
2460 boundary.SetSize(NumOfBdrElements);
2461 for (
int el = 0; el < NumOfBdrElements; ++el)
2463 boundary[el] = elements_0D[el];
2465 for (
size_t el = 0; el < ho_el_order_1D.size(); el++)
2467 mesh_order = max(mesh_order, ho_el_order_1D[el]);
2472 MFEM_ABORT(
"Gmsh file : no elements found");
2486 this->FinalizeTopology();
2492 BasisType::ClosedUniform);
2500 for (
int el = 0; el < NumOfElements; el++)
2502 const int * vm = NULL;
2504 switch (GetElementType(el))
2506 case Element::SEGMENT:
2507 ho_verts = ho_verts_1D[el];
2508 el_order = ho_el_order_1D[el];
2509 if (!ho_lin[el_order])
2511 ho_lin[el_order] =
new int[ho_verts->
Size()];
2514 vm = ho_lin[el_order];
2516 case Element::TRIANGLE:
2517 ho_verts = ho_verts_2D[el];
2518 el_order = ho_el_order_2D[el];
2519 if (!ho_tri[el_order])
2521 ho_tri[el_order] =
new int[ho_verts->
Size()];
2524 vm = ho_tri[el_order];
2526 case Element::QUADRILATERAL:
2527 ho_verts = ho_verts_2D[el];
2528 el_order = ho_el_order_2D[el];
2529 if (!ho_sqr[el_order])
2531 ho_sqr[el_order] =
new int[ho_verts->
Size()];
2534 vm = ho_sqr[el_order];
2536 case Element::TETRAHEDRON:
2537 ho_verts = ho_verts_3D[el];
2538 el_order = ho_el_order_3D[el];
2539 if (!ho_tet[el_order])
2541 ho_tet[el_order] =
new int[ho_verts->
Size()];
2544 vm = ho_tet[el_order];
2546 case Element::HEXAHEDRON:
2547 ho_verts = ho_verts_3D[el];
2548 el_order = ho_el_order_3D[el];
2549 if (!ho_hex[el_order])
2551 ho_hex[el_order] =
new int[ho_verts->
Size()];
2554 vm = ho_hex[el_order];
2556 case Element::WEDGE:
2557 ho_verts = ho_verts_3D[el];
2558 el_order = ho_el_order_3D[el];
2559 if (!ho_wdg[el_order])
2561 ho_wdg[el_order] =
new int[ho_verts->
Size()];
2564 vm = ho_wdg[el_order];
2577 MFEM_WARNING(
"Unsupported Gmsh element type.");
2580 int nv = (ho_verts) ? ho_verts->
Size() : 0;
2582 for (
int v = 0; v<nv; v++)
2584 double * c = GetVertex((*ho_verts)[vm[v]]);
2585 for (
int d=0; d<spaceDim; d++)
2587 Nodes_gf(spaceDim * (o + v) + d) = c[d];
2595 for (
size_t el=0; el<ho_verts_1D.size(); el++)
2597 delete ho_verts_1D[el];
2599 for (
size_t el=0; el<ho_verts_2D.size(); el++)
2601 delete ho_verts_2D[el];
2603 for (
size_t el=0; el<ho_verts_3D.size(); el++)
2605 delete ho_verts_3D[el];
2609 for (
int ord=4; ord<ho_lin.
Size(); ord++)
2611 if (ho_lin[ord] != NULL) {
delete [] ho_lin[ord]; }
2613 for (
int ord=4; ord<ho_tri.
Size(); ord++)
2615 if (ho_tri[ord] != NULL) {
delete [] ho_tri[ord]; }
2617 for (
int ord=4; ord<ho_sqr.
Size(); ord++)
2619 if (ho_sqr[ord] != NULL) {
delete [] ho_sqr[ord]; }
2621 for (
int ord=4; ord<ho_tet.
Size(); ord++)
2623 if (ho_tet[ord] != NULL) {
delete [] ho_tet[ord]; }
2625 for (
int ord=4; ord<ho_hex.
Size(); ord++)
2627 if (ho_hex[ord] != NULL) {
delete [] ho_hex[ord]; }
2629 for (
int ord=4; ord<ho_wdg.
Size(); ord++)
2631 if (ho_wdg[ord] != NULL) {
delete [] ho_wdg[ord]; }
2633 for (
int ord=4; ord<ho_pyr.
Size(); ord++)
2635 if (ho_pyr[ord] != NULL) {
delete [] ho_pyr[ord]; }
2643 else if (buff ==
"$Periodic")
2650 for (
int i = 0; i < v2v.
Size(); i++)
2656 input >> num_per_ent;
2657 getline(input, buff);
2658 for (
int i = 0; i < num_per_ent; i++)
2660 getline(input, buff);
2661 getline(input, buff);
2662 if (!strncmp(buff.c_str(),
"Affine", 6))
2668 num_nodes = atoi(buff.c_str());
2670 for (
int j=0; j<num_nodes; j++)
2673 input >> slave >> master;
2674 v2v[slave - 1] = master - 1;
2676 getline(input, buff);
2683 for (
int slave = 0; slave < v2v.
Size(); slave++)
2685 int master = v2v[slave];
2686 if (master != slave)
2689 while (v2v[master] != master && master != slave)
2691 master = v2v[master];
2693 if (master == slave)
2702 v2v[slave] = master;
2708 if (mesh_order == 1)
2710 this->FinalizeTopology();
2712 this->SetCurvature(1,
true, spaceDim, Ordering::byVDIM);
2717 for (
int i = 0; i < this->GetNE(); i++)
2719 Element *el = this->GetElement(i);
2722 for (
int j = 0; j < nv; j++)
2729 for (
int i = 0; i < this->GetNBE(); i++)
2731 Element *el = this->GetBdrElement(i);
2734 for (
int j = 0; j < nv; j++)
2742 this->RemoveUnusedVertices();
2745 this->RemoveInternalBoundaries();
2747 this->FinalizeTopology();
2752 SetCurvature(mesh_order, periodic, spaceDim, Ordering::byVDIM);
2755 Nodes->ProjectCoefficient(NodesCoef);
2760 #ifdef MFEM_USE_NETCDF
2761 void Mesh::ReadCubit(
const char *filename,
int &curved,
int &read_gf)
2768 const int sideMapTri3[3][2] =
2775 const int sideMapQuad4[4][2] =
2783 const int sideMapTri6[3][3] =
2790 const int sideMapQuad9[4][3] =
2798 const int sideMapTet4[4][3] =
2806 const int sideMapTet10[4][6] =
2814 const int sideMapHex8[6][4] =
2824 const int sideMapHex27[6][9] =
2826 {1,2,6,5,9,14,17,13,26},
2827 {2,3,7,6,10,15,18,14,25},
2828 {4,3,7,8,11,15,19,16,27},
2829 {1,4,8,5,12,16,20,13,24},
2830 {1,4,3,2,12,11,10,9,22},
2831 {5,8,7,6,20,19,18,17,23}
2836 const int mfemToGenesisTet10[10] = {1,2,3,4,5,7,8,6,9,10};
2839 const int mfemToGenesisHex27[27] = {1,2,3,4,5,6,7,8,9,10,11,
2841 12,17,18,19,20,13,14,15,
2843 16,22,26,25,27,24,23,21
2846 const int mfemToGenesisTri6[6] = {1,2,3,4,5,6};
2847 const int mfemToGenesisQuad9[9] = {1,2,3,4,5,6,7,8,9};
2854 char str_dummy[256];
2861 if ((retval = nc_open(filename, NC_NOWRITE, &ncid)))
2863 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
2869 size_t num_dim=0, num_nodes=0, num_elem=0, num_el_blk=0, num_side_sets=0;
2871 if ((retval = nc_inq_dimid(ncid,
"num_dim", &
id)) ||
2872 (retval = nc_inq_dim(ncid,
id, str_dummy, &num_dim)) ||
2874 (retval = nc_inq_dimid(ncid,
"num_nodes", &
id)) ||
2875 (retval = nc_inq_dim(ncid,
id, str_dummy, &num_nodes)) ||
2877 (retval = nc_inq_dimid(ncid,
"num_elem", &
id)) ||
2878 (retval = nc_inq_dim(ncid,
id, str_dummy, &num_elem)) ||
2880 (retval = nc_inq_dimid(ncid,
"num_el_blk", &
id)) ||
2881 (retval = nc_inq_dim(ncid,
id, str_dummy, &num_el_blk)))
2883 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
2885 if ((retval = nc_inq_dimid(ncid,
"num_side_sets", &
id)) ||
2886 (retval = nc_inq_dim(ncid,
id, str_dummy, &num_side_sets)))
2894 size_t *num_el_in_blk =
new size_t[num_el_blk];
2895 size_t num_node_per_el;
2897 int previous_num_node_per_el = 0;
2898 for (
int i = 0; i < (int) num_el_blk; i++)
2900 sprintf(temp_str,
"num_el_in_blk%d", i+1);
2901 if ((retval = nc_inq_dimid(ncid, temp_str, &temp_id)) ||
2902 (retval = nc_inq_dim(ncid, temp_id, str_dummy, &num_el_in_blk[i])))
2904 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
2907 sprintf(temp_str,
"num_nod_per_el%d", i+1);
2908 if ((retval = nc_inq_dimid(ncid, temp_str, &temp_id)) ||
2909 (retval = nc_inq_dim(ncid, temp_id, str_dummy, &num_node_per_el)))
2911 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
2918 if ((
int) num_node_per_el != previous_num_node_per_el)
2920 MFEM_ABORT(
"Element blocks of different element types not supported");
2923 previous_num_node_per_el = num_node_per_el;
2927 enum CubitElementType
2949 CubitElementType cubit_element_type = ELEMENT_TRI3;
2950 CubitFaceType cubit_face_type = FACE_EDGE2;
2951 int num_element_linear_nodes = 0;
2955 switch (num_node_per_el)
2959 cubit_element_type = ELEMENT_TRI3;
2960 cubit_face_type = FACE_EDGE2;
2961 num_element_linear_nodes = 3;
2966 cubit_element_type = ELEMENT_TRI6;
2967 cubit_face_type = FACE_EDGE3;
2968 num_element_linear_nodes = 3;
2973 cubit_element_type = ELEMENT_QUAD4;
2974 cubit_face_type = FACE_EDGE2;
2975 num_element_linear_nodes = 4;
2980 cubit_element_type = ELEMENT_QUAD9;
2981 cubit_face_type = FACE_EDGE3;
2982 num_element_linear_nodes = 4;
2987 MFEM_ABORT(
"Don't know what to do with a " << num_node_per_el <<
2988 " node 2D element\n");
2992 else if (num_dim == 3)
2994 switch (num_node_per_el)
2998 cubit_element_type = ELEMENT_TET4;
2999 cubit_face_type = FACE_TRI3;
3000 num_element_linear_nodes = 4;
3005 cubit_element_type = ELEMENT_TET10;
3006 cubit_face_type = FACE_TRI6;
3007 num_element_linear_nodes = 4;
3012 cubit_element_type = ELEMENT_HEX8;
3013 cubit_face_type = FACE_QUAD4;
3014 num_element_linear_nodes = 8;
3019 cubit_element_type = ELEMENT_HEX27;
3020 cubit_face_type = FACE_QUAD9;
3021 num_element_linear_nodes = 8;
3026 MFEM_ABORT(
"Don't know what to do with a " << num_node_per_el <<
3027 " node 3D element\n");
3033 MFEM_ABORT(
"Invalid dimension: num_dim = " << num_dim);
3038 if (cubit_element_type == ELEMENT_TRI3 || cubit_element_type == ELEMENT_QUAD4 ||
3039 cubit_element_type == ELEMENT_TET4 || cubit_element_type == ELEMENT_HEX8)
3043 else if (cubit_element_type == ELEMENT_TRI6 ||
3044 cubit_element_type == ELEMENT_QUAD9 ||
3045 cubit_element_type == ELEMENT_TET10 || cubit_element_type == ELEMENT_HEX27)
3051 size_t *num_side_in_ss =
new size_t[num_side_sets];
3052 for (
int i = 0; i < (int) num_side_sets; i++)
3054 sprintf(temp_str,
"num_side_ss%d", i+1);
3055 if ((retval = nc_inq_dimid(ncid, temp_str, &temp_id)) ||
3056 (retval = nc_inq_dim(ncid, temp_id, str_dummy, &num_side_in_ss[i])))
3058 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
3063 double *coordx =
new double[num_nodes];
3064 double *coordy =
new double[num_nodes];
3065 double *coordz =
new double[num_nodes];
3067 if ((retval = nc_inq_varid(ncid,
"coordx", &
id)) ||
3068 (retval = nc_get_var_double(ncid,
id, coordx)) ||
3069 (retval = nc_inq_varid(ncid,
"coordy", &
id)) ||
3070 (retval = nc_get_var_double(ncid,
id, coordy)))
3072 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
3077 if ((retval = nc_inq_varid(ncid,
"coordz", &
id)) ||
3078 (retval = nc_get_var_double(ncid,
id, coordz)))
3080 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
3085 int **elem_blk =
new int*[num_el_blk];
3086 for (
int i = 0; i < (int) num_el_blk; i++)
3088 elem_blk[i] =
new int[num_el_in_blk[i] * num_node_per_el];
3089 sprintf(temp_str,
"connect%d", i+1);
3090 if ((retval = nc_inq_varid(ncid, temp_str, &temp_id)) ||
3091 (retval = nc_get_var_int(ncid, temp_id, elem_blk[i])))
3093 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
3096 int *ebprop =
new int[num_el_blk];
3097 if ((retval = nc_inq_varid(ncid,
"eb_prop1", &
id)) ||
3098 (retval = nc_get_var_int(ncid,
id, ebprop)))
3100 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
3105 int **elem_ss =
new int*[num_side_sets];
3106 int **side_ss =
new int*[num_side_sets];
3108 for (
int i = 0; i < (int) num_side_sets; i++)
3110 elem_ss[i] =
new int[num_side_in_ss[i]];
3111 side_ss[i] =
new int[num_side_in_ss[i]];
3113 sprintf(temp_str,
"elem_ss%d", i+1);
3114 if ((retval = nc_inq_varid(ncid, temp_str, &temp_id)) ||
3115 (retval = nc_get_var_int(ncid, temp_id, elem_ss[i])))
3117 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
3120 sprintf(temp_str,
"side_ss%d",i+1);
3121 if ((retval = nc_inq_varid(ncid, temp_str, &temp_id)) ||
3122 (retval = nc_get_var_int(ncid, temp_id, side_ss[i])))
3124 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
3128 int *ssprop =
new int[num_side_sets];
3129 if ((num_side_sets > 0) &&
3130 ((retval = nc_inq_varid(ncid,
"ss_prop1", &
id)) ||
3131 (retval = nc_get_var_int(ncid,
id, ssprop))))
3133 MFEM_ABORT(
"Fatal NetCDF error: " << nc_strerror(retval));
3139 int num_face_nodes = 0;
3140 int num_face_linear_nodes = 0;
3142 switch (cubit_face_type)
3147 num_face_linear_nodes = 2;
3153 num_face_linear_nodes = 2;
3159 num_face_linear_nodes = 3;
3165 num_face_linear_nodes = 3;
3171 num_face_linear_nodes = 4;
3177 num_face_linear_nodes = 4;
3184 int *start_of_block =
new int[num_el_blk+1];
3185 start_of_block[0] = 0;
3186 for (
int i = 1; i < (int) num_el_blk+1; i++)
3188 start_of_block[i] = start_of_block[i-1] + num_el_in_blk[i-1];
3191 int **ss_node_id =
new int*[num_side_sets];
3193 for (
int i = 0; i < (int) num_side_sets; i++)
3195 ss_node_id[i] =
new int[num_side_in_ss[i]*num_face_nodes];
3196 for (
int j = 0; j < (int) num_side_in_ss[i]; j++)
3198 int glob_ind = elem_ss[i][j]-1;
3201 while (iblk < (
int) num_el_blk && glob_ind >= start_of_block[iblk+1])
3205 if (iblk >= (
int) num_el_blk)
3207 MFEM_ABORT(
"Sideset element does not exist");
3209 loc_ind = glob_ind - start_of_block[iblk];
3210 int this_side = side_ss[i][j];
3211 int ielem = loc_ind*num_node_per_el;
3213 for (
int k = 0; k < num_face_nodes; k++)
3216 switch (cubit_element_type)
3218 case (ELEMENT_TRI3):
3220 inode = sideMapTri3[this_side-1][k];
3223 case (ELEMENT_TRI6):
3225 inode = sideMapTri6[this_side-1][k];
3228 case (ELEMENT_QUAD4):
3230 inode = sideMapQuad4[this_side-1][k];
3233 case (ELEMENT_QUAD9):
3235 inode = sideMapQuad9[this_side-1][k];
3238 case (ELEMENT_TET4):
3240 inode = sideMapTet4[this_side-1][k];
3243 case (ELEMENT_TET10):
3245 inode = sideMapTet10[this_side-1][k];
3248 case (ELEMENT_HEX8):
3250 inode = sideMapHex8[this_side-1][k];
3253 case (ELEMENT_HEX27):
3255 inode = sideMapHex27[this_side-1][k];
3259 ss_node_id[i][j*num_face_nodes+k] =
3260 elem_blk[iblk][ielem + inode - 1];
3266 std::vector<int> uniqueVertexID;
3268 for (
int iblk = 0; iblk < (int) num_el_blk; iblk++)
3270 for (
int i = 0; i < (int) num_el_in_blk[iblk]; i++)
3272 for (
int j = 0; j < num_element_linear_nodes; j++)
3274 uniqueVertexID.push_back(elem_blk[iblk][i*num_node_per_el + j]);
3278 std::sort(uniqueVertexID.begin(), uniqueVertexID.end());
3279 std::vector<int>::iterator newEnd;
3280 newEnd = std::unique(uniqueVertexID.begin(), uniqueVertexID.end());
3281 uniqueVertexID.resize(std::distance(uniqueVertexID.begin(), newEnd));
3287 std::map<int,int> cubitToMFEMVertMap;
3288 for (
int i = 0; i < (int) uniqueVertexID.size(); i++)
3290 cubitToMFEMVertMap[uniqueVertexID[i]] = i+1;
3292 MFEM_ASSERT(cubitToMFEMVertMap.size() == uniqueVertexID.size(),
3293 "This should never happen\n");
3299 NumOfVertices = uniqueVertexID.size();
3300 vertices.SetSize(NumOfVertices);
3301 for (
int i = 0; i < (int) uniqueVertexID.size(); i++)
3303 vertices[i](0) = coordx[uniqueVertexID[i] - 1];
3304 vertices[i](1) = coordy[uniqueVertexID[i] - 1];
3307 vertices[i](2) = coordz[uniqueVertexID[i] - 1];
3311 NumOfElements = num_elem;
3312 elements.SetSize(num_elem);
3314 int renumberedVertID[8];
3315 for (
int iblk = 0; iblk < (int) num_el_blk; iblk++)
3317 int NumNodePerEl = num_node_per_el;
3318 for (
int i = 0; i < (int) num_el_in_blk[iblk]; i++)
3320 for (
int j = 0; j < num_element_linear_nodes; j++)
3322 renumberedVertID[j] =
3323 cubitToMFEMVertMap[elem_blk[iblk][i*NumNodePerEl+j]]-1;
3326 switch (cubit_element_type)
3328 case (ELEMENT_TRI3):
3329 case (ELEMENT_TRI6):
3331 elements[elcount] =
new Triangle(renumberedVertID,ebprop[iblk]);
3334 case (ELEMENT_QUAD4):
3335 case (ELEMENT_QUAD9):
3337 elements[elcount] =
new Quadrilateral(renumberedVertID,ebprop[iblk]);
3340 case (ELEMENT_TET4):
3341 case (ELEMENT_TET10):
3343 #ifdef MFEM_USE_MEMALLOC
3344 elements[elcount] = TetMemory.Alloc();
3345 elements[elcount]->SetVertices(renumberedVertID);
3346 elements[elcount]->SetAttribute(ebprop[iblk]);
3348 elements[elcount] =
new Tetrahedron(renumberedVertID,
3353 case (ELEMENT_HEX8):
3354 case (ELEMENT_HEX27):
3356 elements[elcount] =
new Hexahedron(renumberedVertID,ebprop[iblk]);
3366 NumOfBdrElements = 0;
3367 for (
int iss = 0; iss < (int) num_side_sets; iss++)
3369 NumOfBdrElements += num_side_in_ss[iss];
3371 boundary.SetSize(NumOfBdrElements);
3373 for (
int iss = 0; iss < (int) num_side_sets; iss++)
3375 for (
int i = 0; i < (int) num_side_in_ss[iss]; i++)
3377 for (
int j = 0; j < num_face_linear_nodes; j++)
3379 renumberedVertID[j] =
3380 cubitToMFEMVertMap[ss_node_id[iss][i*num_face_nodes+j]] - 1;
3382 switch (cubit_face_type)
3387 boundary[sidecount] =
new Segment(renumberedVertID,ssprop[iss]);
3393 boundary[sidecount] =
new Triangle(renumberedVertID,ssprop[iss]);
3399 boundary[sidecount] =
new Quadrilateral(renumberedVertID,ssprop[iss]);
3412 switch (cubit_element_type)
3414 case (ELEMENT_TRI6):
3416 mymap = (
int *) mfemToGenesisTri6;
3419 case (ELEMENT_QUAD9):
3421 mymap = (
int *) mfemToGenesisQuad9;
3424 case (ELEMENT_TET10):
3426 mymap = (
int *) mfemToGenesisTet10;
3429 case (ELEMENT_HEX27):
3431 mymap = (
int *) mfemToGenesisHex27;
3434 case (ELEMENT_TRI3):
3435 case (ELEMENT_QUAD4):
3436 case (ELEMENT_TET4):
3437 case (ELEMENT_HEX8):
3439 MFEM_ABORT(
"Something went wrong. Linear elements detected when order is 2.");
3451 Nodes->MakeOwner(fec);
3459 for (
int i = 0; i < NumOfElements; i++)
3466 for (
int l = 0; l < dofs.
Size(); l++) { vdofs[l] = dofs[l]; }
3470 while (iblk < (
int) num_el_blk && i >= start_of_block[iblk+1]) { iblk++; }
3471 loc_ind = i - start_of_block[iblk];
3472 for (
int j = 0; j < dofs.
Size(); j++)
3474 int point_id = elem_blk[iblk][loc_ind*num_node_per_el + mymap[j] - 1] - 1;
3475 (*Nodes)(vdofs[j]) = coordx[point_id];
3476 (*Nodes)(vdofs[j]+1) = coordy[point_id];
3479 (*Nodes)(vdofs[j]+2) = coordz[point_id];
3489 for (
int i = 0; i < (int) num_side_sets; i++)
3491 delete [] elem_ss[i];
3492 delete [] side_ss[i];
3497 delete [] num_el_in_blk;
3498 delete [] num_side_in_ss;
3503 for (
int i = 0; i < (int) num_el_blk; i++)
3505 delete [] elem_blk[i];
3509 delete [] start_of_block;
3511 for (
int i = 0; i < (int) num_side_sets; i++)
3513 delete [] ss_node_id[i];
3515 delete [] ss_node_id;
3520 #endif // #ifdef MFEM_USE_NETCDF
Abstract class for all finite elements.
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Mesh * Make3D(int nsteps, double rstep, double aspect, int order, bool sfc)
int Size() const
Return the logical size of the array.
void Load(std::istream &in, int fmt=0)
Read an Array from the stream in using format fmt. The format fmt can be:
Class for grid function - Vector with associated FE space.
virtual void GetVertices(Array< int > &v) const =0
Returns element's vertices.
void GmshHOSegmentMapping(int order, int *map)
Generate Gmsh vertex mapping for a Segment.
T * GetData()
Returns the data.
int Size() const
Returns the size of the vector.
void MakeOwner(FiniteElementCollection *fec_)
Make the GridFunction the owner of fec and fes.
virtual void SetVertices(const int *ind)
Set the vertices according to the given input.
Data type quadrilateral element.
double * GetData() const
Return a pointer to the beginning of the Vector data.
void skip_comment_lines(std::istream &is, const char comment_char)
Check if the stream starts with comment_char. If so skip it.
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...
Mesh * Make2D(int nsteps, double rstep, double phi, double aspect, int order, bool sfc)
Class for standard nodal finite elements.
void Load(std::istream **in, int np, int *dim)
Reads a vector from multiple files.
Data type hexahedron element.
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Type
Constants for the classes derived from Element.
Data type triangle element.
size_t NumBase64Chars(size_t nbytes)
Return the number of characters needed to encode nbytes in base-64.
const T * Read(const Memory< T > &mem, int size, bool on_dev=true)
Get a pointer for read access to mem with the mfem::Device's DeviceMemoryClass, if on_dev = true...
void GmshHOQuadrilateralMapping(int order, int *map)
Generate Gmsh vertex mapping for a Quadrilateral.
virtual DofTransformation * GetElementDofs(int elem, Array< int > &dofs) const
Returns indices of degrees of freedom of element 'elem'.
Data type tetrahedron element.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
void GmshHOTriangleMapping(int order, int *map)
Generate Gmsh vertex mapping for a Triangle.
void GmshHOTetrahedronMapping(int order, int *map)
Generate Gmsh vertex mapping for a Tetrahedron.
const Array< int > & GetLexicographicOrdering() const
Get an Array<int> that maps lexicographically ordered indices to the indices of the respective nodes/...
bool StringCompare(const char *s1, const char *s2)
const char * VTKByteOrder()
Determine the byte order and return either "BigEndian" or "LittleEndian".
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
int index(int i, int j, int nx, int ny)
virtual int GetNVertices() const =0
void SetAttribute(const int attr)
Set element's attribute.
Piecewise-(bi)quadratic continuous finite elements.
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the GridFunction.
void filter_dos(std::string &line)
Check for, and remove, a trailing '\r' from and std::string.
Vector coefficient defined by a vector GridFunction.
void GmshHOHexahedronMapping(int order, int *map)
Generate Gmsh vertex mapping for a Hexahedron.
void CreateVTKElementConnectivity(Array< int > &con, Geometry::Type geom, int ref)
Create the VTK element connectivity array for a given element geometry and refinement level...
Arbitrary order H1-conforming (continuous) finite elements.
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...
Abstract data type element.
Data type line segment element.
Arbitrary order "L2-conforming" discontinuous finite elements.
void DofsToVDofs(Array< int > &dofs, int ndofs=-1) const