16 #include <deal.II/grid/grid_out.h> 17 #include <deal.II/base/parameter_handler.h> 18 #include <deal.II/base/exceptions.h> 19 #include <deal.II/base/point.h> 20 #include <deal.II/base/quadrature.h> 21 #include <deal.II/base/qprojector.h> 22 #include <deal.II/base/geometry_info.h> 23 #include <deal.II/grid/tria.h> 24 #include <deal.II/grid/tria_accessor.h> 25 #include <deal.II/grid/tria_iterator.h> 26 #include <deal.II/fe/mapping.h> 27 #include <deal.II/numerics/data_out.h> 39 DEAL_II_NAMESPACE_OPEN
45 const bool write_faces,
46 const bool write_diameter,
47 const bool write_measure,
48 const bool write_all_faces) :
49 write_cells (write_cells),
50 write_faces (write_faces),
51 write_diameter (write_diameter),
52 write_measure (write_measure),
53 write_all_faces (write_all_faces)
59 "Write the mesh connectivity as DX grid cells");
61 "Write faces of cells. These may be boundary faces " 62 "or all faces between mesh cells, according to " 63 "\"Write all faces\"");
65 "If cells are written, additionally write their" 66 " diameter as data for visualization");
68 "Write the volume of each cell as data");
70 "Write all faces, not only boundary");
84 const bool write_lines) :
85 write_faces (write_faces),
86 write_lines (write_lines)
104 const bool write_faces,
105 const bool write_lines) :
106 write_preamble (write_preamble),
107 write_faces (write_faces),
108 write_lines (write_lines)
130 const unsigned int n_boundary_face_points,
131 const bool curved_inner_cells) :
132 write_cell_numbers (write_cell_numbers),
133 n_boundary_face_points(n_boundary_face_points),
134 curved_inner_cells(curved_inner_cells)
154 const unsigned int size,
155 const double line_width,
156 const bool color_lines_on_user_flag,
157 const unsigned int n_boundary_face_points,
158 const bool color_lines_level) :
159 size_type (size_type),
161 line_width (line_width),
162 color_lines_on_user_flag(color_lines_on_user_flag),
163 n_boundary_face_points(n_boundary_face_points),
164 color_lines_level(color_lines_level)
172 "Depending on this parameter, either the" 174 "of the eps is scaled to \"Size\"");
176 "Size of the output in points");
178 "Width of the lines drawn in points");
180 "Draw lines with user flag set in different color");
182 "Number of points on boundary edges. " 183 "Increase this beyond 2 to see curved boundaries.");
185 "Draw different colors according to grid level.");
191 if (param.
get(
"Size by") == std::string(
"width"))
193 else if (param.
get(
"Size by") == std::string(
"height"))
205 const unsigned int size,
206 const double line_width,
207 const bool color_lines_on_user_flag,
208 const unsigned int n_boundary_face_points)
211 color_lines_on_user_flag,
212 n_boundary_face_points)
228 const unsigned int size,
229 const double line_width,
230 const bool color_lines_on_user_flag,
231 const unsigned int n_boundary_face_points,
232 const bool write_cell_numbers,
233 const bool write_cell_number_level,
234 const bool write_vertex_numbers,
235 const bool color_lines_level
239 color_lines_on_user_flag,
240 n_boundary_face_points,
242 write_cell_numbers (write_cell_numbers),
243 write_cell_number_level (write_cell_number_level),
244 write_vertex_numbers (write_vertex_numbers)
251 "(2D only) Write cell numbers" 252 " into the centers of cells");
254 "(2D only) if \"Cell number\" is true, write" 255 "numbers in the form level.number");
257 "Write numbers for each vertex");
264 write_cell_numbers = param.
get_bool(
"Cell number");
265 write_cell_number_level = param.
get_bool(
"Level number");
266 write_vertex_numbers = param.
get_bool(
"Vertex number");
272 const unsigned int size,
273 const double line_width,
274 const bool color_lines_on_user_flag,
275 const unsigned int n_boundary_face_points,
276 const double azimut_angle,
277 const double turn_angle)
280 color_lines_on_user_flag,
281 n_boundary_face_points),
282 azimut_angle (azimut_angle),
283 turn_angle (turn_angle)
290 "Azimuth of the viw point, that is, the angle " 291 "in the plane from the x-axis.");
293 "Elevation of the view point above the xy-plane.");
300 azimut_angle = 90- param.
get_double(
"Elevation");
309 color_by(material_id),
311 n_boundary_face_points(0),
317 boundary_thickness(3)
349 const unsigned int boundary_line_thickness,
352 const int azimuth_angle,
353 const int polar_angle,
355 const bool convert_level_number_to_height,
356 const bool label_level_number,
357 const bool label_cell_index,
358 const bool label_material_id,
359 const bool label_subdomain_id,
360 const bool draw_colorbar,
361 const bool draw_legend)
365 line_thickness(line_thickness),
366 boundary_line_thickness(boundary_line_thickness),
368 background(background),
369 azimuth_angle(azimuth_angle),
370 polar_angle(polar_angle),
372 convert_level_number_to_height(convert_level_number_to_height),
373 level_height_factor(0.3f),
374 cell_font_scaling(1.f),
375 label_level_number(label_level_number),
376 label_cell_index(label_cell_index),
377 label_material_id(label_material_id),
378 label_subdomain_id(label_subdomain_id),
379 label_level_subdomain_id(false),
380 draw_colorbar(draw_colorbar),
381 draw_legend(draw_legend)
386 draw_bounding_box (false)
404 default_format (none)
486 switch (output_format)
529 if (format_name ==
"none" || format_name ==
"false")
532 if (format_name ==
"dx")
535 if (format_name ==
"ucd")
538 if (format_name ==
"gnuplot")
541 if (format_name ==
"eps")
544 if (format_name ==
"xfig")
547 if (format_name ==
"msh")
550 if (format_name ==
"svg")
553 if (format_name ==
"mathgl")
556 if (format_name ==
"vtk")
559 if (format_name ==
"vtu")
571 return "none|dx|gnuplot|eps|ucd|xfig|msh|svg|mathgl|vtk|vtu";
690 std::ostream &)
const 697 std::ostream &)
const 704 std::ostream &)
const 711 template <
int dim,
int spacedim>
713 std::ostream &out)
const 719 const std::vector<Point<spacedim> > &vertices = tria.
get_vertices();
729 std::vector<unsigned int> renumber(vertices.size());
732 unsigned int new_number=0;
733 for (
unsigned int i=0; i<vertices.size(); ++i)
735 renumber[i]=new_number++;
743 out <<
"object \"vertices\" class array type float rank 1 shape " << dim
744 <<
" items " << n_vertices <<
" data follows" 747 for (
unsigned int i=0; i<vertices.size(); ++i)
749 out <<
'\t' << vertices[i] <<
'\n';
764 out <<
"object \"cells\" class array type int rank 1 shape " 765 << n_vertices_per_cell
766 <<
" items " << n_cells <<
" data follows" <<
'\n';
770 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
774 out <<
"attribute \"element type\" string \"";
775 if (dim==1) out <<
"lines";
776 if (dim==2) out <<
"quads";
777 if (dim==3) out <<
"cubes";
779 <<
"attribute \"ref\" string \"positions\"" <<
'\n' <<
'\n';
783 out <<
"object \"material\" class array type int rank 0 items " 784 << n_cells <<
" data follows" <<
'\n';
786 out <<
' ' << (
unsigned int)cell->material_id();
788 <<
"attribute \"dep\" string \"connections\"" <<
'\n' <<
'\n';
790 out <<
"object \"level\" class array type int rank 0 items " 791 << n_cells <<
" data follows" <<
'\n';
793 out <<
' ' << cell->level();
795 <<
"attribute \"dep\" string \"connections\"" <<
'\n' <<
'\n';
799 out <<
"object \"measure\" class array type float rank 0 items " 800 << n_cells <<
" data follows" <<
'\n';
802 out <<
'\t' << cell->measure();
804 <<
"attribute \"dep\" string \"connections\"" <<
'\n' <<
'\n';
809 out <<
"object \"diameter\" class array type float rank 0 items " 810 << n_cells <<
" data follows" <<
'\n';
812 out <<
'\t' << cell->diameter();
814 <<
"attribute \"dep\" string \"connections\"" <<
'\n' <<
'\n';
820 out <<
"object \"faces\" class array type int rank 1 shape " 821 << n_vertices_per_face
822 <<
" items " << n_faces <<
" data follows" 827 for (
unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
831 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_face; ++v)
836 out <<
"attribute \"element type\" string \"";
837 if (dim==2) out <<
"lines";
838 if (dim==3) out <<
"quads";
840 <<
"attribute \"ref\" string \"positions\"" <<
'\n' <<
'\n';
845 out <<
"object \"boundary\" class array type int rank 0 items " 846 << n_faces <<
" data follows" <<
'\n';
851 for (
unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
852 out <<
' ' << (
int)(
signed char)cell->face(f)->boundary_id();
855 out <<
"attribute \"dep\" string \"connections\"" <<
'\n' <<
'\n';
859 out <<
"object \"face measure\" class array type float rank 0 items " 860 << n_faces <<
" data follows" <<
'\n';
863 for (
unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
864 out <<
' ' << cell->face(f)->measure();
867 out <<
"attribute \"dep\" string \"connections\"" <<
'\n' <<
'\n';
872 out <<
"object \"face diameter\" class array type float rank 0 items " 873 << n_faces <<
" data follows" <<
'\n';
876 for (
unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
877 out <<
' ' << cell->face(f)->diameter();
880 out <<
"attribute \"dep\" string \"connections\"" <<
'\n' <<
'\n';
897 out <<
"object \"deal data\" class field" <<
'\n' 898 <<
"component \"positions\" value \"vertices\"" <<
'\n' 899 <<
"component \"connections\" value \"cells\"" <<
'\n';
903 out <<
"object \"cell data\" class field" <<
'\n' 904 <<
"component \"positions\" value \"vertices\"" <<
'\n' 905 <<
"component \"connections\" value \"cells\"" <<
'\n';
906 out <<
"component \"material\" value \"material\"" <<
'\n';
907 out <<
"component \"level\" value \"level\"" <<
'\n';
909 out <<
"component \"measure\" value \"measure\"" <<
'\n';
911 out <<
"component \"diameter\" value \"diameter\"" <<
'\n';
916 out <<
"object \"face data\" class field" <<
'\n' 917 <<
"component \"positions\" value \"vertices\"" <<
'\n' 918 <<
"component \"connections\" value \"faces\"" <<
'\n';
919 out <<
"component \"boundary\" value \"boundary\"" <<
'\n';
921 out <<
"component \"measure\" value \"face measure\"" <<
'\n';
923 out <<
"component \"diameter\" value \"face diameter\"" <<
'\n';
927 <<
"object \"grid data\" class group" <<
'\n';
929 out <<
"member \"cells\" value \"cell data\"" <<
'\n';
931 out <<
"member \"faces\" value \"face data\"" <<
'\n';
932 out <<
"end" <<
'\n';
944 template <
int dim,
int spacedim>
946 std::ostream &out)
const 953 const std::vector<Point<spacedim> > &vertices = tria.
get_vertices();
977 out <<
"@f$NOD" <<
'\n' 978 << n_vertices <<
'\n';
983 for (
unsigned int i=0; i<vertices.size(); ++i)
989 for (
unsigned int d=spacedim+1; d<=3; ++d)
995 out <<
"@f$ENDNOD" <<
'\n' 1036 unsigned int elm_type;
1056 out << cell->active_cell_index()+1 <<
' ' << elm_type <<
' ' 1057 <<
static_cast<unsigned int>(cell->material_id()) <<
' ' 1058 << cell->subdomain_id() <<
' ' 1063 for (
unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_cell;
1080 out <<
"@f$ENDELM\n";
1090 template <
int dim,
int spacedim>
1092 std::ostream &out)
const 1099 const std::vector<Point<spacedim> > &vertices = tria.
get_vertices();
1113 std::time_t time1= std::time (
nullptr);
1114 std::tm *time = std::localtime(&time1);
1115 out <<
"# This file was generated by the deal.II library." <<
'\n' 1117 << time->tm_year+1900 <<
"/" 1118 << time->tm_mon+1 <<
"/" 1119 << time->tm_mday <<
'\n' 1121 << time->tm_hour <<
":" 1122 << std::setw(2) << time->tm_min <<
":" 1123 << std::setw(2) << time->tm_sec <<
'\n' 1125 <<
"# For a description of the UCD format see the AVS Developer's guide." 1131 out << n_vertices <<
' ' 1142 for (
unsigned int i=0; i<vertices.size(); ++i)
1148 for (
unsigned int d=spacedim+1; d<=3; ++d)
1157 out << cell->active_cell_index()+1 <<
' ' 1158 <<
static_cast<unsigned int>(cell->material_id())
1190 for (
unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_cell;
1216 template <
int dim,
int spacedim>
1234 const int spacedim = 2;
1242 out <<
"#FIG 3.2\nLandscape\nCenter\nInches" << std::endl
1243 <<
"A4\n100.00\nSingle" << std::endl
1245 <<
"-3" << std::endl
1246 <<
"# generated by deal.II GridOut class" << std::endl
1247 <<
"# reduce first number to scale up image" << std::endl
1248 <<
"1200 2" << std::endl;
1251 unsigned int colno = 32;
1252 out <<
"0 " << colno++ <<
" #ff0000" << std::endl;
1253 out <<
"0 " << colno++ <<
" #ff8000" << std::endl;
1254 out <<
"0 " << colno++ <<
" #ffd000" << std::endl;
1255 out <<
"0 " << colno++ <<
" #ffff00" << std::endl;
1256 out <<
"0 " << colno++ <<
" #c0ff00" << std::endl;
1257 out <<
"0 " << colno++ <<
" #80ff00" << std::endl;
1258 out <<
"0 " << colno++ <<
" #00f000" << std::endl;
1259 out <<
"0 " << colno++ <<
" #00f0c0" << std::endl;
1260 out <<
"0 " << colno++ <<
" #00f0ff" << std::endl;
1261 out <<
"0 " << colno++ <<
" #00c0ff" << std::endl;
1262 out <<
"0 " << colno++ <<
" #0080ff" << std::endl;
1263 out <<
"0 " << colno++ <<
" #0040ff" << std::endl;
1264 out <<
"0 " << colno++ <<
" #0000c0" << std::endl;
1265 out <<
"0 " << colno++ <<
" #5000ff" << std::endl;
1266 out <<
"0 " << colno++ <<
" #8000ff" << std::endl;
1267 out <<
"0 " << colno++ <<
" #b000ff" << std::endl;
1268 out <<
"0 " << colno++ <<
" #ff00ff" << std::endl;
1269 out <<
"0 " << colno++ <<
" #ff80ff" << std::endl;
1271 for (
unsigned int i=0; i<8; ++i)
1272 out <<
"0 " << colno++ <<
" #" << std::hex << 32*i+31 << 32*i+31 << 32*i+31 << std::dec << std::endl;
1274 for (
unsigned int i=1; i<16; ++i)
1275 out <<
"0 " << colno++ <<
" #00" << std::hex << 16*i+15 << std::dec <<
"00" << std::endl;
1277 for (
unsigned int i=1; i<16; ++i)
1278 out <<
"0 " << colno++ <<
" #" << std::hex << 16*i+15 << 16*i+15 << std::dec <<
"00" << std::endl;
1280 for (
unsigned int i=1; i<16; ++i)
1281 out <<
"0 " << colno++ <<
" #" << std::hex << 16*i+15 << std::dec <<
"0000" << std::endl;
1283 for (
unsigned int i=1; i<16; ++i)
1284 out <<
"0 " << colno++ <<
" #" << std::hex << 16*i+15 <<
"00" << 16*i+15 << std::dec << std::endl;
1286 for (
unsigned int i=1; i<16; ++i)
1287 out <<
"0 " << colno++ <<
" #0000" << std::hex << 16*i+15 << std::dec << std::endl;
1289 for (
unsigned int i=1; i<16; ++i)
1290 out <<
"0 " << colno++ <<
" #00" << std::hex << 16*i+15 << 16*i+15 << std::dec << std::endl;
1300 for (; cell != end; ++cell)
1316 out << static_cast<unsigned int>(cell->material_id()) + 32;
1319 out << cell->level() + 8;
1322 out << cell->subdomain_id() + 32;
1325 out << cell->level_subdomain_id() + 32;
1334 ? (900-cell->level())
1335 : (900+cell->material_id()))
1341 << nv+1 << std::endl;
1347 for (
unsigned int k=0; k<=nv; ++k)
1351 for (
unsigned int d=0; d<static_cast<unsigned int>(dim); ++
d)
1355 out <<
'\t' << ((
d==0) ? val : -val);
1360 static const unsigned int face_reorder[4]= {2,1,3,0};
1362 for (
unsigned int f=0; f<nf; ++f)
1365 face = cell->face(face_reorder[f]);
1374 << (1 + (
unsigned int) bi);
1379 ? (800-cell->level())
1386 << nvf << std::endl;
1393 for (
unsigned int k=0; k<nvf; ++k)
1396 for (
unsigned int d=0; d<static_cast<unsigned int>(dim); ++
d)
1400 out <<
'\t' << ((
d==0) ? val : -val);
1417 template <
int dim,
int spacedim>
1419 std::ostream &)
const 1428 unsigned int n_materials = 0;
1429 unsigned int n_levels = 0;
1430 unsigned int n_subdomains = 0;
1431 unsigned int n_level_subdomains = 0;
1435 unsigned int min_level, max_level;
1443 Assert (height != 0 || width != 0,
ExcMessage(
"You have to set at least one of width and height"));
1445 unsigned int margin_in_percent = 0;
1447 margin_in_percent = 8;
1450 unsigned int cell_label_font_size;
1467 float x_max_perspective, x_min_perspective;
1468 float y_max_perspective, y_min_perspective;
1470 float x_dimension_perspective, y_dimension_perspective;
1474 double x_min = tria.
begin()->vertex(0)[0];
1475 double x_max = x_min;
1476 double y_min = tria.
begin()->vertex(0)[1];
1477 double y_max = y_min;
1479 double x_dimension, y_dimension;
1481 min_level = max_level = tria.
begin()->level();
1484 unsigned int materials[256];
1485 for (
unsigned int material_index = 0; material_index < 256; material_index++)
1486 materials[material_index] = 0;
1489 unsigned int levels[256];
1490 for (
unsigned int level_index = 0; level_index < 256; level_index++)
1491 levels[level_index] = 0;
1494 unsigned int subdomains[256];
1495 for (
unsigned int subdomain_index = 0; subdomain_index < 256; subdomain_index++)
1496 subdomains[subdomain_index] = 0;
1499 int level_subdomains[256];
1500 for (
int level_subdomain_index = 0; level_subdomain_index < 256; level_subdomain_index++)
1501 level_subdomains[level_subdomain_index] = 0;
1509 for (
unsigned int vertex_index = 0; vertex_index < 4; vertex_index++)
1511 if (cell->vertex(vertex_index)[0] < x_min) x_min = cell->vertex(vertex_index)[0];
1512 if (cell->vertex(vertex_index)[0] > x_max) x_max = cell->vertex(vertex_index)[0];
1514 if (cell->vertex(vertex_index)[1] < y_min) y_min = cell->vertex(vertex_index)[1];
1515 if (cell->vertex(vertex_index)[1] > y_max) y_max = cell->vertex(vertex_index)[1];
1518 if ((
unsigned int)cell->level() < min_level) min_level = cell->level();
1519 if ((
unsigned int)cell->level() > max_level) max_level = cell->level();
1521 materials[(
unsigned int)cell->material_id()] = 1;
1522 levels[(
unsigned int)cell->level()] = 1;
1524 subdomains[cell->subdomain_id()+2] = 1;
1525 level_subdomains[cell->level_subdomain_id()+2] = 1;
1528 x_dimension = x_max - x_min;
1529 y_dimension = y_max - y_min;
1532 for (
unsigned int material_index = 0; material_index < 256; material_index++)
1534 if (materials[material_index]) n_materials++;
1538 for (
unsigned int level_index = 0; level_index < 256; level_index++)
1540 if (levels[level_index]) n_levels++;
1544 for (
unsigned int subdomain_index = 0; subdomain_index < 256; subdomain_index++)
1546 if (subdomains[subdomain_index]) n_subdomains++;
1550 for (
int level_subdomain_index = 0; level_subdomain_index < 256; level_subdomain_index++)
1552 if (level_subdomains[level_subdomain_index]) n_level_subdomains++;
1567 n = n_level_subdomains;
1574 camera_position[0] = 0;
1575 camera_position[1] = 0;
1576 camera_position[2] = 2. * std::max(x_dimension, y_dimension);
1578 camera_direction[0] = 0;
1579 camera_direction[1] = 0;
1580 camera_direction[2] = -1;
1582 camera_horizontal[0] = 1;
1583 camera_horizontal[1] = 0;
1584 camera_horizontal[2] = 0;
1586 camera_focus = .5 * std::max(x_dimension, y_dimension);
1592 const double angle_factor = 3.14159265 / 180.;
1604 camera_position[1] = camera_position_temp[1];
1605 camera_position[2] = camera_position_temp[2];
1607 camera_direction[1] = camera_direction_temp[1];
1608 camera_direction[2] = camera_direction_temp[2];
1610 camera_horizontal[1] = camera_horizontal_temp[1];
1611 camera_horizontal[2] = camera_horizontal_temp[2];
1623 camera_position[0] = camera_position_temp[0];
1624 camera_position[1] = camera_position_temp[1];
1626 camera_direction[0] = camera_direction_temp[0];
1627 camera_direction[1] = camera_direction_temp[1];
1629 camera_horizontal[0] = camera_horizontal_temp[0];
1630 camera_horizontal[1] = camera_horizontal_temp[1];
1633 camera_position[0] = x_min + .5 * x_dimension;
1634 camera_position[1] = y_min + .5 * y_dimension;
1641 point[0] = tria.
begin()->vertex(0)[0];
1642 point[1] = tria.
begin()->vertex(0)[1];
1645 float min_level_min_vertex_distance = 0;
1652 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
1654 x_max_perspective = projection_decomposition[0];
1655 x_min_perspective = projection_decomposition[0];
1657 y_max_perspective = projection_decomposition[1];
1658 y_min_perspective = projection_decomposition[1];
1662 point[0] = cell->vertex(0)[0];
1663 point[1] = cell->vertex(0)[1];
1671 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
1673 if (x_max_perspective < projection_decomposition[0]) x_max_perspective = projection_decomposition[0];
1674 if (x_min_perspective > projection_decomposition[0]) x_min_perspective = projection_decomposition[0];
1676 if (y_max_perspective < projection_decomposition[1]) y_max_perspective = projection_decomposition[1];
1677 if (y_min_perspective > projection_decomposition[1]) y_min_perspective = projection_decomposition[1];
1679 point[0] = cell->vertex(1)[0];
1680 point[1] = cell->vertex(1)[1];
1682 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
1684 if (x_max_perspective < projection_decomposition[0]) x_max_perspective = projection_decomposition[0];
1685 if (x_min_perspective > projection_decomposition[0]) x_min_perspective = projection_decomposition[0];
1687 if (y_max_perspective < projection_decomposition[1]) y_max_perspective = projection_decomposition[1];
1688 if (y_min_perspective > projection_decomposition[1]) y_min_perspective = projection_decomposition[1];
1690 point[0] = cell->vertex(2)[0];
1691 point[1] = cell->vertex(2)[1];
1693 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
1695 if (x_max_perspective < projection_decomposition[0]) x_max_perspective = projection_decomposition[0];
1696 if (x_min_perspective > projection_decomposition[0]) x_min_perspective = projection_decomposition[0];
1698 if (y_max_perspective < projection_decomposition[1]) y_max_perspective = projection_decomposition[1];
1699 if (y_min_perspective > projection_decomposition[1]) y_min_perspective = projection_decomposition[1];
1701 point[0] = cell->vertex(3)[0];
1702 point[1] = cell->vertex(3)[1];
1704 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
1706 if (x_max_perspective < projection_decomposition[0]) x_max_perspective = projection_decomposition[0];
1707 if (x_min_perspective > projection_decomposition[0]) x_min_perspective = projection_decomposition[0];
1709 if (y_max_perspective < projection_decomposition[1]) y_max_perspective = projection_decomposition[1];
1710 if (y_min_perspective > projection_decomposition[1]) y_min_perspective = projection_decomposition[1];
1712 if ((
unsigned int)cell->level() == min_level) min_level_min_vertex_distance = cell->minimum_vertex_distance();
1715 x_dimension_perspective = x_max_perspective - x_min_perspective;
1716 y_dimension_perspective = y_max_perspective - y_min_perspective;
1720 width =
static_cast<unsigned int>(.5 + height * (x_dimension_perspective / y_dimension_perspective));
1721 else if (height == 0)
1722 height =
static_cast<unsigned int>(.5 + width * (y_dimension_perspective / x_dimension_perspective));
1723 unsigned int additional_width = 0;
1725 unsigned int font_size =
static_cast<unsigned int>(.5 + (height/100.) * 1.75);
1726 cell_label_font_size =
static_cast<unsigned int>(.5 +
1729 * min_level_min_vertex_distance
1730 / std::min(x_dimension, y_dimension)));
1734 additional_width =
static_cast<unsigned int>(.5 + height * .4);
1738 additional_width =
static_cast<unsigned int>(.5 + height * .175);
1749 out <<
"<svg width=\"" << width + additional_width <<
"\" height=\"" << height <<
"\" xmlns=\"http://www.w3.org/2000/svg\" version=\"1.1\">" 1755 out <<
" <linearGradient id=\"background_gradient\" gradientUnits=\"userSpaceOnUse\" x1=\"0\" y1=\"0\" x2=\"0\" y2=\"" << height <<
"\">" <<
'\n' 1756 <<
" <stop offset=\"0\" style=\"stop-color:white\"/>" <<
'\n' 1757 <<
" <stop offset=\"1\" style=\"stop-color:lightsteelblue\"/>" <<
'\n' 1758 <<
" </linearGradient>" <<
'\n';
1764 out <<
"<!-- internal style sheet -->" <<
'\n' 1765 <<
"<style type=\"text/css\"><![CDATA[" <<
'\n';
1770 else out <<
" rect.background{fill:none}" <<
'\n';
1774 <<
" text{font-family:Helvetica; text-anchor:middle; fill:rgb(25,25,25)}" <<
'\n' 1782 unsigned int labeling_index = 0;
1784 for (
unsigned int index = 0; index < n; index++)
1788 if (n != 1) h = .6 - (index / (n-1.)) * .6;
1795 unsigned int i =
static_cast<unsigned int>(h * 6);
1797 double f = h * 6 - i;
1804 r = 255, g =
static_cast<unsigned int>(.5 + 255*t);
1807 r =
static_cast<unsigned int>(.5 + 255*q), g = 255;
1810 g = 255, b =
static_cast<unsigned int>(.5 + 255*t);
1813 g =
static_cast<unsigned int>(.5 + 255*q), b = 255;
1816 r =
static_cast<unsigned int>(.5 + 255*t), b = 255;
1819 r = 255, b =
static_cast<unsigned int>(.5 + 255*q);
1828 while (!materials[labeling_index]) labeling_index++;
1831 while (!levels[labeling_index]) labeling_index++;
1834 while (!subdomains[labeling_index]) labeling_index++;
1837 while (!level_subdomains[labeling_index]) labeling_index++;
1843 out <<
" path.p" << labeling_index
1844 <<
"{fill:rgb(" << r <<
',' << g <<
',' << b <<
"); " 1847 out <<
" path.ps" << labeling_index
1848 <<
"{fill:rgb(" <<
static_cast<unsigned int>(.5 + .75 * r) <<
',' << static_cast<unsigned int>(.5 + .75 * g) <<
',' <<
static_cast<unsigned int>(.5 + .75 * b) <<
"); " 1851 out <<
" rect.r" << labeling_index
1852 <<
"{fill:rgb(" << r <<
',' << g <<
',' << b <<
"); " 1859 out <<
"]]></style>" <<
'\n' <<
'\n';
1862 out <<
" <rect class=\"background\" width=\"" << width <<
"\" height=\"" << height <<
"\"/>" <<
'\n';
1866 unsigned int x_offset = 0;
1868 if (
svg_flags.
margin) x_offset =
static_cast<unsigned int>(.5 + (height/100.) * (margin_in_percent/2.));
1869 else x_offset =
static_cast<unsigned int>(.5 + height * .025);
1871 out <<
" <text x=\"" << x_offset <<
"\" y=\"" <<
static_cast<unsigned int>(.5 + height * .0525) <<
'\"' 1872 <<
" style=\"font-weight:100; fill:lightsteelblue; text-anchor:start; font-family:Courier; font-size:" <<
static_cast<unsigned int>(.5 + height * .045) <<
"px\">" 1873 <<
"deal.II" <<
"</text>" <<
'\n';
1887 out <<
" <!-- cells -->" <<
'\n';
1889 for (
unsigned int level_index = min_level; level_index <= max_level; level_index++)
1892 cell = tria.
begin(level_index),
1893 endc = tria.
end(level_index);
1895 for (; cell != endc; ++cell)
1904 out <<
" class=\"p";
1911 out << (
unsigned int)cell->material_id();
1914 out << (
unsigned int)cell->level();
1918 out << cell->subdomain_id() + 2;
1923 out << cell->level_subdomain_id() + 2;
1934 point[0] = cell->vertex(0)[0];
1935 point[1] = cell->vertex(0)[1];
1943 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
1945 out << static_cast<unsigned int>(.5 + ((projection_decomposition[0] - x_min_perspective) / x_dimension_perspective) * (width - (width/100.) * 2. * margin_in_percent) + ((width/100.) * margin_in_percent)) <<
' ' 1946 <<
static_cast<unsigned int>(.5 + height - (height/100.) * margin_in_percent - ((projection_decomposition[1] - y_min_perspective) / y_dimension_perspective) * (height - (height/100.) * 2. * margin_in_percent));
1950 point[0] = cell->vertex(1)[0];
1951 point[1] = cell->vertex(1)[1];
1953 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
1955 out << static_cast<unsigned int>(.5 + ((projection_decomposition[0] - x_min_perspective) / x_dimension_perspective) * (width - (width/100.) * 2. * margin_in_percent) + ((width/100.) * margin_in_percent)) <<
' ' 1956 <<
static_cast<unsigned int>(.5 + height - (height/100.) * margin_in_percent - ((projection_decomposition[1] - y_min_perspective) / y_dimension_perspective) * (height - (height/100.) * 2. * margin_in_percent));
1960 point[0] = cell->vertex(3)[0];
1961 point[1] = cell->vertex(3)[1];
1963 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
1965 out << static_cast<unsigned int>(.5 + ((projection_decomposition[0] - x_min_perspective) / x_dimension_perspective) * (width - (width/100.) * 2. * margin_in_percent) + ((width/100.) * margin_in_percent)) <<
' ' 1966 <<
static_cast<unsigned int>(.5 + height - (height/100.) * margin_in_percent - ((projection_decomposition[1] - y_min_perspective) / y_dimension_perspective) * (height - (height/100.) * 2. * margin_in_percent));
1970 point[0] = cell->vertex(2)[0];
1971 point[1] = cell->vertex(2)[1];
1973 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
1975 out << static_cast<unsigned int>(.5 + ((projection_decomposition[0] - x_min_perspective) / x_dimension_perspective) * (width - (width/100.) * 2. * margin_in_percent) + ((width/100.) * margin_in_percent)) <<
' ' 1976 <<
static_cast<unsigned int>(.5 + height - (height/100.) * margin_in_percent - ((projection_decomposition[1] - y_min_perspective) / y_dimension_perspective) * (height - (height/100.) * 2. * margin_in_percent));
1980 point[0] = cell->vertex(0)[0];
1981 point[1] = cell->vertex(0)[1];
1983 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
1985 out << static_cast<unsigned int>(.5 + ((projection_decomposition[0] - x_min_perspective) / x_dimension_perspective) * (width - (width/100.) * 2. * margin_in_percent) + ((width/100.) * margin_in_percent)) <<
' ' 1986 <<
static_cast<unsigned int>(.5 + height - (height/100.) * margin_in_percent - ((projection_decomposition[1] - y_min_perspective) / y_dimension_perspective) * (height - (height/100.) * 2. * margin_in_percent));
1988 out <<
"\"/>" <<
'\n';
1993 point[0] = cell->center()[0];
1994 point[1] = cell->center()[1];
2002 float distance_to_camera = sqrt(pow(point[0] - camera_position[0], 2.) + pow(point[1] - camera_position[1], 2.) + pow(point[2] - camera_position[2], 2.));
2003 float distance_factor = distance_to_camera / (2. * std::max(x_dimension, y_dimension));
2005 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
2007 const unsigned int font_size_this_cell =
static_cast<unsigned int>(.5 + cell_label_font_size * pow(.5, (
float)cell->level() - 4. + 3.5 * distance_factor));
2010 <<
" x=\"" <<
static_cast<unsigned int>(.5 + ((projection_decomposition[0] - x_min_perspective) / x_dimension_perspective) * (width - (width/100.) * 2. * margin_in_percent) + ((width/100.) * margin_in_percent))
2011 <<
"\" y=\"" <<
static_cast<unsigned int>(.5 + height - (height/100.) * margin_in_percent - ((projection_decomposition[1] - y_min_perspective) / y_dimension_perspective) * (height - (height/100.) * 2. * margin_in_percent) + 0.5 * font_size_this_cell)
2012 <<
"\" style=\"font-size:" << font_size_this_cell
2017 out << cell->level();
2023 out << cell->index();
2029 out << (int)cell->material_id();
2039 out << static_cast<int>(cell->subdomain_id());
2051 out << static_cast<int>(cell->level_subdomain_id());
2054 out <<
"</text>" <<
'\n';
2060 for (
unsigned int faceIndex = 0; faceIndex < 4; faceIndex++)
2062 if (cell->at_boundary(faceIndex))
2065 point[0] = cell->face(faceIndex)->vertex(0)[0];
2066 point[1] = cell->face(faceIndex)->vertex(0)[1];
2074 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
2076 out <<
" <line x1=\"" 2077 <<
static_cast<unsigned int>(.5 + ((projection_decomposition[0] - x_min_perspective) / x_dimension_perspective) * (width - (width/100.) * 2. * margin_in_percent) + ((width/100.) * margin_in_percent))
2079 <<
static_cast<unsigned int>(.5 + height - (height/100.) * margin_in_percent - ((projection_decomposition[1] - y_min_perspective) / y_dimension_perspective) * (height - (height/100.) * 2. * margin_in_percent));
2081 point[0] = cell->face(faceIndex)->vertex(1)[0];
2082 point[1] = cell->face(faceIndex)->vertex(1)[1];
2090 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
2093 <<
static_cast<unsigned int>(.5 + ((projection_decomposition[0] - x_min_perspective) / x_dimension_perspective) * (width - (width/100.) * 2. * margin_in_percent) + ((width/100.) * margin_in_percent))
2095 <<
static_cast<unsigned int>(.5 + height - (height/100.) * margin_in_percent - ((projection_decomposition[1] - y_min_perspective) / y_dimension_perspective) * (height - (height/100.) * 2. * margin_in_percent))
2107 additional_width = 0;
2108 if (!
svg_flags.
margin) additional_width =
static_cast<unsigned int>(.5 + (height/100.) * 2.5);
2113 unsigned int line_offset = 0;
2114 out <<
" <rect x=\"" << width + additional_width <<
"\" y=\"" <<
static_cast<unsigned int>(.5 + (height/100.) * margin_in_percent)
2115 <<
"\" width=\"" << static_cast<unsigned int>(.5 + (height/100.) * (40. - margin_in_percent)) <<
"\" height=\"" <<
static_cast<unsigned int>(.5 + height * .165) <<
"\"/>" <<
'\n';
2117 out <<
" <text x=\"" << width + additional_width +
static_cast<unsigned int>(.5 + (height/100.) * 1.25)
2118 <<
"\" y=\"" <<
static_cast<unsigned int>(.5 + (height/100.) * margin_in_percent + (++line_offset) * 1.5 * font_size)
2119 <<
"\" style=\"text-anchor:start; font-weight:bold; font-size:" << font_size
2120 <<
"px\">" <<
"cell label" 2121 <<
"</text>" <<
'\n';
2125 out <<
" <text x=\"" << width + additional_width +
static_cast<unsigned int>(.5 + (height/100.) * 2.)
2126 <<
"\" y=\"" <<
static_cast<unsigned int>(.5 + (height/100.) * margin_in_percent + (++line_offset) * 1.5 * font_size)
2127 <<
"\" style=\"text-anchor:start; font-style:oblique; font-size:" << font_size
2128 <<
"px\">" <<
"level_number";
2133 out <<
"</text>" <<
'\n';
2138 out <<
" <text x=\"" << width + additional_width +
static_cast<unsigned int>(.5 + (height/100.) * 2.)
2139 <<
"\" y=\"" <<
static_cast<unsigned int>(.5 + (height/100.) * margin_in_percent + (++line_offset) * 1.5 * font_size )
2140 <<
"\" style=\"text-anchor:start; font-style:oblique; font-size:" << font_size
2147 out <<
"</text>" <<
'\n';
2152 out <<
" <text x=\"" << width + additional_width +
static_cast<unsigned int>(.5 + (height/100.) * 2.)
2153 <<
"\" y=\"" <<
static_cast<unsigned int>(.5 + (height/100.) * margin_in_percent + (++line_offset) * 1.5 * font_size )
2154 <<
"\" style=\"text-anchor:start; font-style:oblique; font-size:" << font_size
2161 out <<
"</text>" <<
'\n';
2166 out <<
" <text x= \"" << width + additional_width +
static_cast<unsigned int>(.5 + (height/100.) * 2.)
2167 <<
"\" y=\"" <<
static_cast<unsigned int>(.5 + (height/100.) * margin_in_percent + (++line_offset) * 1.5 * font_size )
2168 <<
"\" style=\"text-anchor:start; font-style:oblique; font-size:" << font_size
2175 out <<
"</text>" <<
'\n';
2180 out <<
" <text x= \"" << width + additional_width +
static_cast<unsigned int>(.5 + (height/100.) * 2.)
2181 <<
"\" y=\"" <<
static_cast<unsigned int>(.5 + (height/100.) * margin_in_percent + (++line_offset) * 1.5 * font_size )
2182 <<
"\" style=\"text-anchor:start; font-style:oblique; font-size:" << font_size
2184 <<
"level_subdomain_id" 2185 <<
"</text>" <<
'\n';
2192 out <<
" <text x=\"" << width + additional_width
2193 <<
"\" y=\"" <<
static_cast<unsigned int>(.5 + (height/100.) * margin_in_percent + 10.75 * font_size)
2194 <<
"\" style=\"text-anchor:start; font-size:" << font_size <<
"px\">" 2202 out <<
'\n' <<
" <!-- colorbar -->" <<
'\n';
2204 out <<
" <text x=\"" << width + additional_width
2205 <<
"\" y=\"" <<
static_cast<unsigned int>(.5 + (height/100.) * (margin_in_percent + 29.) - (font_size/1.25))
2206 <<
"\" style=\"text-anchor:start; font-weight:bold; font-size:" << font_size <<
"px\">";
2211 out <<
"material_id";
2214 out <<
"level_number";
2217 out <<
"subdomain_id";
2220 out <<
"level_subdomain_id";
2226 out <<
"</text>" <<
'\n';
2228 unsigned int element_height =
static_cast<unsigned int>(((height/100.) * (71. - 2.*margin_in_percent)) / n);
2229 unsigned int element_width =
static_cast<unsigned int>(.5 + (height/100.) * 2.5);
2231 int labeling_index = 0;
2233 for (
unsigned int index = 0; index < n; index++)
2238 while (!materials[labeling_index]) labeling_index++;
2241 while (!levels[labeling_index]) labeling_index++;
2244 while (!subdomains[labeling_index]) labeling_index++;
2247 while (!level_subdomains[labeling_index]) labeling_index++;
2253 out <<
" <rect class=\"r" << labeling_index
2254 <<
"\" x=\"" << width + additional_width
2255 <<
"\" y=\"" <<
static_cast<unsigned int>(.5 + (height/100.) * (margin_in_percent + 29)) + (n-index-1) * element_height
2256 <<
"\" width=\"" << element_width
2257 <<
"\" height=\"" << element_height
2260 out <<
" <text x=\"" << width + additional_width + 1.5 * element_width
2261 <<
"\" y=\"" <<
static_cast<unsigned int>(.5 + (height/100.) * (margin_in_percent + 29)) + (n-index-1 + .5) * element_height +
static_cast<unsigned int>(.5 + font_size * .35) <<
"\"" 2262 <<
" style=\"text-anchor:start; font-size:" <<
static_cast<unsigned int>(.5 + font_size) <<
"px";
2264 if (index == 0 || index == n-1) out <<
"; font-weight:bold";
2266 out <<
"\">" << labeling_index;
2268 if (index == n-1) out <<
" max";
2269 if (index == 0) out <<
" min";
2271 out <<
"</text>" <<
'\n';
2279 out <<
'\n' <<
"</svg>";
2287 std::ostream &)
const 2294 template <
int dim,
int spacedim>
2296 std::ostream &out)
const 2304 const std::time_t time1 = std::time (
nullptr);
2305 const std::tm *time = std::localtime (&time1);
2308 <<
"\n# This file was generated by the deal.II library." 2310 << time->tm_year+1900 <<
"/" 2311 << std::setfill(
'0') << std::setw (2) << time->tm_mon+1 <<
"/" 2312 << std::setfill(
'0') << std::setw (2) << time->tm_mday
2314 << std::setfill(
'0') << std::setw (2) << time->tm_hour <<
":" 2315 << std::setfill(
'0') << std::setw (2) << time->tm_min <<
":" 2316 << std::setfill(
'0') << std::setw (2) << time->tm_sec
2318 <<
"\n# For a description of the MathGL script format see the MathGL manual. " 2320 <<
"\n# Note: This file is understood by MathGL v2.1 and higher only, and can " 2321 <<
"\n# be quickly viewed in a graphical environment using \'mglview\'. " 2328 const std::string axes =
"xyz";
2343 out <<
"\nsetsize 800 800";
2344 out <<
"\nrotate 0 0";
2347 out <<
"\nsetsize 800 800";
2348 out <<
"\nrotate 60 40";
2357 <<
"\n# Vertex ordering." 2358 <<
"\n# list <vertex order> <vertex indices>" 2367 out <<
"\nlist f 0 1 2 3" 2371 out <<
"\nlist f 0 2 4 6 | 1 3 5 7 | 0 4 1 5 | 2 6 3 7 | 0 1 2 3 | 4 5 6 7" 2380 <<
"\n# List of vertices." 2381 <<
"\n# list <id> <vertices>" 2386 typename ::Triangulation<dim, spacedim>::active_cell_iterator
2391 for (; cell!=endc; ++cell)
2393 for (
unsigned int i=0; i<dim; ++i)
2400 out <<
"\nlist " << axes[i] << cell->active_cell_index() <<
" ";
2401 for (
unsigned int j=0; j<GeometryInfo<dim>::vertices_per_cell; ++j)
2402 out << cell->vertex(j)[i] <<
" ";
2409 <<
"\n# List of cells to quadplot." 2410 <<
"\n# quadplot <vertex order> <id> <style>" 2414 out <<
"\nquadplot f ";
2415 for (
unsigned int j=0; j<dim; ++j)
2416 out << axes[j] << i <<
" ";
2441 template <
int dim,
int spacedim,
typename ITERATOR,
typename END>
2444 ITERATOR cell, END end)
2447 for (; cell != end; ++cell)
2453 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
2455 patch.
vertices[v] = cell->vertex(v);
2456 patch.
data(0,v) = cell->level();
2457 patch.
data(1,v) =
static_cast<int>(cell->manifold_id());
2458 patch.
data(2,v) = cell->material_id();
2459 if (!cell->has_children())
2460 patch.
data(3,v) =
static_cast<int>(cell->subdomain_id());
2462 patch.
data(3,v) = -1;
2463 patch.
data(4,v) =
static_cast<int>(cell->level_subdomain_id());
2465 patches.push_back (patch);
2469 std::vector<std::string> triangulation_patch_data_names ()
2471 std::vector<std::string> v(5);
2476 v[4] =
"level_subdomain";
2483 template <
int dim,
int spacedim>
2485 std::ostream &out)
const 2493 std::vector<DataOutBase::Patch<dim,spacedim> > patches;
2495 generate_triangulation_patches(patches, tria.
begin_active(), tria.
end());
2497 triangulation_patch_data_names(),
2498 std::vector<std::tuple<unsigned int, unsigned int, std::string> >(),
2507 template <
int dim,
int spacedim>
2509 std::ostream &out)
const 2517 std::vector<DataOutBase::Patch<dim,spacedim> > patches;
2519 generate_triangulation_patches(patches, tria.
begin_active(), tria.
end());
2521 triangulation_patch_data_names(),
2522 std::vector<std::tuple<unsigned int, unsigned int, std::string> >(),
2531 template <
int dim,
int spacedim>
2533 const std::string &filename_without_extension,
2534 const bool view_levels,
2535 const bool include_artificial)
const 2537 std::vector<DataOutBase::Patch<dim,spacedim> > patches;
2538 const unsigned int n_datasets=4;
2539 std::vector<std::string> data_names;
2540 data_names.emplace_back(
"level");
2541 data_names.emplace_back(
"subdomain");
2542 data_names.emplace_back(
"level_subdomain");
2543 data_names.emplace_back(
"proc_writing");
2548 for (cell=tria.
begin(), endc=tria.
end();
2549 cell != endc; ++cell)
2553 if (cell->has_children())
2555 if (!include_artificial &&
2559 else if (!include_artificial)
2561 if (cell->has_children() &&
2564 else if (!cell->has_children() &&
2574 for (
unsigned int vertex=0; vertex<n_q_points; ++vertex)
2576 patch.
vertices[vertex] = cell->vertex(vertex);
2577 patch.
data(0,vertex) = cell->level();
2578 if (!cell->has_children())
2579 patch.
data(1,vertex) =
2580 (double)static_cast<int>(cell->subdomain_id());
2582 patch.
data(1,vertex) = -1.0;
2583 patch.
data(2,vertex) =
2584 (double)static_cast<int>(cell->level_subdomain_id());
2588 for (
unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
2590 patches.push_back(patch);
2596 std::string new_file = filename_without_extension +
".vtu";
2600 new_file = filename_without_extension +
".proc" +
2605 if (tr->locally_owned_subdomain() == 0)
2607 std::vector<std::string> filenames;
2611 std::size_t pos = filename_without_extension.find_last_of(
'/');
2612 if (pos == std::string::npos)
2617 for (
unsigned int i=0; i<n_procs; ++i)
2618 filenames.push_back (filename_without_extension.substr(pos) +
2622 const std::string pvtu_master_filename = (filename_without_extension +
".pvtu");
2623 std::ofstream pvtu_master (pvtu_master_filename.c_str());
2630 Vector<float> dummy_vector (tr->n_active_cells());
2631 data_out.add_data_vector (dummy_vector,
"level");
2632 data_out.add_data_vector (dummy_vector,
"subdomain");
2633 data_out.add_data_vector (dummy_vector,
"level_subdomain");
2634 data_out.add_data_vector (dummy_vector,
"proc_writing");
2636 data_out.build_patches ();
2638 data_out.write_pvtu_record (pvtu_master, filenames);
2642 std::ofstream out(new_file.c_str());
2643 std::vector<std::tuple<unsigned int, unsigned int, std::string> > vector_data_ranges;
2697 template <
int dim,
int spacedim>
2701 unsigned int n_faces = 0;
2704 face != endf; ++face)
2705 if ((face->at_boundary()) &&
2706 (face->boundary_id() != 0))
2714 template <
int dim,
int spacedim>
2720 std::vector<bool> line_flags;
2722 .save_user_flags_line (line_flags);
2724 .clear_user_flags_line ();
2726 unsigned int n_lines = 0;
2731 cell != endc; ++cell)
2732 for (
unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++l)
2733 if (cell->line(l)->at_boundary()
2735 (cell->line(l)->boundary_id() != 0)
2737 (cell->line(l)->user_flag_set() ==
false))
2740 cell->line(l)->set_user_flag();
2746 .load_user_flags_line (line_flags);
2756 const unsigned int next_element_index,
2757 std::ostream &)
const 2759 return next_element_index;
2765 const unsigned int next_element_index,
2766 std::ostream &)
const 2768 return next_element_index;
2773 const unsigned int next_element_index,
2774 std::ostream &)
const 2776 return next_element_index;
2782 const unsigned int next_element_index,
2783 std::ostream &)
const 2785 return next_element_index;
2790 const unsigned int next_element_index,
2791 std::ostream &)
const 2793 return next_element_index;
2799 const unsigned int next_element_index,
2800 std::ostream &)
const 2802 return next_element_index;
2808 const unsigned int next_element_index,
2809 std::ostream &)
const 2811 return next_element_index;
2816 const unsigned int next_element_index,
2817 std::ostream &)
const 2819 return next_element_index;
2825 template <
int dim,
int spacedim>
2828 const unsigned int next_element_index,
2829 std::ostream &out)
const 2831 unsigned int current_element_index = next_element_index;
2835 face != endf; ++face)
2836 if (face->at_boundary() &&
2837 (face->boundary_id() != 0))
2839 out << current_element_index <<
' ';
2851 out << static_cast<unsigned int>(face->boundary_id())
2853 << static_cast<unsigned int>(face->boundary_id())
2856 for (
unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_face; ++vertex)
2861 ++current_element_index;
2863 return current_element_index;
2867 template <
int dim,
int spacedim>
2870 const unsigned int next_element_index,
2871 std::ostream &out)
const 2873 unsigned int current_element_index = next_element_index;
2878 std::vector<bool> line_flags;
2880 .save_user_flags_line (line_flags);
2882 .clear_user_flags_line ();
2887 cell != endc; ++cell)
2888 for (
unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++l)
2889 if (cell->line(l)->at_boundary()
2891 (cell->line(l)->boundary_id() != 0)
2893 (cell->line(l)->user_flag_set() ==
false))
2895 out << next_element_index <<
" 1 ";
2896 out << static_cast<unsigned int>(cell->line(l)->boundary_id())
2898 << static_cast<unsigned int>(cell->line(l)->boundary_id())
2901 for (
unsigned int vertex=0; vertex<2; ++vertex)
2909 ++current_element_index;
2910 cell->line(l)->set_user_flag();
2916 .load_user_flags_line (line_flags);
2918 return current_element_index;
2926 const unsigned int next_element_index,
2927 std::ostream &)
const 2929 return next_element_index;
2934 const unsigned int next_element_index,
2935 std::ostream &)
const 2937 return next_element_index;
2942 const unsigned int next_element_index,
2943 std::ostream &)
const 2945 return next_element_index;
2950 const unsigned int next_element_index,
2951 std::ostream &)
const 2953 return next_element_index;
2958 const unsigned int next_element_index,
2959 std::ostream &)
const 2961 return next_element_index;
2967 const unsigned int next_element_index,
2968 std::ostream &)
const 2970 return next_element_index;
2976 const unsigned int next_element_index,
2977 std::ostream &)
const 2979 return next_element_index;
2984 const unsigned int next_element_index,
2985 std::ostream &)
const 2987 return next_element_index;
2992 template <
int dim,
int spacedim>
2995 const unsigned int next_element_index,
2996 std::ostream &out)
const 2998 unsigned int current_element_index = next_element_index;
3002 face != endf; ++face)
3003 if (face->at_boundary() &&
3004 (face->boundary_id() != 0))
3006 out << current_element_index <<
" " 3007 <<
static_cast<unsigned int>(face->boundary_id())
3021 for (
unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_face; ++vertex)
3025 ++current_element_index;
3027 return current_element_index;
3032 template <
int dim,
int spacedim>
3035 const unsigned int next_element_index,
3036 std::ostream &out)
const 3038 unsigned int current_element_index = next_element_index;
3043 std::vector<bool> line_flags;
3045 .save_user_flags_line (line_flags);
3047 .clear_user_flags_line ();
3052 cell != endc; ++cell)
3053 for (
unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++l)
3054 if (cell->line(l)->at_boundary()
3056 (cell->line(l)->boundary_id() != 0)
3058 (cell->line(l)->user_flag_set() ==
false))
3060 out << current_element_index <<
" " 3061 <<
static_cast<unsigned int>(cell->line(l)->boundary_id())
3064 for (
unsigned int vertex=0; vertex<2; ++vertex)
3072 ++current_element_index;
3073 cell->line(l)->set_user_flag();
3079 .load_user_flags_line (line_flags);
3080 return current_element_index;
3088 camera_vertical[0] = camera_horizontal[1] * camera_direction[2] - camera_horizontal[2] * camera_direction[1];
3089 camera_vertical[1] = camera_horizontal[2] * camera_direction[0] - camera_horizontal[0] * camera_direction[2];
3090 camera_vertical[2] = camera_horizontal[0] * camera_direction[1] - camera_horizontal[1] * camera_direction[0];
3094 phi /= (point[0] - camera_position[0]) * camera_direction[0] + (point[1] - camera_position[1]) * camera_direction[1] + (point[2] - camera_position[2]) * camera_direction[2];
3097 projection[0] = camera_position[0] + phi * (point[0] - camera_position[0]);
3098 projection[1] = camera_position[1] + phi * (point[1] - camera_position[1]);
3099 projection[2] = camera_position[2] + phi * (point[2] - camera_position[2]);
3102 projection_decomposition[0] = (projection[0] - camera_position[0] - camera_focus * camera_direction[0]) * camera_horizontal[0];
3103 projection_decomposition[0] += (projection[1] - camera_position[1] - camera_focus * camera_direction[1]) * camera_horizontal[1];
3104 projection_decomposition[0] += (projection[2] - camera_position[2] - camera_focus * camera_direction[2]) * camera_horizontal[2];
3106 projection_decomposition[1] = (projection[0] - camera_position[0] - camera_focus * camera_direction[0]) * camera_vertical[0];
3107 projection_decomposition[1] += (projection[1] - camera_position[1] - camera_focus * camera_direction[1]) * camera_vertical[1];
3108 projection_decomposition[1] += (projection[2] - camera_position[2] - camera_focus * camera_direction[2]) * camera_vertical[2];
3110 return projection_decomposition;
3119 template <
int spacedim>
3120 void write_gnuplot (const ::Triangulation<1,spacedim> &tria,
3129 typename ::Triangulation<dim,spacedim>::active_cell_iterator
3130 cell=tria.begin_active();
3131 const typename ::Triangulation<dim,spacedim>::active_cell_iterator
3133 for (; cell!=endc; ++cell)
3136 out <<
"# cell " << cell <<
'\n';
3138 out << cell->vertex(0)
3139 <<
' ' << cell->level()
3140 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3142 <<
' ' << cell->level()
3143 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3156 template <
int spacedim>
3157 void write_gnuplot (const ::Triangulation<2,spacedim> &tria,
3166 const unsigned int n_additional_points=
3168 const unsigned int n_points=2+n_additional_points;
3170 typename ::Triangulation<dim,spacedim>::active_cell_iterator
3171 cell=tria.begin_active();
3172 const typename ::Triangulation<dim,spacedim>::active_cell_iterator
3181 std::vector<
Point<dim-1> > boundary_points;
3182 if (mapping!=
nullptr)
3184 boundary_points.resize(n_points);
3185 boundary_points[0][0]=0;
3186 boundary_points[n_points-1][0]=1;
3187 for (
unsigned int i=1; i<n_points-1; ++i)
3188 boundary_points[i](0)= 1.*i/(n_points-1);
3190 std::vector<double> dummy_weights(n_points, 1./n_points);
3191 Quadrature<dim-1> quadrature(boundary_points, dummy_weights);
3196 for (; cell!=endc; ++cell)
3199 out <<
"# cell " << cell <<
'\n';
3201 if (mapping==
nullptr ||
3211 for (
unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
3213 <<
' ' << cell->level()
3214 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n';
3215 out << cell->vertex(0)
3216 <<
' ' << cell->level()
3217 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3228 for (
unsigned int face_no=0;
3229 face_no<GeometryInfo<dim>::faces_per_cell; ++face_no)
3231 const typename ::Triangulation<dim,spacedim>::face_iterator
3232 face = cell->face(face_no);
3240 const unsigned int offset=face_no*n_points;
3241 for (
unsigned int i=0; i<n_points; ++i)
3243 (cell, q_projector->
point(offset+i)))
3244 <<
' ' << cell->level()
3245 <<
' ' <<
static_cast<unsigned int>(cell->material_id())
3258 out << face->vertex(0)
3259 <<
' ' << cell->level()
3260 <<
' ' <<
static_cast<unsigned int>(cell->material_id())
3263 <<
' ' << cell->level()
3264 <<
' ' <<
static_cast<unsigned int>(cell->material_id())
3273 if (q_projector !=
nullptr)
3285 template <
int spacedim>
3286 void write_gnuplot (const ::Triangulation<3,spacedim> &tria,
3295 const unsigned int n_additional_points=
3297 const unsigned int n_points=2+n_additional_points;
3299 typename ::Triangulation<dim,spacedim>::active_cell_iterator
3300 cell=tria.begin_active();
3301 const typename ::Triangulation<dim,spacedim>::active_cell_iterator
3310 std::vector<Point<1> > boundary_points;
3311 if (mapping!=
nullptr)
3313 boundary_points.resize(n_points);
3314 boundary_points[0][0]=0;
3315 boundary_points[n_points-1][0]=1;
3316 for (
unsigned int i=1; i<n_points-1; ++i)
3317 boundary_points[i](0)= 1.*i/(n_points-1);
3319 std::vector<double> dummy_weights(n_points, 1./n_points);
3324 QIterated<dim-1> quadrature(quadrature1d, 1);
3328 for (; cell!=endc; ++cell)
3331 out <<
"# cell " << cell <<
'\n';
3333 if (mapping==
nullptr || n_points==2 ||
3337 out << cell->vertex(0)
3338 <<
' ' << cell->level()
3339 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3341 <<
' ' << cell->level()
3342 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3344 <<
' ' << cell->level()
3345 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3347 <<
' ' << cell->level()
3348 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3350 <<
' ' << cell->level()
3351 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3354 out << cell->vertex(2)
3355 <<
' ' << cell->level()
3356 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3358 <<
' ' << cell->level()
3359 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3361 <<
' ' << cell->level()
3362 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3364 <<
' ' << cell->level()
3365 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3367 <<
' ' << cell->level()
3368 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3372 out << cell->vertex(0)
3373 <<
' ' << cell->level()
3374 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3376 <<
' ' << cell->level()
3377 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3379 out << cell->vertex(1)
3380 <<
' ' << cell->level()
3381 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3383 <<
' ' << cell->level()
3384 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3386 out << cell->vertex(5)
3387 <<
' ' << cell->level()
3388 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3390 <<
' ' << cell->level()
3391 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3393 out << cell->vertex(4)
3394 <<
' ' << cell->level()
3395 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3397 <<
' ' << cell->level()
3398 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3403 for (
unsigned int face_no=0; face_no<GeometryInfo<dim>::faces_per_cell; ++face_no)
3405 const typename ::Triangulation<dim,spacedim>::face_iterator
3406 face = cell->face(face_no);
3408 if (face->at_boundary())
3410 const unsigned int offset=face_no*n_points*n_points;
3411 for (
unsigned int i=0; i<n_points-1; ++i)
3412 for (
unsigned int j=0; j<n_points-1; ++j)
3415 cell, q_projector->
point(offset+i*n_points+j));
3417 <<
' ' << cell->level()
3418 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n';
3420 cell, q_projector->
point(offset+(i+1)*n_points+j)))
3421 <<
' ' << cell->level()
3422 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n';
3424 cell, q_projector->
point(offset+(i+1)*n_points+j+1)))
3425 <<
' ' << cell->level()
3426 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n';
3428 cell, q_projector->
point(offset+i*n_points+j+1)))
3429 <<
' ' << cell->level()
3430 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n';
3436 <<
' ' << cell->level()
3437 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n';
3438 out <<
'\n' <<
'\n';
3443 for (
unsigned int l=0; l<GeometryInfo<dim>::lines_per_face; ++
l)
3445 const typename ::Triangulation<dim,spacedim>::line_iterator
3449 &v1=line->vertex(1);
3463 for (
unsigned int i=0; i<n_points; ++i)
3465 (cell, (1-boundary_points[i][0])*u0+boundary_points[i][0]*u1))
3466 <<
' ' << cell->level()
3467 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n';
3471 <<
' ' << cell->level()
3472 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3474 <<
' ' << cell->level()
3475 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n';
3477 out <<
'\n' <<
'\n';
3484 if (q_projector !=
nullptr)
3499 template <
int dim,
int spacedim>
3505 internal::write_gnuplot (tria, out, mapping,
gnuplot_flags);
3523 const unsigned int l)
3525 first(f), second(s),
3526 colorize(c), level(l)
3531 void write_eps (const ::Triangulation<1> &,
3540 void write_eps (const ::Triangulation<1,2> &,
3549 void write_eps (const ::Triangulation<1,3> &,
3558 void write_eps (const ::Triangulation<2,3> &,
3569 template <
int dim,
int spacedim>
3570 void write_eps (const ::Triangulation<dim, spacedim> &tria,
3576 typedef std::list<LineEntry> LineList;
3588 static_cast<const GridOutFlags::EpsFlagsBase &>(eps_flags_3);
3616 for (typename ::Triangulation<dim, spacedim>::active_cell_iterator
3617 cell=tria.begin_active();
3618 cell!=tria.end(); ++cell)
3619 for (
unsigned int line_no=0;
3620 line_no<GeometryInfo<dim>::lines_per_cell; ++line_no)
3622 typename ::Triangulation<dim, spacedim>::line_iterator
3623 line=cell->line(line_no);
3635 if (!line->has_children() &&
3636 (mapping==
nullptr || !line->at_boundary()))
3656 line_list.emplace_back(
Point<2>(line->vertex(0)(0),
3657 line->vertex(0)(1)),
3659 line->vertex(1)(1)),
3660 line->user_flag_set(),
3670 if (mapping!=
nullptr)
3677 std::vector<
Point<dim-1> > boundary_points (n_points);
3679 for (
unsigned int i=0; i<n_points; ++i)
3680 boundary_points[i](0) = 1.*(i+1)/(n_points+1);
3682 Quadrature<dim-1> quadrature (boundary_points);
3689 for (typename ::Triangulation<dim, spacedim>::active_cell_iterator
3690 cell=tria.begin_active();
3691 cell!=tria.end(); ++cell)
3692 for (
unsigned int face_no=0; face_no<GeometryInfo<dim>::faces_per_cell; ++face_no)
3694 const typename ::Triangulation<dim, spacedim>::face_iterator
3695 face = cell->face(face_no);
3697 if (face->at_boundary())
3700 Point<2> p0 (p0_dim(0), p0_dim(1));
3707 const unsigned int offset=face_no*n_points;
3708 for (
unsigned int i=0; i<n_points; ++i)
3711 (cell, q_projector.
point(offset+i)));
3712 const Point<2> p1 (p1_dim(0), p1_dim(1));
3714 line_list.emplace_back (p0,
3716 face->user_flag_set(),
3723 const Point<2> p1 (p1_dim(0), p1_dim(1));
3724 line_list.emplace_back (p0,
3726 face->user_flag_set(),
3741 typename ::Triangulation<dim, spacedim>::active_cell_iterator
3742 cell=tria.begin_active(),
3764 const double turn_angle = eps_flags_3.
turn_angle;
3765 const Point<dim> view_direction(-std::sin(z_angle * 2.*pi / 360.) * std::sin(turn_angle * 2.*pi / 360.),
3766 +std::sin(z_angle * 2.*pi / 360.) * std::cos(turn_angle * 2.*pi / 360.),
3767 -std::cos(z_angle * 2.*pi / 360.));
3783 - ((
Point<dim>(1,0,0) * view_direction) * view_direction)
3784 - ((
Point<dim>(1,0,0) * unit_vector1) * unit_vector1));
3788 for (; cell!=endc; ++cell)
3789 for (
unsigned int line_no=0;
3790 line_no<GeometryInfo<dim>::lines_per_cell; ++line_no)
3792 typename ::Triangulation<dim, spacedim>::line_iterator
3793 line=cell->line(line_no);
3794 line_list.emplace_back (
Point<2>(line->vertex(0) * unit_vector2,
3795 line->vertex(0) * unit_vector1),
3796 Point<2>(line->vertex(1) * unit_vector2,
3797 line->vertex(1) * unit_vector1),
3798 line->user_flag_set(),
3814 double x_min = tria.begin_active()->vertex(0)(0);
3815 double x_max = x_min;
3816 double y_min = tria.begin_active()->vertex(0)(1);
3817 double y_max = y_min;
3818 unsigned int max_level = line_list.begin()->level;
3820 for (LineList::const_iterator line=line_list.begin();
3821 line!=line_list.end(); ++line)
3823 x_min = std::min (x_min, line->first(0));
3824 x_min = std::min (x_min, line->second(0));
3826 x_max = std::max (x_max, line->first(0));
3827 x_max = std::max (x_max, line->second(0));
3829 y_min = std::min (y_min, line->first(1));
3830 y_min = std::min (y_min, line->second(1));
3832 y_max = std::max (y_max, line->first(1));
3833 y_max = std::max (y_max, line->second(1));
3835 max_level = std::max (max_level, line->level);
3843 const double scale = (eps_flags_base.
size /
3855 std::time_t time1= std::time (
nullptr);
3856 std::tm *time = std::localtime(&time1);
3857 out <<
"%!PS-Adobe-2.0 EPSF-1.2" <<
'\n' 3858 <<
"%%Title: deal.II Output" <<
'\n' 3859 <<
"%%Creator: the deal.II library" <<
'\n' 3860 <<
"%%Creation Date: " 3861 << time->tm_year+1900 <<
"/" 3862 << time->tm_mon+1 <<
"/" 3863 << time->tm_mday <<
" - " 3864 << time->tm_hour <<
":" 3865 << std::setw(2) << time->tm_min <<
":" 3866 << std::setw(2) << time->tm_sec <<
'\n' 3867 <<
"%%BoundingBox: " 3871 <<
static_cast<unsigned int>(std::floor(( (x_max-x_min) * scale )+1))
3873 << static_cast<unsigned int>(std::floor(( (y_max-y_min) * scale )+1))
3882 out <<
"/m {moveto} bind def" <<
'\n' 3883 <<
"/x {lineto stroke} bind def" <<
'\n' 3884 <<
"/b {0 0 0 setrgbcolor} def" <<
'\n' 3885 <<
"/r {1 0 0 setrgbcolor} def" <<
'\n';
3895 << (0.66666/std::max(1U,(max_level-1)))
3896 <<
" mul 1 0.8 sethsbcolor} def" <<
'\n';
3909 out << (
"/R {rmoveto} bind def\n" 3910 "/Symbol-Oblique /Symbol findfont [1 0 .167 1 0 0] makefont\n" 3911 "dup length dict begin {1 index /FID eq {pop pop} {def} ifelse} forall\n" 3912 "currentdict end definefont\n" 3913 "/MFshow {{dup dup 0 get findfont exch 1 get scalefont setfont\n" 3914 "[ currentpoint ] exch dup 2 get 0 exch rmoveto dup dup 5 get exch 4 get\n" 3915 "{show} {stringwidth pop 0 rmoveto}ifelse dup 3 get\n" 3916 "{2 get neg 0 exch rmoveto pop} {pop aload pop moveto}ifelse} forall} bind def\n" 3917 "/MFwidth {0 exch {dup 3 get{dup dup 0 get findfont exch 1 get scalefont setfont\n" 3918 "5 get stringwidth pop add}\n" 3919 "{pop} ifelse} forall} bind def\n" 3920 "/MCshow { currentpoint stroke m\n" 3921 "exch dup MFwidth -2 div 3 -1 roll R MFshow } def\n")
3925 out <<
"%%EndProlog" <<
'\n' 3929 out << eps_flags_base.
line_width <<
" setlinewidth" <<
'\n';
3933 const Point<2> offset(x_min, y_min);
3935 for (LineList::const_iterator line=line_list.begin();
3936 line!=line_list.end(); ++line)
3944 << (line->first - offset) *
scale <<
" m " 3945 << (line->second - offset) * scale <<
" x" <<
'\n';
3948 << (line->first - offset) * scale <<
" m " 3949 << (line->second - offset) *
scale <<
" x" <<
'\n';
3955 out <<
"(Helvetica) findfont 140 scalefont setfont" 3958 typename ::Triangulation<dim, spacedim>::active_cell_iterator
3959 cell = tria.begin_active (),
3961 for (; cell!=endc; ++cell)
3963 out << (cell->center()(0)-offset(0))*scale <<
' ' 3964 << (cell->center()(1)-offset(1))*scale
3966 <<
"[ [(Helvetica) 12.0 0.0 true true (";
3970 out << cell->index();
3981 out <<
"(Helvetica) findfont 140 scalefont setfont" 3988 std::set<unsigned int> treated_vertices;
3989 typename ::Triangulation<dim, spacedim>::active_cell_iterator
3990 cell = tria.begin_active (),
3992 for (; cell!=endc; ++cell)
3993 for (
unsigned int vertex=0;
3994 vertex<GeometryInfo<dim>::vertices_per_cell;
3996 if (treated_vertices.find(cell->vertex_index(vertex))
3998 treated_vertices.end())
4000 treated_vertices.insert (cell->vertex_index(vertex));
4002 out << (cell->vertex(vertex)(0)-offset(0))*scale <<
' ' 4003 << (cell->vertex(vertex)(1)-offset(1))*scale
4005 <<
"[ [(Helvetica) 10.0 0.0 true true (" 4006 << cell->vertex_index(vertex)
4013 out <<
"showpage" <<
'\n';
4025 template <
int dim,
int spacedim>
4030 internal::write_eps (tria, out, mapping,
4035 template <
int dim,
int spacedim>
4041 switch (output_format)
4091 template <
int dim,
int spacedim>
4101 #include "grid_out.inst" 4104 DEAL_II_NAMESPACE_CLOSE
void parse_parameters(ParameterHandler ¶m)
void parse_parameters(ParameterHandler ¶m)
unsigned int n_boundary_faces(const Triangulation< dim, spacedim > &tria) const
unsigned int n_active_cells() const
long int get_integer(const std::string &entry_string) const
unsigned int n_used_vertices() const
static const unsigned int invalid_unsigned_int
OutputFormat default_format
static void declare_parameters(ParameterHandler &prm)
void write_vtk(const std::vector< Patch< dim, spacedim > > &patches, const std::vector< std::string > &data_names, const std::vector< std::tuple< unsigned int, unsigned int, std::string > > &vector_data_ranges, const VtkFlags &flags, std::ostream &out)
DX(const bool write_cells=true, const bool write_faces=false, const bool write_diameter=false, const bool write_measure=false, const bool write_all_faces=true)
active_face_iterator begin_active_face() const
static OutputFormat parse_output_format(const std::string &format_name)
void parse_parameters(ParameterHandler ¶m)
bool margin
Margin around the plotted area.
Gnuplot(const bool write_cell_number=false, const unsigned int n_boundary_face_points=2, const bool curved_inner_cells=false)
void attach_triangulation(const Triangulation< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &)
void write_xfig(const Triangulation< dim, spacedim > &tria, std::ostream &out, const Mapping< dim, spacedim > *mapping=nullptr) const
static ::ExceptionBase & ExcIO()
unsigned int n_boundary_face_points
static void declare_parameters(ParameterHandler ¶m)
GridOutFlags::Eps< 2 > eps_flags_2
unsigned int height
Height of the plot in SVG units, computed from width if zero. Defaults to 1000.
bool draw_legend
Draw a legend next to the plotted grid, explaining the label of the cells.
void write_vtu(const Triangulation< dim, spacedim > &tria, std::ostream &out) const
std::string get(const std::string &entry_string) const
unsigned int neighbors[dim > 0 ? GeometryInfo< dim >::faces_per_cell :1]
void parse_parameters(ParameterHandler ¶m)
void write_mesh_per_processor_as_vtu(const Triangulation< dim, spacedim > &tria, const std::string &filename_without_extension, const bool view_levels=false, const bool include_artificial=false) const
const Point< dim > & point(const unsigned int i) const
static void declare_parameters(ParameterHandler ¶m)
Convert the level number into the cell color.
write() calls write_eps()
active_cell_iterator begin_active(const unsigned int level=0) const
unsigned int line_thickness
Thickness of the lines between cells.
unsigned int n_boundary_face_points
bool convert_level_number_to_height
Interpret the level number of the cells as altitude over the x-y-plane (useful in the perspective vie...
#define AssertThrow(cond, exc)
numbers::NumberTraits< Number >::real_type norm() const
static Point< 2 > svg_project_point(Point< 3 > point, Point< 3 > camera_position, Point< 3 > camera_direction, Point< 3 > camera_horizontal, float camera_focus)
void write_mathgl(const Triangulation< dim, spacedim > &tria, std::ostream &out) const
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const =0
unsigned int boundary_line_thickness
Thickness of lines at the boundary.
Point< spacedim > vertices[GeometryInfo< dim >::vertices_per_cell]
cell_iterator begin(const unsigned int level=0) const
unsigned int write_msh_faces(const Triangulation< dim, spacedim > &tria, const unsigned int next_element_index, std::ostream &out) const
void parse_parameters(ParameterHandler ¶m)
void parse_parameters(ParameterHandler ¶m)
write() calls write_ucd()
void write_gnuplot(const Triangulation< dim, spacedim > &tria, std::ostream &out, const Mapping< dim, spacedim > *mapping=nullptr) const
cell_iterator end() const
unsigned int n_boundary_lines(const Triangulation< dim, spacedim > &tria) const
void enter_subsection(const std::string &subsection)
GridOutFlags::Vtk vtk_flags
bool label_cell_index
Write cell index into each cell. Defaults to true.
GridOutFlags::Gnuplot gnuplot_flags
void write_svg(const Triangulation< 2, 2 > &tria, std::ostream &out) const
Convert the global subdomain id into the cell color.
Convert the material id into the cell color.
static ::ExceptionBase & ExcMessage(std::string arg1)
static void declare_parameters(ParameterHandler ¶m)
std::size_t memory_consumption() const
double get_double(const std::string &entry_name) const
GridOutFlags::MathGL mathgl_flags
write() calls write_mathgl()
GridOutFlags::DX dx_flags
GridOutFlags::Msh msh_flags
GridOutFlags::Eps< 1 > eps_flags_1
void reinit(const TableIndices< N > &new_size, const bool omit_default_initialization=false)
write() calls write_gnuplot()
#define Assert(cond, exc)
void write_gnuplot(const std::vector< Patch< dim, spacedim > > &patches, const std::vector< std::string > &data_names, const std::vector< std::tuple< unsigned int, unsigned int, std::string > > &vector_data_ranges, const GnuplotFlags &flags, std::ostream &out)
bool write_cell_number_level
void parse_parameters(const ParameterHandler &prm)
bool label_material_id
Write material id of each cell. Defaults to false.
EpsFlagsBase(const SizeType size_type=width, const unsigned int size=300, const double line_width=0.5, const bool color_lines_on_user_flag=false, const unsigned int n_boundary_face_points=2, const bool color_lines_level=false)
Abstract base class for mapping classes.
GridOutFlags::XFig xfig_flags
bool label_subdomain_id
Write subdomain id of each cell. Defaults to false.
Use a gradient from white (top) to steelblue (bottom), and add date and time plus a deal...
Ucd(const bool write_preamble=false, const bool write_faces=false, const bool write_lines=false)
const std::vector< Point< spacedim > > & get_vertices() const
static void declare_parameters(ParameterHandler ¶m)
bool label_level_subdomain_id
Write level subdomain id of each cell. Defaults to false.
static std::string get_output_format_names()
Convert the subdomain id into the cell color.
Msh(const bool write_faces=false, const bool write_lines=false)
unsigned int n_subdivisions
bool get_bool(const std::string &entry_name) const
void write_vtk(const Triangulation< dim, spacedim > &tria, std::ostream &out) const
GridOutFlags::Ucd ucd_flags
static void declare_parameters(ParameterHandler ¶m)
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)
static void declare_parameters(ParameterHandler ¶m)
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
float cell_font_scaling
Scaling of the font for cell annotations. Defaults to 1.
Convert the level subdomain id into the cell color.
write() calls write_xfig()
bool write_vertex_numbers
Convert the level subdomain id into the cell color.
static Quadrature< dim > project_to_all_faces(const SubQuadrature &quadrature)
const types::subdomain_id artificial_subdomain_id
Convert the level into the cell color.
void parse_parameters(ParameterHandler ¶m)
unsigned int n_boundary_face_points
static void declare_parameters(ParameterHandler ¶m)
unsigned int write_ucd_faces(const Triangulation< dim, spacedim > &tria, const unsigned int next_element_index, std::ostream &out) const
write() calls write_msh()
const std::vector< bool > & get_used_vertices() const
void write_ucd(const Triangulation< dim, spacedim > &tria, std::ostream &out) const
unsigned int write_ucd_lines(const Triangulation< dim, spacedim > &tria, const unsigned int next_element_index, std::ostream &out) const
write() calls write_svg()
unsigned int write_msh_lines(const Triangulation< dim, spacedim > &tria, const unsigned int next_element_index, std::ostream &out) const
GridOutFlags::Svg svg_flags
Svg(const unsigned int line_thickness=2, const unsigned int boundary_line_thickness=4, bool margin=true, const Background background=white, const int azimuth_angle=0, const int polar_angle=0, const Coloring coloring=level_number, const bool convert_level_number_to_height=false, const bool label_level_number=true, const bool label_cell_index=true, const bool label_material_id=false, const bool label_subdomain_id=false, const bool draw_colorbar=true, const bool draw_legend=true)
bool color_lines_on_user_flag
void parse_parameters(ParameterHandler ¶m)
void write_eps(const Triangulation< dim, spacedim > &tria, std::ostream &out, const Mapping< dim, spacedim > *mapping=nullptr) const
void declare_entry(const std::string &entry, const std::string &default_value, const Patterns::PatternBase &pattern=Patterns::Anything(), const std::string &documentation=std::string())
static ::ExceptionBase & ExcNotImplemented()
Convert the material id into the cell color (default)
float level_height_factor
The factor determining the vertical distance between levels (default = 0.3)
void write_eps(const std::vector< Patch< 2, spacedim > > &patches, const std::vector< std::string > &data_names, const std::vector< std::tuple< unsigned int, unsigned int, std::string > > &vector_data_ranges, const EpsFlags &flags, std::ostream &out)
GridOutFlags::Eps< 3 > eps_flags_3
face_iterator end_face() const
static ::ExceptionBase & ExcInvalidState()
const types::boundary_id internal_face_boundary_id
write() calls write_vtu()
void write_vtu(const std::vector< Patch< dim, spacedim > > &patches, const std::vector< std::string > &data_names, const std::vector< std::tuple< unsigned int, unsigned int, std::string > > &vector_data_ranges, const VtkFlags &flags, std::ostream &out)
bool points_are_available
static void declare_parameters(ParameterHandler ¶m)
void set_flags(const GridOutFlags::DX &flags)
bool label_level_number
Write level number into each cell. Defaults to true.
void parse_parameters(ParameterHandler ¶m)
bool draw_colorbar
Draw a colorbar next to the plotted grid with respect to the chosen coloring of the cells...
void write_dx(const Triangulation< dim, spacedim > &tria, std::ostream &out) const
void write(const Triangulation< dim, spacedim > &tria, std::ostream &out, const OutputFormat output_format, const Mapping< dim, spacedim > *mapping=nullptr) const
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const =0
GridOutFlags::Vtu vtu_flags
unsigned int width
The width of the plot. Computed automatically from height if zero (default)
write() calls write_vtk()
virtual types::subdomain_id locally_owned_subdomain() const
void parse_parameters(ParameterHandler ¶m)
void parse_parameters(ParameterHandler ¶m)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
static ::ExceptionBase & ExcInternalError()
void write_msh(const Triangulation< dim, spacedim > &tria, std::ostream &out) const
std::string default_suffix() const