26 #include <boost/archive/binary_iarchive.hpp>
27 #include <boost/io/ios_state.hpp>
28 #include <boost/property_tree/ptree.hpp>
29 #include <boost/property_tree/xml_parser.hpp>
30 #include <boost/serialization/serialization.hpp>
39 #ifdef DEAL_II_WITH_NETCDF
40 # include <netcdfcpp.h>
43 #ifdef DEAL_II_WITH_ASSIMP
44 # include <assimp/Importer.hpp>
45 # include <assimp/postprocess.h>
46 # include <assimp/scene.h>
63 template <
int spacedim>
65 assign_1d_boundary_ids(
66 const std::map<unsigned int, types::boundary_id> &boundary_ids,
69 if (boundary_ids.size() > 0)
70 for (
const auto &cell :
triangulation.active_cell_iterators())
72 if (boundary_ids.find(cell->vertex_index(f)) != boundary_ids.end())
77 "You are trying to prescribe boundary ids on the face "
78 "of a 1d cell (i.e., on a vertex), but this face is not actually at "
79 "the boundary of the mesh. This is not allowed."));
80 cell->face(f)->set_boundary_id(
81 boundary_ids.find(cell->vertex_index(f))->second);
86 template <
int dim,
int spacedim>
88 assign_1d_boundary_ids(
const std::map<unsigned int, types::boundary_id> &,
97 template <
int dim,
int spacedim>
99 : tria(nullptr, typeid(*this).name())
104 template <
int dim,
int spacedim>
113 template <
int dim,
int spacedim>
125 text[0] =
"# vtk DataFile Version 3.0";
128 text[3] =
"DATASET UNSTRUCTURED_GRID";
130 for (
unsigned int i = 0; i < 4; ++i)
135 line.compare(text[i]) == 0,
138 "While reading VTK file, failed to find a header line with text <") +
145 std::vector<Point<spacedim>>
vertices;
146 std::vector<CellData<dim>> cells;
155 if (keyword ==
"POINTS")
157 unsigned int n_vertices;
162 for (
unsigned int vertex = 0; vertex < n_vertices; ++vertex)
166 in >> x(0) >> x(1) >> x(2);
169 for (
unsigned int d = 0;
d < spacedim; ++
d)
177 "While reading VTK file, failed to find POINTS section"));
181 unsigned int n_geometric_objects = 0;
184 if (keyword ==
"CELLS")
186 in >> n_geometric_objects;
191 for (
unsigned int count = 0; count < n_geometric_objects; count++)
204 cells.emplace_back();
206 for (
unsigned int j = 0; j < type; j++)
207 in >> cells.back().vertices[j];
209 cells.back().material_id = 0;
221 for (
unsigned int j = 0; j < type;
231 for (
unsigned int j = 0; j < type;
242 "While reading VTK file, unknown file type encountered"));
248 for (
unsigned int count = 0; count < n_geometric_objects; count++)
260 cells.emplace_back();
262 for (
unsigned int j = 0; j < type; j++)
263 in >> cells.back().vertices[j];
265 cells.back().material_id = 0;
274 for (
unsigned int j = 0; j < type;
287 "While reading VTK file, unknown cell type encountered"));
292 for (
unsigned int count = 0; count < n_geometric_objects; count++)
300 "While reading VTK file, unknown cell type encountered"));
301 cells.emplace_back();
303 for (
unsigned int j = 0; j < type; j++)
304 in >> cells.back().vertices[j];
306 cells.back().material_id = 0;
312 "While reading VTK file, failed to find CELLS section"));
320 keyword ==
"CELL_TYPES",
322 "While reading VTK file, missing CELL_TYPES section. Found <" +
323 keyword +
"> instead.")));
327 n_ints == n_geometric_objects,
328 ExcMessage(
"The VTK reader found a CELL_DATA statement "
329 "that lists a total of " +
331 " cell data objects, but this needs to "
332 "equal the number of cells (which is " +
334 ") plus the number of quads (" +
336 " in 3d or the number of lines (" +
341 for (
unsigned int i = 0; i < n_ints; ++i)
345 while (in >> keyword)
346 if (keyword ==
"CELL_DATA")
352 ExcMessage(
"The VTK reader found a CELL_DATA statement "
353 "that lists a total of " +
355 " cell data objects, but this needs to "
356 "equal the number of cells (which is " +
358 ") plus the number of quads (" +
361 " in 3d or the number of lines (" +
366 const std::vector<std::string> data_sets{
"MaterialID",
369 for (
unsigned int i = 0; i < data_sets.size(); ++i)
372 while (in >> keyword)
373 if (keyword ==
"SCALARS")
378 std::string field_name;
380 if (std::find(data_sets.begin(),
382 field_name) == data_sets.end())
394 std::getline(in, line);
397 std::min(
static_cast<std::size_t
>(3),
398 line.size() - 1)) ==
"int",
400 "While reading VTK file, material- and manifold IDs can only have type 'int'."));
404 keyword ==
"LOOKUP_TABLE",
406 "While reading VTK file, missing keyword 'LOOKUP_TABLE'."));
410 keyword ==
"default",
412 "While reading VTK file, missing keyword 'default'."));
419 for (
unsigned int i = 0; i < cells.size(); i++)
423 if (field_name ==
"MaterialID")
424 cells[i].material_id =
426 else if (field_name ==
"ManifoldID")
427 cells[i].manifold_id =
439 if (field_name ==
"MaterialID")
440 boundary_quad.material_id =
442 else if (field_name ==
"ManifoldID")
443 boundary_quad.manifold_id =
452 if (field_name ==
"MaterialID")
453 boundary_line.material_id =
455 else if (field_name ==
"ManifoldID")
456 boundary_line.manifold_id =
468 if (field_name ==
"MaterialID")
469 boundary_line.material_id =
471 else if (field_name ==
"ManifoldID")
472 boundary_line.manifold_id =
491 tria->create_triangulation_compatibility(
vertices, cells, subcelldata);
498 "While reading VTK file, failed to find CELLS section"));
503 template <
int dim,
int spacedim>
507 namespace pt = boost::property_tree;
509 pt::read_xml(in, tree);
510 auto section = tree.get_optional<std::string>(
"VTKFile.dealiiData");
514 "While reading a VTU file, failed to find dealiiData section. "
515 "Notice that we can only read grid files in .vtu format that "
516 "were created by the deal.II library, using a call to "
517 "GridOut::write_vtu(), where the flag "
518 "GridOutFlags::Vtu::serialize_triangulation is set to true."));
522 const auto string_archive =
524 std::istringstream in_stream(string_archive);
525 boost::archive::binary_iarchive ia(in_stream);
530 template <
int dim,
int spacedim>
538 skip_comment_lines(in,
'#');
549 AssertThrow(tmp == 2411, ExcUnknownSectionType(tmp));
551 std::vector<Point<spacedim>>
vertices;
573 in >> x[0] >> x[1] >> x[2];
577 for (
unsigned int d = 0;
d < spacedim;
d++)
592 AssertThrow(tmp == 2412, ExcUnknownSectionType(tmp));
594 std::vector<CellData<dim>> cells;
623 AssertThrow((type == 11) || (type == 44) || (type == 94) || (type == 115),
624 ExcUnknownElementType(type));
626 if ((((type == 44) || (type == 94)) && (dim == 2)) ||
627 ((type == 115) && (dim == 3)))
629 cells.emplace_back();
633 in >> cells.back().vertices[v];
635 cells.back().material_id = 0;
638 cells.back().vertices[v] =
vertex_indices[cells.back().vertices[v]];
640 cell_indices[no] = no_cell;
644 else if (((type == 11) && (dim == 2)) ||
645 ((type == 11) && (dim == 3)))
653 for (
unsigned int &vertex :
659 for (
unsigned int &vertex :
663 line_indices[no] = no_line;
667 else if (((type == 44) || (type == 94)) && (dim == 3))
672 for (
unsigned int &vertex :
678 for (
unsigned int &vertex :
682 quad_indices[no] = no_quad;
690 "> when running in dim=" +
709 AssertThrow((tmp == 2467) || (tmp == 2477), ExcUnknownSectionType(tmp));
731 const unsigned int n_lines =
732 (n_entities % 2 == 0) ? (n_entities / 2) : ((n_entities + 1) / 2);
734 for (
unsigned int line = 0; line < n_lines; line++)
736 unsigned int n_fragments;
738 if (line == n_lines - 1)
739 n_fragments = (n_entities % 2 == 0) ? (2) : (1);
743 for (
unsigned int no_fragment = 0; no_fragment < n_fragments;
749 if (cell_indices.count(no) > 0)
750 cells[cell_indices[no]].material_id = id;
752 if (line_indices.count(no) > 0)
756 if (quad_indices.count(no) > 0)
774 tria->create_triangulation_compatibility(
vertices, cells, subcelldata);
779 template <
int dim,
int spacedim>
782 const bool apply_all_indicators_to_manifolds)
788 skip_comment_lines(in,
'#');
791 unsigned int n_vertices;
801 std::vector<Point<spacedim>>
vertices(n_vertices);
807 for (
unsigned int vertex = 0; vertex < n_vertices; ++vertex)
814 in >> vertex_number >> x[0] >> x[1] >> x[2];
817 for (
unsigned int d = 0;
d < spacedim; ++
d)
826 std::vector<CellData<dim>> cells;
829 for (
unsigned int cell = 0; cell <
n_cells; ++cell)
838 std::string cell_type;
848 if (((cell_type ==
"line") && (dim == 1)) ||
849 ((cell_type ==
"quad") && (dim == 2)) ||
850 ((cell_type ==
"hex") && (dim == 3)))
854 cells.emplace_back();
856 in >> cells.back().vertices[i];
867 if (apply_all_indicators_to_manifolds)
868 cells.back().manifold_id =
878 cells.back().vertices[i] =
885 cells.back().vertices[i]));
890 else if ((cell_type ==
"line") && ((dim == 2) || (dim == 3)))
910 if (apply_all_indicators_to_manifolds)
927 for (
unsigned int &vertex :
939 else if ((cell_type ==
"quad") && (dim == 3))
961 if (apply_all_indicators_to_manifolds)
978 for (
unsigned int &vertex :
992 AssertThrow(
false, ExcUnknownIdentifier(cell_type));
1004 if (dim == spacedim)
1008 tria->create_triangulation_compatibility(
vertices, cells, subcelldata);
1013 template <
int dim,
int spacedim>
1020 read_in_abaqus(std::istream &in);
1022 write_out_avs_ucd(std::ostream &out)
const;
1025 const double tolerance;
1028 get_global_node_numbers(
const int face_cell_no,
1029 const int face_cell_face_no)
const;
1032 std::vector<std::vector<double>> node_list;
1035 std::vector<std::vector<double>> cell_list;
1037 std::vector<std::vector<double>> face_list;
1040 std::map<std::string, std::vector<int>> elsets_list;
1044 template <
int dim,
int spacedim>
1047 const bool apply_all_indicators_to_manifolds)
1054 Assert((spacedim == 2 && dim == spacedim) ||
1055 (spacedim == 3 && (dim == spacedim || dim == spacedim - 1)),
1061 Abaqus_to_UCD<dim, spacedim> abaqus_to_ucd;
1062 abaqus_to_ucd.read_in_abaqus(in);
1064 std::stringstream in_ucd;
1065 abaqus_to_ucd.write_out_avs_ucd(in_ucd);
1073 read_ucd(in_ucd, apply_all_indicators_to_manifolds);
1075 catch (std::exception &exc)
1077 std::cerr <<
"Exception on processing internal UCD data: " << std::endl
1078 << exc.what() << std::endl;
1083 "Internal conversion from ABAQUS file to UCD format was unsuccessful. "
1084 "More information is provided in an error message printed above. "
1085 "Are you sure that your ABAQUS mesh file conforms with the requirements "
1086 "listed in the documentation?"));
1093 "Internal conversion from ABAQUS file to UCD format was unsuccessful. "
1094 "Are you sure that your ABAQUS mesh file conforms with the requirements "
1095 "listed in the documentation?"));
1100 template <
int dim,
int spacedim>
1110 skip_comment_lines(in,
'#');
1116 AssertThrow(line ==
"MeshVersionFormatted 0", ExcInvalidDBMESHInput(line));
1118 skip_empty_lines(in);
1122 AssertThrow(line ==
"Dimension", ExcInvalidDBMESHInput(line));
1123 unsigned int dimension;
1125 AssertThrow(dimension == dim, ExcDBMESHWrongDimension(dimension));
1126 skip_empty_lines(in);
1139 while (getline(in, line), line.find(
"# END") == std::string::npos)
1141 skip_empty_lines(in);
1146 AssertThrow(line ==
"Vertices", ExcInvalidDBMESHInput(line));
1148 unsigned int n_vertices;
1152 std::vector<Point<spacedim>>
vertices(n_vertices);
1153 for (
unsigned int vertex = 0; vertex < n_vertices; ++vertex)
1156 for (
unsigned int d = 0;
d < dim; ++
d)
1163 skip_empty_lines(in);
1169 AssertThrow(line ==
"Edges", ExcInvalidDBMESHInput(line));
1171 unsigned int n_edges;
1173 for (
unsigned int edge = 0; edge < n_edges; ++edge)
1182 skip_empty_lines(in);
1191 AssertThrow(line ==
"CrackedEdges", ExcInvalidDBMESHInput(line));
1194 for (
unsigned int edge = 0; edge < n_edges; ++edge)
1203 skip_empty_lines(in);
1209 AssertThrow(line ==
"Quadrilaterals", ExcInvalidDBMESHInput(line));
1211 std::vector<CellData<dim>> cells;
1215 for (
unsigned int cell = 0; cell <
n_cells; ++cell)
1219 cells.emplace_back();
1222 in >> cells.back().vertices[i];
1225 (
static_cast<unsigned int>(cells.back().vertices[i]) <=
1229 --cells.back().vertices[i];
1237 skip_empty_lines(in);
1245 while (getline(in, line), ((line.find(
"End") == std::string::npos) && (in)))
1262 tria->create_triangulation_compatibility(
vertices, cells, subcelldata);
1267 template <
int dim,
int spacedim>
1288 unsigned int n_vertices;
1299 for (
unsigned int i = 0; i < 8; ++i)
1303 std::vector<CellData<2>> cells(
n_cells);
1306 for (
unsigned int cell = 0; cell <
n_cells; ++cell)
1316 for (
unsigned int &vertex : cells[cell].
vertices)
1323 std::vector<Point<2>>
vertices(n_vertices);
1324 for (
unsigned int vertex = 0; vertex < n_vertices; ++vertex)
1329 in >> x[0] >> x[1] >> x[2];
1332 for (
unsigned int d = 0;
d < 2; ++
d)
1342 tria->create_triangulation_compatibility(
vertices, cells, subcelldata);
1354 static const unsigned int xda_to_dealII_map[] = {0, 1, 5, 4, 3, 2, 6, 7};
1361 unsigned int n_vertices;
1372 for (
unsigned int i = 0; i < 8; ++i)
1376 std::vector<CellData<3>> cells(
n_cells);
1379 for (
unsigned int cell = 0; cell <
n_cells; ++cell)
1389 unsigned int xda_ordered_nodes[8];
1391 for (
unsigned int &xda_ordered_node : xda_ordered_nodes)
1392 in >> xda_ordered_node;
1394 for (
unsigned int i = 0; i < 8; i++)
1395 cells[cell].
vertices[i] = xda_ordered_nodes[xda_to_dealII_map[i]];
1401 std::vector<Point<3>>
vertices(n_vertices);
1402 for (
unsigned int vertex = 0; vertex < n_vertices; ++vertex)
1407 in >> x[0] >> x[1] >> x[2];
1410 for (
unsigned int d = 0;
d < 3; ++
d)
1420 tria->create_triangulation_compatibility(
vertices, cells, subcelldata);
1425 template <
int dim,
int spacedim>
1432 unsigned int n_vertices;
1439 std::array<std::map<int, int>, 4> tag_maps;
1444 unsigned int gmsh_file_format = 0;
1445 if (line ==
"@f$NOD")
1446 gmsh_file_format = 10;
1447 else if (line ==
"@f$MeshFormat")
1448 gmsh_file_format = 20;
1454 if (gmsh_file_format == 20)
1457 unsigned int file_type, data_size;
1459 in >> version >> file_type >> data_size;
1462 gmsh_file_format =
static_cast<unsigned int>(version * 10);
1470 AssertThrow(line ==
"@f$EndMeshFormat", ExcInvalidGMSHInput(line));
1474 if (line ==
"@f$PhysicalNames")
1480 while (line !=
"@f$EndPhysicalNames");
1485 if (line ==
"@f$Entities")
1487 unsigned long n_points, n_curves, n_surfaces, n_volumes;
1489 in >> n_points >> n_curves >> n_surfaces >> n_volumes;
1490 for (
unsigned int i = 0; i < n_points; ++i)
1494 unsigned int n_physicals;
1495 double box_min_x, box_min_y, box_min_z, box_max_x, box_max_y,
1499 if (gmsh_file_format > 40)
1501 in >> tag >> box_min_x >> box_min_y >> box_min_z >>
1503 box_max_x = box_min_x;
1504 box_max_y = box_min_y;
1505 box_max_z = box_min_z;
1509 in >> tag >> box_min_x >> box_min_y >> box_min_z >>
1510 box_max_x >> box_max_y >> box_max_z >> n_physicals;
1514 ExcMessage(
"More than one tag is not supported!"));
1516 int physical_tag = 0;
1517 for (
unsigned int j = 0; j < n_physicals; ++j)
1519 tag_maps[0][tag] = physical_tag;
1521 for (
unsigned int i = 0; i < n_curves; ++i)
1525 unsigned int n_physicals;
1526 double box_min_x, box_min_y, box_min_z, box_max_x, box_max_y,
1530 in >> tag >> box_min_x >> box_min_y >> box_min_z >> box_max_x >>
1531 box_max_y >> box_max_z >> n_physicals;
1534 ExcMessage(
"More than one tag is not supported!"));
1536 int physical_tag = 0;
1537 for (
unsigned int j = 0; j < n_physicals; ++j)
1539 tag_maps[1][tag] = physical_tag;
1543 for (
unsigned int j = 0; j < n_points; ++j)
1547 for (
unsigned int i = 0; i < n_surfaces; ++i)
1551 unsigned int n_physicals;
1552 double box_min_x, box_min_y, box_min_z, box_max_x, box_max_y,
1556 in >> tag >> box_min_x >> box_min_y >> box_min_z >> box_max_x >>
1557 box_max_y >> box_max_z >> n_physicals;
1560 ExcMessage(
"More than one tag is not supported!"));
1562 int physical_tag = 0;
1563 for (
unsigned int j = 0; j < n_physicals; ++j)
1565 tag_maps[2][tag] = physical_tag;
1569 for (
unsigned int j = 0; j < n_curves; ++j)
1572 for (
unsigned int i = 0; i < n_volumes; ++i)
1576 unsigned int n_physicals;
1577 double box_min_x, box_min_y, box_min_z, box_max_x, box_max_y,
1581 in >> tag >> box_min_x >> box_min_y >> box_min_z >> box_max_x >>
1582 box_max_y >> box_max_z >> n_physicals;
1585 ExcMessage(
"More than one tag is not supported!"));
1587 int physical_tag = 0;
1588 for (
unsigned int j = 0; j < n_physicals; ++j)
1590 tag_maps[3][tag] = physical_tag;
1594 for (
unsigned int j = 0; j < n_surfaces; ++j)
1598 AssertThrow(line ==
"@f$EndEntities", ExcInvalidGMSHInput(line));
1603 if (line ==
"@f$PartitionedEntities")
1609 while (line !=
"@f$EndPartitionedEntities");
1616 AssertThrow(line ==
"@f$Nodes", ExcInvalidGMSHInput(line));
1620 int n_entity_blocks = 1;
1621 if (gmsh_file_format > 40)
1625 in >> n_entity_blocks >> n_vertices >> min_node_tag >> max_node_tag;
1627 else if (gmsh_file_format == 40)
1629 in >> n_entity_blocks >> n_vertices;
1633 std::vector<Point<spacedim>>
vertices(n_vertices);
1640 unsigned int global_vertex = 0;
1641 for (
int entity_block = 0; entity_block < n_entity_blocks; ++entity_block)
1644 unsigned long numNodes;
1646 if (gmsh_file_format < 40)
1648 numNodes = n_vertices;
1655 int tagEntity, dimEntity;
1656 in >> tagEntity >> dimEntity >> parametric >> numNodes;
1659 std::vector<int> vertex_numbers;
1661 if (gmsh_file_format > 40)
1662 for (
unsigned long vertex_per_entity = 0;
1663 vertex_per_entity < numNodes;
1664 ++vertex_per_entity)
1666 in >> vertex_number;
1667 vertex_numbers.push_back(vertex_number);
1670 for (
unsigned long vertex_per_entity = 0; vertex_per_entity < numNodes;
1671 ++vertex_per_entity, ++global_vertex)
1677 if (gmsh_file_format > 40)
1679 vertex_number = vertex_numbers[vertex_per_entity];
1680 in >> x[0] >> x[1] >> x[2];
1683 in >> vertex_number >> x[0] >> x[1] >> x[2];
1685 for (
unsigned int d = 0;
d < spacedim; ++
d)
1691 if (parametric != 0)
1706 static const std::string end_nodes_marker[] = {
"@f$ENDNOD",
"@f$EndNodes"};
1707 AssertThrow(line == end_nodes_marker[gmsh_file_format == 10 ? 0 : 1],
1708 ExcInvalidGMSHInput(line));
1712 static const std::string begin_elements_marker[] = {
"@f$ELM",
"@f$Elements"};
1713 AssertThrow(line == begin_elements_marker[gmsh_file_format == 10 ? 0 : 1],
1714 ExcInvalidGMSHInput(line));
1717 if (gmsh_file_format > 40)
1721 in >> n_entity_blocks >>
n_cells >> min_node_tag >> max_node_tag;
1723 else if (gmsh_file_format == 40)
1725 in >> n_entity_blocks >>
n_cells;
1729 n_entity_blocks = 1;
1736 std::vector<CellData<dim>> cells;
1738 std::map<unsigned int, types::boundary_id> boundary_ids_1d;
1741 unsigned int global_cell = 0;
1742 for (
int entity_block = 0; entity_block < n_entity_blocks; ++entity_block)
1745 unsigned long numElements;
1748 if (gmsh_file_format < 40)
1754 else if (gmsh_file_format == 40)
1756 int tagEntity, dimEntity;
1757 in >> tagEntity >> dimEntity >> cell_type >> numElements;
1763 int tagEntity, dimEntity;
1764 in >> dimEntity >> tagEntity >> cell_type >> numElements;
1768 for (
unsigned int cell_per_entity = 0; cell_per_entity < numElements;
1769 ++cell_per_entity, ++global_cell)
1778 unsigned int nod_num;
1798 unsigned int elm_number = 0;
1799 if (gmsh_file_format < 40)
1805 if (gmsh_file_format < 20)
1811 else if (gmsh_file_format < 40)
1816 unsigned int n_tags;
1823 for (
unsigned int i = 1; i < n_tags; ++i)
1852 if (((cell_type == 1) && (dim == 1)) ||
1853 ((cell_type == 3) && (dim == 2)) ||
1854 ((cell_type == 5) && (dim == 3)))
1859 "Number of nodes does not coincide with the "
1860 "number required for this object"));
1863 cells.emplace_back();
1865 in >> cells.back().vertices[i];
1887 ExcInvalidVertexIndexGmsh(cell_per_entity,
1889 cells.back().vertices[i]));
1892 cells.back().vertices[i] =
1896 else if ((cell_type == 1) && ((dim == 2) || (dim == 3)))
1920 for (
unsigned int &vertex :
1934 else if ((cell_type == 3) && (dim == 3))
1960 for (
unsigned int &vertex :
1973 else if (cell_type == 15)
1976 unsigned int node_index = 0;
1977 if (gmsh_file_format < 20)
1979 for (
unsigned int i = 0; i < nod_num; ++i)
2001 ExcMessage(
"Found triangles while reading a file "
2002 "in gmsh format. deal.II does not "
2003 "support triangles"));
2005 ExcMessage(
"Found tetrahedra while reading a file "
2006 "in gmsh format. deal.II does not "
2007 "support tetrahedra"));
2009 AssertThrow(
false, ExcGmshUnsupportedGeometry(cell_type));
2017 static const std::string end_elements_marker[] = {
"@f$ENDELM",
"@f$EndElements"};
2018 AssertThrow(line == end_elements_marker[gmsh_file_format == 10 ? 0 : 1],
2019 ExcInvalidGMSHInput(line));
2028 AssertThrow(cells.size() > 0, ExcGmshNoCellInformation());
2034 if (dim == spacedim)
2038 tria->create_triangulation_compatibility(
vertices, cells, subcelldata);
2043 assign_1d_boundary_ids(boundary_ids_1d, *tria);
2081 #ifndef DEAL_II_WITH_NETCDF
2085 const unsigned int dim = 2;
2086 const unsigned int spacedim = 2;
2110 const unsigned int coord = 1;
2116 const unsigned int x2d = 0;
2117 const unsigned int y2d = 2;
2126 NcFile nc(filename.c_str());
2130 NcDim *elements_dim = nc.get_dim(
"no_of_elements");
2132 const unsigned int n_cells = elements_dim->size();
2136 NcDim *marker_dim = nc.get_dim(
"no_of_markers");
2138 const unsigned int n_markers = marker_dim->size();
2140 NcVar *marker_var = nc.get_var(
"marker");
2143 AssertThrow(
static_cast<unsigned int>(marker_var->get_dim(0)->size()) ==
2147 std::vector<int> marker(n_markers);
2150 marker_var->get(&*marker.begin(), n_markers);
2155 NcDim *bquads_dim = nc.get_dim(
"no_of_surfacequadrilaterals");
2157 const unsigned int n_bquads = bquads_dim->size();
2159 NcVar *bmarker_var = nc.get_var(
"boundarymarker_of_surfaces");
2162 AssertThrow(
static_cast<unsigned int>(bmarker_var->get_dim(0)->size()) ==
2166 std::vector<int> bmarker(n_bquads);
2167 bmarker_var->get(&*bmarker.begin(), n_bquads);
2172 std::map<int, unsigned int> n_bquads_per_bmarker;
2173 for (
unsigned int i = 0; i < n_markers; ++i)
2177 AssertThrow(n_bquads_per_bmarker.find(marker[i]) ==
2178 n_bquads_per_bmarker.end(),
2181 n_bquads_per_bmarker[marker[i]] =
2182 count(bmarker.begin(), bmarker.end(), marker[i]);
2201 NcDim *quad_vertices_dim = nc.get_dim(
"points_per_surfacequadrilateral");
2203 const unsigned int vertices_per_quad = quad_vertices_dim->size();
2207 NcVar *vertex_indices_var = nc.get_var(
"points_of_surfacequadrilaterals");
2211 vertex_indices_var->get_dim(0)->size()) == n_bquads,
2214 vertex_indices_var->get_dim(1)->size()) == vertices_per_quad,
2229 NcDim *vertices_dim = nc.get_dim(
"no_of_points");
2231 const unsigned int n_vertices = vertices_dim->size();
2233 NcVar *points_xc = nc.get_var(
"points_xc");
2234 NcVar *points_yc = nc.get_var(
"points_yc");
2235 NcVar *points_zc = nc.get_var(
"points_zc");
2242 AssertThrow(points_yc->get_dim(0)->size() ==
static_cast<int>(n_vertices),
2244 AssertThrow(points_zc->get_dim(0)->size() ==
static_cast<int>(n_vertices),
2246 AssertThrow(points_xc->get_dim(0)->size() ==
static_cast<int>(n_vertices),
2248 std::vector<std::vector<double>> point_values(
2249 3, std::vector<double>(n_vertices));
2250 points_xc->get(point_values[0].data(), n_vertices);
2251 points_yc->get(point_values[1].data(), n_vertices);
2252 points_zc->get(point_values[2].data(), n_vertices);
2255 std::vector<Point<spacedim>>
vertices(n_vertices);
2256 for (
unsigned int i = 0; i < n_vertices; ++i)
2258 vertices[i](0) = point_values[x2d][i];
2259 vertices[i](1) = point_values[y2d][i];
2265 std::map<int, bool> zero_plane_markers;
2266 for (
unsigned int quad = 0; quad < n_bquads; ++quad)
2268 bool zero_plane =
true;
2269 for (
unsigned int i = 0; i < vertices_per_quad; ++i)
2270 if (point_values[coord][
vertex_indices[quad * vertices_per_quad + i]] !=
2278 zero_plane_markers[bmarker[quad]] =
true;
2280 unsigned int sum_of_zero_plane_cells = 0;
2281 for (std::map<int, bool>::const_iterator iter = zero_plane_markers.begin();
2282 iter != zero_plane_markers.end();
2284 sum_of_zero_plane_cells += n_bquads_per_bmarker[iter->first];
2290 std::vector<CellData<dim>> cells(
n_cells);
2291 for (
unsigned int quad = 0, cell = 0; quad < n_bquads; ++quad)
2293 bool zero_plane =
false;
2294 for (std::map<int, bool>::const_iterator iter =
2295 zero_plane_markers.begin();
2296 iter != zero_plane_markers.end();
2298 if (bmarker[quad] == iter->first)
2306 for (
unsigned int i = 0; i < vertices_per_quad; ++i)
2312 cells[cell].vertices[i] =
2322 tria->create_triangulation_compatibility(
vertices, cells, subcelldata);
2331 #ifndef DEAL_II_WITH_NETCDF
2338 const unsigned int dim = 3;
2339 const unsigned int spacedim = 3;
2345 NcFile nc(filename.c_str());
2349 NcDim *elements_dim = nc.get_dim(
"no_of_elements");
2351 const unsigned int n_cells = elements_dim->size();
2353 NcDim *hexes_dim = nc.get_dim(
"no_of_hexaeders");
2355 const unsigned int n_hexes = hexes_dim->size();
2357 ExcMessage(
"deal.II can handle purely hexaedral grids, only."));
2363 NcDim *hex_vertices_dim = nc.get_dim(
"points_per_hexaeder");
2365 const unsigned int vertices_per_hex = hex_vertices_dim->size();
2369 NcVar *vertex_indices_var = nc.get_var(
"points_of_hexaeders");
2373 vertex_indices_var->get_dim(0)->size()) ==
n_cells,
2376 vertex_indices_var->get_dim(1)->size()) == vertices_per_hex,
2391 NcDim *vertices_dim = nc.get_dim(
"no_of_points");
2393 const unsigned int n_vertices = vertices_dim->size();
2395 NcVar *points_xc = nc.get_var(
"points_xc");
2396 NcVar *points_yc = nc.get_var(
"points_yc");
2397 NcVar *points_zc = nc.get_var(
"points_zc");
2404 AssertThrow(points_yc->get_dim(0)->size() ==
static_cast<int>(n_vertices),
2406 AssertThrow(points_zc->get_dim(0)->size() ==
static_cast<int>(n_vertices),
2408 AssertThrow(points_xc->get_dim(0)->size() ==
static_cast<int>(n_vertices),
2410 std::vector<std::vector<double>> point_values(
2411 3, std::vector<double>(n_vertices));
2412 points_xc->get(point_values[0].data(), n_vertices);
2413 points_yc->get(point_values[1].data(), n_vertices);
2414 points_zc->get(point_values[2].data(), n_vertices);
2417 std::vector<Point<spacedim>>
vertices(n_vertices);
2418 for (
unsigned int i = 0; i < n_vertices; ++i)
2420 vertices[i](0) = point_values[0][i];
2421 vertices[i](1) = point_values[1][i];
2422 vertices[i](2) = point_values[2][i];
2426 std::vector<CellData<dim>> cells(
n_cells);
2427 for (
unsigned int cell = 0; cell <
n_cells; ++cell)
2428 for (
unsigned int i = 0; i < vertices_per_hex; ++i)
2440 NcDim *quad_vertices_dim = nc.get_dim(
"points_per_surfacequadrilateral");
2442 const unsigned int vertices_per_quad = quad_vertices_dim->size();
2446 NcVar *bvertex_indices_var = nc.get_var(
"points_of_surfacequadrilaterals");
2449 const unsigned int n_bquads = bvertex_indices_var->get_dim(0)->size();
2451 bvertex_indices_var->get_dim(1)->size()) ==
2455 std::vector<int> bvertex_indices(n_bquads * vertices_per_quad);
2456 bvertex_indices_var->get(&*bvertex_indices.begin(),
2463 NcDim *bquads_dim = nc.get_dim(
"no_of_surfacequadrilaterals");
2465 AssertThrow(
static_cast<unsigned int>(bquads_dim->size()) == n_bquads,
2468 NcVar *bmarker_var = nc.get_var(
"boundarymarker_of_surfaces");
2471 AssertThrow(
static_cast<unsigned int>(bmarker_var->get_dim(0)->size()) ==
2475 std::vector<int> bmarker(n_bquads);
2476 bmarker_var->get(&*bmarker.begin(), n_bquads);
2482 for (
const int id : bmarker)
2494 for (
unsigned int i = 0; i < n_bquads; ++i)
2496 for (
unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_face; ++v)
2507 tria->create_triangulation_compatibility(
vertices, cells, subcelldata);
2512 template <
int dim,
int spacedim>
2515 std::string & header,
2516 std::vector<unsigned int> &tecplot2deal,
2517 unsigned int & n_vars,
2518 unsigned int & n_vertices,
2520 std::vector<unsigned int> &IJK,
2545 std::transform(header.begin(), header.end(), header.begin(), ::toupper);
2549 std::replace(header.begin(), header.end(),
'\t',
' ');
2550 std::replace(header.begin(), header.end(),
',',
' ');
2551 std::replace(header.begin(), header.end(),
'\n',
' ');
2558 if (header[pos + 1] ==
' ')
2559 header.erase(pos + 1, 1);
2560 else if (header[pos - 1] ==
' ')
2562 header.erase(pos - 1, 1);
2566 pos = header.find(
'=', ++pos);
2569 std::vector<std::string> entries =
2573 for (
unsigned int i = 0; i < entries.size(); ++i)
2582 tecplot2deal[0] = 0;
2585 while (entries[i][0] ==
'"')
2587 if (entries[i] ==
"\"X\"")
2588 tecplot2deal[0] = n_vars;
2589 else if (entries[i] ==
"\"Y\"")
2595 tecplot2deal[1] = n_vars;
2597 else if (entries[i] ==
"\"Z\"")
2603 tecplot2deal[2] = n_vars;
2615 "Tecplot file must contain at least one variable for each dimension"));
2616 for (
unsigned int d = 1;
d < dim; ++
d)
2618 tecplot2deal[
d] > 0,
2620 "Tecplot file must contain at least one variable for each dimension."));
2625 "ZONETYPE=FELINESEG") &&
2629 "ZONETYPE=FEQUADRILATERAL") &&
2633 "ZONETYPE=FEBRICK") &&
2641 "The tecplot file contains an unsupported ZONETYPE."));
2644 "DATAPACKING=POINT"))
2647 "DATAPACKING=BLOCK"))
2670 "ET=QUADRILATERAL") &&
2682 "The tecplot file contains an unsupported ElementType."));
2690 dim > 1 || IJK[1] == 1,
2692 "Parameter 'J=' found in tecplot, although this is only possible for dimensions greater than 1."));
2698 dim > 2 || IJK[2] == 1,
2700 "Parameter 'K=' found in tecplot, although this is only possible for dimensions greater than 2."));
2715 for (
unsigned int d = 0;
d < dim; ++
d)
2720 "Tecplot file does not contain a complete and consistent set of parameters"));
2721 n_vertices *= IJK[
d];
2730 "Tecplot file does not contain a complete and consistent set of parameters"));
2736 n_cells = *std::max_element(IJK.begin(), IJK.end());
2740 "Tecplot file does not contain a complete and consistent set of parameters"));
2750 const unsigned int dim = 2;
2751 const unsigned int spacedim = 2;
2756 skip_comment_lines(in,
'#');
2759 std::string line, header;
2766 std::string letters =
"abcdfghijklmnopqrstuvwxyzABCDFGHIJKLMNOPQRSTUVWXYZ";
2769 while (line.find_first_of(letters) != std::string::npos)
2771 header +=
" " + line;
2778 std::vector<unsigned int> tecplot2deal(dim);
2779 std::vector<unsigned int> IJK(dim);
2780 unsigned int n_vars, n_vertices,
n_cells;
2781 bool structured, blocked;
2783 parse_tecplot_header(header,
2798 std::vector<Point<spacedim>>
vertices(n_vertices + 1);
2801 std::vector<CellData<dim>> cells(
n_cells);
2817 unsigned int next_index = 0;
2821 if (tecplot2deal[0] == 0)
2825 std::vector<std::string> first_var =
2828 for (
unsigned int i = 1; i < first_var.size() + 1; ++i)
2829 vertices[i](0) = std::strtod(first_var[i - 1].c_str(), &endptr);
2834 for (
unsigned int j = first_var.size() + 1; j < n_vertices + 1; ++j)
2842 for (
unsigned int i = 1; i < n_vars; ++i)
2851 if (next_index == dim && structured)
2854 if ((next_index < dim) && (i == tecplot2deal[next_index]))
2857 for (
unsigned int j = 1; j < n_vertices + 1; ++j)
2865 for (
unsigned int j = 1; j < n_vertices + 1; ++j)
2877 std::vector<double> vars(n_vars);
2882 std::vector<std::string> first_vertex =
2885 for (
unsigned int d = 0;
d < dim; ++
d)
2887 std::strtod(first_vertex[tecplot2deal[
d]].c_str(), &endptr);
2891 for (
unsigned int v = 2; v < n_vertices + 1; ++v)
2893 for (
unsigned int i = 0; i < n_vars; ++i)
2899 for (
unsigned int i = 0; i < dim; ++i)
2900 vertices[v](i) = vars[tecplot2deal[i]];
2908 unsigned int I = IJK[0], J = IJK[1];
2910 unsigned int cell = 0;
2912 for (
unsigned int j = 0; j < J - 1; ++j)
2913 for (
unsigned int i = 1; i < I; ++i)
2915 cells[cell].vertices[0] = i + j * I;
2916 cells[cell].vertices[1] = i + 1 + j * I;
2917 cells[cell].vertices[2] = i + 1 + (j + 1) * I;
2918 cells[cell].vertices[3] = i + (j + 1) * I;
2922 std::vector<unsigned int> boundary_vertices(2 * I + 2 * J - 4);
2924 for (
unsigned int i = 1; i < I + 1; ++i)
2926 boundary_vertices[k] = i;
2928 boundary_vertices[k] = i + (J - 1) * I;
2931 for (
unsigned int j = 1; j < J - 1; ++j)
2933 boundary_vertices[k] = 1 + j * I;
2935 boundary_vertices[k] = I + j * I;
2954 for (
unsigned int i = 0; i <
n_cells; ++i)
2965 for (
unsigned int &vertex : cells[i].
vertices)
2982 tria->create_triangulation_compatibility(
vertices, cells, subcelldata);
2987 template <
int dim,
int spacedim>
2996 template <
int dim,
int spacedim>
2999 const unsigned int mesh_index,
3000 const bool remove_duplicates,
3002 const bool ignore_unsupported_types)
3004 #ifdef DEAL_II_WITH_ASSIMP
3009 Assimp::Importer importer;
3012 const aiScene *scene =
3013 importer.ReadFile(filename.c_str(),
3014 aiProcess_RemoveComponent |
3015 aiProcess_JoinIdenticalVertices |
3016 aiProcess_ImproveCacheLocality | aiProcess_SortByPType |
3017 aiProcess_OptimizeGraph | aiProcess_OptimizeMeshes);
3023 ExcMessage(
"Input file contains no meshes."));
3026 (mesh_index < scene->mNumMeshes),
3029 unsigned int start_mesh =
3031 unsigned int end_mesh =
3036 std::vector<Point<spacedim>>
vertices;
3037 std::vector<CellData<dim>> cells;
3041 unsigned int v_offset = 0;
3042 unsigned int c_offset = 0;
3045 for (
unsigned int m = start_mesh; m < end_mesh; ++m)
3047 const aiMesh *mesh = scene->mMeshes[m];
3051 if ((dim == 2) && mesh->mPrimitiveTypes != aiPrimitiveType_POLYGON)
3058 else if ((dim == 1) && mesh->mPrimitiveTypes != aiPrimitiveType_LINE)
3066 const unsigned int n_vertices = mesh->mNumVertices;
3067 const aiVector3D * mVertices = mesh->mVertices;
3070 const unsigned int n_faces = mesh->mNumFaces;
3071 const aiFace * mFaces = mesh->mFaces;
3073 vertices.resize(v_offset + n_vertices);
3074 cells.resize(c_offset + n_faces);
3076 for (
unsigned int i = 0; i < n_vertices; ++i)
3077 for (
unsigned int d = 0;
d < spacedim; ++
d)
3080 unsigned int valid_cell = c_offset;
3081 for (
unsigned int i = 0; i < n_faces; ++i)
3087 cells[valid_cell].vertices[f] =
3088 mFaces[i].mIndices[f] + v_offset;
3090 cells[valid_cell].material_id = m;
3099 " vertices. We expected only " +
3104 cells.resize(valid_cell);
3109 v_offset += n_vertices;
3110 c_offset = valid_cell;
3114 if (cells.size() == 0)
3117 if (remove_duplicates)
3122 unsigned int n_verts = 0;
3126 std::vector<unsigned int> considered_vertices;
3128 vertices, cells, subcelldata, considered_vertices, tol);
3133 if (dim == spacedim)
3139 tria->create_triangulation_compatibility(
vertices, cells, subcelldata);
3141 tria->create_triangulation(
vertices, cells, subcelldata);
3145 (void)remove_duplicates;
3147 (void)ignore_unsupported_types;
3153 template <
int dim,
int spacedim>
3167 if (std::find_if(line.begin(), line.end(), [](
const char c) {
3172 for (
int i = line.length() - 1; i >= 0; --i)
3173 in.putback(line[i]);
3183 template <
int dim,
int spacedim>
3186 const char comment_start)
3191 while (in.get(c) && c == comment_start)
3194 while (in.get() !=
'\n')
3204 skip_empty_lines(in);
3209 template <
int dim,
int spacedim>
3227 double min_x =
vertices[cells[0].vertices[0]](0),
3229 min_y =
vertices[cells[0].vertices[0]](1),
3232 for (
unsigned int i = 0; i < cells.size(); ++i)
3234 for (
const auto vertex : cells[i].
vertices)
3248 out <<
"# cell " << i << std::endl;
3250 for (
const auto vertex : cells[i].
vertices)
3254 out <<
"set label \"" << i <<
"\" at " <<
center(0) <<
',' <<
center(1)
3255 <<
" center" << std::endl;
3258 for (
unsigned int f = 0; f < 2; ++f)
3260 <<
vertices[cells[i].vertices[f]](1) <<
" to "
3264 for (
unsigned int f = 2; f < 4; ++f)
3268 <<
vertices[cells[i].vertices[f]](1) << std::endl;
3274 <<
"set nokey" << std::endl
3275 <<
"pl [" << min_x <<
':' << max_x <<
"][" << min_y <<
':' << max_y
3276 <<
"] " << min_y << std::endl
3277 <<
"pause -1" << std::endl;
3288 for (
const auto &cell : cells)
3291 out <<
vertices[cell.vertices[0]] << std::endl
3292 <<
vertices[cell.vertices[1]] << std::endl
3296 out <<
vertices[cell.vertices[1]] << std::endl
3297 <<
vertices[cell.vertices[2]] << std::endl
3301 out <<
vertices[cell.vertices[3]] << std::endl
3302 <<
vertices[cell.vertices[2]] << std::endl
3306 out <<
vertices[cell.vertices[0]] << std::endl
3307 <<
vertices[cell.vertices[3]] << std::endl
3311 out <<
vertices[cell.vertices[4]] << std::endl
3312 <<
vertices[cell.vertices[5]] << std::endl
3316 out <<
vertices[cell.vertices[5]] << std::endl
3317 <<
vertices[cell.vertices[6]] << std::endl
3321 out <<
vertices[cell.vertices[7]] << std::endl
3322 <<
vertices[cell.vertices[6]] << std::endl
3326 out <<
vertices[cell.vertices[4]] << std::endl
3327 <<
vertices[cell.vertices[7]] << std::endl
3331 out <<
vertices[cell.vertices[0]] << std::endl
3332 <<
vertices[cell.vertices[4]] << std::endl
3336 out <<
vertices[cell.vertices[1]] << std::endl
3337 <<
vertices[cell.vertices[5]] << std::endl
3341 out <<
vertices[cell.vertices[2]] << std::endl
3342 <<
vertices[cell.vertices[6]] << std::endl
3346 out <<
vertices[cell.vertices[3]] << std::endl
3347 <<
vertices[cell.vertices[7]] << std::endl
3355 template <
int dim,
int spacedim>
3363 if (format == Default)
3364 name = search.
find(filename);
3368 std::ifstream in(name.c_str());
3370 if (format == Default)
3374 if (dotpos < name.length() &&
3375 (dotpos > slashpos || slashpos == std::string::npos))
3377 std::string ext = name.substr(dotpos + 1);
3378 format = parse_format(ext);
3381 if (format == netcdf)
3382 read_netcdf(filename);
3388 template <
int dim,
int spacedim>
3392 if (format == Default)
3431 ExcMessage(
"There is no read_netcdf(istream &) function. "
3432 "Use the read_netcdf(string &filename) "
3433 "functions, instead."));
3442 ExcMessage(
"There is no read_assimp(istream &) function. "
3443 "Use the read_assimp(string &filename, ...) "
3444 "functions, instead."));
3455 template <
int dim,
int spacedim>
3484 return ".unknown_format";
3490 template <
int dim,
int spacedim>
3494 if (format_name ==
"dbmesh")
3497 if (format_name ==
"msh")
3500 if (format_name ==
"unv")
3503 if (format_name ==
"vtk")
3506 if (format_name ==
"vtu")
3510 if (format_name ==
"inp")
3513 if (format_name ==
"ucd")
3516 if (format_name ==
"xda")
3519 if (format_name ==
"netcdf")
3522 if (format_name ==
"nc")
3525 if (format_name ==
"tecplot")
3528 if (format_name ==
"dat")
3531 if (format_name ==
"plt")
3550 template <
int dim,
int spacedim>
3554 return "dbmesh|msh|unv|vtk|vtu|ucd|abaqus|xda|netcdf|tecplot|assimp";
3561 template <
int dim,
int spacedim>
3562 Abaqus_to_UCD<dim, spacedim>::Abaqus_to_UCD()
3575 from_string(
T &t,
const std::string &s, std::ios_base &(*f)(std::ios_base &))
3577 std::istringstream iss(s);
3578 return !(iss >> f >> t).fail();
3585 extract_int(
const std::string &s)
3588 for (
const char c : s)
3597 from_string(number, tmp, std::dec);
3603 template <
int dim,
int spacedim>
3605 Abaqus_to_UCD<dim, spacedim>::read_in_abaqus(std::istream &input_stream)
3614 while (std::getline(input_stream, line))
3617 std::transform(line.begin(), line.end(), line.begin(), ::toupper);
3619 if (line.compare(
"*HEADING") == 0 || line.compare(0, 2,
"**") == 0 ||
3620 line.compare(0, 5,
"*PART") == 0)
3623 while (std::getline(input_stream, line))
3629 else if (line.compare(0, 5,
"*NODE") == 0)
3638 while (std::getline(input_stream, line))
3643 std::vector<double> node(spacedim + 1);
3645 std::istringstream iss(line);
3647 for (
unsigned int i = 0; i < spacedim + 1; ++i)
3648 iss >> node[i] >> comma;
3650 node_list.push_back(node);
3653 else if (line.compare(0, 8,
"*ELEMENT") == 0)
3668 const std::string before_material =
"ELSET=EB";
3669 const std::size_t idx = line.find(before_material);
3670 if (idx != std::string::npos)
3672 from_string(material,
3673 line.substr(idx + before_material.size()),
3679 while (std::getline(input_stream, line))
3684 std::istringstream iss(line);
3690 const unsigned int n_data_per_cell =
3692 std::vector<double> cell(n_data_per_cell);
3693 for (
unsigned int i = 0; i < n_data_per_cell; ++i)
3694 iss >> cell[i] >> comma;
3697 cell[0] =
static_cast<double>(material);
3698 cell_list.push_back(cell);
3701 else if (line.compare(0, 8,
"*SURFACE") == 0)
3712 const std::string name_key =
"NAME=";
3713 const std::size_t name_idx_start =
3714 line.find(name_key) + name_key.size();
3715 std::size_t name_idx_end = line.find(
',', name_idx_start);
3716 if (name_idx_end == std::string::npos)
3718 name_idx_end = line.size();
3720 const int b_indicator = extract_int(
3721 line.substr(name_idx_start, name_idx_end - name_idx_start));
3728 while (std::getline(input_stream, line))
3742 std::istringstream iss(line);
3749 std::vector<double> quad_node_list;
3750 const std::string elset_name = line.substr(0, line.find(
','));
3751 if (elsets_list.count(elset_name) != 0)
3755 iss >> stmp >> temp >> face_number;
3757 const std::vector<int> cells = elsets_list[elset_name];
3758 for (
const int cell : cells)
3762 get_global_node_numbers(el_idx, face_number);
3763 quad_node_list.insert(quad_node_list.begin(),
3766 face_list.push_back(quad_node_list);
3773 iss >> el_idx >> comma >> temp >> face_number;
3775 get_global_node_numbers(el_idx, face_number);
3776 quad_node_list.insert(quad_node_list.begin(), b_indicator);
3778 face_list.push_back(quad_node_list);
3782 else if (line.compare(0, 6,
"*ELSET") == 0)
3786 std::string elset_name;
3788 const std::string elset_key =
"*ELSET, ELSET=";
3789 const std::size_t idx = line.find(elset_key);
3790 if (idx != std::string::npos)
3792 const std::string comma =
",";
3793 const std::size_t first_comma = line.find(comma);
3794 const std::size_t second_comma =
3795 line.find(comma, first_comma + 1);
3796 const std::size_t elset_name_start =
3797 line.find(elset_key) + elset_key.size();
3798 elset_name = line.substr(elset_name_start,
3799 second_comma - elset_name_start);
3809 std::vector<int> elements;
3810 const std::size_t generate_idx = line.find(
"GENERATE");
3811 if (generate_idx != std::string::npos)
3814 std::getline(input_stream, line);
3815 std::istringstream iss(line);
3824 iss >> elid_start >> comma >> elid_end;
3828 "While reading an ABAQUS file, the reader "
3829 "expected a comma but found a <") +
3830 comma +
"> in the line <" + line +
">."));
3832 elid_start <= elid_end,
3835 "While reading an ABAQUS file, the reader encountered "
3836 "a GENERATE statement in which the upper bound <") +
3838 "> for the element numbers is not larger or equal "
3839 "than the lower bound <" +
3843 if (iss.rdbuf()->in_avail() != 0)
3844 iss >> comma >> elis_step;
3848 "While reading an ABAQUS file, the reader "
3849 "expected a comma but found a <") +
3850 comma +
"> in the line <" + line +
">."));
3852 for (
int i = elid_start; i <= elid_end; i += elis_step)
3853 elements.push_back(i);
3854 elsets_list[elset_name] = elements;
3856 std::getline(input_stream, line);
3861 while (std::getline(input_stream, line))
3866 std::istringstream iss(line);
3871 iss >> elid >> comma;
3876 "While reading an ABAQUS file, the reader "
3877 "expected a comma but found a <") +
3878 comma +
"> in the line <" + line +
">."));
3880 elements.push_back(elid);
3884 elsets_list[elset_name] = elements;
3889 else if (line.compare(0, 5,
"*NSET") == 0)
3892 while (std::getline(input_stream, line))
3898 else if (line.compare(0, 14,
"*SOLID SECTION") == 0)
3901 const std::string elset_key =
"ELSET=";
3902 const std::size_t elset_start =
3903 line.find(
"ELSET=") + elset_key.size();
3904 const std::size_t elset_end = line.find(
',', elset_start + 1);
3905 const std::string elset_name =
3906 line.substr(elset_start, elset_end - elset_start);
3911 const std::string material_key =
"MATERIAL=";
3912 const std::size_t last_equal =
3913 line.find(
"MATERIAL=") + material_key.size();
3914 const std::size_t material_id_start = line.find(
'-', last_equal);
3917 line.substr(material_id_start + 1),
3921 const std::vector<int> &elset_cells = elsets_list[elset_name];
3922 for (
const int elset_cell : elset_cells)
3924 const int cell_id = elset_cell - 1;
3932 template <
int dim,
int spacedim>
3934 Abaqus_to_UCD<dim, spacedim>::get_global_node_numbers(
3935 const int face_cell_no,
3936 const int face_cell_face_no)
const
3946 if (face_cell_face_no == 1)
3948 quad_node_list[0] = cell_list[face_cell_no - 1][1];
3949 quad_node_list[1] = cell_list[face_cell_no - 1][2];
3951 else if (face_cell_face_no == 2)
3953 quad_node_list[0] = cell_list[face_cell_no - 1][2];
3954 quad_node_list[1] = cell_list[face_cell_no - 1][3];
3956 else if (face_cell_face_no == 3)
3958 quad_node_list[0] = cell_list[face_cell_no - 1][3];
3959 quad_node_list[1] = cell_list[face_cell_no - 1][4];
3961 else if (face_cell_face_no == 4)
3963 quad_node_list[0] = cell_list[face_cell_no - 1][4];
3964 quad_node_list[1] = cell_list[face_cell_no - 1][1];
3974 if (face_cell_face_no == 1)
3976 quad_node_list[0] = cell_list[face_cell_no - 1][1];
3977 quad_node_list[1] = cell_list[face_cell_no - 1][4];
3978 quad_node_list[2] = cell_list[face_cell_no - 1][3];
3979 quad_node_list[3] = cell_list[face_cell_no - 1][2];
3981 else if (face_cell_face_no == 2)
3983 quad_node_list[0] = cell_list[face_cell_no - 1][5];
3984 quad_node_list[1] = cell_list[face_cell_no - 1][8];
3985 quad_node_list[2] = cell_list[face_cell_no - 1][7];
3986 quad_node_list[3] = cell_list[face_cell_no - 1][6];
3988 else if (face_cell_face_no == 3)
3990 quad_node_list[0] = cell_list[face_cell_no - 1][1];
3991 quad_node_list[1] = cell_list[face_cell_no - 1][2];
3992 quad_node_list[2] = cell_list[face_cell_no - 1][6];
3993 quad_node_list[3] = cell_list[face_cell_no - 1][5];
3995 else if (face_cell_face_no == 4)
3997 quad_node_list[0] = cell_list[face_cell_no - 1][2];
3998 quad_node_list[1] = cell_list[face_cell_no - 1][3];
3999 quad_node_list[2] = cell_list[face_cell_no - 1][7];
4000 quad_node_list[3] = cell_list[face_cell_no - 1][6];
4002 else if (face_cell_face_no == 5)
4004 quad_node_list[0] = cell_list[face_cell_no - 1][3];
4005 quad_node_list[1] = cell_list[face_cell_no - 1][4];
4006 quad_node_list[2] = cell_list[face_cell_no - 1][8];
4007 quad_node_list[3] = cell_list[face_cell_no - 1][7];
4009 else if (face_cell_face_no == 6)
4011 quad_node_list[0] = cell_list[face_cell_no - 1][1];
4012 quad_node_list[1] = cell_list[face_cell_no - 1][5];
4013 quad_node_list[2] = cell_list[face_cell_no - 1][8];
4014 quad_node_list[3] = cell_list[face_cell_no - 1][4];
4027 return quad_node_list;
4030 template <
int dim,
int spacedim>
4032 Abaqus_to_UCD<dim, spacedim>::write_out_avs_ucd(std::ostream &output)
const
4041 const boost::io::ios_base_all_saver formatting_saver(output);
4045 output <<
"# Abaqus to UCD mesh conversion" << std::endl;
4046 output <<
"# Mesh type: AVS UCD" << std::endl;
4071 output << node_list.size() <<
"\t" << (cell_list.size() + face_list.size())
4072 <<
"\t0\t0\t0" << std::endl;
4075 output.precision(8);
4079 for (
const auto &node : node_list)
4082 output << node[0] <<
"\t";
4085 output.setf(std::ios::scientific, std::ios::floatfield);
4086 for (
unsigned int jj = 1; jj < spacedim + 1; ++jj)
4089 if (std::abs(node[jj]) > tolerance)
4090 output <<
static_cast<double>(node[jj]) <<
"\t";
4092 output << 0.0 <<
"\t";
4095 output << 0.0 <<
"\t";
4097 output << std::endl;
4098 output.unsetf(std::ios::floatfield);
4102 for (
unsigned int ii = 0; ii < cell_list.size(); ++ii)
4104 output << ii + 1 <<
"\t" << cell_list[ii][0] <<
"\t"
4105 << (dim == 2 ?
"quad" :
"hex") <<
"\t";
4106 for (
unsigned int jj = 1; jj < GeometryInfo<dim>::vertices_per_cell + 1;
4108 output << cell_list[ii][jj] <<
"\t";
4110 output << std::endl;
4114 for (
unsigned int ii = 0; ii < face_list.size(); ++ii)
4116 output << ii + 1 <<
"\t" << face_list[ii][0] <<
"\t"
4117 << (dim == 2 ?
"line" :
"quad") <<
"\t";
4118 for (
unsigned int jj = 1; jj < GeometryInfo<dim>::vertices_per_face + 1;
4120 output << face_list[ii][jj] <<
"\t";
4122 output << std::endl;
4129 #include "grid_in.inst"