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/grid/tria.h> 23 #include <deal.II/grid/tria_accessor.h> 24 #include <deal.II/grid/tria_iterator.h> 25 #include <deal.II/fe/mapping.h> 36 DEAL_II_NAMESPACE_OPEN
42 const bool write_faces,
43 const bool write_diameter,
44 const bool write_measure,
45 const bool write_all_faces) :
46 write_cells (write_cells),
47 write_faces (write_faces),
48 write_diameter (write_diameter),
49 write_measure (write_measure),
50 write_all_faces (write_all_faces)
56 "Write the mesh connectivity as DX grid cells");
58 "Write faces of cells. These may be boundary faces " 59 "or all faces between mesh cells, according to " 60 "\"Write all faces\"");
62 "If cells are written, additionally write their" 63 " diameter as data for visualization");
65 "Write the volume of each cell as data");
67 "Write all faces, not only boundary");
81 const bool write_lines) :
82 write_faces (write_faces),
83 write_lines (write_lines)
103 write_preamble (write_preamble),
104 write_faces (write_faces),
105 write_lines (write_lines)
127 const unsigned int n_boundary_face_points,
128 const bool curved_inner_cells) :
129 write_cell_numbers (write_cell_numbers),
130 n_boundary_face_points(n_boundary_face_points),
131 curved_inner_cells(curved_inner_cells)
151 const unsigned int size,
152 const double line_width,
153 const bool color_lines_on_user_flag,
155 const bool color_lines_level) :
156 size_type (size_type),
158 line_width (line_width),
159 color_lines_on_user_flag(color_lines_on_user_flag),
160 n_boundary_face_points(n_boundary_face_points),
161 color_lines_level(color_lines_level)
169 "Depending on this parameter, either the" 171 "of the eps is scaled to \"Size\"");
173 "Size of the output in points");
175 "Width of the lines drawn in points");
177 "Draw lines with user flag set in different color");
179 "Number of points on boundary edges. " 180 "Increase this beyond 2 to see curved boundaries.");
182 "Draw different colors according to grid level.");
188 if (param.
get(
"Size by") == std::string(
"width"))
190 else if (param.
get(
"Size by") == std::string(
"height"))
202 const unsigned int size,
208 color_lines_on_user_flag,
209 n_boundary_face_points)
225 const unsigned int size,
226 const double line_width,
227 const bool color_lines_on_user_flag,
228 const unsigned int n_boundary_face_points,
229 const bool write_cell_numbers,
230 const bool write_cell_number_level,
231 const bool write_vertex_numbers,
232 const bool color_lines_level
236 color_lines_on_user_flag,
237 n_boundary_face_points,
239 write_cell_numbers (write_cell_numbers),
240 write_cell_number_level (write_cell_number_level),
241 write_vertex_numbers (write_vertex_numbers)
248 "(2D only) Write cell numbers" 249 " into the centers of cells");
251 "(2D only) if \"Cell number\" is true, write" 252 "numbers in the form level.number");
254 "Write numbers for each vertex");
261 write_cell_numbers = param.
get_bool(
"Cell number");
262 write_cell_number_level = param.
get_bool(
"Level number");
263 write_vertex_numbers = param.
get_bool(
"Vertex number");
269 const unsigned int size,
270 const double line_width,
271 const bool color_lines_on_user_flag,
272 const unsigned int n_boundary_face_points,
273 const double azimut_angle,
274 const double turn_angle)
277 color_lines_on_user_flag,
278 n_boundary_face_points),
279 azimut_angle (azimut_angle),
280 turn_angle (turn_angle)
287 "Azimuth of the viw point, that is, the angle " 288 "in the plane from the x-axis.");
290 "Elevation of the view point above the xy-plane.");
297 azimut_angle = 90- param.
get_double(
"Elevation");
306 color_by(material_id),
308 n_boundary_face_points(0),
314 boundary_thickness(3)
346 const unsigned int boundary_line_thickness,
349 const int azimuth_angle,
350 const int polar_angle,
352 const bool convert_level_number_to_height,
353 const bool label_level_number,
354 const bool label_cell_index,
355 const bool label_material_id,
356 const bool label_subdomain_id,
357 const bool draw_colorbar,
358 const bool draw_legend)
362 line_thickness(line_thickness),
363 boundary_line_thickness(boundary_line_thickness),
365 background(background),
366 azimuth_angle(azimuth_angle),
367 polar_angle(polar_angle),
369 convert_level_number_to_height(convert_level_number_to_height),
370 level_height_factor(0.3f),
371 cell_font_scaling(1.f),
372 label_level_number(label_level_number),
373 label_cell_index(label_cell_index),
374 label_material_id(label_material_id),
375 label_subdomain_id(label_subdomain_id),
376 label_level_subdomain_id(false),
377 draw_colorbar(draw_colorbar),
378 draw_legend(draw_legend)
383 draw_bounding_box (false)
401 default_format (none)
483 switch (output_format)
508 Assert (
false, ExcNotImplemented());
526 if (format_name ==
"none" || format_name ==
"false")
529 if (format_name ==
"dx")
532 if (format_name ==
"ucd")
535 if (format_name ==
"gnuplot")
538 if (format_name ==
"eps")
541 if (format_name ==
"xfig")
544 if (format_name ==
"msh")
547 if (format_name ==
"svg")
550 if (format_name ==
"mathgl")
553 if (format_name ==
"vtk")
556 if (format_name ==
"vtu")
568 return "none|dx|gnuplot|eps|ucd|xfig|msh|svg|mathgl|vtk|vtu";
687 std::ostream &)
const 689 Assert (
false, ExcNotImplemented());
694 std::ostream &)
const 696 Assert (
false, ExcNotImplemented());
701 std::ostream &)
const 703 Assert (
false, ExcNotImplemented());
708 template <
int dim,
int spacedim>
710 std::ostream &out)
const 716 const std::vector<Point<spacedim> > &vertices = tria.
get_vertices();
726 std::vector<unsigned int> renumber(vertices.size());
729 unsigned int new_number=0;
730 for (
unsigned int i=0; i<vertices.size(); ++i)
732 renumber[i]=new_number++;
733 Assert(new_number==n_vertices, ExcInternalError());
740 out <<
"object \"vertices\" class array type float rank 1 shape " << dim
741 <<
" items " << n_vertices <<
" data follows" 744 for (
unsigned int i=0; i<vertices.size(); ++i)
746 out <<
'\t' << vertices[i] <<
'\n';
761 out <<
"object \"cells\" class array type int rank 1 shape " 762 << n_vertices_per_cell
763 <<
" items " << n_cells <<
" data follows" <<
'\n';
767 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
771 out <<
"attribute \"element type\" string \"";
772 if (dim==1) out <<
"lines";
773 if (dim==2) out <<
"quads";
774 if (dim==3) out <<
"cubes";
776 <<
"attribute \"ref\" string \"positions\"" <<
'\n' <<
'\n';
780 out <<
"object \"material\" class array type int rank 0 items " 781 << n_cells <<
" data follows" <<
'\n';
783 out <<
' ' << (
unsigned int)cell->material_id();
785 <<
"attribute \"dep\" string \"connections\"" <<
'\n' <<
'\n';
787 out <<
"object \"level\" class array type int rank 0 items " 788 << n_cells <<
" data follows" <<
'\n';
790 out <<
' ' << cell->level();
792 <<
"attribute \"dep\" string \"connections\"" <<
'\n' <<
'\n';
796 out <<
"object \"measure\" class array type float rank 0 items " 797 << n_cells <<
" data follows" <<
'\n';
799 out <<
'\t' << cell->measure();
801 <<
"attribute \"dep\" string \"connections\"" <<
'\n' <<
'\n';
806 out <<
"object \"diameter\" class array type float rank 0 items " 807 << n_cells <<
" data follows" <<
'\n';
809 out <<
'\t' << cell->diameter();
811 <<
"attribute \"dep\" string \"connections\"" <<
'\n' <<
'\n';
817 out <<
"object \"faces\" class array type int rank 1 shape " 818 << n_vertices_per_face
819 <<
" items " << n_faces <<
" data follows" 824 for (
unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
828 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_face; ++v)
833 out <<
"attribute \"element type\" string \"";
834 if (dim==2) out <<
"lines";
835 if (dim==3) out <<
"quads";
837 <<
"attribute \"ref\" string \"positions\"" <<
'\n' <<
'\n';
842 out <<
"object \"boundary\" class array type int rank 0 items " 843 << n_faces <<
" data follows" <<
'\n';
848 for (
unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
849 out <<
' ' << (
int)(
signed char)cell->face(f)->boundary_id();
852 out <<
"attribute \"dep\" string \"connections\"" <<
'\n' <<
'\n';
856 out <<
"object \"face measure\" class array type float rank 0 items " 857 << n_faces <<
" data follows" <<
'\n';
860 for (
unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
861 out <<
' ' << cell->face(f)->measure();
864 out <<
"attribute \"dep\" string \"connections\"" <<
'\n' <<
'\n';
869 out <<
"object \"face diameter\" class array type float rank 0 items " 870 << n_faces <<
" data follows" <<
'\n';
873 for (
unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
874 out <<
' ' << cell->face(f)->diameter();
877 out <<
"attribute \"dep\" string \"connections\"" <<
'\n' <<
'\n';
894 out <<
"object \"deal data\" class field" <<
'\n' 895 <<
"component \"positions\" value \"vertices\"" <<
'\n' 896 <<
"component \"connections\" value \"cells\"" <<
'\n';
900 out <<
"object \"cell data\" class field" <<
'\n' 901 <<
"component \"positions\" value \"vertices\"" <<
'\n' 902 <<
"component \"connections\" value \"cells\"" <<
'\n';
903 out <<
"component \"material\" value \"material\"" <<
'\n';
904 out <<
"component \"level\" value \"level\"" <<
'\n';
906 out <<
"component \"measure\" value \"measure\"" <<
'\n';
908 out <<
"component \"diameter\" value \"diameter\"" <<
'\n';
913 out <<
"object \"face data\" class field" <<
'\n' 914 <<
"component \"positions\" value \"vertices\"" <<
'\n' 915 <<
"component \"connections\" value \"faces\"" <<
'\n';
916 out <<
"component \"boundary\" value \"boundary\"" <<
'\n';
918 out <<
"component \"measure\" value \"face measure\"" <<
'\n';
920 out <<
"component \"diameter\" value \"face diameter\"" <<
'\n';
924 <<
"object \"grid data\" class group" <<
'\n';
926 out <<
"member \"cells\" value \"cell data\"" <<
'\n';
928 out <<
"member \"faces\" value \"face data\"" <<
'\n';
929 out <<
"end" <<
'\n';
941 template <
int dim,
int spacedim>
943 std::ostream &out)
const 950 const std::vector<Point<spacedim> > &vertices = tria.
get_vertices();
974 out <<
"@f$NOD" <<
'\n' 975 << n_vertices <<
'\n';
980 for (
unsigned int i=0; i<vertices.size(); ++i)
986 for (
unsigned int d=spacedim+1; d<=3; ++d)
992 out <<
"@f$ENDNOD" <<
'\n' 1033 unsigned int elm_type;
1046 Assert(
false, ExcNotImplemented());
1053 out << cell->active_cell_index()+1 <<
' ' << elm_type <<
' ' 1054 <<
static_cast<unsigned int>(cell->material_id()) <<
' ' 1055 << cell->subdomain_id() <<
' ' 1060 for (
unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_cell;
1077 out <<
"@f$ENDELM\n";
1087 template <
int dim,
int spacedim>
1089 std::ostream &out)
const 1096 const std::vector<Point<spacedim> > &vertices = tria.
get_vertices();
1110 std::time_t time1= std::time (0);
1111 std::tm *time = std::localtime(&time1);
1112 out <<
"# This file was generated by the deal.II library." <<
'\n' 1114 << time->tm_year+1900 <<
"/" 1115 << time->tm_mon+1 <<
"/" 1116 << time->tm_mday <<
'\n' 1118 << time->tm_hour <<
":" 1119 << std::setw(2) << time->tm_min <<
":" 1120 << std::setw(2) << time->tm_sec <<
'\n' 1122 <<
"# For a description of the UCD format see the AVS Developer's guide." 1128 out << n_vertices <<
' ' 1139 for (
unsigned int i=0; i<vertices.size(); ++i)
1145 for (
unsigned int d=spacedim+1; d<=3; ++d)
1154 out << cell->active_cell_index()+1 <<
' ' 1155 <<
static_cast<unsigned int>(cell->material_id())
1169 Assert (
false, ExcNotImplemented());
1187 for (
unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_cell;
1213 template <
int dim,
int spacedim>
1219 Assert (
false, ExcNotImplemented());
1231 const int spacedim = 2;
1239 out <<
"#FIG 3.2\nLandscape\nCenter\nInches" << std::endl
1240 <<
"A4\n100.00\nSingle" << std::endl
1242 <<
"-3" << std::endl
1243 <<
"# generated by deal.II GridOut class" << std::endl
1244 <<
"# reduce first number to scale up image" << std::endl
1245 <<
"1200 2" << std::endl;
1248 unsigned int colno = 32;
1249 out <<
"0 " << colno++ <<
" #ff0000" << std::endl;
1250 out <<
"0 " << colno++ <<
" #ff8000" << std::endl;
1251 out <<
"0 " << colno++ <<
" #ffd000" << std::endl;
1252 out <<
"0 " << colno++ <<
" #ffff00" << std::endl;
1253 out <<
"0 " << colno++ <<
" #c0ff00" << std::endl;
1254 out <<
"0 " << colno++ <<
" #80ff00" << std::endl;
1255 out <<
"0 " << colno++ <<
" #00f000" << std::endl;
1256 out <<
"0 " << colno++ <<
" #00f0c0" << std::endl;
1257 out <<
"0 " << colno++ <<
" #00f0ff" << std::endl;
1258 out <<
"0 " << colno++ <<
" #00c0ff" << std::endl;
1259 out <<
"0 " << colno++ <<
" #0080ff" << std::endl;
1260 out <<
"0 " << colno++ <<
" #0040ff" << std::endl;
1261 out <<
"0 " << colno++ <<
" #0000c0" << std::endl;
1262 out <<
"0 " << colno++ <<
" #5000ff" << std::endl;
1263 out <<
"0 " << colno++ <<
" #8000ff" << std::endl;
1264 out <<
"0 " << colno++ <<
" #b000ff" << std::endl;
1265 out <<
"0 " << colno++ <<
" #ff00ff" << std::endl;
1266 out <<
"0 " << colno++ <<
" #ff80ff" << std::endl;
1268 for (
unsigned int i=0; i<8; ++i)
1269 out <<
"0 " << colno++ <<
" #" << std::hex << 32*i+31 << 32*i+31 << 32*i+31 << std::dec << std::endl;
1271 for (
unsigned int i=1; i<16; ++i)
1272 out <<
"0 " << colno++ <<
" #00" << std::hex << 16*i+15 << std::dec <<
"00" << std::endl;
1274 for (
unsigned int i=1; i<16; ++i)
1275 out <<
"0 " << colno++ <<
" #" << std::hex << 16*i+15 << 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 << std::dec <<
"0000" << std::endl;
1280 for (
unsigned int i=1; i<16; ++i)
1281 out <<
"0 " << colno++ <<
" #" << std::hex << 16*i+15 <<
"00" << 16*i+15 << std::dec << std::endl;
1283 for (
unsigned int i=1; i<16; ++i)
1284 out <<
"0 " << colno++ <<
" #0000" << std::hex << 16*i+15 << std::dec << std::endl;
1286 for (
unsigned int i=1; i<16; ++i)
1287 out <<
"0 " << colno++ <<
" #00" << std::hex << 16*i+15 << 16*i+15 << std::dec << std::endl;
1297 for (; cell != end; ++cell)
1313 out << static_cast<unsigned int>(cell->material_id()) + 32;
1316 out << cell->level() + 8;
1319 out << cell->subdomain_id() + 32;
1322 out << cell->level_subdomain_id() + 32;
1325 Assert(
false, ExcInternalError());
1331 ? (900-cell->level())
1332 : (900+cell->material_id()))
1338 << nv+1 << std::endl;
1344 for (
unsigned int k=0; k<=nv; ++k)
1348 for (
unsigned int d=0; d<static_cast<unsigned int>(dim); ++d)
1352 out <<
'\t' << ((d==0) ? val : -val);
1357 static const unsigned int face_reorder[4]= {2,1,3,0};
1359 for (
unsigned int f=0; f<nf; ++f)
1362 face = cell->face(face_reorder[f]);
1371 << (1 + (
unsigned int) bi);
1376 ? (800-cell->level())
1383 << nvf << std::endl;
1390 for (
unsigned int k=0; k<nvf; ++k)
1393 for (
unsigned int d=0; d<static_cast<unsigned int>(dim); ++d)
1397 out <<
'\t' << ((d==0) ? val : -val);
1414 template <
int dim,
int spacedim>
1416 std::ostream &)
const 1418 Assert(
false, ExcNotImplemented());
1425 unsigned int n_materials = 0;
1426 unsigned int n_levels = 0;
1427 unsigned int n_subdomains = 0;
1428 unsigned int n_level_subdomains = 0;
1432 unsigned int min_level, max_level;
1440 Assert (height != 0 || width != 0,
ExcMessage(
"You have to set at least one of width and height"));
1442 unsigned int margin_in_percent = 0;
1444 margin_in_percent = 8;
1447 unsigned int cell_label_font_size;
1464 float x_max_perspective, x_min_perspective;
1465 float y_max_perspective, y_min_perspective;
1467 float x_dimension_perspective, y_dimension_perspective;
1471 double x_min = tria.
begin()->vertex(0)[0];
1472 double x_max = tria.
begin()->vertex(0)[0];
1473 double y_min = tria.
begin()->vertex(0)[1];
1474 double y_max = tria.
begin()->vertex(0)[1];
1476 double x_dimension, y_dimension;
1478 min_level = max_level = tria.
begin()->level();
1481 unsigned int materials[256];
1482 for (
unsigned int material_index = 0; material_index < 256; material_index++)
1483 materials[material_index] = 0;
1486 unsigned int levels[256];
1487 for (
unsigned int level_index = 0; level_index < 256; level_index++)
1488 levels[level_index] = 0;
1491 unsigned int subdomains[256];
1492 for (
unsigned int subdomain_index = 0; subdomain_index < 256; subdomain_index++)
1493 subdomains[subdomain_index] = 0;
1496 int level_subdomains[256];
1497 for (
int level_subdomain_index = 0; level_subdomain_index < 256; level_subdomain_index++)
1498 level_subdomains[level_subdomain_index] = 0;
1506 for (
unsigned int vertex_index = 0; vertex_index < 4; vertex_index++)
1508 if (cell->vertex(vertex_index)[0] < x_min) x_min = cell->vertex(vertex_index)[0];
1509 if (cell->vertex(vertex_index)[0] > x_max) x_max = cell->vertex(vertex_index)[0];
1511 if (cell->vertex(vertex_index)[1] < y_min) y_min = cell->vertex(vertex_index)[1];
1512 if (cell->vertex(vertex_index)[1] > y_max) y_max = cell->vertex(vertex_index)[1];
1515 if ((
unsigned int)cell->level() < min_level) min_level = cell->level();
1516 if ((
unsigned int)cell->level() > max_level) max_level = cell->level();
1518 materials[(
unsigned int)cell->material_id()] = 1;
1519 levels[(
unsigned int)cell->level()] = 1;
1521 subdomains[cell->subdomain_id()+2] = 1;
1522 level_subdomains[cell->level_subdomain_id()+2] = 1;
1525 x_dimension = x_max - x_min;
1526 y_dimension = y_max - y_min;
1529 for (
unsigned int material_index = 0; material_index < 256; material_index++)
1531 if (materials[material_index]) n_materials++;
1535 for (
unsigned int level_index = 0; level_index < 256; level_index++)
1537 if (levels[level_index]) n_levels++;
1541 for (
unsigned int subdomain_index = 0; subdomain_index < 256; subdomain_index++)
1543 if (subdomains[subdomain_index]) n_subdomains++;
1547 for (
int level_subdomain_index = 0; level_subdomain_index < 256; level_subdomain_index++)
1549 if (level_subdomains[level_subdomain_index]) n_level_subdomains++;
1564 n = n_level_subdomains;
1571 camera_position[0] = 0;
1572 camera_position[1] = 0;
1573 camera_position[2] = 2. * std::max(x_dimension, y_dimension);
1575 camera_direction[0] = 0;
1576 camera_direction[1] = 0;
1577 camera_direction[2] = -1;
1579 camera_horizontal[0] = 1;
1580 camera_horizontal[1] = 0;
1581 camera_horizontal[2] = 0;
1583 camera_focus = .5 * std::max(x_dimension, y_dimension);
1589 const double angle_factor = 3.14159265 / 180.;
1601 camera_position[1] = camera_position_temp[1];
1602 camera_position[2] = camera_position_temp[2];
1604 camera_direction[1] = camera_direction_temp[1];
1605 camera_direction[2] = camera_direction_temp[2];
1607 camera_horizontal[1] = camera_horizontal_temp[1];
1608 camera_horizontal[2] = camera_horizontal_temp[2];
1620 camera_position[0] = camera_position_temp[0];
1621 camera_position[1] = camera_position_temp[1];
1623 camera_direction[0] = camera_direction_temp[0];
1624 camera_direction[1] = camera_direction_temp[1];
1626 camera_horizontal[0] = camera_horizontal_temp[0];
1627 camera_horizontal[1] = camera_horizontal_temp[1];
1630 camera_position[0] = x_min + .5 * x_dimension;
1631 camera_position[1] = y_min + .5 * y_dimension;
1638 point[0] = tria.
begin()->vertex(0)[0];
1639 point[1] = tria.
begin()->vertex(0)[1];
1642 float min_level_min_vertex_distance = 0;
1649 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
1651 x_max_perspective = projection_decomposition[0];
1652 x_min_perspective = projection_decomposition[0];
1654 y_max_perspective = projection_decomposition[1];
1655 y_min_perspective = projection_decomposition[1];
1659 point[0] = cell->vertex(0)[0];
1660 point[1] = cell->vertex(0)[1];
1668 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
1670 if (x_max_perspective < projection_decomposition[0]) x_max_perspective = projection_decomposition[0];
1671 if (x_min_perspective > projection_decomposition[0]) x_min_perspective = projection_decomposition[0];
1673 if (y_max_perspective < projection_decomposition[1]) y_max_perspective = projection_decomposition[1];
1674 if (y_min_perspective > projection_decomposition[1]) y_min_perspective = projection_decomposition[1];
1676 point[0] = cell->vertex(1)[0];
1677 point[1] = cell->vertex(1)[1];
1679 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
1681 if (x_max_perspective < projection_decomposition[0]) x_max_perspective = projection_decomposition[0];
1682 if (x_min_perspective > projection_decomposition[0]) x_min_perspective = projection_decomposition[0];
1684 if (y_max_perspective < projection_decomposition[1]) y_max_perspective = projection_decomposition[1];
1685 if (y_min_perspective > projection_decomposition[1]) y_min_perspective = projection_decomposition[1];
1687 point[0] = cell->vertex(2)[0];
1688 point[1] = cell->vertex(2)[1];
1690 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
1692 if (x_max_perspective < projection_decomposition[0]) x_max_perspective = projection_decomposition[0];
1693 if (x_min_perspective > projection_decomposition[0]) x_min_perspective = projection_decomposition[0];
1695 if (y_max_perspective < projection_decomposition[1]) y_max_perspective = projection_decomposition[1];
1696 if (y_min_perspective > projection_decomposition[1]) y_min_perspective = projection_decomposition[1];
1698 point[0] = cell->vertex(3)[0];
1699 point[1] = cell->vertex(3)[1];
1701 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
1703 if (x_max_perspective < projection_decomposition[0]) x_max_perspective = projection_decomposition[0];
1704 if (x_min_perspective > projection_decomposition[0]) x_min_perspective = projection_decomposition[0];
1706 if (y_max_perspective < projection_decomposition[1]) y_max_perspective = projection_decomposition[1];
1707 if (y_min_perspective > projection_decomposition[1]) y_min_perspective = projection_decomposition[1];
1709 if ((
unsigned int)cell->level() == min_level) min_level_min_vertex_distance = cell->minimum_vertex_distance();
1712 x_dimension_perspective = x_max_perspective - x_min_perspective;
1713 y_dimension_perspective = y_max_perspective - y_min_perspective;
1717 width =
static_cast<unsigned int>(.5 + height * (x_dimension_perspective / y_dimension_perspective));
1718 else if (height == 0)
1719 height =
static_cast<unsigned int>(.5 + width * (y_dimension_perspective / x_dimension_perspective));
1720 unsigned int additional_width = 0;
1722 unsigned int font_size =
static_cast<unsigned int>(.5 + (height/100.) * 1.75);
1723 cell_label_font_size =
static_cast<unsigned int>(.5 +
1726 * min_level_min_vertex_distance
1727 / std::min(x_dimension, y_dimension)));
1731 additional_width =
static_cast<unsigned int>(.5 + height * .4);
1735 additional_width =
static_cast<unsigned int>(.5 + height * .175);
1746 out <<
"<svg width=\"" << width + additional_width <<
"\" height=\"" << height <<
"\" xmlns=\"http://www.w3.org/2000/svg\" version=\"1.1\">" 1752 out <<
" <linearGradient id=\"background_gradient\" gradientUnits=\"userSpaceOnUse\" x1=\"0\" y1=\"0\" x2=\"0\" y2=\"" << height <<
"\">" <<
'\n' 1753 <<
" <stop offset=\"0\" style=\"stop-color:white\"/>" <<
'\n' 1754 <<
" <stop offset=\"1\" style=\"stop-color:lightsteelblue\"/>" <<
'\n' 1755 <<
" </linearGradient>" <<
'\n';
1761 out <<
"<!-- internal style sheet -->" <<
'\n' 1762 <<
"<style type=\"text/css\"><![CDATA[" <<
'\n';
1767 else out <<
" rect.background{fill:none}" <<
'\n';
1771 <<
" text{font-family:Helvetica; text-anchor:middle; fill:rgb(25,25,25)}" <<
'\n' 1779 unsigned int labeling_index = 0;
1781 for (
unsigned int index = 0; index < n; index++)
1785 if (n != 1) h = .6 - (index / (n-1.)) * .6;
1792 unsigned int i =
static_cast<unsigned int>(h * 6);
1794 double f = h * 6 - i;
1801 r = 255, g =
static_cast<unsigned int>(.5 + 255*t);
1804 r =
static_cast<unsigned int>(.5 + 255*q), g = 255;
1807 g = 255, b =
static_cast<unsigned int>(.5 + 255*t);
1810 g =
static_cast<unsigned int>(.5 + 255*q), b = 255;
1813 r =
static_cast<unsigned int>(.5 + 255*t), b = 255;
1816 r = 255, b =
static_cast<unsigned int>(.5 + 255*q);
1825 while (!materials[labeling_index]) labeling_index++;
1828 while (!levels[labeling_index]) labeling_index++;
1831 while (!subdomains[labeling_index]) labeling_index++;
1834 while (!level_subdomains[labeling_index]) labeling_index++;
1840 out <<
" path.p" << labeling_index
1841 <<
"{fill:rgb(" << r <<
',' << g <<
',' << b <<
"); " 1844 out <<
" path.ps" << labeling_index
1845 <<
"{fill:rgb(" <<
static_cast<unsigned int>(.5 + .75 * r) <<
',' << static_cast<unsigned int>(.5 + .75 * g) <<
',' <<
static_cast<unsigned int>(.5 + .75 * b) <<
"); " 1848 out <<
" rect.r" << labeling_index
1849 <<
"{fill:rgb(" << r <<
',' << g <<
',' << b <<
"); " 1856 out <<
"]]></style>" <<
'\n' <<
'\n';
1859 out <<
" <rect class=\"background\" width=\"" << width <<
"\" height=\"" << height <<
"\"/>" <<
'\n';
1863 unsigned int x_offset = 0;
1865 if (
svg_flags.
margin) x_offset =
static_cast<unsigned int>(.5 + (height/100.) * (margin_in_percent/2.));
1866 else x_offset =
static_cast<unsigned int>(.5 + height * .025);
1868 out <<
" <text x=\"" << x_offset <<
"\" y=\"" <<
static_cast<unsigned int>(.5 + height * .0525) <<
'\"' 1869 <<
" style=\"font-weight:100; fill:lightsteelblue; text-anchor:start; font-family:Courier; font-size:" <<
static_cast<unsigned int>(.5 + height * .045) <<
"px\">" 1870 <<
"deal.II" <<
"</text>" <<
'\n';
1884 out <<
" <!-- cells -->" <<
'\n';
1886 for (
unsigned int level_index = min_level; level_index <= max_level; level_index++)
1889 cell = tria.
begin(level_index),
1890 endc = tria.
end(level_index);
1892 for (; cell != endc; ++cell)
1901 out <<
" class=\"p";
1908 out << (
unsigned int)cell->material_id();
1911 out << (
unsigned int)cell->level();
1915 out << cell->subdomain_id() + 2;
1920 out << cell->level_subdomain_id() + 2;
1931 point[0] = cell->vertex(0)[0];
1932 point[1] = cell->vertex(0)[1];
1940 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
1942 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)) <<
' ' 1943 <<
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));
1947 point[0] = cell->vertex(1)[0];
1948 point[1] = cell->vertex(1)[1];
1950 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
1952 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)) <<
' ' 1953 <<
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));
1957 point[0] = cell->vertex(3)[0];
1958 point[1] = cell->vertex(3)[1];
1960 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
1962 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)) <<
' ' 1963 <<
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));
1967 point[0] = cell->vertex(2)[0];
1968 point[1] = cell->vertex(2)[1];
1970 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
1972 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)) <<
' ' 1973 <<
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));
1977 point[0] = cell->vertex(0)[0];
1978 point[1] = cell->vertex(0)[1];
1980 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
1982 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)) <<
' ' 1983 <<
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));
1985 out <<
"\"/>" <<
'\n';
1990 point[0] = cell->center()[0];
1991 point[1] = cell->center()[1];
1999 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.));
2000 float distance_factor = distance_to_camera / (2. * std::max(x_dimension, y_dimension));
2002 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
2004 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));
2007 <<
" 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))
2008 <<
"\" 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)
2009 <<
"\" style=\"font-size:" << font_size_this_cell
2014 out << cell->level();
2020 out << cell->index();
2026 out << (int)cell->material_id();
2036 out << static_cast<int>(cell->subdomain_id());
2048 out << static_cast<int>(cell->level_subdomain_id());
2051 out <<
"</text>" <<
'\n';
2057 for (
unsigned int faceIndex = 0; faceIndex < 4; faceIndex++)
2059 if (cell->at_boundary(faceIndex))
2062 point[0] = cell->face(faceIndex)->vertex(0)[0];
2063 point[1] = cell->face(faceIndex)->vertex(0)[1];
2071 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
2073 out <<
" <line x1=\"" 2074 <<
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))
2076 <<
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));
2078 point[0] = cell->face(faceIndex)->vertex(1)[0];
2079 point[1] = cell->face(faceIndex)->vertex(1)[1];
2087 projection_decomposition =
GridOut::svg_project_point(point, camera_position, camera_direction, camera_horizontal, camera_focus);
2090 <<
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))
2092 <<
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))
2104 unsigned int line_offset = 0;
2106 additional_width = 0;
2107 if (!
svg_flags.
margin) additional_width =
static_cast<unsigned int>(.5 + (height/100.) * 2.5);
2112 out <<
" <rect x=\"" << width + additional_width <<
"\" y=\"" <<
static_cast<unsigned int>(.5 + (height/100.) * margin_in_percent)
2113 <<
"\" width=\"" << static_cast<unsigned int>(.5 + (height/100.) * (40. - margin_in_percent)) <<
"\" height=\"" <<
static_cast<unsigned int>(.5 + height * .165) <<
"\"/>" <<
'\n';
2115 out <<
" <text x=\"" << width + additional_width +
static_cast<unsigned int>(.5 + (height/100.) * 1.25)
2116 <<
"\" y=\"" <<
static_cast<unsigned int>(.5 + (height/100.) * margin_in_percent + (++line_offset) * 1.5 * font_size)
2117 <<
"\" style=\"text-anchor:start; font-weight:bold; font-size:" << font_size
2118 <<
"px\">" <<
"cell label" 2119 <<
"</text>" <<
'\n';
2123 out <<
" <text x=\"" << width + additional_width +
static_cast<unsigned int>(.5 + (height/100.) * 2.)
2124 <<
"\" y=\"" <<
static_cast<unsigned int>(.5 + (height/100.) * margin_in_percent + (++line_offset) * 1.5 * font_size)
2125 <<
"\" style=\"text-anchor:start; font-style:oblique; font-size:" << font_size
2126 <<
"px\">" <<
"level_number";
2131 out <<
"</text>" <<
'\n';
2136 out <<
" <text x=\"" << width + additional_width +
static_cast<unsigned int>(.5 + (height/100.) * 2.)
2137 <<
"\" y=\"" <<
static_cast<unsigned int>(.5 + (height/100.) * margin_in_percent + (++line_offset) * 1.5 * font_size )
2138 <<
"\" style=\"text-anchor:start; font-style:oblique; font-size:" << font_size
2145 out <<
"</text>" <<
'\n';
2150 out <<
" <text x=\"" << width + additional_width +
static_cast<unsigned int>(.5 + (height/100.) * 2.)
2151 <<
"\" y=\"" <<
static_cast<unsigned int>(.5 + (height/100.) * margin_in_percent + (++line_offset) * 1.5 * font_size )
2152 <<
"\" style=\"text-anchor:start; font-style:oblique; font-size:" << font_size
2159 out <<
"</text>" <<
'\n';
2164 out <<
" <text x= \"" << width + additional_width +
static_cast<unsigned int>(.5 + (height/100.) * 2.)
2165 <<
"\" y=\"" <<
static_cast<unsigned int>(.5 + (height/100.) * margin_in_percent + (++line_offset) * 1.5 * font_size )
2166 <<
"\" style=\"text-anchor:start; font-style:oblique; font-size:" << font_size
2173 out <<
"</text>" <<
'\n';
2178 out <<
" <text x= \"" << width + additional_width +
static_cast<unsigned int>(.5 + (height/100.) * 2.)
2179 <<
"\" y=\"" <<
static_cast<unsigned int>(.5 + (height/100.) * margin_in_percent + (++line_offset) * 1.5 * font_size )
2180 <<
"\" style=\"text-anchor:start; font-style:oblique; font-size:" << font_size
2182 <<
"level_subdomain_id" 2183 <<
"</text>" <<
'\n';
2190 out <<
" <text x=\"" << width + additional_width
2191 <<
"\" y=\"" <<
static_cast<unsigned int>(.5 + (height/100.) * margin_in_percent + 10.75 * font_size)
2192 <<
"\" style=\"text-anchor:start; font-size:" << font_size <<
"px\">" 2200 out <<
'\n' <<
" <!-- colorbar -->" <<
'\n';
2202 out <<
" <text x=\"" << width + additional_width
2203 <<
"\" y=\"" <<
static_cast<unsigned int>(.5 + (height/100.) * (margin_in_percent + 29.) - (font_size/1.25))
2204 <<
"\" style=\"text-anchor:start; font-weight:bold; font-size:" << font_size <<
"px\">";
2209 out <<
"material_id";
2212 out <<
"level_number";
2215 out <<
"subdomain_id";
2218 out <<
"level_subdomain_id";
2224 out <<
"</text>" <<
'\n';
2226 unsigned int element_height =
static_cast<unsigned int>(((height/100.) * (71. - 2.*margin_in_percent)) / n);
2227 unsigned int element_width =
static_cast<unsigned int>(.5 + (height/100.) * 2.5);
2229 int labeling_index = 0;
2231 for (
unsigned int index = 0; index < n; index++)
2236 while (!materials[labeling_index]) labeling_index++;
2239 while (!levels[labeling_index]) labeling_index++;
2242 while (!subdomains[labeling_index]) labeling_index++;
2245 while (!level_subdomains[labeling_index]) labeling_index++;
2251 out <<
" <rect class=\"r" << labeling_index
2252 <<
"\" x=\"" << width + additional_width
2253 <<
"\" y=\"" <<
static_cast<unsigned int>(.5 + (height/100.) * (margin_in_percent + 29)) + (n-index-1) * element_height
2254 <<
"\" width=\"" << element_width
2255 <<
"\" height=\"" << element_height
2258 out <<
" <text x=\"" << width + additional_width + 1.5 * element_width
2259 <<
"\" 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) <<
"\"" 2260 <<
" style=\"text-anchor:start; font-size:" <<
static_cast<unsigned int>(.5 + font_size) <<
"px";
2262 if (index == 0 || index == n-1) out <<
"; font-weight:bold";
2264 out <<
"\">" << labeling_index;
2266 if (index == n-1) out <<
" max";
2267 if (index == 0) out <<
" min";
2269 out <<
"</text>" <<
'\n';
2277 out <<
'\n' <<
"</svg>";
2285 std::ostream &)
const 2288 Assert (
false, ExcNotImplemented());
2292 template <
int dim,
int spacedim>
2294 std::ostream &out)
const 2302 const std::time_t time1 = std::time (0);
2303 const std::tm *time = std::localtime (&time1);
2306 <<
"\n# This file was generated by the deal.II library." 2308 << time->tm_year+1900 <<
"/" 2309 << std::setfill(
'0') << std::setw (2) << time->tm_mon+1 <<
"/" 2310 << std::setfill(
'0') << std::setw (2) << time->tm_mday
2312 << std::setfill(
'0') << std::setw (2) << time->tm_hour <<
":" 2313 << std::setfill(
'0') << std::setw (2) << time->tm_min <<
":" 2314 << std::setfill(
'0') << std::setw (2) << time->tm_sec
2316 <<
"\n# For a description of the MathGL script format see the MathGL manual. " 2318 <<
"\n# Note: This file is understood by MathGL v2.1 and higher only, and can " 2319 <<
"\n# be quickly viewed in a graphical environment using \'mglview\'. " 2326 const std::string axes =
"xyz";
2341 out <<
"\nsetsize 800 800";
2342 out <<
"\nrotate 0 0";
2345 out <<
"\nsetsize 800 800";
2346 out <<
"\nrotate 60 40";
2349 Assert (
false, ExcNotImplemented ());
2355 <<
"\n# Vertex ordering." 2356 <<
"\n# list <vertex order> <vertex indices>" 2365 out <<
"\nlist f 0 1 2 3" 2369 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" 2373 Assert (
false, ExcNotImplemented ());
2378 <<
"\n# List of vertices." 2379 <<
"\n# list <id> <vertices>" 2384 typename ::Triangulation<dim, spacedim>::active_cell_iterator
2389 for (; cell!=endc; ++cell)
2391 for (
unsigned int i=0; i<dim; ++i)
2398 out <<
"\nlist " << axes[i] << cell->active_cell_index() <<
" ";
2399 for (
unsigned int j=0; j<GeometryInfo<dim>::vertices_per_cell; ++j)
2400 out << cell->vertex(j)[i] <<
" ";
2407 <<
"\n# List of cells to quadplot." 2408 <<
"\n# quadplot <vertex order> <id> <style>" 2412 out <<
"\nquadplot f ";
2413 for (
unsigned int j=0; j<dim; ++j)
2414 out << axes[j] << i <<
" ";
2439 template <
int dim,
int spacedim,
typename ITERATOR,
typename END>
2442 ITERATOR cell, END end)
2445 for (; cell != end; ++cell)
2451 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
2453 patch.
vertices[v] = cell->vertex(v);
2454 patch.
data(0,v) = cell->level();
2455 patch.
data(1,v) =
static_cast<int>(cell->manifold_id());
2456 patch.
data(2,v) = cell->material_id();
2457 if (!cell->has_children())
2458 patch.
data(3,v) =
static_cast<int>(cell->subdomain_id());
2460 patch.
data(3,v) = -1;
2461 patch.
data(4,v) =
static_cast<int>(cell->level_subdomain_id());
2463 patches.push_back (patch);
2467 std::vector<std::string> triangulation_patch_data_names ()
2469 std::vector<std::string> v(5);
2474 v[4] =
"level_subdomain";
2481 template <
int dim,
int spacedim>
2483 std::ostream &out)
const 2491 std::vector<DataOutBase::Patch<dim,spacedim> > patches;
2493 generate_triangulation_patches(patches, tria.
begin_active(), tria.
end());
2495 triangulation_patch_data_names(),
2496 std::vector<std_cxx11::tuple<unsigned int, unsigned int, std::string> >(),
2505 template <
int dim,
int spacedim>
2507 std::ostream &out)
const 2515 std::vector<DataOutBase::Patch<dim,spacedim> > patches;
2517 generate_triangulation_patches(patches, tria.
begin_active(), tria.
end());
2519 triangulation_patch_data_names(),
2520 std::vector<std_cxx11::tuple<unsigned int, unsigned int, std::string> >(),
2572 template <
int dim,
int spacedim>
2576 unsigned int n_faces = 0;
2579 face != endf; ++face)
2580 if ((face->at_boundary()) &&
2581 (face->boundary_id() != 0))
2589 template <
int dim,
int spacedim>
2595 std::vector<bool> line_flags;
2597 .save_user_flags_line (line_flags);
2599 .clear_user_flags_line ();
2601 unsigned int n_lines = 0;
2606 cell != endc; ++cell)
2607 for (
unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++l)
2608 if (cell->line(l)->at_boundary()
2610 (cell->line(l)->boundary_id() != 0)
2612 (cell->line(l)->user_flag_set() ==
false))
2615 cell->line(l)->set_user_flag();
2621 .load_user_flags_line (line_flags);
2631 const unsigned int next_element_index,
2632 std::ostream &)
const 2634 return next_element_index;
2640 const unsigned int next_element_index,
2641 std::ostream &)
const 2643 return next_element_index;
2648 const unsigned int next_element_index,
2649 std::ostream &)
const 2651 return next_element_index;
2657 const unsigned int next_element_index,
2658 std::ostream &)
const 2660 return next_element_index;
2665 const unsigned int next_element_index,
2666 std::ostream &)
const 2668 return next_element_index;
2674 const unsigned int next_element_index,
2675 std::ostream &)
const 2677 return next_element_index;
2683 const unsigned int next_element_index,
2684 std::ostream &)
const 2686 return next_element_index;
2691 const unsigned int next_element_index,
2692 std::ostream &)
const 2694 return next_element_index;
2700 template <
int dim,
int spacedim>
2703 const unsigned int next_element_index,
2704 std::ostream &out)
const 2706 unsigned int current_element_index = next_element_index;
2710 face != endf; ++face)
2711 if (face->at_boundary() &&
2712 (face->boundary_id() != 0))
2714 out << current_element_index <<
' ';
2724 Assert (
false, ExcNotImplemented());
2726 out << static_cast<unsigned int>(face->boundary_id())
2728 << static_cast<unsigned int>(face->boundary_id())
2731 for (
unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_face; ++vertex)
2736 ++current_element_index;
2738 return current_element_index;
2742 template <
int dim,
int spacedim>
2745 const unsigned int next_element_index,
2746 std::ostream &out)
const 2748 unsigned int current_element_index = next_element_index;
2753 std::vector<bool> line_flags;
2755 .save_user_flags_line (line_flags);
2757 .clear_user_flags_line ();
2762 cell != endc; ++cell)
2763 for (
unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++l)
2764 if (cell->line(l)->at_boundary()
2766 (cell->line(l)->boundary_id() != 0)
2768 (cell->line(l)->user_flag_set() ==
false))
2770 out << next_element_index <<
" 1 ";
2771 out << static_cast<unsigned int>(cell->line(l)->boundary_id())
2773 << static_cast<unsigned int>(cell->line(l)->boundary_id())
2776 for (
unsigned int vertex=0; vertex<2; ++vertex)
2784 ++current_element_index;
2785 cell->line(l)->set_user_flag();
2791 .load_user_flags_line (line_flags);
2793 return current_element_index;
2801 const unsigned int next_element_index,
2802 std::ostream &)
const 2804 return next_element_index;
2809 const unsigned int next_element_index,
2810 std::ostream &)
const 2812 return next_element_index;
2817 const unsigned int next_element_index,
2818 std::ostream &)
const 2820 return next_element_index;
2825 const unsigned int next_element_index,
2826 std::ostream &)
const 2828 return next_element_index;
2833 const unsigned int next_element_index,
2834 std::ostream &)
const 2836 return next_element_index;
2842 const unsigned int next_element_index,
2843 std::ostream &)
const 2845 return next_element_index;
2851 const unsigned int next_element_index,
2852 std::ostream &)
const 2854 return next_element_index;
2859 const unsigned int next_element_index,
2860 std::ostream &)
const 2862 return next_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;
2877 face != endf; ++face)
2878 if (face->at_boundary() &&
2879 (face->boundary_id() != 0))
2881 out << current_element_index <<
" " 2882 <<
static_cast<unsigned int>(face->boundary_id())
2893 Assert (
false, ExcNotImplemented());
2896 for (
unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_face; ++vertex)
2900 ++current_element_index;
2902 return current_element_index;
2907 template <
int dim,
int spacedim>
2910 const unsigned int next_element_index,
2911 std::ostream &out)
const 2913 unsigned int current_element_index = next_element_index;
2918 std::vector<bool> line_flags;
2920 .save_user_flags_line (line_flags);
2922 .clear_user_flags_line ();
2927 cell != endc; ++cell)
2928 for (
unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++l)
2929 if (cell->line(l)->at_boundary()
2931 (cell->line(l)->boundary_id() != 0)
2933 (cell->line(l)->user_flag_set() ==
false))
2935 out << current_element_index <<
" " 2936 <<
static_cast<unsigned int>(cell->line(l)->boundary_id())
2939 for (
unsigned int vertex=0; vertex<2; ++vertex)
2947 ++current_element_index;
2948 cell->line(l)->set_user_flag();
2954 .load_user_flags_line (line_flags);
2955 return current_element_index;
2963 camera_vertical[0] = camera_horizontal[1] * camera_direction[2] - camera_horizontal[2] * camera_direction[1];
2964 camera_vertical[1] = camera_horizontal[2] * camera_direction[0] - camera_horizontal[0] * camera_direction[2];
2965 camera_vertical[2] = camera_horizontal[0] * camera_direction[1] - camera_horizontal[1] * camera_direction[0];
2969 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];
2972 projection[0] = camera_position[0] + phi * (point[0] - camera_position[0]);
2973 projection[1] = camera_position[1] + phi * (point[1] - camera_position[1]);
2974 projection[2] = camera_position[2] + phi * (point[2] - camera_position[2]);
2977 projection_decomposition[0] = (projection[0] - camera_position[0] - camera_focus * camera_direction[0]) * camera_horizontal[0];
2978 projection_decomposition[0] += (projection[1] - camera_position[1] - camera_focus * camera_direction[1]) * camera_horizontal[1];
2979 projection_decomposition[0] += (projection[2] - camera_position[2] - camera_focus * camera_direction[2]) * camera_horizontal[2];
2981 projection_decomposition[1] = (projection[0] - camera_position[0] - camera_focus * camera_direction[0]) * camera_vertical[0];
2982 projection_decomposition[1] += (projection[1] - camera_position[1] - camera_focus * camera_direction[1]) * camera_vertical[1];
2983 projection_decomposition[1] += (projection[2] - camera_position[2] - camera_focus * camera_direction[2]) * camera_vertical[2];
2985 return projection_decomposition;
2994 template <
int spacedim>
2995 void write_gnuplot (const ::Triangulation<1,spacedim> &tria,
3004 typename ::Triangulation<dim,spacedim>::active_cell_iterator
3005 cell=tria.begin_active();
3006 const typename ::Triangulation<dim,spacedim>::active_cell_iterator
3008 for (; cell!=endc; ++cell)
3011 out <<
"# cell " << cell <<
'\n';
3013 out << cell->vertex(0)
3014 <<
' ' << cell->level()
3015 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3017 <<
' ' << cell->level()
3018 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3031 template <
int spacedim>
3032 void write_gnuplot (const ::Triangulation<2,spacedim> &tria,
3041 const unsigned int n_additional_points=
3043 const unsigned int n_points=2+n_additional_points;
3045 typename ::Triangulation<dim,spacedim>::active_cell_iterator
3046 cell=tria.begin_active();
3047 const typename ::Triangulation<dim,spacedim>::active_cell_iterator
3056 std::vector<
Point<dim-1> > boundary_points;
3059 boundary_points.resize(n_points);
3060 boundary_points[0][0]=0;
3061 boundary_points[n_points-1][0]=1;
3062 for (
unsigned int i=1; i<n_points-1; ++i)
3063 boundary_points[i](0)= 1.*i/(n_points-1);
3065 std::vector<double> dummy_weights(n_points, 1./n_points);
3066 Quadrature<dim-1> quadrature(boundary_points, dummy_weights);
3071 for (; cell!=endc; ++cell)
3074 out <<
"# cell " << cell <<
'\n';
3086 for (
unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
3088 <<
' ' << cell->level()
3089 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n';
3090 out << cell->vertex(0)
3091 <<
' ' << cell->level()
3092 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3103 for (
unsigned int face_no=0;
3104 face_no<GeometryInfo<dim>::faces_per_cell; ++face_no)
3106 const typename ::Triangulation<dim,spacedim>::face_iterator
3107 face = cell->face(face_no);
3115 const unsigned int offset=face_no*n_points;
3116 for (
unsigned int i=0; i<n_points; ++i)
3118 (cell, q_projector->
point(offset+i)))
3119 <<
' ' << cell->level()
3120 <<
' ' <<
static_cast<unsigned int>(cell->material_id())
3133 out << face->vertex(0)
3134 <<
' ' << cell->level()
3135 <<
' ' <<
static_cast<unsigned int>(cell->material_id())
3138 <<
' ' << cell->level()
3139 <<
' ' <<
static_cast<unsigned int>(cell->material_id())
3148 if (q_projector != 0)
3160 template <
int spacedim>
3161 void write_gnuplot (const ::Triangulation<3,spacedim> &tria,
3170 const unsigned int n_additional_points=
3172 const unsigned int n_points=2+n_additional_points;
3174 typename ::Triangulation<dim,spacedim>::active_cell_iterator
3175 cell=tria.begin_active();
3176 const typename ::Triangulation<dim,spacedim>::active_cell_iterator
3185 std::vector<Point<1> > boundary_points;
3188 boundary_points.resize(n_points);
3189 boundary_points[0][0]=0;
3190 boundary_points[n_points-1][0]=1;
3191 for (
unsigned int i=1; i<n_points-1; ++i)
3192 boundary_points[i](0)= 1.*i/(n_points-1);
3194 std::vector<double> dummy_weights(n_points, 1./n_points);
3199 QIterated<dim-1> quadrature(quadrature1d, 1);
3203 for (; cell!=endc; ++cell)
3206 out <<
"# cell " << cell <<
'\n';
3208 if (mapping==0 || n_points==2 ||
3212 out << cell->vertex(0)
3213 <<
' ' << cell->level()
3214 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3216 <<
' ' << cell->level()
3217 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3219 <<
' ' << cell->level()
3220 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3222 <<
' ' << cell->level()
3223 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3225 <<
' ' << cell->level()
3226 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3229 out << cell->vertex(2)
3230 <<
' ' << cell->level()
3231 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3233 <<
' ' << cell->level()
3234 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3236 <<
' ' << cell->level()
3237 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3239 <<
' ' << cell->level()
3240 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3242 <<
' ' << cell->level()
3243 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3247 out << cell->vertex(0)
3248 <<
' ' << cell->level()
3249 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3251 <<
' ' << cell->level()
3252 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3254 out << cell->vertex(1)
3255 <<
' ' << cell->level()
3256 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3258 <<
' ' << cell->level()
3259 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3261 out << cell->vertex(5)
3262 <<
' ' << cell->level()
3263 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3265 <<
' ' << cell->level()
3266 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3268 out << cell->vertex(4)
3269 <<
' ' << cell->level()
3270 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3272 <<
' ' << cell->level()
3273 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3278 for (
unsigned int face_no=0; face_no<GeometryInfo<dim>::faces_per_cell; ++face_no)
3280 const typename ::Triangulation<dim,spacedim>::face_iterator
3281 face = cell->face(face_no);
3283 if (face->at_boundary())
3285 const unsigned int offset=face_no*n_points*n_points;
3286 for (
unsigned int i=0; i<n_points-1; ++i)
3287 for (
unsigned int j=0; j<n_points-1; ++j)
3290 cell, q_projector->
point(offset+i*n_points+j));
3292 <<
' ' << cell->level()
3293 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n';
3295 cell, q_projector->
point(offset+(i+1)*n_points+j)))
3296 <<
' ' << cell->level()
3297 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n';
3299 cell, q_projector->
point(offset+(i+1)*n_points+j+1)))
3300 <<
' ' << cell->level()
3301 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n';
3303 cell, q_projector->
point(offset+i*n_points+j+1)))
3304 <<
' ' << cell->level()
3305 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n';
3311 <<
' ' << cell->level()
3312 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n';
3313 out <<
'\n' <<
'\n';
3318 for (
unsigned int l=0; l<GeometryInfo<dim>::lines_per_face; ++l)
3320 const typename ::Triangulation<dim,spacedim>::line_iterator
3324 &v1=line->vertex(1);
3338 for (
unsigned int i=0; i<n_points; ++i)
3340 (cell, (1-boundary_points[i][0])*u0+boundary_points[i][0]*u1))
3341 <<
' ' << cell->level()
3342 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n';
3346 <<
' ' << cell->level()
3347 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n' 3349 <<
' ' << cell->level()
3350 <<
' ' <<
static_cast<unsigned int>(cell->material_id()) <<
'\n';
3352 out <<
'\n' <<
'\n';
3359 if (q_projector != 0)
3374 template <
int dim,
int spacedim>
3380 internal::write_gnuplot (tria, out, mapping,
gnuplot_flags);
3398 const unsigned int l)
3400 first(f), second(s),
3401 colorize(c), level(l)
3406 void write_eps (const ::Triangulation<1> &,
3412 Assert(
false, ExcNotImplemented());
3415 void write_eps (const ::Triangulation<1,2> &,
3421 Assert(
false, ExcNotImplemented());
3424 void write_eps (const ::Triangulation<1,3> &,
3430 Assert(
false, ExcNotImplemented());
3433 void write_eps (const ::Triangulation<2,3> &,
3439 Assert(
false, ExcNotImplemented());
3444 template <
int dim,
int spacedim>
3445 void write_eps (const ::Triangulation<dim, spacedim> &tria,
3451 typedef std::list<LineEntry> LineList;
3459 &eps_flags_base = (dim==2 ?
3462 static_cast<const GridOutFlags::EpsFlagsBase &>(eps_flags_3) :
3485 Assert(
false, ExcInternalError());
3491 for (typename ::Triangulation<dim, spacedim>::active_cell_iterator
3492 cell=tria.begin_active();
3493 cell!=tria.end(); ++cell)
3494 for (
unsigned int line_no=0;
3495 line_no<GeometryInfo<dim>::lines_per_cell; ++line_no)
3497 typename ::Triangulation<dim, spacedim>::line_iterator
3498 line=cell->line(line_no);
3510 if (!line->has_children() &&
3511 (mapping==0 || !line->at_boundary()))
3531 line_list.push_back (LineEntry(
Point<2>(line->vertex(0)(0),
3532 line->vertex(0)(1)),
3534 line->vertex(1)(1)),
3535 line->user_flag_set(),
3552 std::vector<
Point<dim-1> > boundary_points (n_points);
3554 for (
unsigned int i=0; i<n_points; ++i)
3555 boundary_points[i](0) = 1.*(i+1)/(n_points+1);
3557 Quadrature<dim-1> quadrature (boundary_points);
3564 for (typename ::Triangulation<dim, spacedim>::active_cell_iterator
3565 cell=tria.begin_active();
3566 cell!=tria.end(); ++cell)
3567 for (
unsigned int face_no=0; face_no<GeometryInfo<dim>::faces_per_cell; ++face_no)
3569 const typename ::Triangulation<dim, spacedim>::face_iterator
3570 face = cell->face(face_no);
3572 if (face->at_boundary())
3575 Point<2> p0 (p0_dim(0), p0_dim(1));
3582 const unsigned int offset=face_no*n_points;
3583 for (
unsigned int i=0; i<n_points; ++i)
3586 (cell, q_projector.
point(offset+i)));
3587 const Point<2> p1 (p1_dim(0), p1_dim(1));
3589 line_list.push_back (LineEntry(p0, p1,
3590 face->user_flag_set(),
3597 const Point<2> p1 (p1_dim(0), p1_dim(1));
3598 line_list.push_back (LineEntry(p0, p1,
3599 face->user_flag_set(),
3612 Assert (mapping == 0, ExcNotImplemented());
3614 typename ::Triangulation<dim, spacedim>::active_cell_iterator
3615 cell=tria.begin_active(),
3637 const double turn_angle = eps_flags_3.
turn_angle;
3638 const Point<dim> view_direction(-std::sin(z_angle * 2.*pi / 360.) * std::sin(turn_angle * 2.*pi / 360.),
3639 +std::sin(z_angle * 2.*pi / 360.) * std::cos(turn_angle * 2.*pi / 360.),
3640 -std::cos(z_angle * 2.*pi / 360.));
3656 - ((
Point<dim>(1,0,0) * view_direction) * view_direction)
3657 - ((
Point<dim>(1,0,0) * unit_vector1) * unit_vector1));
3661 for (; cell!=endc; ++cell)
3662 for (
unsigned int line_no=0;
3663 line_no<GeometryInfo<dim>::lines_per_cell; ++line_no)
3665 typename ::Triangulation<dim, spacedim>::line_iterator
3666 line=cell->line(line_no);
3667 line_list.push_back (LineEntry(
Point<2>(line->vertex(0) * unit_vector2,
3668 line->vertex(0) * unit_vector1),
3669 Point<2>(line->vertex(1) * unit_vector2,
3670 line->vertex(1) * unit_vector1),
3671 line->user_flag_set(),
3679 Assert (
false, ExcNotImplemented());
3687 double x_min = tria.begin_active()->vertex(0)(0);
3688 double x_max = x_min;
3689 double y_min = tria.begin_active()->vertex(0)(1);
3690 double y_max = y_min;
3691 unsigned int max_level = line_list.begin()->level;
3693 for (LineList::const_iterator line=line_list.begin();
3694 line!=line_list.end(); ++line)
3696 x_min = std::min (x_min, line->first(0));
3697 x_min = std::min (x_min, line->second(0));
3699 x_max = std::max (x_max, line->first(0));
3700 x_max = std::max (x_max, line->second(0));
3702 y_min = std::min (y_min, line->first(1));
3703 y_min = std::min (y_min, line->second(1));
3705 y_max = std::max (y_max, line->first(1));
3706 y_max = std::max (y_max, line->second(1));
3708 max_level = std::max (max_level, line->level);
3716 const double scale = (eps_flags_base.
size /
3717 (eps_flags_base.
size_type==GridOutFlags::EpsFlagsBase::width ?
3728 std::time_t time1= std::time (0);
3729 std::tm *time = std::localtime(&time1);
3730 out <<
"%!PS-Adobe-2.0 EPSF-1.2" <<
'\n' 3731 <<
"%%Title: deal.II Output" <<
'\n' 3732 <<
"%%Creator: the deal.II library" <<
'\n' 3733 <<
"%%Creation Date: " 3734 << time->tm_year+1900 <<
"/" 3735 << time->tm_mon+1 <<
"/" 3736 << time->tm_mday <<
" - " 3737 << time->tm_hour <<
":" 3738 << std::setw(2) << time->tm_min <<
":" 3739 << std::setw(2) << time->tm_sec <<
'\n' 3740 <<
"%%BoundingBox: " 3744 <<
static_cast<unsigned int>(std::floor(( (x_max-x_min) * scale )+1))
3746 << static_cast<unsigned int>(std::floor(( (y_max-y_min) * scale )+1))
3755 out <<
"/m {moveto} bind def" <<
'\n' 3756 <<
"/x {lineto stroke} bind def" <<
'\n' 3757 <<
"/b {0 0 0 setrgbcolor} def" <<
'\n' 3758 <<
"/r {1 0 0 setrgbcolor} def" <<
'\n';
3768 << (0.66666/std::max(1U,(max_level-1)))
3769 <<
" mul 1 0.8 sethsbcolor} def" <<
'\n';
3782 out << (
"/R {rmoveto} bind def\n" 3783 "/Symbol-Oblique /Symbol findfont [1 0 .167 1 0 0] makefont\n" 3784 "dup length dict begin {1 index /FID eq {pop pop} {def} ifelse} forall\n" 3785 "currentdict end definefont\n" 3786 "/MFshow {{dup dup 0 get findfont exch 1 get scalefont setfont\n" 3787 "[ currentpoint ] exch dup 2 get 0 exch rmoveto dup dup 5 get exch 4 get\n" 3788 "{show} {stringwidth pop 0 rmoveto}ifelse dup 3 get\n" 3789 "{2 get neg 0 exch rmoveto pop} {pop aload pop moveto}ifelse} forall} bind def\n" 3790 "/MFwidth {0 exch {dup 3 get{dup dup 0 get findfont exch 1 get scalefont setfont\n" 3791 "5 get stringwidth pop add}\n" 3792 "{pop} ifelse} forall} bind def\n" 3793 "/MCshow { currentpoint stroke m\n" 3794 "exch dup MFwidth -2 div 3 -1 roll R MFshow } def\n")
3798 out <<
"%%EndProlog" <<
'\n' 3802 out << eps_flags_base.
line_width <<
" setlinewidth" <<
'\n';
3806 const Point<2> offset(x_min, y_min);
3808 for (LineList::const_iterator line=line_list.begin();
3809 line!=line_list.end(); ++line)
3817 << (line->first - offset) * scale <<
" m " 3818 << (line->second - offset) * scale <<
" x" <<
'\n';
3821 << (line->first - offset) * scale <<
" m " 3822 << (line->second - offset) * scale <<
" x" <<
'\n';
3828 out <<
"(Helvetica) findfont 140 scalefont setfont" 3831 typename ::Triangulation<dim, spacedim>::active_cell_iterator
3832 cell = tria.begin_active (),
3834 for (; cell!=endc; ++cell)
3836 out << (cell->center()(0)-offset(0))*scale <<
' ' 3837 << (cell->center()(1)-offset(1))*scale
3839 <<
"[ [(Helvetica) 12.0 0.0 true true (";
3843 out << cell->index();
3854 out <<
"(Helvetica) findfont 140 scalefont setfont" 3861 std::set<unsigned int> treated_vertices;
3862 typename ::Triangulation<dim, spacedim>::active_cell_iterator
3863 cell = tria.begin_active (),
3865 for (; cell!=endc; ++cell)
3866 for (
unsigned int vertex=0;
3867 vertex<GeometryInfo<dim>::vertices_per_cell;
3869 if (treated_vertices.find(cell->vertex_index(vertex))
3871 treated_vertices.end())
3873 treated_vertices.insert (cell->vertex_index(vertex));
3875 out << (cell->vertex(vertex)(0)-offset(0))*scale <<
' ' 3876 << (cell->vertex(vertex)(1)-offset(1))*scale
3878 <<
"[ [(Helvetica) 10.0 0.0 true true (" 3879 << cell->vertex_index(vertex)
3886 out <<
"showpage" <<
'\n';
3898 template <
int dim,
int spacedim>
3903 internal::write_eps (tria, out, mapping,
3908 template <
int dim,
int spacedim>
3914 switch (output_format)
3960 Assert (
false, ExcInternalError());
3964 template <
int dim,
int spacedim>
3974 #include "grid_out.inst" 3977 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
OutputFormat default_format
static void declare_parameters(ParameterHandler &prm)
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)
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.
::ExceptionBase & ExcMessage(std::string arg1)
void write_vtu(const Triangulation< dim, spacedim > &tria, std::ostream &out) const
std::string get(const std::string &entry_string) const
void parse_parameters(ParameterHandler ¶m)
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.
void write_eps(const Triangulation< dim, spacedim > &tria, std::ostream &out, const Mapping< dim, spacedim > *mapping=0) const
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()
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.
void write_gnuplot(const Triangulation< dim, spacedim > &tria, std::ostream &out, const Mapping< dim, spacedim > *mapping=0) const
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)
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)
static void declare_parameters(ParameterHandler ¶m)
float cell_font_scaling
Scaling of the font for cell annotations. Defaults to 1.
Convert the level subdomain id into the cell color.
void write_vtk(const std::vector< Patch< dim, spacedim > > &patches, const std::vector< std::string > &data_names, const std::vector< std_cxx11::tuple< unsigned int, unsigned int, std::string > > &vector_data_ranges, const VtkFlags &flags, std::ostream &out)
write() calls write_xfig()
void write_vtu(const std::vector< Patch< dim, spacedim > > &patches, const std::vector< std::string > &data_names, const std::vector< std_cxx11::tuple< unsigned int, unsigned int, std::string > > &vector_data_ranges, const VtkFlags &flags, std::ostream &out)
bool write_vertex_numbers
Convert the level subdomain id into the cell color.
static Quadrature< dim > project_to_all_faces(const SubQuadrature &quadrature)
void write(const Triangulation< dim, spacedim > &tria, std::ostream &out, const OutputFormat output_format, const Mapping< dim, spacedim > *mapping=0) const
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 declare_entry(const std::string &entry, const std::string &default_value, const Patterns::PatternBase &pattern=Patterns::Anything(), const std::string &documentation=std::string())
Convert the material id into the cell color (default)
float level_height_factor
The factor determining the vertical distance between levels (default = 0.3)
unsigned char boundary_id
GridOutFlags::Eps< 3 > eps_flags_3
face_iterator end_face() const
const types::boundary_id internal_face_boundary_id
write() calls write_vtu()
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
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)
void write_xfig(const Triangulation< dim, spacedim > &tria, std::ostream &out, const Mapping< dim, spacedim > *mapping=0) const
write() calls write_vtk()
void parse_parameters(ParameterHandler ¶m)
void parse_parameters(ParameterHandler ¶m)
void write_msh(const Triangulation< dim, spacedim > &tria, std::ostream &out) const
std::string default_suffix() const