17 #include <deal.II/base/path_search.h> 18 #include <deal.II/base/utilities.h> 19 #include <deal.II/base/exceptions.h> 21 #include <deal.II/grid/grid_in.h> 22 #include <deal.II/grid/tria.h> 23 #include <deal.II/grid/grid_reordering.h> 24 #include <deal.II/grid/grid_tools.h> 33 #ifdef DEAL_II_WITH_NETCDF 34 #include <netcdfcpp.h> 38 DEAL_II_NAMESPACE_OPEN
51 template <
int spacedim>
53 assign_1d_boundary_ids (
const std::map<unsigned int, types::boundary_id> &boundary_ids,
56 if (boundary_ids.size() > 0)
59 cell != triangulation.
end(); ++cell)
60 for (
unsigned int f=0; f<GeometryInfo<1>::faces_per_cell; ++f)
61 if (boundary_ids.find(cell->vertex_index(f)) != boundary_ids.end())
64 ExcMessage (
"You are trying to prescribe boundary ids on the face " 65 "of a 1d cell (i.e., on a vertex), but this face is not actually at " 66 "the boundary of the mesh. This is not allowed."));
67 cell->face(f)->set_boundary_id (boundary_ids.find(cell->vertex_index(f))->second);
72 template <
int dim,
int spacedim>
74 assign_1d_boundary_ids (
const std::map<unsigned int, types::boundary_id> &,
83 template <
int dim,
int spacedim>
85 tria(0, typeid(*this).name()), default_format(ucd)
89 template <
int dim,
int spacedim>
97 template<
int dim,
int spacedim>
109 text[0] =
"# vtk DataFile Version 3.0";
112 text[3] =
"DATASET UNSTRUCTURED_GRID";
114 for (
unsigned int i = 0; i < 4; ++i)
119 ExcMessage(std::string(
"While reading VTK file, failed to find a header line with text <") +
126 std::vector< Point<spacedim> > vertices;
127 std::vector< CellData<dim> > cells;
136 if (keyword ==
"POINTS")
138 unsigned int n_vertices;
141 in.ignore(256,
'\n');
143 for (
unsigned int vertex = 0; vertex < n_vertices; ++vertex)
147 in >> x(0) >> x(1) >> x(2);
150 for (
unsigned int d=0; d<spacedim; ++d)
151 vertices.back()(d) = x(d);
157 ExcMessage (
"While reading VTK file, failed to find POINTS section"));
161 std::string checkline;
163 in.ignore(256,
'\n');
165 getline(in,checkline);
166 if (checkline.compare(
"") != 0)
175 if (keyword ==
"CELLS")
177 unsigned int n_geometric_objects;
178 in >> n_geometric_objects;
183 for (
unsigned int count = 0; count < n_geometric_objects; count++)
199 for (
unsigned int j = 0; j < type; j++)
200 in >> cells.back().vertices[j];
202 cells.back().material_id = 0;
209 for (
unsigned int j = 0; j < type; j++)
217 ExcMessage (
"While reading VTK file, unknown file type encountered"));
223 for (
unsigned int count = 0; count < n_geometric_objects; count++)
237 for (
unsigned int j = 0; j < type; j++)
238 in >> cells.back().vertices[j];
240 cells.back().material_id = 0;
249 for (
unsigned int j = 0; j < type; j++)
259 ExcMessage (
"While reading VTK file, unknown file type encountered"));
264 ExcMessage (
"While reading VTK file, failed to find CELLS section"));
270 if (keyword ==
"CELL_TYPES")
272 in.ignore(256,
'\n');
274 while (!in.fail() && !in.eof())
277 if (keyword !=
"12" && keyword !=
"9")
286 if (keyword ==
"CELL_DATA")
297 ExcMessage (
"The VTK reader found a CELL_DATA statement " 298 "that lists a total of " 300 " cell data objects, but this needs to " 301 "equal the number of cells (which is " 303 ") plus the number of quads (" 305 " in 3d or the number of lines (" 311 std::string textnew[2];
312 textnew[0] =
"SCALARS MaterialID double";
313 textnew[1] =
"LOOKUP_TABLE default";
315 in.ignore(256,
'\n');
317 for (
unsigned int i = 0; i < 2; i++)
319 getline(in, linenew);
321 if (linenew.size() > textnew[0].size())
322 linenew.resize(textnew[0].size());
325 ExcMessage (std::string(
"While reading VTK file, failed to find <") +
326 textnew[i] +
"> section"));
333 for (
unsigned int i = 0; i < cells.size(); i++)
337 cells[i].material_id = id;
342 for (
unsigned int i = 0; i < subcelldata.
boundary_quads.size(); i++)
351 for (
unsigned int i = 0; i < subcelldata.
boundary_lines.size(); i++)
371 tria->create_triangulation_compatibility(vertices,
379 ExcMessage (
"While reading VTK file, failed to find CELLS section"));
384 template<
int dim,
int spacedim>
391 skip_comment_lines(in,
'#');
402 AssertThrow(tmp == 2411, ExcUnknownSectionType(tmp));
404 std::vector< Point<spacedim> > vertices;
405 std::map<int, int> vertex_indices;
422 in >> dummy >> dummy >> dummy;
425 in >> x[0] >> x[1] >> x[2];
429 for (
unsigned int d = 0; d < spacedim; d++)
430 vertices.back()(d) = x[d];
432 vertex_indices[no] = no_vertex;
444 AssertThrow(tmp == 2412, ExcUnknownSectionType(tmp));
446 std::vector< CellData<dim> > cells;
449 std::map<int, int> cell_indices;
450 std::map<int, int> line_indices;
451 std::map<int, int> quad_indices;
470 in >> type >> dummy >> dummy >> dummy >> dummy;
472 AssertThrow((type == 11)||(type == 44)||(type == 94)||(type == 115), ExcUnknownElementType(type));
474 if ( (((type == 44)||(type == 94))&&(dim == 2)) || ((type == 115)&&(dim == 3)) )
479 for (
unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_cell; v++)
480 in >> cells.back().vertices[v];
482 cells.back().material_id = 0;
484 for (
unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_cell; v++)
485 cells.back().vertices[v] = vertex_indices[cells.back().vertices[v]];
487 cell_indices[no] = no_cell;
491 else if ( ((type == 11)&&(dim == 2)) || ((type == 11)&&(dim == 3)) )
494 in >> dummy >> dummy >> dummy;
499 for (
unsigned int v = 0; v < 2; v++)
504 for (
unsigned int v = 0; v < 2; v++)
507 line_indices[no] = no_line;
511 else if ( ((type == 44)||(type == 94)) && (dim == 3) )
516 for (
unsigned int v = 0; v < 4; v++)
521 for (
unsigned int v = 0; v < 4; v++)
524 quad_indices[no] = no_quad;
532 +
"> when running in dim=" 551 AssertThrow((tmp == 2467)||(tmp == 2477), ExcUnknownSectionType(tmp));
567 in >> dummy >> dummy >> dummy >> dummy >> dummy >> dummy >> n_entities;
572 const unsigned int n_lines = (n_entities%2 == 0)?(n_entities/2):((n_entities+1)/2);
574 for (
unsigned int line = 0; line < n_lines; line++)
576 unsigned int n_fragments;
578 if (line == n_lines-1)
579 n_fragments = (n_entities%2 == 0)?(2):(1);
583 for (
unsigned int no_fragment = 0; no_fragment < n_fragments; no_fragment++)
586 in >> dummy >> no >> dummy >> dummy;
588 if ( cell_indices.count(no) > 0 )
589 cells[cell_indices[no]].material_id =
id;
591 if ( line_indices.count(no) > 0 )
594 if ( quad_indices.count(no) > 0 )
613 tria->create_triangulation_compatibility(vertices,
620 template <
int dim,
int spacedim>
622 const bool apply_all_indicators_to_manifolds)
628 skip_comment_lines (in,
'#');
631 unsigned int n_vertices;
632 unsigned int n_cells;
643 std::vector<Point<spacedim> > vertices (n_vertices);
647 std::map<int,int> vertex_indices;
649 for (
unsigned int vertex=0; vertex<n_vertices; ++vertex)
657 >> x[0] >> x[1] >> x[2];
660 for (
unsigned int d=0; d<spacedim; ++d)
661 vertices[vertex](d) = x[d];
665 vertex_indices[vertex_number] = vertex;
669 std::vector<CellData<dim> > cells;
672 for (
unsigned int cell=0; cell<n_cells; ++cell)
681 std::string cell_type;
685 unsigned int material_id;
691 if (((cell_type ==
"line") && (dim == 1)) ||
692 ((cell_type ==
"quad") && (dim == 2)) ||
693 ((cell_type ==
"hex" ) && (dim == 3)))
698 for (
unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
699 in >> cells.back().vertices[i];
702 Assert(material_id<= std::numeric_limits<types::material_id>::max(),
703 ExcIndexRange(material_id,0,std::numeric_limits<types::material_id>::max()));
712 for (
unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
713 if (vertex_indices.find (cells.back().vertices[i]) != vertex_indices.end())
715 cells.back().vertices[i] = vertex_indices[cells.back().vertices[i]];
725 else if ((cell_type ==
"line") && ((dim == 2) || (dim == 3)))
733 Assert(material_id<= std::numeric_limits<types::boundary_id>::max(),
734 ExcIndexRange(material_id,0,std::numeric_limits<types::boundary_id>::max()));
739 if (apply_all_indicators_to_manifolds)
748 for (
unsigned int i=0; i<2; ++i)
749 if (vertex_indices.find (subcelldata.
boundary_lines.back().vertices[i]) !=
750 vertex_indices.end())
764 else if ((cell_type ==
"quad") && (dim == 3))
774 Assert(material_id<= std::numeric_limits<types::boundary_id>::max(),
775 ExcIndexRange(material_id,0,std::numeric_limits<types::boundary_id>::max()));
780 if (apply_all_indicators_to_manifolds)
789 for (
unsigned int i=0; i<4; ++i)
790 if (vertex_indices.find (subcelldata.
boundary_quads.back().vertices[i]) !=
791 vertex_indices.end())
808 AssertThrow (
false, ExcUnknownIdentifier(cell_type));
823 tria->create_triangulation_compatibility (vertices, cells, subcelldata);
834 void read_in_abaqus (std::istream &in);
835 void write_out_avs_ucd (std::ostream &out)
const;
838 const double tolerance;
840 std::vector<double> get_global_node_numbers (
const int face_cell_no,
841 const int face_cell_face_no)
const;
844 std::vector< std::vector<double> > node_list;
846 std::vector< std::vector<double> > cell_list;
848 std::vector< std::vector<double> > face_list;
850 std::map< std::string, std::vector<int> > elsets_list;
854 template <
int dim,
int spacedim>
856 const bool apply_all_indicators_to_manifolds)
864 Abaqus_to_UCD<dim> abaqus_to_ucd;
865 abaqus_to_ucd.read_in_abaqus(in);
867 std::stringstream in_ucd;
868 abaqus_to_ucd.write_out_avs_ucd(in_ucd);
876 read_ucd(in_ucd, apply_all_indicators_to_manifolds);
878 catch (std::exception &exc)
881 <<
"Exception on processing internal UCD data: " << std::endl
885 AssertThrow(
false,
ExcMessage(
"Internal conversion from ABAQUS file to UCD format was unsuccessful. " 886 "More information is provided in an error message printed above. " 887 "Are you sure that your ABAQUS mesh file conforms with the requirements " 888 "listed in the documentation?"));
892 AssertThrow(
false,
ExcMessage(
"Internal conversion from ABAQUS file to UCD format was unsuccessful. " 893 "Are you sure that your ABAQUS mesh file conforms with the requirements " 894 "listed in the documentation?"));
899 template <
int dim,
int spacedim>
908 skip_comment_lines (in,
'#');
915 ExcInvalidDBMESHInput(line));
917 skip_empty_lines (in);
921 AssertThrow (line==
"Dimension", ExcInvalidDBMESHInput(line));
922 unsigned int dimension;
924 AssertThrow (dimension == dim, ExcDBMESHWrongDimension(dimension));
925 skip_empty_lines (in);
938 while (getline(in,line), line.find(
"# END")==std::string::npos)
940 skip_empty_lines (in);
945 AssertThrow (line==
"Vertices", ExcInvalidDBMESHInput(line));
947 unsigned int n_vertices;
951 std::vector<Point<spacedim> > vertices (n_vertices);
952 for (
unsigned int vertex=0; vertex<n_vertices; ++vertex)
955 for (
unsigned int d=0; d<dim; ++d)
956 in >> vertices[vertex][d];
962 skip_empty_lines(in);
968 AssertThrow (line==
"Edges", ExcInvalidDBMESHInput(line));
970 unsigned int n_edges;
972 for (
unsigned int edge=0; edge<n_edges; ++edge)
975 in >> dummy >> dummy;
981 skip_empty_lines(in);
990 AssertThrow (line==
"CrackedEdges", ExcInvalidDBMESHInput(line));
993 for (
unsigned int edge=0; edge<n_edges; ++edge)
996 in >> dummy >> dummy;
1002 skip_empty_lines(in);
1008 AssertThrow (line==
"Quadrilaterals", ExcInvalidDBMESHInput(line));
1010 std::vector<CellData<dim> > cells;
1012 unsigned int n_cells;
1014 for (
unsigned int cell=0; cell<n_cells; ++cell)
1019 for (
unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
1021 in >> cells.back().vertices[i];
1025 (static_cast<unsigned int>(cells.back().vertices[i]) <= vertices.size()),
1028 --cells.back().vertices[i];
1036 skip_empty_lines(in);
1044 while (getline(in,line), ((line.find(
"End")==std::string::npos) && (in)))
1060 tria->create_triangulation_compatibility (vertices, cells, subcelldata);
1065 template <
int dim,
int spacedim>
1084 unsigned int n_vertices;
1085 unsigned int n_cells;
1095 for (
unsigned int i=0; i<8; ++i)
1099 std::vector<CellData<2> > cells (n_cells);
1102 for (
unsigned int cell=0; cell<n_cells; ++cell)
1113 for (
unsigned int i=0; i<4; ++i)
1114 in >> cells[cell].vertices[i];
1120 std::vector<Point<2> > vertices (n_vertices);
1121 for (
unsigned int vertex=0; vertex<n_vertices; ++vertex)
1126 in >> x[0] >> x[1] >> x[2];
1129 for (
unsigned int d=0;
d<2; ++
d)
1130 vertices[vertex](d) = x[
d];
1139 tria->create_triangulation_compatibility (vertices, cells, subcelldata);
1150 static const unsigned int xda_to_dealII_map[] = {0,1,5,4,3,2,6,7};
1157 unsigned int n_vertices;
1158 unsigned int n_cells;
1168 for (
unsigned int i=0; i<8; ++i)
1172 std::vector<CellData<3> > cells (n_cells);
1175 for (
unsigned int cell=0; cell<n_cells; ++cell)
1186 unsigned int xda_ordered_nodes[8];
1188 for (
unsigned int i=0; i<8; ++i)
1189 in >> xda_ordered_nodes[i];
1191 for (
unsigned int i=0; i<8; i++)
1192 cells[cell].vertices[i] = xda_ordered_nodes[xda_to_dealII_map[i]];
1198 std::vector<Point<3> > vertices (n_vertices);
1199 for (
unsigned int vertex=0; vertex<n_vertices; ++vertex)
1204 in >> x[0] >> x[1] >> x[2];
1207 for (
unsigned int d=0;
d<3; ++
d)
1208 vertices[vertex](d) = x[
d];
1217 tria->create_triangulation_compatibility (vertices, cells, subcelldata);
1223 template <
int dim,
int spacedim>
1229 unsigned int n_vertices;
1230 unsigned int n_cells;
1237 unsigned int gmsh_file_format = 0;
1238 if (line ==
"@f$NOD")
1239 gmsh_file_format = 1;
1240 else if (line ==
"@f$MeshFormat")
1241 gmsh_file_format = 2;
1248 if (gmsh_file_format == 2)
1251 unsigned int file_type, data_size;
1253 in >> version >> file_type >> data_size;
1255 Assert ( (version >= 2.0) &&
1267 ExcInvalidGMSHInput(line));
1272 if (line ==
"@f$PhysicalNames")
1278 while (line !=
"@f$EndPhysicalNames");
1286 ExcInvalidGMSHInput(line));
1291 std::vector<Point<spacedim> > vertices (n_vertices);
1295 std::map<int,int> vertex_indices;
1297 for (
unsigned int vertex=0; vertex<n_vertices; ++vertex)
1304 >> x[0] >> x[1] >> x[2];
1306 for (
unsigned int d=0; d<spacedim; ++d)
1307 vertices[vertex](d) = x[d];
1309 vertex_indices[vertex_number] = vertex;
1314 static const std::string end_nodes_marker[] = {
"@f$ENDNOD",
"@f$EndNodes" };
1315 AssertThrow (line==end_nodes_marker[gmsh_file_format-1],
1316 ExcInvalidGMSHInput(line));
1320 static const std::string begin_elements_marker[] = {
"@f$ELM",
"@f$Elements" };
1321 AssertThrow (line==begin_elements_marker[gmsh_file_format-1],
1322 ExcInvalidGMSHInput(line));
1329 std::vector<CellData<dim> > cells;
1331 std::map<unsigned int, types::boundary_id> boundary_ids_1d;
1333 for (
unsigned int cell=0; cell<n_cells; ++cell)
1342 unsigned int cell_type;
1343 unsigned int material_id;
1344 unsigned int nod_num;
1359 unsigned int elm_number;
1363 switch (gmsh_file_format)
1378 unsigned int n_tags;
1385 for (
unsigned int i=1; i<n_tags; ++i)
1413 if (((cell_type == 1) && (dim == 1)) ||
1414 ((cell_type == 3) && (dim == 2)) ||
1415 ((cell_type == 5) && (dim == 3)))
1419 ExcMessage (
"Number of nodes does not coincide with the " 1420 "number required for this object"));
1424 for (
unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
1425 in >> cells.back().vertices[i];
1428 Assert(material_id<= std::numeric_limits<types::material_id>::max(),
1429 ExcIndexRange(material_id,0,std::numeric_limits<types::material_id>::max()));
1438 for (
unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
1440 AssertThrow (vertex_indices.find (cells.back().vertices[i]) !=
1441 vertex_indices.end(),
1442 ExcInvalidVertexIndexGmsh(cell, elm_number,
1443 cells.back().vertices[i]));
1446 cells.back().vertices[i] = vertex_indices[cells.back().vertices[i]];
1449 else if ((cell_type == 1) && ((dim == 2) || (dim == 3)))
1457 Assert(material_id<= std::numeric_limits<types::boundary_id>::max(),
1458 ExcIndexRange(material_id,0,std::numeric_limits<types::boundary_id>::max()));
1468 for (
unsigned int i=0; i<2; ++i)
1469 if (vertex_indices.find (subcelldata.
boundary_lines.back().vertices[i]) !=
1470 vertex_indices.end())
1484 else if ((cell_type == 3) && (dim == 3))
1494 Assert(material_id<= std::numeric_limits<types::boundary_id>::max(),
1495 ExcIndexRange(material_id,0,std::numeric_limits<types::boundary_id>::max()));
1505 for (
unsigned int i=0; i<4; ++i)
1506 if (vertex_indices.find (subcelldata.
boundary_quads.back().vertices[i]) !=
1507 vertex_indices.end())
1522 else if (cell_type == 15)
1525 unsigned int node_index;
1526 switch (gmsh_file_format)
1530 for (
unsigned int i=0; i<nod_num; ++i)
1544 boundary_ids_1d[vertex_indices[node_index]] = material_id;
1555 ExcMessage(
"Found triangles while reading a file " 1556 "in gmsh format. deal.II does not " 1557 "support triangles"));
1559 ExcMessage(
"Found tetrahedra while reading a file " 1560 "in gmsh format. deal.II does not " 1561 "support tetrahedra"));
1563 AssertThrow (
false, ExcGmshUnsupportedGeometry(cell_type));
1569 static const std::string end_elements_marker[] = {
"@f$ENDELM",
"@f$EndElements" };
1570 AssertThrow (line==end_elements_marker[gmsh_file_format-1],
1571 ExcInvalidGMSHInput(line));
1580 AssertThrow(cells.size() > 0, ExcGmshNoCellInformation());
1589 tria->create_triangulation_compatibility (vertices, cells, subcelldata);
1594 assign_1d_boundary_ids (boundary_ids_1d, *tria);
1627 #ifndef DEAL_II_WITH_NETCDF 1631 const unsigned int dim=2;
1632 const unsigned int spacedim=2;
1633 const bool output=
false;
1657 const unsigned int coord=1;
1663 const unsigned int x2d=0;
1664 const unsigned int y2d=2;
1673 NcFile nc (filename.c_str());
1677 NcDim *elements_dim=nc.get_dim(
"no_of_elements");
1679 const unsigned int n_cells=elements_dim->size();
1683 NcDim *marker_dim=nc.get_dim(
"no_of_markers");
1685 const unsigned int n_markers=marker_dim->size();
1687 NcVar *marker_var=nc.get_var(
"marker");
1691 marker_var->get_dim(0)->size())==n_markers,
ExcIO());
1693 std::vector<int> marker(n_markers);
1696 marker_var->get(&*marker.begin(), n_markers);
1700 std::cout <<
"n_cell=" << n_cells << std::endl;
1701 std::cout <<
"marker: ";
1702 for (
unsigned int i=0; i<n_markers; ++i)
1703 std::cout << marker[i] <<
" ";
1704 std::cout << std::endl;
1710 NcDim *bquads_dim=nc.get_dim(
"no_of_surfacequadrilaterals");
1712 const unsigned int n_bquads=bquads_dim->size();
1714 NcVar *bmarker_var=nc.get_var(
"boundarymarker_of_surfaces");
1718 bmarker_var->get_dim(0)->size())==n_bquads,
ExcIO());
1720 std::vector<int> bmarker(n_bquads);
1721 bmarker_var->get(&*bmarker.begin(), n_bquads);
1726 std::map<int, unsigned int> n_bquads_per_bmarker;
1727 for (
unsigned int i=0; i<n_markers; ++i)
1731 AssertThrow(n_bquads_per_bmarker.find(marker[i])==
1732 n_bquads_per_bmarker.end(),
ExcIO());
1734 n_bquads_per_bmarker[marker[i]]=
1735 count(bmarker.begin(), bmarker.end(), marker[i]);
1751 std::cout <<
"n_bquads_per_bmarker: " << std::endl;
1752 std::map<int, unsigned int>::const_iterator
1753 iter=n_bquads_per_bmarker.begin();
1754 for (; iter!=n_bquads_per_bmarker.end(); ++iter)
1755 std::cout <<
" n_bquads_per_bmarker[" << iter->first
1756 <<
"]=" << iter->second << std::endl;
1763 NcDim *quad_vertices_dim=nc.get_dim(
"points_per_surfacequadrilateral");
1765 const unsigned int vertices_per_quad=quad_vertices_dim->size();
1768 NcVar *vertex_indices_var=nc.get_var(
"points_of_surfacequadrilaterals");
1772 vertex_indices_var->get_dim(0)->size())==n_bquads,
ExcIO());
1774 vertex_indices_var->get_dim(1)->size())==vertices_per_quad,
ExcIO());
1776 std::vector<int> vertex_indices(n_bquads*vertices_per_quad);
1777 vertex_indices_var->get(&*vertex_indices.begin(), n_bquads, vertices_per_quad);
1779 for (
unsigned int i=0; i<vertex_indices.size(); ++i)
1784 std::cout <<
"vertex_indices:" << std::endl;
1785 for (
unsigned int i=0, v=0; i<n_bquads; ++i)
1787 for (
unsigned int j=0; j<vertices_per_quad; ++j)
1788 std::cout << vertex_indices[v++] <<
" ";
1789 std::cout << std::endl;
1797 NcDim *vertices_dim=nc.get_dim(
"no_of_points");
1799 const unsigned int n_vertices=vertices_dim->size();
1801 std::cout <<
"n_vertices=" << n_vertices << std::endl;
1803 NcVar *points_xc=nc.get_var(
"points_xc");
1804 NcVar *points_yc=nc.get_var(
"points_yc");
1805 NcVar *points_zc=nc.get_var(
"points_zc");
1813 static_cast<int>(n_vertices),
ExcIO());
1815 static_cast<int>(n_vertices),
ExcIO());
1817 static_cast<int>(n_vertices),
ExcIO());
1818 std::vector<std::vector<double> > point_values(
1819 3, std::vector<double> (n_vertices));
1820 points_xc->get(&*point_values[0].begin(), n_vertices);
1821 points_yc->get(&*point_values[1].begin(), n_vertices);
1822 points_zc->get(&*point_values[2].begin(), n_vertices);
1825 std::vector<Point<spacedim> > vertices (n_vertices);
1826 for (
unsigned int i=0; i<n_vertices; ++i)
1828 vertices[i](0)=point_values[x2d][i];
1829 vertices[i](1)=point_values[y2d][i];
1835 std::map<int, bool> zero_plane_markers;
1836 for (
unsigned int quad=0; quad<n_bquads; ++quad)
1838 bool zero_plane=
true;
1839 for (
unsigned int i=0; i<vertices_per_quad; ++i)
1840 if (point_values[coord][vertex_indices[quad*vertices_per_quad+i]]!=0)
1847 zero_plane_markers[bmarker[quad]]=
true;
1849 unsigned int sum_of_zero_plane_cells=0;
1850 for (std::map<int, bool>::const_iterator iter=zero_plane_markers.begin();
1851 iter != zero_plane_markers.end(); ++iter)
1853 sum_of_zero_plane_cells+=n_bquads_per_bmarker[iter->first];
1855 std::cout <<
"bmarker=" << iter->first << std::endl;
1862 std::vector<CellData<dim> > cells(n_cells);
1863 for (
unsigned int quad=0, cell=0; quad<n_bquads; ++quad)
1865 bool zero_plane=
false;
1866 for (std::map<int, bool>::const_iterator iter=zero_plane_markers.begin();
1867 iter != zero_plane_markers.end(); ++iter)
1868 if (bmarker[quad]==iter->first)
1876 for (
unsigned int i=0; i<vertices_per_quad; ++i)
1878 Assert(point_values[coord][vertex_indices[
1880 cells[cell].vertices[i]=vertex_indices[quad*vertices_per_quad+i];
1889 tria->create_triangulation_compatibility (vertices, cells, subcelldata);
1897 #ifndef DEAL_II_WITH_NETCDF 1904 const unsigned int dim=3;
1905 const unsigned int spacedim=3;
1906 const bool output=
false;
1912 NcFile nc (filename.c_str());
1916 NcDim *elements_dim=nc.get_dim(
"no_of_elements");
1918 const unsigned int n_cells=elements_dim->size();
1920 std::cout <<
"n_cell=" << n_cells << std::endl;
1922 NcDim *hexes_dim=nc.get_dim(
"no_of_hexaeders");
1924 const unsigned int n_hexes=hexes_dim->size();
1926 ExcMessage(
"deal.II can handle purely hexaedral grids, only."));
1932 NcDim *hex_vertices_dim=nc.get_dim(
"points_per_hexaeder");
1934 const unsigned int vertices_per_hex=hex_vertices_dim->size();
1937 NcVar *vertex_indices_var=nc.get_var(
"points_of_hexaeders");
1941 vertex_indices_var->get_dim(0)->size())==n_cells,
ExcIO());
1943 vertex_indices_var->get_dim(1)->size())==vertices_per_hex,
ExcIO());
1945 std::vector<int> vertex_indices(n_cells*vertices_per_hex);
1948 vertex_indices_var->get(&*vertex_indices.begin(), n_cells, vertices_per_hex);
1950 for (
unsigned int i=0; i<vertex_indices.size(); ++i)
1955 std::cout <<
"vertex_indices:" << std::endl;
1956 for (
unsigned int cell=0, v=0; cell<n_cells; ++cell)
1958 for (
unsigned int i=0; i<vertices_per_hex; ++i)
1959 std::cout << vertex_indices[v++] <<
" ";
1960 std::cout << std::endl;
1968 NcDim *vertices_dim=nc.get_dim(
"no_of_points");
1970 const unsigned int n_vertices=vertices_dim->size();
1972 std::cout <<
"n_vertices=" << n_vertices << std::endl;
1974 NcVar *points_xc=nc.get_var(
"points_xc");
1975 NcVar *points_yc=nc.get_var(
"points_yc");
1976 NcVar *points_zc=nc.get_var(
"points_zc");
1984 static_cast<int>(n_vertices),
ExcIO());
1986 static_cast<int>(n_vertices),
ExcIO());
1988 static_cast<int>(n_vertices),
ExcIO());
1989 std::vector<std::vector<double> > point_values(
1990 3, std::vector<double> (n_vertices));
1992 const bool switch_y_z=
false;
1993 points_xc->get(&*point_values[0].begin(), n_vertices);
1996 points_yc->get(&*point_values[2].begin(), n_vertices);
1997 points_zc->get(&*point_values[1].begin(), n_vertices);
2001 points_yc->get(&*point_values[1].begin(), n_vertices);
2002 points_zc->get(&*point_values[2].begin(), n_vertices);
2006 std::vector<Point<spacedim> > vertices (n_vertices);
2007 for (
unsigned int i=0; i<n_vertices; ++i)
2009 vertices[i](0)=point_values[0][i];
2010 vertices[i](1)=point_values[1][i];
2011 vertices[i](2)=point_values[2][i];
2015 std::vector<CellData<dim> > cells(n_cells);
2016 for (
unsigned int cell=0; cell<n_cells; ++cell)
2017 for (
unsigned int i=0; i<vertices_per_hex; ++i)
2018 cells[cell].vertices[i]=vertex_indices[cell*vertices_per_hex+i];
2029 NcDim *quad_vertices_dim=nc.get_dim(
"points_per_surfacequadrilateral");
2031 const unsigned int vertices_per_quad=quad_vertices_dim->size();
2034 NcVar *bvertex_indices_var=nc.get_var(
"points_of_surfacequadrilaterals");
2037 const unsigned int n_bquads=bvertex_indices_var->get_dim(0)->size();
2039 bvertex_indices_var->get_dim(1)->size())==
2042 std::vector<int> bvertex_indices(n_bquads*vertices_per_quad);
2043 bvertex_indices_var->get(&*bvertex_indices.begin(), n_bquads, vertices_per_quad);
2047 std::cout <<
"bquads: ";
2048 for (
unsigned int i=0; i<n_bquads; ++i)
2050 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_face; ++v)
2051 std::cout << bvertex_indices[
2053 std::cout << std::endl;
2060 NcDim *bquads_dim=nc.get_dim(
"no_of_surfacequadrilaterals");
2063 bquads_dim->size())==n_bquads,
ExcIO());
2065 NcVar *bmarker_var=nc.get_var(
"boundarymarker_of_surfaces");
2069 bmarker_var->get_dim(0)->size())==n_bquads,
ExcIO());
2071 std::vector<int> bmarker(n_bquads);
2072 bmarker_var->get(&*bmarker.begin(), n_bquads);
2078 for (
unsigned int i=0; i<bmarker.size(); ++i)
2085 for (
unsigned int i=0; i<n_bquads; ++i)
2087 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_face; ++v)
2091 = static_cast<types::boundary_id>(bmarker[i]);
2097 tria->create_triangulation_compatibility (vertices, cells, subcelldata);
2102 template <
int dim,
int spacedim>
2104 std::vector<unsigned int> &tecplot2deal,
2105 unsigned int &n_vars,
2106 unsigned int &n_vertices,
2107 unsigned int &n_cells,
2108 std::vector<unsigned int> &IJK,
2131 std::transform(header.begin(),header.end(),header.begin(),::toupper);
2135 std::replace(header.begin(),header.end(),
'\t',
' ');
2136 std::replace(header.begin(),header.end(),
',',
' ');
2137 std::replace(header.begin(),header.end(),
'\n',
' ');
2141 std::string::size_type pos=header.find(
"=");
2143 while (pos!=static_cast<std::string::size_type>(std::string::npos))
2144 if (header[pos+1]==
' ')
2145 header.erase(pos+1,1);
2146 else if (header[pos-1]==
' ')
2148 header.erase(pos-1,1);
2152 pos=header.find(
"=",++pos);
2158 for (
unsigned int i=0; i<entries.size(); ++i)
2170 while (entries[i][0]==
'"')
2172 if (entries[i]==
"\"X\"")
2173 tecplot2deal[0]=n_vars;
2174 else if (entries[i]==
"\"Y\"")
2180 tecplot2deal[1]=n_vars;
2182 else if (entries[i]==
"\"Z\"")
2188 tecplot2deal[2]=n_vars;
2198 ExcMessage(
"Tecplot file must contain at least one variable for each dimension"));
2199 for (
unsigned int d=1; d<dim; ++d)
2201 ExcMessage(
"Tecplot file must contain at least one variable for each dimension."));
2255 ExcMessage(
"Parameter 'J=' found in tecplot, although this is only possible for dimensions greater than 1."));
2261 ExcMessage(
"Parameter 'K=' found in tecplot, although this is only possible for dimensions greater than 2."));
2276 for (
unsigned int d=0; d<dim; ++d)
2279 ExcMessage(
"Tecplot file does not contain a complete and consistent set of parameters"));
2281 n_cells*=(IJK[d]-1);
2287 ExcMessage(
"Tecplot file does not contain a complete and consistent set of parameters"));
2293 n_cells=*std::max_element(IJK.begin(),IJK.end());
2295 ExcMessage(
"Tecplot file does not contain a complete and consistent set of parameters"));
2305 const unsigned int dim=2;
2306 const unsigned int spacedim=2;
2311 skip_comment_lines (in,
'#');
2314 std::string line, header;
2321 std::string letters =
"abcdfghijklmnopqrstuvwxyzABCDFGHIJKLMNOPQRSTUVWXYZ";
2324 while (line.find_first_of(letters)!=std::string::npos)
2333 std::vector<unsigned int> tecplot2deal(dim);
2334 std::vector<unsigned int> IJK(dim);
2335 unsigned int n_vars,
2341 parse_tecplot_header(header,
2342 tecplot2deal,n_vars,n_vertices,n_cells,IJK,
2343 structured,blocked);
2351 std::vector<Point<spacedim> > vertices(n_vertices+1);
2354 std::vector<CellData<dim> > cells(n_cells);
2370 unsigned int next_index=0;
2374 if (tecplot2deal[0]==0)
2380 for (
unsigned int i=1; i<first_var.size()+1; ++i)
2381 vertices[i](0) = std::strtod (first_var[i-1].c_str(), &endptr);
2386 for (
unsigned int j=first_var.size()+1; j<n_vertices+1; ++j)
2387 in>>vertices[j](next_index);
2394 for (
unsigned int i=1; i<n_vars; ++i)
2403 if (next_index==dim && structured)
2406 if ((next_index<dim) && (i==tecplot2deal[next_index]))
2409 for (
unsigned int j=1; j<n_vertices+1; ++j)
2410 in>>vertices[j](next_index);
2417 for (
unsigned int j=1; j<n_vertices+1; ++j)
2429 std::vector<double> vars(n_vars);
2436 for (
unsigned int d=0;
d<dim; ++
d)
2437 vertices[1](d) = std::strtod (first_vertex[tecplot2deal[d]].c_str(), &endptr);
2441 for (
unsigned int v=2; v<n_vertices+1; ++v)
2443 for (
unsigned int i=0; i<n_vars; ++i)
2449 for (
unsigned int i=0; i<dim; ++i)
2450 vertices[v](i)=vars[tecplot2deal[i]];
2458 unsigned int I=IJK[0],
2461 unsigned int cell=0;
2463 for (
unsigned int j=0; j<J-1; ++j)
2464 for (
unsigned int i=1; i<I; ++i)
2466 cells[cell].vertices[0]=i+ j *I;
2467 cells[cell].vertices[1]=i+1+j *I;
2468 cells[cell].vertices[2]=i+1+(j+1)*I;
2469 cells[cell].vertices[3]=i +(j+1)*I;
2473 std::vector<unsigned int> boundary_vertices(2*I+2*J-4);
2475 for (
unsigned int i=1; i<I+1; ++i)
2477 boundary_vertices[k]=i;
2479 boundary_vertices[k]=i+(J-1)*I;
2482 for (
unsigned int j=1; j<J-1; ++j)
2484 boundary_vertices[k]=1+j*I;
2486 boundary_vertices[k]=I+j*I;
2502 for (
unsigned int i=0; i<n_cells; ++i)
2513 for (
unsigned int j=0; j<GeometryInfo<dim>::vertices_per_cell; ++j)
2514 in>>cells[i].vertices[j];
2529 tria->create_triangulation_compatibility (vertices, cells, subcelldata);
2534 template <
int dim,
int spacedim>
2541 template <
int dim,
int spacedim>
2554 if (std::find_if (line.begin(), line.end(),
2555 std_cxx11::bind (std::not_equal_to<char>(),std_cxx11::_1,
' '))
2559 for (
int i=line.length()-1; i>=0; --i)
2560 in.putback (line[i]);
2570 template <
int dim,
int spacedim>
2572 const char comment_start)
2577 while ((c=in.get()) == comment_start)
2580 while (in.get() !=
'\n')
2589 skip_empty_lines(in);
2594 template <
int dim,
int spacedim>
2607 const std::vector<
Point<2> > &vertices,
2610 double min_x = vertices[cells[0].vertices[0]](0),
2611 max_x = vertices[cells[0].vertices[0]](0),
2612 min_y = vertices[cells[0].vertices[0]](1),
2613 max_y = vertices[cells[0].vertices[0]](1);
2615 for (
unsigned int i=0; i<cells.size(); ++i)
2617 for (
unsigned int v=0; v<4; ++v)
2619 const Point<2> &p = vertices[cells[i].vertices[v]];
2631 out <<
"# cell " << i << std::endl;
2633 for (
unsigned int f=0; f<4; ++f)
2634 center += vertices[cells[i].vertices[f]];
2637 out <<
"set label \"" << i <<
"\" at " 2638 << center(0) <<
',' << center(1)
2643 for (
unsigned int f=0; f<2; ++f)
2644 out <<
"set arrow from " 2645 << vertices[cells[i].vertices[f]](0) <<
',' 2646 << vertices[cells[i].vertices[f]](1)
2648 << vertices[cells[i].vertices[(f+1)%4]](0) <<
',' 2649 << vertices[cells[i].vertices[(f+1)%4]](1)
2652 for (
unsigned int f=2; f<4; ++f)
2653 out <<
"set arrow from " 2654 << vertices[cells[i].vertices[(f+1)%4]](0) <<
',' 2655 << vertices[cells[i].vertices[(f+1)%4]](1)
2657 << vertices[cells[i].vertices[f]](0) <<
',' 2658 << vertices[cells[i].vertices[f]](1)
2665 <<
"set nokey" << std::endl
2666 <<
"pl [" << min_x <<
':' << max_x <<
"][" 2667 << min_y <<
':' << max_y <<
"] " 2668 << min_y << std::endl
2669 <<
"pause -1" << std::endl;
2677 const std::vector<
Point<3> > &vertices,
2680 for (
unsigned int cell=0; cell<cells.size(); ++cell)
2683 out << vertices[cells[cell].vertices[0]]
2685 << vertices[cells[cell].vertices[1]]
2686 << std::endl << std::endl << std::endl;
2688 out << vertices[cells[cell].vertices[1]]
2690 << vertices[cells[cell].vertices[2]]
2691 << std::endl << std::endl << std::endl;
2693 out << vertices[cells[cell].vertices[3]]
2695 << vertices[cells[cell].vertices[2]]
2696 << std::endl << std::endl << std::endl;
2698 out << vertices[cells[cell].vertices[0]]
2700 << vertices[cells[cell].vertices[3]]
2701 << std::endl << std::endl << std::endl;
2703 out << vertices[cells[cell].vertices[4]]
2705 << vertices[cells[cell].vertices[5]]
2706 << std::endl << std::endl << std::endl;
2708 out << vertices[cells[cell].vertices[5]]
2710 << vertices[cells[cell].vertices[6]]
2711 << std::endl << std::endl << std::endl;
2713 out << vertices[cells[cell].vertices[7]]
2715 << vertices[cells[cell].vertices[6]]
2716 << std::endl << std::endl << std::endl;
2718 out << vertices[cells[cell].vertices[4]]
2720 << vertices[cells[cell].vertices[7]]
2721 << std::endl << std::endl << std::endl;
2723 out << vertices[cells[cell].vertices[0]]
2725 << vertices[cells[cell].vertices[4]]
2726 << std::endl << std::endl << std::endl;
2728 out << vertices[cells[cell].vertices[1]]
2730 << vertices[cells[cell].vertices[5]]
2731 << std::endl << std::endl << std::endl;
2733 out << vertices[cells[cell].vertices[2]]
2735 << vertices[cells[cell].vertices[6]]
2736 << std::endl << std::endl << std::endl;
2738 out << vertices[cells[cell].vertices[3]]
2740 << vertices[cells[cell].vertices[7]]
2741 << std::endl << std::endl << std::endl;
2747 template <
int dim,
int spacedim>
2755 if (format == Default)
2756 name = search.
find(filename);
2758 name = search.
find(filename, default_suffix(format));
2760 std::ifstream in(name.c_str());
2762 if (format == Default)
2764 const std::string::size_type slashpos = name.find_last_of(
'/');
2765 const std::string::size_type dotpos = name.find_last_of(
'.');
2766 if (dotpos < name.length()
2767 && (dotpos > slashpos || slashpos == std::string::npos))
2769 std::string ext = name.substr(dotpos+1);
2770 format = parse_format(ext);
2773 if (format == netcdf)
2774 read_netcdf(filename);
2780 template <
int dim,
int spacedim>
2784 if (format == Default)
2785 format = default_format;
2819 "Use the read(_netcdf)(string &filename) " 2820 "functions, instead."));
2835 template <
int dim,
int spacedim>
2861 return ".unknown_format";
2867 template <
int dim,
int spacedim>
2871 if (format_name ==
"dbmesh")
2874 if (format_name ==
"msh")
2877 if (format_name ==
"unv")
2880 if (format_name ==
"vtk")
2884 if (format_name ==
"inp")
2887 if (format_name ==
"ucd")
2890 if (format_name ==
"xda")
2893 if (format_name ==
"netcdf")
2896 if (format_name ==
"nc")
2899 if (format_name ==
"tecplot")
2902 if (format_name ==
"dat")
2905 if (format_name ==
"plt")
2924 template <
int dim,
int spacedim>
2927 return "dbmesh|msh|unv|vtk|ucd|abaqus|xda|netcdf|tecplot";
2933 Abaqus_to_UCD<dim>::Abaqus_to_UCD ()
2941 template <
class T>
bool 2943 const std::string &s,
2944 std::ios_base& (*f) (std::ios_base &))
2946 std::istringstream iss (s);
2947 return ! (iss >> f >> t).fail();
2952 extract_int (
const std::string &s)
2955 for (
unsigned int i = 0; i<s.size(); ++i)
2964 from_string(number, tmp, std::dec);
2970 Abaqus_to_UCD<dim>::read_in_abaqus (std::istream &input_stream)
2978 std::getline (input_stream, line);
2980 while (!input_stream.fail() && !input_stream.eof())
2982 std::transform(line.begin(), line.end(), line.begin(), ::toupper);
2984 if (line.compare (
"*HEADING") == 0 ||
2985 line.compare (0, 2,
"**") == 0 ||
2986 line.compare (0, 5,
"*PART") == 0)
2989 while (!input_stream.fail() && !input_stream.eof())
2991 std::getline (input_stream, line);
2996 else if (line.compare (0, 5,
"*NODE") == 0)
3005 while (!input_stream.fail() && !input_stream.eof())
3007 std::getline (input_stream, line);
3011 std::vector <double> node (dim+1);
3013 std::istringstream iss (line);
3015 for (
unsigned int i = 0; i < dim+1; ++i)
3016 iss >> node[i] >> comma;
3018 node_list.push_back (node);
3021 else if (line.compare (0, 8,
"*ELEMENT") == 0)
3036 const std::string before_material =
"ELSET=EB";
3037 const std::size_t idx = line.find (before_material);
3038 if (idx != std::string::npos)
3040 from_string (material, line.substr (idx + before_material.size()), std::dec);
3045 std::getline (input_stream, line);
3046 while (!input_stream.fail() && !input_stream.eof())
3051 std::istringstream iss (line);
3058 std::vector <double> cell (n_data_per_cell);
3059 for (
unsigned int i = 0; i < n_data_per_cell; ++i)
3060 iss >> cell[i] >> comma;
3063 cell[0] =
static_cast<double> (material);
3064 cell_list.push_back (cell);
3066 std::getline (input_stream, line);
3069 else if (line.compare (0, 8,
"*SURFACE") == 0)
3080 const std::string name_key =
"NAME=";
3081 const std::size_t name_idx_start = line.find(name_key) + name_key.size();
3082 std::size_t name_idx_end = line.find(
',', name_idx_start);
3083 if (name_idx_end == std::string::npos)
3085 name_idx_end = line.size();
3087 const int b_indicator = extract_int(line.substr(name_idx_start, name_idx_end - name_idx_start));
3093 std::getline (input_stream, line);
3094 while (!input_stream.fail() && !input_stream.eof())
3100 std::transform(line.begin(), line.end(), line.begin(), ::toupper);
3104 std::istringstream iss (line);
3110 std::vector <double> quad_node_list;
3111 const std::string elset_name = line.substr(0, line.find(
','));
3112 if (elsets_list.count(elset_name) != 0)
3116 iss >> stmp >> temp >> face_number;
3118 const std::vector<int> cells = elsets_list[elset_name];
3119 for (
unsigned int i = 0; i <cells.size(); ++i)
3122 quad_node_list = get_global_node_numbers (el_idx, face_number);
3123 quad_node_list.insert (quad_node_list.begin(), b_indicator);
3125 face_list.push_back (quad_node_list);
3132 iss >> el_idx >> comma >> temp >> face_number;
3133 quad_node_list = get_global_node_numbers (el_idx, face_number);
3134 quad_node_list.insert (quad_node_list.begin(), b_indicator);
3136 face_list.push_back (quad_node_list);
3139 std::getline (input_stream, line);
3142 else if (line.compare (0, 6,
"*ELSET") == 0)
3146 std::string elset_name;
3148 const std::string elset_key =
"*ELSET, ELSET=";
3149 const std::size_t idx = line.find(elset_key);
3150 if (idx != std::string::npos)
3152 const std::string comma =
",";
3153 const std::size_t first_comma = line.find(comma);
3154 const std::size_t second_comma = line.find(comma, first_comma+1);
3155 const std::size_t elset_name_start = line.find(elset_key) + elset_key.size();
3156 elset_name = line.substr(elset_name_start, second_comma-elset_name_start);
3165 std::vector<int> elements;
3166 const std::size_t generate_idx = line.find(
"GENERATE");
3167 if (generate_idx != std::string::npos)
3170 std::getline (input_stream, line);
3171 std::istringstream iss (line);
3178 iss >> elid_start >> comma >> elid_end;
3180 if (iss.rdbuf()->in_avail() != 0)
3181 iss >> comma >> elis_step;
3182 for (
int i = elid_start; i <= elid_end; i+= elis_step)
3184 elements.push_back(i);
3186 elsets_list[elset_name] = elements;
3188 std::getline (input_stream, line);
3193 std::getline (input_stream, line);
3194 while (!input_stream.fail() && !input_stream.eof())
3199 std::istringstream iss (line);
3204 iss >> elid >> comma;
3205 elements.push_back (elid);
3208 std::getline (input_stream, line);
3211 elsets_list[elset_name] = elements;
3216 else if (line.compare (0, 5,
"*NSET") == 0)
3219 while (!input_stream.fail() && !input_stream.eof())
3221 std::getline (input_stream, line);
3226 else if (line.compare(0, 14,
"*SOLID SECTION") == 0)
3229 const std::string elset_key =
"ELSET=";
3230 const std::size_t elset_start = line.find(
"ELSET=") + elset_key.size();
3231 const std::size_t elset_end = line.find(
',', elset_start+1);
3232 const std::string elset_name = line.substr(elset_start, elset_end-elset_start);
3237 const std::string material_key =
"MATERIAL=";
3238 const std::size_t last_equal = line.find(
"MATERIAL=") + material_key.size();
3239 const std::size_t material_id_start = line.find(
'-', last_equal);
3241 from_string(material_id, line.substr(material_id_start+1), std::dec);
3244 const std::vector<int> &elset_cells = elsets_list[elset_name];
3245 for (
unsigned int i = 0; i < elset_cells.size(); ++i)
3247 const int cell_id = elset_cells[i] - 1;
3253 std::getline (input_stream, line);
3262 Abaqus_to_UCD<dim>::get_global_node_numbers (
const int face_cell_no,
3263 const int face_cell_face_no)
const 3273 if (face_cell_face_no == 1)
3275 quad_node_list[0] = cell_list[face_cell_no - 1][1];
3276 quad_node_list[1] = cell_list[face_cell_no - 1][2];
3278 else if (face_cell_face_no == 2)
3280 quad_node_list[0] = cell_list[face_cell_no - 1][2];
3281 quad_node_list[1] = cell_list[face_cell_no - 1][3];
3283 else if (face_cell_face_no == 3)
3285 quad_node_list[0] = cell_list[face_cell_no - 1][3];
3286 quad_node_list[1] = cell_list[face_cell_no - 1][4];
3288 else if (face_cell_face_no == 4)
3290 quad_node_list[0] = cell_list[face_cell_no - 1][4];
3291 quad_node_list[1] = cell_list[face_cell_no - 1][1];
3300 if (face_cell_face_no == 1)
3302 quad_node_list[0] = cell_list[face_cell_no - 1][1];
3303 quad_node_list[1] = cell_list[face_cell_no - 1][4];
3304 quad_node_list[2] = cell_list[face_cell_no - 1][3];
3305 quad_node_list[3] = cell_list[face_cell_no - 1][2];
3307 else if (face_cell_face_no == 2)
3309 quad_node_list[0] = cell_list[face_cell_no - 1][5];
3310 quad_node_list[1] = cell_list[face_cell_no - 1][8];
3311 quad_node_list[2] = cell_list[face_cell_no - 1][7];
3312 quad_node_list[3] = cell_list[face_cell_no - 1][6];
3314 else if (face_cell_face_no == 3)
3316 quad_node_list[0] = cell_list[face_cell_no - 1][1];
3317 quad_node_list[1] = cell_list[face_cell_no - 1][2];
3318 quad_node_list[2] = cell_list[face_cell_no - 1][6];
3319 quad_node_list[3] = cell_list[face_cell_no - 1][5];
3321 else if (face_cell_face_no == 4)
3323 quad_node_list[0] = cell_list[face_cell_no - 1][2];
3324 quad_node_list[1] = cell_list[face_cell_no - 1][3];
3325 quad_node_list[2] = cell_list[face_cell_no - 1][7];
3326 quad_node_list[3] = cell_list[face_cell_no - 1][6];
3328 else if (face_cell_face_no == 5)
3330 quad_node_list[0] = cell_list[face_cell_no - 1][3];
3331 quad_node_list[1] = cell_list[face_cell_no - 1][4];
3332 quad_node_list[2] = cell_list[face_cell_no - 1][8];
3333 quad_node_list[3] = cell_list[face_cell_no - 1][7];
3335 else if (face_cell_face_no == 6)
3337 quad_node_list[0] = cell_list[face_cell_no - 1][1];
3338 quad_node_list[1] = cell_list[face_cell_no - 1][5];
3339 quad_node_list[2] = cell_list[face_cell_no - 1][8];
3340 quad_node_list[3] = cell_list[face_cell_no - 1][4];
3352 return quad_node_list;
3357 Abaqus_to_UCD<dim>::write_out_avs_ucd (std::ostream &output)
const 3367 output <<
"# Abaqus to UCD mesh conversion" << std::endl;
3368 output <<
"# Mesh type: AVS UCD" << std::endl;
3394 << node_list.size() <<
"\t" 3395 << (cell_list.size() + face_list.size()) <<
"\t0\t0\t0" 3399 for (
unsigned int ii = 0; ii < node_list.size(); ++ii)
3401 for (
unsigned int jj = 0; jj < dim + 1; ++jj)
3406 output << node_list[ii][jj] <<
"\t";
3411 output.setf (std::ios::scientific,
3412 std::ios::floatfield);
3413 output.precision (8);
3414 if (std::abs (node_list[ii][jj]) > tolerance)
3415 output << static_cast<double> (node_list[ii][jj]) <<
"\t";
3417 output << 0.0 <<
"\t";
3421 output << 0.0 <<
"\t";
3425 output.unsetf (std::ios::floatfield);
3429 for (
unsigned int ii = 0; ii < cell_list.size(); ++ii)
3433 << cell_list[ii][0] <<
"\t" 3434 << (dim == 2 ?
"quad" :
"hex") <<
"\t";
3435 for (
unsigned int jj = 1; jj < GeometryInfo<dim>::vertices_per_cell + 1; ++jj)
3437 << cell_list[ii][jj] <<
"\t";
3444 for (
unsigned int ii = 0; ii < face_list.size(); ++ii)
3448 << face_list[ii][0] <<
"\t" 3449 << (dim == 2 ?
"line" :
"quad") <<
"\t";
3450 for (
unsigned int jj = 1; jj < GeometryInfo<dim>::vertices_per_face + 1; ++jj)
3452 << face_list[ii][jj] <<
"\t";
3462 #include "grid_in.inst" 3464 DEAL_II_NAMESPACE_CLOSE
std::vector< CellData< 1 > > boundary_lines
static std::string get_format_names()
static const unsigned int invalid_unsigned_int
static ::ExceptionBase & ExcIO()
unsigned char material_id
std::pair< int, unsigned int > get_integer_at_position(const std::string &name, const unsigned int position)
std::vector< std::string > break_text_into_lines(const std::string &original_text, const unsigned int width, const char delimiter=' ')
active_cell_iterator begin_active(const unsigned int level=0) const
#define AssertThrow(cond, exc)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
static Format parse_format(const std::string &format_name)
static void skip_empty_lines(std::istream &in)
void read_vtk(std::istream &in)
cell_iterator end() const
void read_dbmesh(std::istream &in)
void read_tecplot(std::istream &in)
static ::ExceptionBase & ExcInvalidState()
static ::ExceptionBase & ExcMessage(std::string arg1)
bool check_consistency(const unsigned int dim) const
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
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)
static ::ExceptionBase & ExcNeedsNetCDF()
#define Assert(cond, exc)
static void invert_all_cells_of_negative_grid(const std::vector< Point< spacedim > > &all_vertices, std::vector< CellData< dim > > &original_cells)
static void skip_comment_lines(std::istream &in, const char comment_start)
static ::ExceptionBase & ExcInvalidVertexIndex(int arg1, int arg2, int arg3)
bool match_at_string_start(const std::string &name, const std::string &pattern)
void read_ucd(std::istream &in, const bool apply_all_indicators_to_manifolds=false)
static void reorder_cells(std::vector< CellData< dim > > &original_cells, const bool use_new_style_ordering=false)
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
std::string find(const std::string &filename, const char *open_mode="r")
void read_xda(std::istream &in)
void read_msh(std::istream &in)
std::vector< CellData< 2 > > boundary_quads
void read_netcdf(const std::string &filename)
static std::string default_suffix(const Format format)
void read_unv(std::istream &in)
void read(std::istream &in, Format format=Default)
static ::ExceptionBase & ExcNotImplemented()
unsigned char boundary_id
const types::boundary_id internal_face_boundary_id
void read_abaqus(std::istream &in, const bool apply_all_indicators_to_manifolds=false)
static void debug_output_grid(const std::vector< CellData< dim > > &cells, const std::vector< Point< spacedim > > &vertices, std::ostream &out)
const types::material_id invalid_material_id
void attach_triangulation(Triangulation< dim, spacedim > &tria)
static ::ExceptionBase & ExcNoTriangulationSelected()
static ::ExceptionBase & ExcInternalError()