27#include <boost/algorithm/string.hpp>
28#include <boost/archive/binary_iarchive.hpp>
29#include <boost/io/ios_state.hpp>
30#include <boost/property_tree/ptree.hpp>
31#include <boost/property_tree/xml_parser.hpp>
32#include <boost/serialization/serialization.hpp>
34#ifdef DEAL_II_GMSH_WITH_API
45#ifdef DEAL_II_WITH_ASSIMP
46# include <assimp/Importer.hpp>
47# include <assimp/postprocess.h>
48# include <assimp/scene.h>
51#ifdef DEAL_II_TRILINOS_WITH_SEACAS
78 template <
int spacedim>
80 assign_1d_boundary_ids(
87 if (cell->face(face_no)->at_boundary())
88 for (const auto &pair : boundary_ids)
89 if (cell->face(face_no)->vertex(0) == pair.
first)
91 cell->face(face_no)->set_boundary_id(pair.second);
97 template <
int dim,
int spacedim>
99 assign_1d_boundary_ids(
111 template <
int dim,
int spacedim>
119 const auto n_hypercube_vertices =
120 ReferenceCells::get_hypercube<dim>().n_vertices();
121 bool is_only_hypercube =
true;
123 if (cell.
vertices.size() != n_hypercube_vertices)
125 is_only_hypercube =
false;
133 if (is_only_hypercube)
138template <
int dim,
int spacedim>
140 :
tria(nullptr, typeid(*this).name())
141 , default_format(ucd)
146template <
int dim,
int spacedim>
148 :
tria(&t, typeid(*this).name())
149 , default_format(ucd)
154template <
int dim,
int spacedim>
163template <
int dim,
int spacedim>
175 text[0] =
"# vtk DataFile Version 3.0";
178 text[3] =
"DATASET UNSTRUCTURED_GRID";
180 for (
unsigned int i = 0; i < 4; ++i)
185 line.compare(text[i]) == 0,
188 "While reading VTK file, failed to find a header line with text <") +
195 std::vector<Point<spacedim>>
vertices;
196 std::vector<CellData<dim>> cells;
205 if (keyword ==
"POINTS")
207 unsigned int n_vertices;
212 for (
unsigned int vertex = 0; vertex < n_vertices; ++vertex)
216 in >> x(0) >> x(1) >> x(2);
219 for (
unsigned int d = 0; d < spacedim; ++d)
227 "While reading VTK file, failed to find POINTS section"));
231 unsigned int n_geometric_objects = 0;
234 if (keyword ==
"CELLS")
237 std::vector<unsigned int> cell_types;
239 std::streampos oldpos = in.tellg();
242 while (in >> keyword)
243 if (keyword ==
"CELL_TYPES")
247 cell_types.resize(n_ints);
249 for (
unsigned int i = 0; i < n_ints; ++i)
258 in >> n_geometric_objects;
263 for (
unsigned int count = 0; count < n_geometric_objects; ++count)
265 unsigned int n_vertices;
269 if (cell_types[count] == 10 || cell_types[count] == 12)
277 cells.emplace_back(n_vertices);
279 for (
unsigned int j = 0; j < n_vertices;
281 in >> cells.back().vertices[j];
285 if (cell_types[count] == 12)
287 std::swap(cells.back().vertices[2],
288 cells.back().vertices[3]);
289 std::swap(cells.back().vertices[6],
290 cells.back().vertices[7]);
293 cells.back().material_id = 0;
296 else if (cell_types[count] == 5 || cell_types[count] == 9)
305 for (
unsigned int j = 0; j < n_vertices;
312 else if (cell_types[count] == 3)
316 for (
unsigned int j = 0; j < n_vertices;
327 "While reading VTK file, unknown cell type encountered"));
332 for (
unsigned int count = 0; count < n_geometric_objects; ++count)
334 unsigned int n_vertices;
338 if (cell_types[count] == 5 || cell_types[count] == 9)
345 cells.emplace_back(n_vertices);
347 for (
unsigned int j = 0; j < n_vertices;
349 in >> cells.back().vertices[j];
353 if (cell_types[count] == 9)
357 std::swap(cells.back().vertices[2],
358 cells.back().vertices[3]);
361 cells.back().material_id = 0;
364 else if (cell_types[count] == 3)
370 for (
unsigned int j = 0; j < n_vertices;
383 "While reading VTK file, unknown cell type encountered"));
388 for (
unsigned int count = 0; count < n_geometric_objects; ++count)
394 cell_types[count] == 3 && type == 2,
396 "While reading VTK file, unknown cell type encountered"));
397 cells.emplace_back(type);
399 for (
unsigned int j = 0; j < type; ++j)
400 in >> cells.back().vertices[j];
402 cells.back().material_id = 0;
408 "While reading VTK file, failed to find CELLS section"));
415 keyword ==
"CELL_TYPES",
417 "While reading VTK file, missing CELL_TYPES section. Found <" +
418 keyword +
"> instead.")));
422 n_ints == n_geometric_objects,
423 ExcMessage(
"The VTK reader found a CELL_DATA statement "
424 "that lists a total of " +
426 " cell data objects, but this needs to "
427 "equal the number of cells (which is " +
429 ") plus the number of quads (" +
431 " in 3d or the number of lines (" +
436 for (
unsigned int i = 0; i < n_ints; ++i)
440 while (in >> keyword)
441 if (keyword ==
"CELL_DATA")
447 ExcMessage(
"The VTK reader found a CELL_DATA statement "
448 "that lists a total of " +
450 " cell data objects, but this needs to "
451 "equal the number of cells (which is " +
453 ") plus the number of quads (" +
456 " in 3d or the number of lines (" +
461 const std::vector<std::string> data_sets{
"MaterialID",
464 for (
unsigned int i = 0; i < data_sets.size(); ++i)
467 while (in >> keyword)
468 if (keyword ==
"SCALARS")
473 std::string field_name;
475 if (std::find(data_sets.begin(),
477 field_name) == data_sets.end())
489 std::getline(in, line);
492 std::min(
static_cast<std::size_t
>(3),
493 line.size() - 1)) ==
"int",
495 "While reading VTK file, material- and manifold IDs can only have type 'int'."));
499 keyword ==
"LOOKUP_TABLE",
501 "While reading VTK file, missing keyword 'LOOKUP_TABLE'."));
505 keyword ==
"default",
507 "While reading VTK file, missing keyword 'default'."));
514 for (
unsigned int i = 0; i < cells.size(); ++i)
518 if (field_name ==
"MaterialID")
519 cells[i].material_id =
521 else if (field_name ==
"ManifoldID")
522 cells[i].manifold_id =
534 if (field_name ==
"MaterialID")
535 boundary_quad.material_id =
537 else if (field_name ==
"ManifoldID")
538 boundary_quad.manifold_id =
547 if (field_name ==
"MaterialID")
548 boundary_line.material_id =
550 else if (field_name ==
"ManifoldID")
551 boundary_line.manifold_id =
563 if (field_name ==
"MaterialID")
564 boundary_line.material_id =
566 else if (field_name ==
"ManifoldID")
567 boundary_line.manifold_id =
577 apply_grid_fixup_functions(
vertices, cells, subcelldata);
583 "While reading VTK file, failed to find CELLS section"));
588template <
int dim,
int spacedim>
592 namespace pt = boost::property_tree;
594 pt::read_xml(in, tree);
595 auto section = tree.get_optional<std::string>(
"VTKFile.dealiiData");
599 "While reading a VTU file, failed to find dealiiData section. "
600 "Notice that we can only read grid files in .vtu format that "
601 "were created by the deal.II library, using a call to "
602 "GridOut::write_vtu(), where the flag "
603 "GridOutFlags::Vtu::serialize_triangulation is set to true."));
607 const auto string_archive =
609 std::istringstream in_stream(string_archive);
610 boost::archive::binary_iarchive ia(in_stream);
615template <
int dim,
int spacedim>
619 Assert(
tria !=
nullptr, ExcNoTriangulationSelected());
623 skip_comment_lines(in,
'#');
634 "Expected '-1' before and after a section."));
646 std::getline(in, line);
648 boost::algorithm::trim(line);
649 if (line.compare(
"-1") == 0)
660 std::vector<Point<spacedim>>
vertices;
679 in >> dummy >> dummy >> dummy;
682 in >> x[0] >> x[1] >> x[2];
686 for (
unsigned int d = 0; d < spacedim; ++d)
701 AssertThrow(tmp == 2412, ExcUnknownSectionType(tmp));
703 std::vector<CellData<dim>> cells;
730 in >> type >> dummy >> dummy >> dummy >> dummy;
732 AssertThrow((type == 11) || (type == 44) || (type == 94) || (type == 115),
733 ExcUnknownElementType(type));
735 if ((((type == 44) || (type == 94)) && (dim == 2)) ||
736 ((type == 115) && (dim == 3)))
738 const auto reference_cell = ReferenceCells::get_hypercube<dim>();
739 cells.emplace_back();
744 .vertices[reference_cell.unv_vertex_to_deal_vertex(v)];
746 cells.back().material_id = 0;
749 cells.back().vertices[v] =
vertex_indices[cells.back().vertices[v]];
751 cell_indices[no] = no_cell;
755 else if (((type == 11) && (dim == 2)) ||
756 ((type == 11) && (dim == 3)))
759 in >> dummy >> dummy >> dummy;
764 for (
unsigned int &vertex :
770 for (
unsigned int &vertex :
774 line_indices[no] = no_line;
778 else if (((type == 44) || (type == 94)) && (dim == 3))
789 .vertices[reference_cell.unv_vertex_to_deal_vertex(v)];
793 for (
unsigned int &vertex :
797 quad_indices[no] = no_quad;
805 "> when running in dim=" +
824 AssertThrow((tmp == 2467) || (tmp == 2477), ExcUnknownSectionType(tmp));
840 in >> dummy >> dummy >> dummy >> dummy >> dummy >> dummy >>
851 std::getline(in, line);
854 "The line before the line containing an ID has too "
855 "many entries. This is not a valid UNV file."));
857 std::getline(in, line);
859 std::istringstream id_stream(line);
862 !id_stream.fail() && id_stream.eof(),
864 "The given UNV file contains a boundary or material id set to '" +
866 "', which cannot be parsed as a fixed-width integer, whereas "
867 "deal.II only supports integer boundary and material ids. To fix "
868 "this, ensure that all such ids are given integer values."));
871 id <=
long(std::numeric_limits<types::material_id>::max()),
872 ExcMessage(
"The provided integer id '" + std::to_string(
id) +
873 "' is not convertible to either types::material_id nor "
874 "types::boundary_id."));
876 const unsigned int n_lines =
877 (n_entities % 2 == 0) ? (n_entities / 2) : ((n_entities + 1) / 2);
879 for (
unsigned int line = 0; line < n_lines; ++line)
881 unsigned int n_fragments;
883 if (line == n_lines - 1)
884 n_fragments = (n_entities % 2 == 0) ? (2) : (1);
888 for (
unsigned int no_fragment = 0; no_fragment < n_fragments;
892 in >> dummy >> no >> dummy >> dummy;
894 if (cell_indices.count(no) > 0)
895 cells[cell_indices[no]].material_id = id;
897 if (line_indices.count(no) > 0)
901 if (quad_indices.count(no) > 0)
909 apply_grid_fixup_functions(
vertices, cells, subcelldata);
915template <
int dim,
int spacedim>
918 const bool apply_all_indicators_to_manifolds)
920 Assert(
tria !=
nullptr, ExcNoTriangulationSelected());
924 skip_comment_lines(in,
'#');
927 unsigned int n_vertices;
928 unsigned int n_cells;
931 in >> n_vertices >> n_cells >> dummy
937 std::vector<Point<spacedim>>
vertices(n_vertices);
943 for (
unsigned int vertex = 0; vertex < n_vertices; ++vertex)
950 in >> vertex_number >> x[0] >> x[1] >> x[2];
953 for (
unsigned int d = 0; d < spacedim; ++d)
962 std::vector<CellData<dim>> cells;
965 for (
unsigned int cell = 0; cell < n_cells; ++cell)
974 std::string cell_type;
978 unsigned int material_id;
984 if (((cell_type ==
"line") && (dim == 1)) ||
985 ((cell_type ==
"quad") && (dim == 2)) ||
986 ((cell_type ==
"hex") && (dim == 3)))
990 cells.emplace_back();
995 Assert(material_id <= std::numeric_limits<types::material_id>::max(),
998 std::numeric_limits<types::material_id>::max()));
1003 if (apply_all_indicators_to_manifolds)
1004 cells.back().manifold_id =
1006 cells.back().material_id = material_id;
1014 cells.back().vertices[i] =
1020 ExcInvalidVertexIndex(cell,
1021 cells.back().vertices[i]));
1026 else if ((cell_type ==
"line") && ((dim == 2) || (dim == 3)))
1034 Assert(material_id <= std::numeric_limits<types::boundary_id>::max(),
1037 std::numeric_limits<types::boundary_id>::max()));
1046 if (apply_all_indicators_to_manifolds)
1063 for (
unsigned int &vertex :
1071 AssertThrow(
false, ExcInvalidVertexIndex(cell, vertex));
1075 else if ((cell_type ==
"quad") && (dim == 3))
1084 Assert(material_id <= std::numeric_limits<types::boundary_id>::max(),
1087 std::numeric_limits<types::boundary_id>::max()));
1096 if (apply_all_indicators_to_manifolds)
1113 for (
unsigned int &vertex :
1121 Assert(
false, ExcInvalidVertexIndex(cell, vertex));
1127 AssertThrow(
false, ExcUnknownIdentifier(cell_type));
1132 apply_grid_fixup_functions(
vertices, cells, subcelldata);
1138 template <
int dim,
int spacedim>
1145 read_in_abaqus(std::istream &in);
1147 write_out_avs_ucd(std::ostream &out)
const;
1150 const double tolerance;
1153 get_global_node_numbers(
const int face_cell_no,
1154 const int face_cell_face_no)
const;
1157 std::vector<std::vector<double>> node_list;
1160 std::vector<std::vector<double>> cell_list;
1162 std::vector<std::vector<double>> face_list;
1165 std::map<std::string, std::vector<int>> elsets_list;
1169template <
int dim,
int spacedim>
1172 const bool apply_all_indicators_to_manifolds)
1174 Assert(
tria !=
nullptr, ExcNoTriangulationSelected());
1179 Assert((spacedim == 2 && dim == spacedim) ||
1180 (spacedim == 3 && (dim == spacedim || dim == spacedim - 1)),
1186 Abaqus_to_UCD<dim, spacedim> abaqus_to_ucd;
1187 abaqus_to_ucd.read_in_abaqus(in);
1189 std::stringstream in_ucd;
1190 abaqus_to_ucd.write_out_avs_ucd(in_ucd);
1198 read_ucd(in_ucd, apply_all_indicators_to_manifolds);
1200 catch (std::exception &exc)
1202 std::cerr <<
"Exception on processing internal UCD data: " << std::endl
1203 << exc.what() << std::endl;
1208 "Internal conversion from ABAQUS file to UCD format was unsuccessful. "
1209 "More information is provided in an error message printed above. "
1210 "Are you sure that your ABAQUS mesh file conforms with the requirements "
1211 "listed in the documentation?"));
1218 "Internal conversion from ABAQUS file to UCD format was unsuccessful. "
1219 "Are you sure that your ABAQUS mesh file conforms with the requirements "
1220 "listed in the documentation?"));
1225template <
int dim,
int spacedim>
1229 Assert(
tria !=
nullptr, ExcNoTriangulationSelected());
1235 skip_comment_lines(in,
'#');
1241 AssertThrow(line ==
"MeshVersionFormatted 0", ExcInvalidDBMESHInput(line));
1243 skip_empty_lines(in);
1247 AssertThrow(line ==
"Dimension", ExcInvalidDBMESHInput(line));
1248 unsigned int dimension;
1250 AssertThrow(dimension == dim, ExcDBMESHWrongDimension(dimension));
1251 skip_empty_lines(in);
1264 while (getline(in, line), line.find(
"# END") == std::string::npos)
1266 skip_empty_lines(in);
1271 AssertThrow(line ==
"Vertices", ExcInvalidDBMESHInput(line));
1273 unsigned int n_vertices;
1277 std::vector<Point<spacedim>>
vertices(n_vertices);
1278 for (
unsigned int vertex = 0; vertex < n_vertices; ++vertex)
1281 for (
unsigned int d = 0; d < dim; ++d)
1288 skip_empty_lines(in);
1294 AssertThrow(line ==
"Edges", ExcInvalidDBMESHInput(line));
1296 unsigned int n_edges;
1298 for (
unsigned int edge = 0; edge < n_edges; ++edge)
1301 in >> dummy >> dummy;
1307 skip_empty_lines(in);
1316 AssertThrow(line ==
"CrackedEdges", ExcInvalidDBMESHInput(line));
1319 for (
unsigned int edge = 0; edge < n_edges; ++edge)
1322 in >> dummy >> dummy;
1328 skip_empty_lines(in);
1334 AssertThrow(line ==
"Quadrilaterals", ExcInvalidDBMESHInput(line));
1336 static constexpr std::array<unsigned int, 8> local_vertex_numbering = {
1337 {0, 1, 5, 4, 2, 3, 7, 6}};
1338 std::vector<CellData<dim>> cells;
1340 unsigned int n_cells;
1342 for (
unsigned int cell = 0; cell < n_cells; ++cell)
1346 cells.emplace_back();
1350 cells.back().vertices[dim == 3 ? local_vertex_numbering[i] :
1354 (
static_cast<unsigned int>(cells.back().vertices[i]) <=
1356 ExcInvalidVertexIndex(cell, cells.back().vertices[i]));
1358 --cells.back().vertices[i];
1366 skip_empty_lines(in);
1374 while (getline(in, line), ((line.find(
"End") == std::string::npos) && (in)))
1380 apply_grid_fixup_functions(
vertices, cells, subcelldata);
1386template <
int dim,
int spacedim>
1390 Assert(
tria !=
nullptr, ExcNoTriangulationSelected());
1393 const auto reference_cell = ReferenceCells::get_hypercube<dim>();
1397 std::getline(in, line);
1399 unsigned int n_vertices;
1400 unsigned int n_cells;
1404 std::getline(in, line);
1407 std::getline(in, line);
1410 for (
unsigned int i = 0; i < 8; ++i)
1411 std::getline(in, line);
1414 std::vector<CellData<dim>> cells(n_cells);
1425 for (
unsigned int i = 0; i < GeometryInfo<dim>::vertices_per_cell; ++i)
1426 in >> cell.vertices[reference_cell.exodusii_vertex_to_deal_vertex(i)];
1430 std::vector<Point<spacedim>>
vertices(n_vertices);
1433 for (
unsigned int d = 0; d < spacedim; ++d)
1435 for (
unsigned int d = spacedim; d < 3; ++d)
1444 apply_grid_fixup_functions(
vertices, cells, subcelldata);
1450template <
int dim,
int spacedim>
1454 Assert(
tria !=
nullptr, ExcNoTriangulationSelected());
1467 std::stringstream whole_file;
1472 std::getline(in, line);
1480 if ((line.size() > 0) && (line.back() ==
'\r'))
1481 line.erase(line.size() - 1);
1486 if (line.find(
'#') != std::string::npos)
1487 line.erase(line.find(
'#'), std::string::npos);
1488 while ((line.size() > 0) && (line.back() ==
' '))
1489 line.erase(line.size() - 1);
1491 if (line.size() > 0)
1492 whole_file <<
'\n' << line;
1511 unsigned int version_major, version_minor;
1512 whole_file >> version_major >> version_minor;
1513 AssertThrow((version_major == 0) && (version_minor == 1),
1514 ExcMessage(
"deal.II can currently only read version 0.1 "
1515 "of the mphtxt file format."));
1520 unsigned int n_tags;
1521 whole_file >> n_tags;
1522 for (
unsigned int i = 0; i < n_tags; ++i)
1525 while (whole_file.peek() ==
'\n')
1527 std::getline(whole_file, dummy);
1533 unsigned int n_types;
1534 whole_file >> n_types;
1535 for (
unsigned int i = 0; i < n_types; ++i)
1538 while (whole_file.peek() ==
'\n')
1540 std::getline(whole_file, dummy);
1560 whole_file >> dummy >> dummy >> dummy;
1564 while (whole_file.peek() ==
'\n')
1566 std::getline(whole_file, s);
1568 ExcMessage(
"Expected '4 Mesh', but got '" + s +
"'."));
1571 unsigned int version;
1572 whole_file >> version;
1576 unsigned int file_space_dim;
1577 whole_file >> file_space_dim;
1581 "The mesh file uses a different number of space dimensions "
1582 "than the triangulation you want to read it into."));
1584 unsigned int n_vertices;
1585 whole_file >> n_vertices;
1587 unsigned int starting_vertex_index;
1588 whole_file >> starting_vertex_index;
1590 std::vector<Point<spacedim>>
vertices(n_vertices);
1591 for (
unsigned int v = 0; v < n_vertices; ++v)
1622 std::vector<CellData<dim>> cells;
1625 unsigned int n_types;
1626 whole_file >> n_types;
1627 for (
unsigned int type = 0; type < n_types; ++type)
1634 whole_file >> dummy;
1638 std::string object_name;
1639 whole_file >> object_name;
1641 static const std::map<std::string, ReferenceCell> name_to_type = {
1651 AssertThrow(name_to_type.find(object_name) != name_to_type.end(),
1652 ExcMessage(
"The input file contains a cell type <" +
1654 "> that the reader does not "
1655 "current support"));
1656 const ReferenceCell object_type = name_to_type.at(object_name);
1658 unsigned int n_vertices_per_element;
1659 whole_file >> n_vertices_per_element;
1661 unsigned int n_elements;
1662 whole_file >> n_elements;
1676 ExcMessage(
"Triangles should not appear in input files "
1684 "Quadrilaterals should not appear in input files "
1691 ExcMessage(
"Tetrahedra should not appear in input files "
1692 "for 1d or 2d meshes."));
1699 "Prisms (wedges) should not appear in input files "
1700 "for 1d or 2d meshes."));
1719 std::vector<unsigned int> vertices_for_this_element(
1720 n_vertices_per_element);
1721 for (
unsigned int e = 0; e < n_elements; ++e)
1724 for (
unsigned int v = 0; v < n_vertices_per_element; ++v)
1726 whole_file >> vertices_for_this_element[v];
1727 vertices_for_this_element[v] -= starting_vertex_index;
1736 cells.emplace_back();
1737 cells.back().vertices = vertices_for_this_element;
1743 vertices_for_this_element;
1751 cells.emplace_back();
1752 cells.back().vertices = vertices_for_this_element;
1758 vertices_for_this_element;
1766 cells.emplace_back();
1767 cells.back().vertices = vertices_for_this_element;
1780 unsigned int n_geom_entity_indices;
1781 whole_file >> n_geom_entity_indices;
1783 (n_geom_entity_indices == n_elements),
1790 if (n_geom_entity_indices != 0)
1792 for (
unsigned int e = 0; e < n_geom_entity_indices; ++e)
1795 unsigned int geometric_entity_index;
1796 whole_file >> geometric_entity_index;
1802 cells[cells.size() - n_elements + e].material_id =
1803 geometric_entity_index;
1808 .boundary_id = geometric_entity_index;
1814 cells[cells.size() - n_elements + e].material_id =
1815 geometric_entity_index;
1820 .boundary_id = geometric_entity_index;
1826 cells[cells.size() - n_elements + e].material_id =
1827 geometric_entity_index;
1853 if (line.vertices[1] < line.vertices[0])
1854 std::swap(line.vertices[0], line.vertices[1]);
1859 return std::lexicographical_compare(a.vertices.begin(),
1881 Assert((face.vertices.size() == 3) || (face.vertices.size() == 4),
1883 std::sort(face.vertices.begin(), face.vertices.end());
1888 return std::lexicographical_compare(a.vertices.begin(),
1899 for (
const auto &face : cell->face_iterators())
1900 if (face->at_boundary())
1906 std::array<unsigned int, 2> face_vertex_indices = {
1907 {face->vertex_index(0), face->vertex_index(1)}};
1908 if (face_vertex_indices[0] > face_vertex_indices[1])
1909 std::swap(face_vertex_indices[0], face_vertex_indices[1]);
1915 face_vertex_indices,
1917 const std::array<unsigned int, 2>
1918 &face_vertex_indices) ->
bool {
1919 return std::lexicographical_compare(
1922 face_vertex_indices.begin(),
1923 face_vertex_indices.end());
1927 (p->vertices[0] == face_vertex_indices[0]) &&
1928 (p->vertices[1] == face_vertex_indices[1]))
1930 face->set_boundary_id(p->boundary_id);
1938 std::vector<unsigned int> face_vertex_indices(
1939 face->n_vertices());
1940 for (
unsigned int v = 0; v < face->n_vertices(); ++v)
1941 face_vertex_indices[v] = face->vertex_index(v);
1942 std::sort(face_vertex_indices.begin(),
1943 face_vertex_indices.end());
1949 face_vertex_indices,
1951 const std::vector<unsigned int>
1952 &face_vertex_indices) ->
bool {
1953 return std::lexicographical_compare(
1956 face_vertex_indices.begin(),
1957 face_vertex_indices.end());
1961 (p->vertices == face_vertex_indices))
1963 face->set_boundary_id(p->boundary_id);
1968 for (
unsigned int e = 0; e < face->n_lines(); ++e)
1970 const auto edge = face->line(e);
1972 std::array<unsigned int, 2> edge_vertex_indices = {
1973 {edge->vertex_index(0), edge->vertex_index(1)}};
1974 if (edge_vertex_indices[0] > edge_vertex_indices[1])
1975 std::swap(edge_vertex_indices[0],
1976 edge_vertex_indices[1]);
1982 edge_vertex_indices,
1984 const std::array<unsigned int, 2>
1985 &edge_vertex_indices) ->
bool {
1986 return std::lexicographical_compare(
1989 edge_vertex_indices.begin(),
1990 edge_vertex_indices.end());
1994 (p->vertices[0] == edge_vertex_indices[0]) &&
1995 (p->vertices[1] == edge_vertex_indices[1]))
1997 edge->set_boundary_id(p->boundary_id);
2006template <
int dim,
int spacedim>
2010 Assert(
tria !=
nullptr, ExcNoTriangulationSelected());
2013 unsigned int n_vertices;
2014 unsigned int n_cells;
2020 std::array<std::map<int, int>, 4> tag_maps;
2025 unsigned int gmsh_file_format = 0;
2026 if (line ==
"@f$NOD")
2027 gmsh_file_format = 10;
2028 else if (line ==
"@f$MeshFormat")
2029 gmsh_file_format = 20;
2035 if (gmsh_file_format == 20)
2038 unsigned int file_type, data_size;
2040 in >> version >> file_type >> data_size;
2043 gmsh_file_format =
static_cast<unsigned int>(version * 10);
2051 AssertThrow(line ==
"@f$EndMeshFormat", ExcInvalidGMSHInput(line));
2055 if (line ==
"@f$PhysicalNames")
2061 while (line !=
"@f$EndPhysicalNames");
2066 if (line ==
"@f$Entities")
2068 unsigned long n_points, n_curves, n_surfaces, n_volumes;
2070 in >> n_points >> n_curves >> n_surfaces >> n_volumes;
2071 for (
unsigned int i = 0; i < n_points; ++i)
2075 unsigned int n_physicals;
2076 double box_min_x, box_min_y, box_min_z, box_max_x, box_max_y,
2080 if (gmsh_file_format > 40)
2082 in >> tag >> box_min_x >> box_min_y >> box_min_z >>
2084 box_max_x = box_min_x;
2085 box_max_y = box_min_y;
2086 box_max_z = box_min_z;
2090 in >> tag >> box_min_x >> box_min_y >> box_min_z >>
2091 box_max_x >> box_max_y >> box_max_z >> n_physicals;
2096 ExcMessage(
"More than one tag is not supported!"));
2098 int physical_tag = 0;
2099 for (
unsigned int j = 0; j < n_physicals; ++j)
2101 tag_maps[0][tag] = physical_tag;
2103 for (
unsigned int i = 0; i < n_curves; ++i)
2107 unsigned int n_physicals;
2108 double box_min_x, box_min_y, box_min_z, box_max_x, box_max_y,
2112 in >> tag >> box_min_x >> box_min_y >> box_min_z >> box_max_x >>
2113 box_max_y >> box_max_z >> n_physicals;
2117 ExcMessage(
"More than one tag is not supported!"));
2119 int physical_tag = 0;
2120 for (
unsigned int j = 0; j < n_physicals; ++j)
2122 tag_maps[1][tag] = physical_tag;
2127 for (
unsigned int j = 0; j < n_points; ++j)
2131 for (
unsigned int i = 0; i < n_surfaces; ++i)
2135 unsigned int n_physicals;
2136 double box_min_x, box_min_y, box_min_z, box_max_x, box_max_y,
2140 in >> tag >> box_min_x >> box_min_y >> box_min_z >> box_max_x >>
2141 box_max_y >> box_max_z >> n_physicals;
2145 ExcMessage(
"More than one tag is not supported!"));
2147 int physical_tag = 0;
2148 for (
unsigned int j = 0; j < n_physicals; ++j)
2150 tag_maps[2][tag] = physical_tag;
2155 for (
unsigned int j = 0; j < n_curves; ++j)
2158 for (
unsigned int i = 0; i < n_volumes; ++i)
2162 unsigned int n_physicals;
2163 double box_min_x, box_min_y, box_min_z, box_max_x, box_max_y,
2167 in >> tag >> box_min_x >> box_min_y >> box_min_z >> box_max_x >>
2168 box_max_y >> box_max_z >> n_physicals;
2172 ExcMessage(
"More than one tag is not supported!"));
2174 int physical_tag = 0;
2175 for (
unsigned int j = 0; j < n_physicals; ++j)
2177 tag_maps[3][tag] = physical_tag;
2182 for (
unsigned int j = 0; j < n_surfaces; ++j)
2186 AssertThrow(line ==
"@f$EndEntities", ExcInvalidGMSHInput(line));
2191 if (line ==
"@f$PartitionedEntities")
2197 while (line !=
"@f$EndPartitionedEntities");
2204 AssertThrow(line ==
"@f$Nodes", ExcInvalidGMSHInput(line));
2208 int n_entity_blocks = 1;
2209 if (gmsh_file_format > 40)
2213 in >> n_entity_blocks >> n_vertices >> min_node_tag >> max_node_tag;
2215 else if (gmsh_file_format == 40)
2217 in >> n_entity_blocks >> n_vertices;
2221 std::vector<Point<spacedim>>
vertices(n_vertices);
2228 unsigned int global_vertex = 0;
2229 for (
int entity_block = 0; entity_block < n_entity_blocks; ++entity_block)
2232 unsigned long numNodes;
2234 if (gmsh_file_format < 40)
2236 numNodes = n_vertices;
2243 int tagEntity, dimEntity;
2244 in >> tagEntity >> dimEntity >> parametric >> numNodes;
2247 std::vector<int> vertex_numbers;
2249 if (gmsh_file_format > 40)
2250 for (
unsigned long vertex_per_entity = 0;
2251 vertex_per_entity < numNodes;
2252 ++vertex_per_entity)
2254 in >> vertex_number;
2255 vertex_numbers.push_back(vertex_number);
2258 for (
unsigned long vertex_per_entity = 0; vertex_per_entity < numNodes;
2259 ++vertex_per_entity, ++global_vertex)
2265 if (gmsh_file_format > 40)
2267 vertex_number = vertex_numbers[vertex_per_entity];
2268 in >> x[0] >> x[1] >> x[2];
2271 in >> vertex_number >> x[0] >> x[1] >> x[2];
2273 for (
unsigned int d = 0; d < spacedim; ++d)
2279 if (parametric != 0)
2294 static const std::string end_nodes_marker[] = {
"@f$ENDNOD",
"@f$EndNodes"};
2295 AssertThrow(line == end_nodes_marker[gmsh_file_format == 10 ? 0 : 1],
2296 ExcInvalidGMSHInput(line));
2300 static const std::string begin_elements_marker[] = {
"@f$ELM",
"@f$Elements"};
2301 AssertThrow(line == begin_elements_marker[gmsh_file_format == 10 ? 0 : 1],
2302 ExcInvalidGMSHInput(line));
2305 if (gmsh_file_format > 40)
2309 in >> n_entity_blocks >> n_cells >> min_node_tag >> max_node_tag;
2311 else if (gmsh_file_format == 40)
2313 in >> n_entity_blocks >> n_cells;
2317 n_entity_blocks = 1;
2324 std::vector<CellData<dim>> cells;
2326 std::map<unsigned int, types::boundary_id> boundary_ids_1d;
2332 std::map<unsigned int, unsigned int> vertex_counts;
2335 static constexpr std::array<unsigned int, 8> local_vertex_numbering = {
2336 {0, 1, 5, 4, 2, 3, 7, 6}};
2337 unsigned int global_cell = 0;
2338 for (
int entity_block = 0; entity_block < n_entity_blocks; ++entity_block)
2340 unsigned int material_id;
2341 unsigned long numElements;
2344 if (gmsh_file_format < 40)
2348 numElements = n_cells;
2350 else if (gmsh_file_format == 40)
2352 int tagEntity, dimEntity;
2353 in >> tagEntity >> dimEntity >> cell_type >> numElements;
2354 material_id = tag_maps[dimEntity][tagEntity];
2359 int tagEntity, dimEntity;
2360 in >> dimEntity >> tagEntity >> cell_type >> numElements;
2361 material_id = tag_maps[dimEntity][tagEntity];
2364 for (
unsigned int cell_per_entity = 0; cell_per_entity < numElements;
2365 ++cell_per_entity, ++global_cell)
2374 unsigned int nod_num;
2395 unsigned int elm_number = 0;
2396 if (gmsh_file_format < 40)
2402 if (gmsh_file_format < 20)
2408 else if (gmsh_file_format < 40)
2413 unsigned int n_tags;
2420 for (
unsigned int i = 1; i < n_tags; ++i)
2425 else if (cell_type == 2)
2427 else if (cell_type == 3)
2429 else if (cell_type == 4)
2431 else if (cell_type == 5)
2442 else if (cell_type == 2)
2444 else if (cell_type == 3)
2446 else if (cell_type == 4)
2448 else if (cell_type == 5)
2474 if (((cell_type == 1) && (dim == 1)) ||
2475 ((cell_type == 2) && (dim == 2)) ||
2476 ((cell_type == 3) && (dim == 2)) ||
2477 ((cell_type == 4) && (dim == 3)) ||
2478 ((cell_type == 5) && (dim == 3)))
2481 unsigned int vertices_per_cell = 0;
2483 vertices_per_cell = 2;
2484 else if (cell_type == 2)
2485 vertices_per_cell = 3;
2486 else if (cell_type == 3)
2487 vertices_per_cell = 4;
2488 else if (cell_type == 4)
2489 vertices_per_cell = 4;
2490 else if (cell_type == 5)
2491 vertices_per_cell = 8;
2495 "Number of nodes does not coincide with the "
2496 "number required for this object"));
2499 cells.emplace_back();
2501 cell.
vertices.resize(vertices_per_cell);
2502 for (
unsigned int i = 0; i < vertices_per_cell; ++i)
2505 if (vertices_per_cell ==
2508 local_vertex_numbering[i] :
2516 std::numeric_limits<types::material_id>::max(),
2520 std::numeric_limits<types::material_id>::max()));
2528 for (
unsigned int i = 0; i < vertices_per_cell; ++i)
2532 ExcInvalidVertexIndexGmsh(cell_per_entity,
2538 vertex_counts[vertex] += 1u;
2542 else if ((cell_type == 1) &&
2543 ((dim == 2) || (dim == 3)))
2552 std::numeric_limits<types::boundary_id>::max(),
2556 std::numeric_limits<types::boundary_id>::max()));
2567 for (
unsigned int &vertex :
2576 ExcInvalidVertexIndex(cell_per_entity,
2581 else if ((cell_type == 2 || cell_type == 3) &&
2585 unsigned int vertices_per_cell = 0;
2588 vertices_per_cell = 3;
2589 else if (cell_type == 3)
2590 vertices_per_cell = 4;
2598 for (
unsigned int i = 0; i < vertices_per_cell; ++i)
2603 std::numeric_limits<types::boundary_id>::max(),
2607 std::numeric_limits<types::boundary_id>::max()));
2618 for (
unsigned int &vertex :
2627 ExcInvalidVertexIndex(cell_per_entity, vertex));
2631 else if (cell_type == 15)
2634 unsigned int node_index = 0;
2635 if (gmsh_file_format < 20)
2655 AssertThrow(
false, ExcGmshUnsupportedGeometry(cell_type));
2663 static const std::string end_elements_marker[] = {
"@f$ENDELM",
"@f$EndElements"};
2664 AssertThrow(line == end_elements_marker[gmsh_file_format == 10 ? 0 : 1],
2665 ExcInvalidGMSHInput(line));
2678 for (
const auto &pair : vertex_counts)
2679 if (pair.second == 1u)
2680 boundary_id_pairs.emplace_back(
vertices[pair.first],
2681 boundary_ids_1d[pair.first]);
2683 apply_grid_fixup_functions(
vertices, cells, subcelldata);
2689 assign_1d_boundary_ids(boundary_id_pairs, *
tria);
2694#ifdef DEAL_II_GMSH_WITH_API
2695template <
int dim,
int spacedim>
2699 Assert(
tria !=
nullptr, ExcNoTriangulationSelected());
2701 const std::map<int, std::uint8_t> gmsh_to_dealii_type = {
2702 {{15, 0}, {1, 1}, {2, 2}, {3, 3}, {4, 4}, {7, 5}, {6, 6}, {5, 7}}};
2705 const std::array<std::vector<unsigned int>, 8> gmsh_to_dealii = {
2712 {{0, 1, 2, 3, 4, 5}},
2713 {{0, 1, 3, 2, 4, 5, 7, 6}}}};
2715 std::vector<Point<spacedim>>
vertices;
2716 std::vector<CellData<dim>> cells;
2718 std::map<unsigned int, types::boundary_id> boundary_ids_1d;
2724 std::map<unsigned int, unsigned int> vertex_counts;
2727 gmsh::option::setNumber(
"General.Verbosity", 0);
2731 ExcMessage(
"You are trying to read a gmsh file with dimension " +
2732 std::to_string(gmsh::model::getDimension()) +
2733 " into a grid of dimension " + std::to_string(dim)));
2738 gmsh::model::mesh::removeDuplicateNodes();
2739 gmsh::model::mesh::renumberNodes();
2740 std::vector<std::size_t> node_tags;
2741 std::vector<double> coord;
2742 std::vector<double> parametricCoord;
2743 gmsh::model::mesh::getNodes(
2744 node_tags, coord, parametricCoord, -1, -1,
false,
false);
2746 for (
unsigned int i = 0; i < node_tags.size(); ++i)
2750 for (
unsigned int d = 0; d < spacedim; ++d)
2754 for (
unsigned int d = spacedim; d < 3; ++d)
2756 ExcMessage(
"The grid you are reading contains nodes that are "
2757 "nonzero in the coordinate with index " +
2759 ", but you are trying to save "
2760 "it on a grid embedded in a " +
2761 std::to_string(spacedim) +
" dimensional space."));
2768 std::vector<std::pair<int, int>> entities;
2769 gmsh::model::getEntities(entities);
2771 for (
const auto &e : entities)
2774 const int &entity_dim = e.first;
2775 const int &entity_tag = e.second;
2781 std::vector<int> physical_tags;
2782 gmsh::model::getPhysicalGroupsForEntity(entity_dim,
2787 if (physical_tags.size())
2788 for (
auto physical_tag : physical_tags)
2791 gmsh::model::getPhysicalName(entity_dim, physical_tag, name);
2795 std::map<std::string, int> id_names;
2797 bool throw_anyway =
false;
2798 bool found_boundary_id =
false;
2801 for (
const auto &it : id_names)
2803 const auto &name = it.first;
2804 const auto &
id = it.second;
2805 if (entity_dim == dim && name ==
"MaterialID")
2808 found_boundary_id =
true;
2810 else if (entity_dim < dim && name ==
"BoundaryID")
2813 found_boundary_id =
true;
2815 else if (name ==
"ManifoldID")
2821 throw_anyway =
true;
2827 if (throw_anyway && !found_boundary_id)
2836 boundary_id = physical_tag;
2841 std::vector<int> element_types;
2842 std::vector<std::vector<std::size_t>> element_ids, element_nodes;
2843 gmsh::model::mesh::getElements(
2844 element_types, element_ids, element_nodes, entity_dim, entity_tag);
2846 for (
unsigned int i = 0; i < element_types.size(); ++i)
2848 const auto &type = gmsh_to_dealii_type.at(element_types[i]);
2849 const auto n_vertices = gmsh_to_dealii[type].size();
2850 const auto &elements = element_ids[i];
2851 const auto &nodes = element_nodes[i];
2852 for (
unsigned int j = 0; j < elements.size(); ++j)
2854 if (entity_dim == dim)
2856 cells.emplace_back(n_vertices);
2857 auto &cell = cells.back();
2858 for (
unsigned int v = 0; v < n_vertices; ++v)
2861 nodes[n_vertices * j + gmsh_to_dealii[type][v]] - 1;
2863 vertex_counts[cell.vertices[v]] += 1u;
2865 cell.manifold_id = manifold_id;
2866 cell.material_id = boundary_id;
2868 else if (entity_dim == 2)
2872 for (
unsigned int v = 0; v < n_vertices; ++v)
2874 nodes[n_vertices * j + gmsh_to_dealii[type][v]] - 1;
2876 face.manifold_id = manifold_id;
2877 face.boundary_id = boundary_id;
2879 else if (entity_dim == 1)
2883 for (
unsigned int v = 0; v < n_vertices; ++v)
2885 nodes[n_vertices * j + gmsh_to_dealii[type][v]] - 1;
2887 line.manifold_id = manifold_id;
2888 line.boundary_id = boundary_id;
2890 else if (entity_dim == 0)
2894 for (
unsigned int j = 0; j < elements.size(); ++j)
2895 boundary_ids_1d[nodes[j] - 1] = boundary_id;
2906 for (
const auto &pair : vertex_counts)
2907 if (pair.second == 1u)
2908 boundary_id_pairs.emplace_back(
vertices[pair.first],
2909 boundary_ids_1d[pair.first]);
2911 apply_grid_fixup_functions(
vertices, cells, subcelldata);
2917 assign_1d_boundary_ids(boundary_id_pairs, *
tria);
2926template <
int dim,
int spacedim>
2929 std::string & header,
2930 std::vector<unsigned int> &tecplot2deal,
2931 unsigned int & n_vars,
2932 unsigned int & n_vertices,
2933 unsigned int & n_cells,
2934 std::vector<unsigned int> &IJK,
2959 std::transform(header.begin(), header.end(), header.begin(), ::toupper);
2963 std::replace(header.begin(), header.end(),
'\t',
' ');
2964 std::replace(header.begin(), header.end(),
',',
' ');
2965 std::replace(header.begin(), header.end(),
'\n',
' ');
2969 std::string::size_type pos = header.find(
'=');
2971 while (pos !=
static_cast<std::string::size_type
>(std::string::npos))
2972 if (header[pos + 1] ==
' ')
2973 header.erase(pos + 1, 1);
2974 else if (header[pos - 1] ==
' ')
2976 header.erase(pos - 1, 1);
2980 pos = header.find(
'=', ++pos);
2983 std::vector<std::string> entries =
2987 for (
unsigned int i = 0; i < entries.size(); ++i)
2996 tecplot2deal[0] = 0;
2999 while (entries[i][0] ==
'"')
3001 if (entries[i] ==
"\"X\"")
3002 tecplot2deal[0] = n_vars;
3003 else if (entries[i] ==
"\"Y\"")
3009 tecplot2deal[1] = n_vars;
3011 else if (entries[i] ==
"\"Z\"")
3017 tecplot2deal[2] = n_vars;
3029 "Tecplot file must contain at least one variable for each dimension"));
3030 for (
unsigned int d = 1; d < dim; ++d)
3032 tecplot2deal[d] > 0,
3034 "Tecplot file must contain at least one variable for each dimension."));
3039 "ZONETYPE=FELINESEG") &&
3043 "ZONETYPE=FEQUADRILATERAL") &&
3047 "ZONETYPE=FEBRICK") &&
3055 "The tecplot file contains an unsupported ZONETYPE."));
3058 "DATAPACKING=POINT"))
3061 "DATAPACKING=BLOCK"))
3084 "ET=QUADRILATERAL") &&
3096 "The tecplot file contains an unsupported ElementType."));
3104 dim > 1 || IJK[1] == 1,
3106 "Parameter 'J=' found in tecplot, although this is only possible for dimensions greater than 1."));
3112 dim > 2 || IJK[2] == 1,
3114 "Parameter 'K=' found in tecplot, although this is only possible for dimensions greater than 2."));
3129 for (
unsigned int d = 0; d < dim; ++d)
3134 "Tecplot file does not contain a complete and consistent set of parameters"));
3135 n_vertices *= IJK[d];
3136 n_cells *= (IJK[d] - 1);
3144 "Tecplot file does not contain a complete and consistent set of parameters"));
3150 n_cells = *std::max_element(IJK.begin(), IJK.end());
3154 "Tecplot file does not contain a complete and consistent set of parameters"));
3164 const unsigned int dim = 2;
3165 const unsigned int spacedim = 2;
3166 Assert(
tria !=
nullptr, ExcNoTriangulationSelected());
3170 skip_comment_lines(in,
'#');
3173 std::string line, header;
3180 std::string letters =
"abcdfghijklmnopqrstuvwxyzABCDFGHIJKLMNOPQRSTUVWXYZ";
3183 while (line.find_first_of(letters) != std::string::npos)
3185 header +=
" " + line;
3192 std::vector<unsigned int> tecplot2deal(dim);
3193 std::vector<unsigned int> IJK(dim);
3194 unsigned int n_vars, n_vertices, n_cells;
3195 bool structured, blocked;
3197 parse_tecplot_header(header,
3212 std::vector<Point<spacedim>>
vertices(n_vertices + 1);
3215 std::vector<CellData<dim>> cells(n_cells);
3231 unsigned int next_index = 0;
3235 if (tecplot2deal[0] == 0)
3239 std::vector<std::string> first_var =
3242 for (
unsigned int i = 1; i < first_var.size() + 1; ++i)
3243 vertices[i](0) = std::strtod(first_var[i - 1].c_str(), &endptr);
3248 for (
unsigned int j = first_var.size() + 1; j < n_vertices + 1; ++j)
3256 for (
unsigned int i = 1; i < n_vars; ++i)
3265 if (next_index == dim && structured)
3268 if ((next_index < dim) && (i == tecplot2deal[next_index]))
3271 for (
unsigned int j = 1; j < n_vertices + 1; ++j)
3279 for (
unsigned int j = 1; j < n_vertices + 1; ++j)
3291 std::vector<double> vars(n_vars);
3296 std::vector<std::string> first_vertex =
3299 for (
unsigned int d = 0; d < dim; ++d)
3301 std::strtod(first_vertex[tecplot2deal[d]].c_str(), &endptr);
3305 for (
unsigned int v = 2; v < n_vertices + 1; ++v)
3307 for (
unsigned int i = 0; i < n_vars; ++i)
3313 for (
unsigned int i = 0; i < dim; ++i)
3314 vertices[v](i) = vars[tecplot2deal[i]];
3322 unsigned int I = IJK[0], J = IJK[1];
3324 unsigned int cell = 0;
3326 for (
unsigned int j = 0; j < J - 1; ++j)
3327 for (
unsigned int i = 1; i < I; ++i)
3329 cells[cell].vertices[0] = i + j * I;
3330 cells[cell].vertices[1] = i + 1 + j * I;
3331 cells[cell].vertices[2] = i + (j + 1) * I;
3332 cells[cell].vertices[3] = i + 1 + (j + 1) * I;
3336 std::vector<unsigned int> boundary_vertices(2 * I + 2 * J - 4);
3338 for (
unsigned int i = 1; i < I + 1; ++i)
3340 boundary_vertices[k] = i;
3342 boundary_vertices[k] = i + (J - 1) * I;
3345 for (
unsigned int j = 1; j < J - 1; ++j)
3347 boundary_vertices[k] = 1 + j * I;
3349 boundary_vertices[k] = I + j * I;
3368 for (
unsigned int i = 0; i < n_cells; ++i)
3385 apply_grid_fixup_functions(
vertices, cells, subcelldata);
3391template <
int dim,
int spacedim>
3400template <
int dim,
int spacedim>
3403 const unsigned int mesh_index,
3404 const bool remove_duplicates,
3406 const bool ignore_unsupported_types)
3408#ifdef DEAL_II_WITH_ASSIMP
3413 Assimp::Importer importer;
3416 const aiScene *scene =
3417 importer.ReadFile(filename.c_str(),
3418 aiProcess_RemoveComponent |
3419 aiProcess_JoinIdenticalVertices |
3420 aiProcess_ImproveCacheLocality | aiProcess_SortByPType |
3421 aiProcess_OptimizeGraph | aiProcess_OptimizeMeshes);
3427 ExcMessage(
"Input file contains no meshes."));
3430 (mesh_index < scene->mNumMeshes),
3433 unsigned int start_mesh =
3435 unsigned int end_mesh =
3440 std::vector<Point<spacedim>>
vertices;
3441 std::vector<CellData<dim>> cells;
3445 unsigned int v_offset = 0;
3446 unsigned int c_offset = 0;
3448 static constexpr std::array<unsigned int, 8> local_vertex_numbering = {
3449 {0, 1, 5, 4, 2, 3, 7, 6}};
3451 for (
unsigned int m = start_mesh; m < end_mesh; ++m)
3453 const aiMesh *mesh = scene->mMeshes[m];
3457 if ((dim == 2) && mesh->mPrimitiveTypes != aiPrimitiveType_POLYGON)
3460 ExcMessage(
"Incompatible mesh " + std::to_string(m) +
3461 "/" + std::to_string(scene->mNumMeshes)));
3464 else if ((dim == 1) && mesh->mPrimitiveTypes != aiPrimitiveType_LINE)
3467 ExcMessage(
"Incompatible mesh " + std::to_string(m) +
3468 "/" + std::to_string(scene->mNumMeshes)));
3472 const unsigned int n_vertices = mesh->mNumVertices;
3473 const aiVector3D * mVertices = mesh->mVertices;
3476 const unsigned int n_faces = mesh->mNumFaces;
3477 const aiFace * mFaces = mesh->mFaces;
3479 vertices.resize(v_offset + n_vertices);
3480 cells.resize(c_offset + n_faces);
3482 for (
unsigned int i = 0; i < n_vertices; ++i)
3483 for (
unsigned int d = 0; d < spacedim; ++d)
3484 vertices[i + v_offset][d] = mVertices[i][d];
3486 unsigned int valid_cell = c_offset;
3487 for (
unsigned int i = 0; i < n_faces; ++i)
3494 .vertices[dim == 3 ? local_vertex_numbering[f] :
3496 mFaces[i].mIndices[f] + v_offset;
3498 cells[valid_cell].material_id = m;
3504 ExcMessage(
"Face " + std::to_string(i) +
" of mesh " +
3505 std::to_string(m) +
" has " +
3506 std::to_string(mFaces[i].mNumIndices) +
3507 " vertices. We expected only " +
3512 cells.resize(valid_cell);
3517 v_offset += n_vertices;
3518 c_offset = valid_cell;
3522 if (cells.size() == 0)
3525 if (remove_duplicates)
3530 unsigned int n_verts = 0;
3534 std::vector<unsigned int> considered_vertices;
3536 vertices, cells, subcelldata, considered_vertices, tol);
3540 apply_grid_fixup_functions(
vertices, cells, subcelldata);
3546 (void)remove_duplicates;
3548 (void)ignore_unsupported_types;
3553#ifdef DEAL_II_TRILINOS_WITH_SEACAS
3560 exodusii_name_to_type(
const std::string &type_name,
3561 const int n_nodes_per_element)
3567 std::string type_name_2 = type_name;
3568 std::transform(type_name_2.begin(),
3570 type_name_2.begin(),
3571 [](
unsigned char c) { return std::toupper(c); });
3572 const std::string
numbers =
"0123456789";
3573 type_name_2.erase(std::find_first_of(type_name_2.begin(),
3580 if (type_name_2 ==
"BAR" || type_name_2 ==
"BEAM" ||
3581 type_name_2 ==
"EDGE" || type_name_2 ==
"TRUSS")
3583 else if (type_name_2 ==
"TRI" || type_name_2 ==
"TRIANGLE")
3585 else if (type_name_2 ==
"QUAD" || type_name_2 ==
"QUADRILATERAL")
3587 else if (type_name_2 ==
"SHELL")
3589 if (n_nodes_per_element == 3)
3594 else if (type_name_2 ==
"TET" || type_name_2 ==
"TETRA" ||
3595 type_name_2 ==
"TETRAHEDRON")
3597 else if (type_name_2 ==
"PYRA" || type_name_2 ==
"PYRAMID")
3599 else if (type_name_2 ==
"WEDGE")
3601 else if (type_name_2 ==
"HEX" || type_name_2 ==
"HEXAHEDRON")
3611 template <
int dim,
int spacedim = dim>
3612 std::pair<SubCellData, std::vector<std::vector<int>>>
3613 read_exodusii_sidesets(
const int ex_id,
3614 const int n_side_sets,
3616 const bool apply_all_indicators_to_manifolds)
3619 std::vector<std::vector<int>> b_or_m_id_to_sideset_ids;
3621 b_or_m_id_to_sideset_ids.emplace_back();
3627 if (dim == spacedim && n_side_sets > 0)
3629 std::vector<int> side_set_ids(n_side_sets);
3630 int ierr = ex_get_ids(ex_id, EX_SIDE_SET, side_set_ids.data());
3631 AssertThrowExodusII(ierr);
3639 std::map<std::size_t, std::vector<int>> face_side_sets;
3640 for (
const int side_set_id : side_set_ids)
3643 int n_distribution_factors = -1;
3645 ierr = ex_get_set_param(ex_id,
3649 &n_distribution_factors);
3650 AssertThrowExodusII(ierr);
3653 std::vector<int> elements(n_sides);
3654 std::vector<int> faces(n_sides);
3655 ierr = ex_get_set(ex_id,
3660 AssertThrowExodusII(ierr);
3668 for (
int side_n = 0; side_n < n_sides; ++side_n)
3670 const long element_n = elements[side_n] - 1;
3671 const long face_n = faces[side_n] - 1;
3672 const std::size_t face_id =
3673 element_n * max_faces_per_cell + face_n;
3674 face_side_sets[face_id].push_back(side_set_id);
3680 std::vector<std::pair<std::size_t, std::vector<int>>>
3681 face_id_to_side_sets;
3682 face_id_to_side_sets.reserve(face_side_sets.size());
3683 for (
auto &pair : face_side_sets)
3686 face_id_to_side_sets.push_back(std::move(pair));
3690 std::sort(face_id_to_side_sets.begin(),
3691 face_id_to_side_sets.end(),
3692 [](
const auto &a,
const auto &b) {
3693 return std::lexicographical_compare(a.second.begin(),
3704 for (
const auto &pair : face_id_to_side_sets)
3706 const std::size_t face_id = pair.first;
3707 const std::vector<int> &face_sideset_ids = pair.second;
3708 if (face_sideset_ids != b_or_m_id_to_sideset_ids.back())
3713 ++current_b_or_m_id;
3714 b_or_m_id_to_sideset_ids.push_back(face_sideset_ids);
3715 Assert(current_b_or_m_id == b_or_m_id_to_sideset_ids.size() - 1,
3719 const unsigned int local_face_n = face_id % max_faces_per_cell;
3720 const CellData<dim> &cell = cells[face_id / max_faces_per_cell];
3723 const unsigned int deal_face_n =
3734 if (apply_all_indicators_to_manifolds)
3735 boundary_line.manifold_id = current_b_or_m_id;
3737 boundary_line.boundary_id = current_b_or_m_id;
3738 for (
unsigned int j = 0; j < face_reference_cell.
n_vertices();
3740 boundary_line.vertices[j] =
3742 deal_face_n, j, 0)];
3749 if (apply_all_indicators_to_manifolds)
3750 boundary_quad.manifold_id = current_b_or_m_id;
3752 boundary_quad.boundary_id = current_b_or_m_id;
3753 for (
unsigned int j = 0; j < face_reference_cell.
n_vertices();
3755 boundary_quad.vertices[j] =
3757 deal_face_n, j, 0)];
3764 return std::make_pair(std::move(subcelldata),
3765 std::move(b_or_m_id_to_sideset_ids));
3770template <
int dim,
int spacedim>
3773 const std::string &filename,
3774 const bool apply_all_indicators_to_manifolds)
3776#ifdef DEAL_II_TRILINOS_WITH_SEACAS
3778 int component_word_size =
sizeof(double);
3780 int floating_point_word_size = 0;
3781 float ex_version = 0.0;
3783 const int ex_id = ex_open(filename.c_str(),
3785 &component_word_size,
3786 &floating_point_word_size,
3789 ExcMessage(
"ExodusII failed to open the specified input file."));
3792 std::vector<char> string_temp(MAX_LINE_LENGTH + 1,
'\0');
3793 int mesh_dimension = 0;
3796 int n_element_blocks = 0;
3797 int n_node_sets = 0;
3798 int n_side_sets = 0;
3800 int ierr = ex_get_init(ex_id,
3808 AssertThrowExodusII(ierr);
3819 std::vector<Point<spacedim>>
vertices;
3822 std::vector<double> xs(n_nodes);
3823 std::vector<double> ys(n_nodes);
3824 std::vector<double> zs(n_nodes);
3826 ierr = ex_get_coord(ex_id, xs.data(), ys.data(), zs.data());
3827 AssertThrowExodusII(ierr);
3829 for (
int vertex_n = 0; vertex_n < n_nodes; ++vertex_n)
3834 vertices.emplace_back(xs[vertex_n]);
3837 vertices.emplace_back(xs[vertex_n], ys[vertex_n]);
3840 vertices.emplace_back(xs[vertex_n], ys[vertex_n], zs[vertex_n]);
3848 std::vector<int> element_block_ids(n_element_blocks);
3849 ierr = ex_get_ids(ex_id, EX_ELEM_BLOCK, element_block_ids.data());
3850 AssertThrowExodusII(ierr);
3852 std::vector<CellData<dim>> cells;
3853 cells.reserve(n_elements);
3857 for (
const int element_block_id : element_block_ids)
3859 std::fill(string_temp.begin(), string_temp.end(),
'\0');
3860 int n_block_elements = 0;
3861 int n_nodes_per_element = 0;
3862 int n_edges_per_element = 0;
3863 int n_faces_per_element = 0;
3864 int n_attributes_per_element = 0;
3867 ierr = ex_get_block(ex_id,
3872 &n_nodes_per_element,
3873 &n_edges_per_element,
3874 &n_faces_per_element,
3875 &n_attributes_per_element);
3876 AssertThrowExodusII(ierr);
3878 exodusii_name_to_type(string_temp.data(), n_nodes_per_element);
3881 "The ExodusII block " + std::to_string(element_block_id) +
3882 " with element type " + std::string(string_temp.data()) +
3884 ", which does not match the topological mesh dimension " +
3885 std::to_string(dim) +
"."));
3892 std::vector<int> connection(n_nodes_per_element * n_block_elements);
3893 ierr = ex_get_conn(ex_id,
3899 AssertThrowExodusII(ierr);
3901 for (
unsigned int elem_n = 0; elem_n < connection.size();
3902 elem_n += n_nodes_per_element)
3908 connection[elem_n + i] - 1;
3911 cells.push_back(std::move(cell));
3916 auto pair = read_exodusii_sidesets<dim, spacedim>(
3917 ex_id, n_side_sets, cells, apply_all_indicators_to_manifolds);
3918 ierr = ex_close(ex_id);
3919 AssertThrowExodusII(ierr);
3921 apply_grid_fixup_functions(
vertices, cells, pair.first);
3928 (void)apply_all_indicators_to_manifolds;
3935template <
int dim,
int spacedim>
3949 if (std::find_if(line.begin(), line.end(), [](
const char c) {
3954 for (
int i = line.size() - 1; i >= 0; --i)
3955 in.putback(line[i]);
3965template <
int dim,
int spacedim>
3968 const char comment_start)
3973 while (in.get(c) && c == comment_start)
3976 while (in.get() !=
'\n')
3986 skip_empty_lines(in);
3991template <
int dim,
int spacedim>
4009 double min_x =
vertices[cells[0].vertices[0]](0),
4011 min_y =
vertices[cells[0].vertices[0]](1),
4014 for (
unsigned int i = 0; i < cells.size(); ++i)
4016 for (
const auto vertex : cells[i].vertices)
4030 out <<
"# cell " << i << std::endl;
4032 for (
const auto vertex : cells[i].vertices)
4036 out <<
"set label \"" << i <<
"\" at " <<
center(0) <<
',' <<
center(1)
4037 <<
" center" << std::endl;
4040 for (
unsigned int f = 0; f < 2; ++f)
4042 <<
vertices[cells[i].vertices[f]](1) <<
" to "
4046 for (
unsigned int f = 2; f < 4; ++f)
4050 <<
vertices[cells[i].vertices[f]](1) << std::endl;
4056 <<
"set nokey" << std::endl
4057 <<
"pl [" << min_x <<
':' << max_x <<
"][" << min_y <<
':' << max_y
4058 <<
"] " << min_y << std::endl
4059 <<
"pause -1" << std::endl;
4070 for (
const auto &cell : cells)
4137template <
int dim,
int spacedim>
4145 if (format == Default)
4146 name = search.
find(filename);
4148 name = search.
find(filename, default_suffix(format));
4151 if (format == Default)
4153 const std::string::size_type slashpos = name.find_last_of(
'/');
4154 const std::string::size_type dotpos = name.find_last_of(
'.');
4155 if (dotpos < name.size() &&
4156 (dotpos > slashpos || slashpos == std::string::npos))
4158 std::string ext = name.substr(dotpos + 1);
4159 format = parse_format(ext);
4163 if (format == assimp)
4167 else if (format == exodusii)
4169 read_exodusii(name);
4173 std::ifstream in(name.c_str());
4179template <
int dim,
int spacedim>
4183 if (format == Default)
4184 format = default_format;
4226 ExcMessage(
"There is no read_assimp(istream &) function. "
4227 "Use the read_assimp(string &filename, ...) "
4228 "functions, instead."));
4233 ExcMessage(
"There is no read_exodusii(istream &) function. "
4234 "Use the read_exodusii(string &filename, ...) "
4235 "function, instead."));
4246template <
int dim,
int spacedim>
4275 return ".unknown_format";
4281template <
int dim,
int spacedim>
4285 if (format_name ==
"dbmesh")
4288 if (format_name ==
"exodusii")
4291 if (format_name ==
"msh")
4294 if (format_name ==
"unv")
4297 if (format_name ==
"vtk")
4300 if (format_name ==
"vtu")
4304 if (format_name ==
"inp")
4307 if (format_name ==
"ucd")
4310 if (format_name ==
"xda")
4313 if (format_name ==
"tecplot")
4316 if (format_name ==
"dat")
4319 if (format_name ==
"plt")
4338template <
int dim,
int spacedim>
4342 return "dbmesh|exodusii|msh|unv|vtk|vtu|ucd|abaqus|xda|tecplot|assimp";
4349 template <
int dim,
int spacedim>
4350 Abaqus_to_UCD<dim, spacedim>::Abaqus_to_UCD()
4363 from_string(T &t,
const std::string &s, std::ios_base &(*f)(std::ios_base &))
4365 std::istringstream iss(s);
4366 return !(iss >> f >> t).fail();
4373 extract_int(
const std::string &s)
4376 for (
const char c : s)
4378 if (isdigit(c) != 0)
4385 from_string(number, tmp, std::dec);
4391 template <
int dim,
int spacedim>
4393 Abaqus_to_UCD<dim, spacedim>::read_in_abaqus(std::istream &input_stream)
4402 while (std::getline(input_stream, line))
4405 std::transform(line.begin(), line.end(), line.begin(), ::toupper);
4407 if (line.compare(
"*HEADING") == 0 || line.compare(0, 2,
"**") == 0 ||
4408 line.compare(0, 5,
"*PART") == 0)
4411 while (std::getline(input_stream, line))
4417 else if (line.compare(0, 5,
"*NODE") == 0)
4426 while (std::getline(input_stream, line))
4431 std::vector<double> node(spacedim + 1);
4433 std::istringstream iss(line);
4435 for (
unsigned int i = 0; i < spacedim + 1; ++i)
4436 iss >> node[i] >> comma;
4438 node_list.push_back(node);
4441 else if (line.compare(0, 8,
"*ELEMENT") == 0)
4456 const std::string before_material =
"ELSET=EB";
4457 const std::size_t idx = line.find(before_material);
4458 if (idx != std::string::npos)
4460 from_string(material,
4461 line.substr(idx + before_material.size()),
4467 while (std::getline(input_stream, line))
4472 std::istringstream iss(line);
4478 const unsigned int n_data_per_cell =
4480 std::vector<double> cell(n_data_per_cell);
4481 for (
unsigned int i = 0; i < n_data_per_cell; ++i)
4482 iss >> cell[i] >> comma;
4485 cell[0] =
static_cast<double>(material);
4486 cell_list.push_back(cell);
4489 else if (line.compare(0, 8,
"*SURFACE") == 0)
4500 const std::string name_key =
"NAME=";
4501 const std::size_t name_idx_start =
4502 line.find(name_key) + name_key.size();
4503 std::size_t name_idx_end = line.find(
',', name_idx_start);
4504 if (name_idx_end == std::string::npos)
4506 name_idx_end = line.size();
4508 const int b_indicator = extract_int(
4509 line.substr(name_idx_start, name_idx_end - name_idx_start));
4516 while (std::getline(input_stream, line))
4522 std::transform(line.begin(),
4530 std::istringstream iss(line);
4537 std::vector<double> quad_node_list;
4538 const std::string elset_name = line.substr(0, line.find(
','));
4539 if (elsets_list.count(elset_name) != 0)
4543 iss >> stmp >> temp >> face_number;
4545 const std::vector<int> cells = elsets_list[elset_name];
4546 for (
const int cell : cells)
4550 get_global_node_numbers(el_idx, face_number);
4551 quad_node_list.insert(quad_node_list.begin(),
4554 face_list.push_back(quad_node_list);
4561 iss >> el_idx >> comma >> temp >> face_number;
4563 get_global_node_numbers(el_idx, face_number);
4564 quad_node_list.insert(quad_node_list.begin(), b_indicator);
4566 face_list.push_back(quad_node_list);
4570 else if (line.compare(0, 6,
"*ELSET") == 0)
4574 std::string elset_name;
4576 const std::string elset_key =
"*ELSET, ELSET=";
4577 const std::size_t idx = line.find(elset_key);
4578 if (idx != std::string::npos)
4580 const std::string comma =
",";
4581 const std::size_t first_comma = line.find(comma);
4582 const std::size_t second_comma =
4583 line.find(comma, first_comma + 1);
4584 const std::size_t elset_name_start =
4585 line.find(elset_key) + elset_key.size();
4586 elset_name = line.substr(elset_name_start,
4587 second_comma - elset_name_start);
4597 std::vector<int> elements;
4598 const std::size_t generate_idx = line.find(
"GENERATE");
4599 if (generate_idx != std::string::npos)
4602 std::getline(input_stream, line);
4603 std::istringstream iss(line);
4613 iss >> elid_start >> comma >> elid_end;
4617 "While reading an ABAQUS file, the reader "
4618 "expected a comma but found a <") +
4619 comma +
"> in the line <" + line +
">."));
4621 elid_start <= elid_end,
4624 "While reading an ABAQUS file, the reader encountered "
4625 "a GENERATE statement in which the upper bound <") +
4627 "> for the element numbers is not larger or equal "
4628 "than the lower bound <" +
4632 if (iss.rdbuf()->in_avail() != 0)
4633 iss >> comma >> elis_step;
4637 "While reading an ABAQUS file, the reader "
4638 "expected a comma but found a <") +
4639 comma +
"> in the line <" + line +
">."));
4641 for (
int i = elid_start; i <= elid_end; i += elis_step)
4642 elements.push_back(i);
4643 elsets_list[elset_name] = elements;
4645 std::getline(input_stream, line);
4650 while (std::getline(input_stream, line))
4655 std::istringstream iss(line);
4660 iss >> elid >> comma;
4665 "While reading an ABAQUS file, the reader "
4666 "expected a comma but found a <") +
4667 comma +
"> in the line <" + line +
">."));
4669 elements.push_back(elid);
4673 elsets_list[elset_name] = elements;
4678 else if (line.compare(0, 5,
"*NSET") == 0)
4681 while (std::getline(input_stream, line))
4687 else if (line.compare(0, 14,
"*SOLID SECTION") == 0)
4691 const std::string elset_key =
"ELSET=";
4692 const std::size_t elset_start =
4693 line.find(
"ELSET=") + elset_key.size();
4694 const std::size_t elset_end = line.find(
',', elset_start + 1);
4695 const std::string elset_name =
4696 line.substr(elset_start, elset_end - elset_start);
4701 const std::string material_key =
"MATERIAL=";
4702 const std::size_t last_equal =
4703 line.find(
"MATERIAL=") + material_key.size();
4704 const std::size_t material_id_start = line.find(
'-', last_equal);
4706 from_string(material_id,
4707 line.substr(material_id_start + 1),
4711 const std::vector<int> &elset_cells = elsets_list[elset_name];
4712 for (
const int elset_cell : elset_cells)
4714 const int cell_id = elset_cell - 1;
4722 template <
int dim,
int spacedim>
4724 Abaqus_to_UCD<dim, spacedim>::get_global_node_numbers(
4725 const int face_cell_no,
4726 const int face_cell_face_no)
const
4736 if (face_cell_face_no == 1)
4738 quad_node_list[0] = cell_list[face_cell_no - 1][1];
4739 quad_node_list[1] = cell_list[face_cell_no - 1][2];
4741 else if (face_cell_face_no == 2)
4743 quad_node_list[0] = cell_list[face_cell_no - 1][2];
4744 quad_node_list[1] = cell_list[face_cell_no - 1][3];
4746 else if (face_cell_face_no == 3)
4748 quad_node_list[0] = cell_list[face_cell_no - 1][3];
4749 quad_node_list[1] = cell_list[face_cell_no - 1][4];
4751 else if (face_cell_face_no == 4)
4753 quad_node_list[0] = cell_list[face_cell_no - 1][4];
4754 quad_node_list[1] = cell_list[face_cell_no - 1][1];
4764 if (face_cell_face_no == 1)
4766 quad_node_list[0] = cell_list[face_cell_no - 1][1];
4767 quad_node_list[1] = cell_list[face_cell_no - 1][4];
4768 quad_node_list[2] = cell_list[face_cell_no - 1][3];
4769 quad_node_list[3] = cell_list[face_cell_no - 1][2];
4771 else if (face_cell_face_no == 2)
4773 quad_node_list[0] = cell_list[face_cell_no - 1][5];
4774 quad_node_list[1] = cell_list[face_cell_no - 1][8];
4775 quad_node_list[2] = cell_list[face_cell_no - 1][7];
4776 quad_node_list[3] = cell_list[face_cell_no - 1][6];
4778 else if (face_cell_face_no == 3)
4780 quad_node_list[0] = cell_list[face_cell_no - 1][1];
4781 quad_node_list[1] = cell_list[face_cell_no - 1][2];
4782 quad_node_list[2] = cell_list[face_cell_no - 1][6];
4783 quad_node_list[3] = cell_list[face_cell_no - 1][5];
4785 else if (face_cell_face_no == 4)
4787 quad_node_list[0] = cell_list[face_cell_no - 1][2];
4788 quad_node_list[1] = cell_list[face_cell_no - 1][3];
4789 quad_node_list[2] = cell_list[face_cell_no - 1][7];
4790 quad_node_list[3] = cell_list[face_cell_no - 1][6];
4792 else if (face_cell_face_no == 5)
4794 quad_node_list[0] = cell_list[face_cell_no - 1][3];
4795 quad_node_list[1] = cell_list[face_cell_no - 1][4];
4796 quad_node_list[2] = cell_list[face_cell_no - 1][8];
4797 quad_node_list[3] = cell_list[face_cell_no - 1][7];
4799 else if (face_cell_face_no == 6)
4801 quad_node_list[0] = cell_list[face_cell_no - 1][1];
4802 quad_node_list[1] = cell_list[face_cell_no - 1][5];
4803 quad_node_list[2] = cell_list[face_cell_no - 1][8];
4804 quad_node_list[3] = cell_list[face_cell_no - 1][4];
4817 return quad_node_list;
4820 template <
int dim,
int spacedim>
4822 Abaqus_to_UCD<dim, spacedim>::write_out_avs_ucd(std::ostream &output)
const
4831 const boost::io::ios_base_all_saver formatting_saver(output);
4835 output <<
"# Abaqus to UCD mesh conversion" << std::endl;
4836 output <<
"# Mesh type: AVS UCD" << std::endl;
4865 output << node_list.size() <<
"\t" << (cell_list.size() + face_list.size())
4866 <<
"\t0\t0\t0" << std::endl;
4869 output.precision(8);
4873 for (
const auto &node : node_list)
4876 output << node[0] <<
"\t";
4879 output.setf(std::ios::scientific, std::ios::floatfield);
4880 for (
unsigned int jj = 1; jj < spacedim + 1; ++jj)
4883 if (
std::abs(node[jj]) > tolerance)
4884 output << static_cast<double>(node[jj]) <<
"\t";
4886 output << 0.0 <<
"\t";
4889 output << 0.0 <<
"\t";
4891 output << std::endl;
4892 output.unsetf(std::ios::floatfield);
4896 for (
unsigned int ii = 0; ii < cell_list.size(); ++ii)
4898 output << ii + 1 <<
"\t" << cell_list[ii][0] <<
"\t"
4899 << (dim == 2 ?
"quad" :
"hex") <<
"\t";
4900 for (
unsigned int jj = 1; jj < GeometryInfo<dim>::vertices_per_cell + 1;
4902 output << cell_list[ii][jj] <<
"\t";
4904 output << std::endl;
4908 for (
unsigned int ii = 0; ii < face_list.size(); ++ii)
4910 output << ii + 1 <<
"\t" << face_list[ii][0] <<
"\t"
4911 << (dim == 2 ?
"line" :
"quad") <<
"\t";
4912 for (
unsigned int jj = 1; jj < GeometryInfo<dim>::vertices_per_face + 1;
4914 output << face_list[ii][jj] <<
"\t";
4916 output << std::endl;
4923#include "grid_in.inst"
void read_vtk(std::istream &in)
static void skip_empty_lines(std::istream &in)
void read_assimp(const std::string &filename, const unsigned int mesh_index=numbers::invalid_unsigned_int, const bool remove_duplicates=true, const double tol=1e-12, const bool ignore_unsupported_element_types=true)
void read_abaqus(std::istream &in, const bool apply_all_indicators_to_manifolds=false)
static std::string default_suffix(const Format format)
void read_xda(std::istream &in)
void read_comsol_mphtxt(std::istream &in)
static void skip_comment_lines(std::istream &in, const char comment_start)
void attach_triangulation(Triangulation< dim, spacedim > &tria)
void read_msh(std::istream &in)
static Format parse_format(const std::string &format_name)
void read_vtu(std::istream &in)
void read_tecplot(std::istream &in)
ExodusIIData read_exodusii(const std::string &filename, const bool apply_all_indicators_to_manifolds=false)
void read_dbmesh(std::istream &in)
void read_ucd(std::istream &in, const bool apply_all_indicators_to_manifolds=false)
void read(std::istream &in, Format format=Default)
static void debug_output_grid(const std::vector< CellData< dim > > &cells, const std::vector< Point< spacedim > > &vertices, std::ostream &out)
static std::string get_format_names()
void read_unv(std::istream &in)
static void parse_tecplot_header(std::string &header, std::vector< unsigned int > &tecplot2deal, unsigned int &n_vars, unsigned int &n_vertices, unsigned int &n_cells, std::vector< unsigned int > &IJK, bool &structured, bool &blocked)
std::string find(const std::string &filename, const char *open_mode="r")
std_cxx20::ranges::iota_view< unsigned int, unsigned int > vertex_indices() const
unsigned int n_vertices() const
unsigned int face_to_cell_vertices(const unsigned int face, const unsigned int vertex, const unsigned char face_orientation) const
unsigned int exodusii_vertex_to_deal_vertex(const unsigned int vertex_n) const
unsigned int get_dimension() const
ReferenceCell face_reference_cell(const unsigned int face_no) const
unsigned int exodusii_face_to_deal_face(const unsigned int face_n) const
static ReferenceCell n_vertices_to_type(const int dim, const unsigned int n_vertices)
virtual void create_triangulation(const std::vector< Point< spacedim > > &vertices, const std::vector< CellData< dim > > &cells, const SubCellData &subcelldata)
void load(Archive &ar, const unsigned int version)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_FALLTHROUGH
IteratorRange< active_cell_iterator > active_cell_iterators() const
static ::ExceptionBase & ExcIO()
static ::ExceptionBase & ExcNeedsAssimp()
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
static ::ExceptionBase & ExcNeedsExodusII()
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcIndexRange(std::size_t arg1, std::size_t arg2, std::size_t arg3)
static ::ExceptionBase & ExcInvalidState()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
constexpr const ReferenceCell Tetrahedron
constexpr const ReferenceCell Quadrilateral
constexpr const ReferenceCell Wedge
constexpr const ReferenceCell Pyramid
constexpr const ReferenceCell Invalid
constexpr const ReferenceCell Triangle
constexpr const ReferenceCell Hexahedron
constexpr const ReferenceCell Vertex
constexpr const ReferenceCell Line
std::pair< int, unsigned int > get_integer_at_position(const std::string &name, const unsigned int position)
std::vector< unsigned char > decode_base64(const std::string &base64_input)
std::vector< std::string > break_text_into_lines(const std::string &original_text, const unsigned int width, const char delimiter=' ')
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
bool match_at_string_start(const std::string &name, const std::string &pattern)
std::string decompress(const std::string &compressed_input)
const types::material_id invalid_material_id
const types::boundary_id internal_face_boundary_id
static const unsigned int invalid_unsigned_int
const types::manifold_id flat_manifold_id
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
std::vector< unsigned int > vertices
types::material_id material_id
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > vertex_indices()
std::vector< std::vector< int > > id_to_sideset_ids
std::vector< CellData< 2 > > boundary_quads
bool check_consistency(const unsigned int dim) const
std::vector< CellData< 1 > > boundary_lines
const ::Triangulation< dim, spacedim > & tria