32#include <boost/algorithm/string.hpp>
33#include <boost/archive/binary_oarchive.hpp>
35#ifdef DEAL_II_GMSH_WITH_API
54 const bool write_faces,
55 const bool write_diameter,
56 const bool write_measure,
57 const bool write_all_faces)
58 : write_cells(write_cells)
59 , write_faces(write_faces)
60 , write_diameter(write_diameter)
61 , write_measure(write_measure)
62 , write_all_faces(write_all_faces)
71 "Write the mesh connectivity as DX grid cells");
75 "Write faces of cells. These may be boundary faces "
76 "or all faces between mesh cells, according to "
77 "\"Write all faces\"");
81 "If cells are written, additionally write their"
82 " diameter as data for visualization");
86 "Write the volume of each cell as data");
90 "Write all faces, not only boundary");
104 Msh::Msh(
const bool write_faces,
const bool write_lines)
105 : write_faces(write_faces)
106 , write_lines(write_lines)
126 const bool write_faces,
127 const bool write_lines)
128 : write_preamble(write_preamble)
129 , write_faces(write_faces)
130 , write_lines(write_lines)
155 const unsigned int n_extra_curved_line_points,
156 const bool curved_inner_cells,
157 const bool write_additional_boundary_lines)
158 : write_cell_numbers(write_cell_numbers)
159 , n_extra_curved_line_points(n_extra_curved_line_points)
160 , curved_inner_cells(curved_inner_cells)
161 , write_additional_boundary_lines(write_additional_boundary_lines)
183 const unsigned int size,
184 const double line_width,
185 const bool color_lines_on_user_flag,
186 const unsigned int n_boundary_face_points,
187 const bool color_lines_level)
188 : size_type(size_type)
190 , line_width(line_width)
191 , color_lines_on_user_flag(color_lines_on_user_flag)
192 , n_boundary_face_points(n_boundary_face_points)
193 , color_lines_level(color_lines_level)
203 "Depending on this parameter, either the "
205 "of the eps is scaled to \"Size\"");
209 "Size of the output in points");
213 "Width of the lines drawn in points");
217 "Draw lines with user flag set in different color");
221 "Number of points on boundary edges. "
222 "Increase this beyond 2 to see curved boundaries.");
226 "Draw different colors according to grid level.");
233 if (param.
get(
"Size by") == std::string(
"width"))
235 else if (param.
get(
"Size by") == std::string(
"height"))
247 const unsigned int size,
248 const double line_width,
249 const bool color_lines_on_user_flag,
250 const unsigned int n_boundary_face_points)
254 color_lines_on_user_flag,
255 n_boundary_face_points)
273 const unsigned int size,
274 const double line_width,
275 const bool color_lines_on_user_flag,
276 const unsigned int n_boundary_face_points,
277 const bool write_cell_numbers,
278 const bool write_cell_number_level,
279 const bool write_vertex_numbers,
280 const bool color_lines_level)
284 color_lines_on_user_flag,
285 n_boundary_face_points,
287 , write_cell_numbers(write_cell_numbers)
288 , write_cell_number_level(write_cell_number_level)
289 , write_vertex_numbers(write_vertex_numbers)
299 "(2D only) Write cell numbers"
300 " into the centers of cells");
304 "(2D only) if \"Cell number\" is true, write "
305 "numbers in the form level.number");
309 "Write numbers for each vertex");
317 write_cell_numbers = param.
get_bool(
"Cell number");
318 write_cell_number_level = param.
get_bool(
"Level number");
319 write_vertex_numbers = param.
get_bool(
"Vertex number");
325 const unsigned int size,
326 const double line_width,
327 const bool color_lines_on_user_flag,
328 const unsigned int n_boundary_face_points,
329 const double azimut_angle,
330 const double turn_angle)
334 color_lines_on_user_flag,
335 n_boundary_face_points)
336 , azimut_angle(azimut_angle)
337 , turn_angle(turn_angle)
347 "Azimuth of the viw point, that is, the angle "
348 "in the plane from the x-axis.");
352 "Elevation of the view point above the xy-plane.");
360 azimut_angle = 90 - param.
get_double(
"Elevation");
367 : draw_boundary(true)
368 , color_by(material_id)
370 , n_boundary_face_points(0)
376 , boundary_thickness(3)
410 const unsigned int boundary_line_thickness,
413 const int azimuth_angle,
414 const int polar_angle,
416 const bool convert_level_number_to_height,
417 const bool label_level_number,
418 const bool label_cell_index,
419 const bool label_material_id,
420 const bool label_subdomain_id,
421 const bool draw_colorbar,
422 const bool draw_legend,
423 const bool label_boundary_id)
426 , line_thickness(line_thickness)
427 , boundary_line_thickness(boundary_line_thickness)
429 , background(background)
430 , azimuth_angle(azimuth_angle)
431 , polar_angle(polar_angle)
433 , convert_level_number_to_height(convert_level_number_to_height)
434 , level_height_factor(0.3f)
435 , cell_font_scaling(1.f)
436 , label_level_number(label_level_number)
437 , label_cell_index(label_cell_index)
438 , label_material_id(label_material_id)
439 , label_subdomain_id(label_subdomain_id)
440 , label_level_subdomain_id(false)
441 , label_boundary_id(label_boundary_id)
442 , draw_colorbar(draw_colorbar)
443 , draw_legend(draw_legend)
447 : draw_bounding_box(false)
466 : default_format(none)
560 switch (output_format)
603 if (format_name ==
"none" || format_name ==
"false")
606 if (format_name ==
"dx")
609 if (format_name ==
"ucd")
612 if (format_name ==
"gnuplot")
615 if (format_name ==
"eps")
618 if (format_name ==
"xfig")
621 if (format_name ==
"msh")
624 if (format_name ==
"svg")
627 if (format_name ==
"mathgl")
630 if (format_name ==
"vtk")
633 if (format_name ==
"vtu")
646 return "none|dx|gnuplot|eps|ucd|xfig|msh|svg|mathgl|vtk|vtu";
779template <
int dim,
int spacedim>
782 std::ostream & out)
const
798 std::vector<unsigned int> renumber(
vertices.size());
801 unsigned int new_number = 0;
802 for (
unsigned int i = 0; i <
vertices.size(); ++i)
804 renumber[i] = new_number++;
808 out <<
"object \"vertices\" class array type float rank 1 shape " << dim
809 <<
" items " << n_vertices <<
" data follows" <<
'\n';
811 for (
unsigned int i = 0; i <
vertices.size(); ++i)
820 const unsigned int n_faces =
828 out <<
"object \"cells\" class array type int rank 1 shape "
829 << n_vertices_per_cell <<
" items " << n_cells <<
" data follows"
840 out <<
"attribute \"element type\" string \"";
848 <<
"attribute \"ref\" string \"positions\"" <<
'\n'
853 out <<
"object \"material\" class array type int rank 0 items " << n_cells
854 <<
" data follows" <<
'\n';
856 out <<
' ' << cell->material_id();
857 out <<
'\n' <<
"attribute \"dep\" string \"connections\"" <<
'\n' <<
'\n';
859 out <<
"object \"level\" class array type int rank 0 items " << n_cells
860 <<
" data follows" <<
'\n';
862 out <<
' ' << cell->level();
863 out <<
'\n' <<
"attribute \"dep\" string \"connections\"" <<
'\n' <<
'\n';
867 out <<
"object \"measure\" class array type float rank 0 items "
868 << n_cells <<
" data follows" <<
'\n';
870 out <<
'\t' << cell->measure();
872 <<
"attribute \"dep\" string \"connections\"" <<
'\n'
878 out <<
"object \"diameter\" class array type float rank 0 items "
879 << n_cells <<
" data follows" <<
'\n';
881 out <<
'\t' << cell->diameter();
883 <<
"attribute \"dep\" string \"connections\"" <<
'\n'
890 out <<
"object \"faces\" class array type int rank 1 shape "
891 << n_vertices_per_face <<
" items " << n_faces <<
" data follows"
896 for (
const unsigned int f : cell->face_indices())
901 for (
unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_face;
904 << renumber[face->vertex_index(
909 out <<
"attribute \"element type\" string \"";
915 <<
"attribute \"ref\" string \"positions\"" <<
'\n'
921 out <<
"object \"boundary\" class array type int rank 0 items " << n_faces
922 <<
" data follows" <<
'\n';
929 <<
static_cast<std::make_signed<types::boundary_id>::type
>(
930 cell->face(f)->boundary_id());
934 out <<
"attribute \"dep\" string \"connections\"" <<
'\n' <<
'\n';
938 out <<
"object \"face measure\" class array type float rank 0 items "
939 << n_faces <<
" data follows" <<
'\n';
943 out <<
' ' << cell->face(f)->measure();
946 out <<
"attribute \"dep\" string \"connections\"" <<
'\n' <<
'\n';
951 out <<
"object \"face diameter\" class array type float rank 0 items "
952 << n_faces <<
" data follows" <<
'\n';
956 out <<
' ' << cell->face(f)->diameter();
959 out <<
"attribute \"dep\" string \"connections\"" <<
'\n' <<
'\n';
972 out <<
"object \"deal data\" class field" <<
'\n'
973 <<
"component \"positions\" value \"vertices\"" <<
'\n'
974 <<
"component \"connections\" value \"cells\"" <<
'\n';
978 out <<
"object \"cell data\" class field" <<
'\n'
979 <<
"component \"positions\" value \"vertices\"" <<
'\n'
980 <<
"component \"connections\" value \"cells\"" <<
'\n';
981 out <<
"component \"material\" value \"material\"" <<
'\n';
982 out <<
"component \"level\" value \"level\"" <<
'\n';
984 out <<
"component \"measure\" value \"measure\"" <<
'\n';
986 out <<
"component \"diameter\" value \"diameter\"" <<
'\n';
991 out <<
"object \"face data\" class field" <<
'\n'
992 <<
"component \"positions\" value \"vertices\"" <<
'\n'
993 <<
"component \"connections\" value \"faces\"" <<
'\n';
994 out <<
"component \"boundary\" value \"boundary\"" <<
'\n';
996 out <<
"component \"measure\" value \"face measure\"" <<
'\n';
998 out <<
"component \"diameter\" value \"face diameter\"" <<
'\n';
1001 out <<
'\n' <<
"object \"grid data\" class group" <<
'\n';
1003 out <<
"member \"cells\" value \"cell data\"" <<
'\n';
1005 out <<
"member \"faces\" value \"face data\"" <<
'\n';
1006 out <<
"end" <<
'\n';
1017template <
int dim,
int spacedim>
1020 std::ostream & out)
const
1048 out <<
"@f$NOD" <<
'\n' << n_vertices <<
'\n';
1053 for (
unsigned int i = 0; i <
vertices.size(); ++i)
1058 for (
unsigned int d = spacedim + 1; d <= 3; ++d)
1064 out <<
"@f$ENDNOD" <<
'\n'
1075 out << cell->active_cell_index() + 1 <<
' '
1076 << cell->reference_cell().gmsh_element_type() <<
' '
1077 << cell->material_id() <<
' ' << cell->subdomain_id() <<
' '
1078 << cell->n_vertices() <<
' ';
1082 for (
const unsigned int vertex : cell->vertex_indices())
1084 if (cell->reference_cell() == ReferenceCells::get_hypercube<dim>())
1088 else if (cell->reference_cell() == ReferenceCells::get_simplex<dim>())
1089 out << cell->vertex_index(vertex) + 1 <<
' ';
1107 out <<
"@f$ENDELM\n";
1117template <
int dim,
int spacedim>
1120 std::ostream & out)
const
1138 std::time_t time1 = std::time(
nullptr);
1139 std::tm * time = std::localtime(&time1);
1141 <<
"# This file was generated by the deal.II library." <<
'\n'
1142 <<
"# Date = " << time->tm_year + 1900 <<
"/" << time->tm_mon + 1
1143 <<
"/" << time->tm_mday <<
'\n'
1144 <<
"# Time = " << time->tm_hour <<
":" << std::setw(2) << time->tm_min
1145 <<
":" << std::setw(2) << time->tm_sec <<
'\n'
1147 <<
"# For a description of the UCD format see the AVS Developer's guide."
1153 out << n_vertices <<
' '
1163 for (
unsigned int i = 0; i <
vertices.size(); ++i)
1168 for (
unsigned int d = spacedim + 1; d <= 3; ++d)
1177 out << cell->active_cell_index() + 1 <<
' ' << cell->material_id() <<
' ';
1234template <
int dim,
int spacedim>
1252 const int spacedim = 2;
1258 out <<
"#FIG 3.2\nLandscape\nCenter\nInches" << std::endl
1259 <<
"A4\n100.00\nSingle"
1262 <<
"-3" << std::endl
1263 <<
"# generated by deal.II GridOut class" << std::endl
1264 <<
"# reduce first number to scale up image" << std::endl
1265 <<
"1200 2" << std::endl;
1268 unsigned int colno = 32;
1269 out <<
"0 " << colno++ <<
" #ff0000" << std::endl;
1270 out <<
"0 " << colno++ <<
" #ff8000" << std::endl;
1271 out <<
"0 " << colno++ <<
" #ffd000" << std::endl;
1272 out <<
"0 " << colno++ <<
" #ffff00" << std::endl;
1273 out <<
"0 " << colno++ <<
" #c0ff00" << std::endl;
1274 out <<
"0 " << colno++ <<
" #80ff00" << std::endl;
1275 out <<
"0 " << colno++ <<
" #00f000" << std::endl;
1276 out <<
"0 " << colno++ <<
" #00f0c0" << std::endl;
1277 out <<
"0 " << colno++ <<
" #00f0ff" << std::endl;
1278 out <<
"0 " << colno++ <<
" #00c0ff" << std::endl;
1279 out <<
"0 " << colno++ <<
" #0080ff" << std::endl;
1280 out <<
"0 " << colno++ <<
" #0040ff" << std::endl;
1281 out <<
"0 " << colno++ <<
" #0000c0" << std::endl;
1282 out <<
"0 " << colno++ <<
" #5000ff" << std::endl;
1283 out <<
"0 " << colno++ <<
" #8000ff" << std::endl;
1284 out <<
"0 " << colno++ <<
" #b000ff" << std::endl;
1285 out <<
"0 " << colno++ <<
" #ff00ff" << std::endl;
1286 out <<
"0 " << colno++ <<
" #ff80ff" << std::endl;
1288 for (
unsigned int i = 0; i < 8; ++i)
1289 out <<
"0 " << colno++ <<
" #" << std::hex << 32 * i + 31 << 32 * i + 31
1290 << 32 * i + 31 << std::dec << std::endl;
1292 for (
unsigned int i = 1; i < 16; ++i)
1293 out <<
"0 " << colno++ <<
" #00" << std::hex << 16 * i + 15 << std::dec
1294 <<
"00" << std::endl;
1296 for (
unsigned int i = 1; i < 16; ++i)
1297 out <<
"0 " << colno++ <<
" #" << std::hex << 16 * i + 15 << 16 * i + 15
1298 << std::dec <<
"00" << std::endl;
1300 for (
unsigned int i = 1; i < 16; ++i)
1301 out <<
"0 " << colno++ <<
" #" << std::hex << 16 * i + 15 << std::dec
1302 <<
"0000" << std::endl;
1304 for (
unsigned int i = 1; i < 16; ++i)
1305 out <<
"0 " << colno++ <<
" #" << std::hex << 16 * i + 15 <<
"00"
1306 << 16 * i + 15 << std::dec << std::endl;
1308 for (
unsigned int i = 1; i < 16; ++i)
1309 out <<
"0 " << colno++ <<
" #0000" << std::hex << 16 * i + 15 << std::dec
1312 for (
unsigned int i = 1; i < 16; ++i)
1313 out <<
"0 " << colno++ <<
" #00" << std::hex << 16 * i + 15 << 16 * i + 15
1314 << std::dec << std::endl;
1336 out << cell->material_id() + 32;
1339 out << cell->level() + 8;
1342 out << cell->subdomain_id() + 32;
1345 out << cell->level_subdomain_id() + 32;
1354 (900 + cell->material_id()))
1360 << nv + 1 << std::endl;
1366 for (
unsigned int k = 0; k <= nv; ++k)
1370 for (
unsigned int d = 0; d < static_cast<unsigned int>(dim); ++d)
1374 out <<
'\t' << ((d == 0) ? val : -val);
1379 static const unsigned int face_reorder[4] = {2, 1, 3, 0};
1381 for (
const unsigned int f : face_reorder)
1409 for (
unsigned int k = 0;
1410 k < GeometryInfo<dim>::vertices_per_face;
1414 for (
unsigned int d = 0; d < static_cast<unsigned int>(dim);
1420 out <<
'\t' << ((d == 0) ? val : -val);
1437#ifdef DEAL_II_GMSH_WITH_API
1438template <
int dim,
int spacedim>
1441 const std::string & filename)
const
1444 const std::array<int, 8> dealii_to_gmsh_type = {{15, 1, 2, 3, 4, 7, 6, 5}};
1447 const std::array<std::vector<unsigned int>, 8> dealii_to_gmsh = {
1454 {{0, 1, 2, 3, 4, 5}},
1455 {{0, 1, 3, 2, 4, 5, 7, 6}}}};
1460 std::vector<double> coords(3 *
vertices.size());
1461 std::vector<std::size_t> nodes(
vertices.size());
1467 for (
unsigned int d = 0; d < spacedim; ++d)
1468 coords[i * 3 + d] = p[d];
1480 using IdPair = std::pair<types::material_id, types::manifold_id>;
1481 std::map<IdPair, int> id_pair_to_entity_tag;
1482 std::vector<IdPair> all_pairs;
1484 std::set<IdPair> set_of_pairs;
1487 set_of_pairs.insert({cell->material_id(), cell->manifold_id()});
1488 for (
const auto &f : cell->face_iterators())
1490 (f->boundary_id() != 0 &&
1492 set_of_pairs.insert({f->boundary_id(), f->manifold_id()});
1494 for (
const auto l : cell->line_indices())
1496 const auto &f = cell->line(l);
1498 (f->boundary_id() != 0 &&
1500 set_of_pairs.insert({f->boundary_id(), f->manifold_id()});
1503 all_pairs = {set_of_pairs.begin(), set_of_pairs.end()};
1506 for (
const auto &p : set_of_pairs)
1507 id_pair_to_entity_tag[p] = entity++;
1510 const auto n_entity_tags = id_pair_to_entity_tag.size();
1513 std::vector<std::vector<std::vector<std::size_t>>> element_ids(
1514 n_entity_tags, std::vector<std::vector<std::size_t>>(8));
1515 std::vector<std::vector<std::vector<std::size_t>>> element_nodes(
1516 n_entity_tags, std::vector<std::vector<std::size_t>>(8));
1519 std::size_t element_id = 1;
1521 const auto add_element = [&](
const auto &element,
const int &entity_tag) {
1522 const auto type = element->reference_cell();
1527 for (
const auto v : element->vertex_indices())
1528 element_nodes[entity_tag - 1][type].emplace_back(
1529 element->vertex_index(dealii_to_gmsh[type][v]) + 1);
1532 element_ids[entity_tag - 1][type].emplace_back(element_id);
1540 std::set<std::pair<int, int>> dim_entity_tag;
1542 auto maybe_add_element =
1543 [&](
const auto & element,
1545 const auto struct_dim = element->structure_dimension;
1546 const auto manifold_id = element->manifold_id();
1549 const bool non_default_boundary_or_material_id =
1550 (boundary_or_material_id != 0 &&
1552 const bool non_default_manifold =
1554 if (struct_dim == dim || non_default_boundary_or_material_id ||
1555 non_default_manifold)
1557 const auto entity_tag =
1558 id_pair_to_entity_tag[{boundary_or_material_id, manifold_id}];
1559 add_element(element, entity_tag);
1560 dim_entity_tag.insert({struct_dim, entity_tag});
1567 maybe_add_element(cell, cell->material_id());
1568 for (
const auto &face : cell->face_iterators())
1569 maybe_add_element(face, face->boundary_id());
1571 for (
const auto l : cell->line_indices())
1572 maybe_add_element(cell->line(l), cell->line(l)->boundary_id());
1577 gmsh::option::setNumber(
"General.Verbosity", 0);
1578 gmsh::model::add(
"Grid generated in deal.II");
1579 for (
const auto &p : dim_entity_tag)
1581 gmsh::model::addDiscreteEntity(p.first, p.second);
1582 gmsh::model::mesh::addNodes(p.first, p.second, nodes, coords);
1585 for (
unsigned int entity_tag = 0; entity_tag < n_entity_tags; ++entity_tag)
1586 for (
unsigned int t = 1; t < 8; ++t)
1588 const auto all_element_ids = element_ids[entity_tag][t];
1589 const auto all_element_nodes = element_nodes[entity_tag][t];
1590 const auto gmsh_t = dealii_to_gmsh_type[t];
1591 if (all_element_ids.size() > 0)
1592 gmsh::model::mesh::addElementsByType(entity_tag + 1,
1600 for (
const auto &it : dim_entity_tag)
1602 const auto &d = it.first;
1603 const auto &entity_tag = it.second;
1604 const auto &boundary_id = all_pairs[entity_tag - 1].first;
1605 const auto &manifold_id = all_pairs[entity_tag - 1].second;
1607 std::string physical_name;
1608 if (d == dim && boundary_id != 0)
1610 static_cast<int>(boundary_id));
1611 else if (d < dim && boundary_id != 0)
1618 std::string sep = physical_name !=
"" ?
", " :
"";
1621 sep +
"ManifoldID:" +
1623 const auto physical_tag =
1624 gmsh::model::addPhysicalGroup(d, {entity_tag}, -1);
1625 if (physical_name !=
"")
1626 gmsh::model::setPhysicalName(d, physical_tag, physical_name);
1630 gmsh::write(filename);
1651 svg_project_point(
const Point<3> & point,
1655 const float camera_focus)
1661 camera_focus / ((point - camera_position) * camera_direction);
1664 camera_position + phi * (point - camera_position);
1666 return {(projection - camera_position - camera_focus * camera_direction) *
1668 (projection - camera_position - camera_focus * camera_direction) *
1675template <
int dim,
int spacedim>
1678 std::ostream & )
const
1681 ExcMessage(
"Mesh output in SVG format is not implemented for anything "
1682 "other than two-dimensional meshes in two-dimensional "
1683 "space. That's because three-dimensional meshes are best "
1684 "viewed in programs that allow changing the viewpoint, "
1685 "but SVG format does not allow this: It is an inherently "
1686 "2d format, and for three-dimensional meshes would "
1687 "require choosing one, fixed viewpoint."
1689 "You probably want to output your mesh in a format such "
1690 "as VTK, VTU, or gnuplot."));
1699 unsigned int min_level, max_level;
1707 Assert(height != 0 || width != 0,
1708 ExcMessage(
"You have to set at least one of width and height"));
1710 unsigned int margin_in_percent = 0;
1712 margin_in_percent = 8;
1715 unsigned int cell_label_font_size;
1728 float x_max_perspective, x_min_perspective;
1729 float y_max_perspective, y_min_perspective;
1731 float x_dimension_perspective, y_dimension_perspective;
1735 double x_min =
tria.
begin()->vertex(0)[0];
1736 double x_max = x_min;
1737 double y_min =
tria.
begin()->vertex(0)[1];
1738 double y_max = y_min;
1740 double x_dimension, y_dimension;
1742 min_level = max_level =
tria.
begin()->level();
1745 std::set<unsigned int> materials;
1748 std::set<unsigned int> levels;
1751 std::set<unsigned int> subdomains;
1754 std::set<int> level_subdomains;
1762 for (
unsigned int vertex_index = 0; vertex_index < cell->n_vertices();
1765 if (cell->vertex(vertex_index)[0] < x_min)
1766 x_min = cell->vertex(vertex_index)[0];
1767 if (cell->vertex(vertex_index)[0] > x_max)
1768 x_max = cell->vertex(vertex_index)[0];
1770 if (cell->vertex(vertex_index)[1] < y_min)
1771 y_min = cell->vertex(vertex_index)[1];
1772 if (cell->vertex(vertex_index)[1] > y_max)
1773 y_max = cell->vertex(vertex_index)[1];
1776 if (
static_cast<unsigned int>(cell->level()) < min_level)
1777 min_level = cell->level();
1778 if (
static_cast<unsigned int>(cell->level()) > max_level)
1779 max_level = cell->level();
1781 materials.insert(cell->material_id());
1782 levels.insert(cell->level());
1783 if (cell->is_active())
1784 subdomains.insert(cell->subdomain_id() + 2);
1785 level_subdomains.insert(cell->level_subdomain_id() + 2);
1788 x_dimension = x_max - x_min;
1789 y_dimension = y_max - y_min;
1792 const unsigned int n_materials = materials.size();
1795 const unsigned int n_levels = levels.size();
1798 const unsigned int n_subdomains = subdomains.size();
1801 const unsigned int n_level_subdomains = level_subdomains.size();
1815 n = n_level_subdomains;
1824 camera_position[0] = 0;
1825 camera_position[1] = 0;
1826 camera_position[2] = 2. *
std::max(x_dimension, y_dimension);
1829 camera_direction[0] = 0;
1830 camera_direction[1] = 0;
1831 camera_direction[2] = -1;
1834 camera_horizontal[0] = 1;
1835 camera_horizontal[1] = 0;
1836 camera_horizontal[2] = 0;
1838 camera_focus = .5 *
std::max(x_dimension, y_dimension);
1844 const double angle_factor = 3.14159265 / 180.;
1847 camera_position_temp[1] =
1850 camera_position_temp[2] =
1854 camera_direction_temp[1] =
1857 camera_direction_temp[2] =
1861 camera_horizontal_temp[1] =
1864 camera_horizontal_temp[2] =
1868 camera_position[1] = camera_position_temp[1];
1869 camera_position[2] = camera_position_temp[2];
1871 camera_direction[1] = camera_direction_temp[1];
1872 camera_direction[2] = camera_direction_temp[2];
1874 camera_horizontal[1] = camera_horizontal_temp[1];
1875 camera_horizontal[2] = camera_horizontal_temp[2];
1878 camera_position_temp[0] =
1881 camera_position_temp[1] =
1885 camera_direction_temp[0] =
1888 camera_direction_temp[1] =
1892 camera_horizontal_temp[0] =
1895 camera_horizontal_temp[1] =
1899 camera_position[0] = camera_position_temp[0];
1900 camera_position[1] = camera_position_temp[1];
1902 camera_direction[0] = camera_direction_temp[0];
1903 camera_direction[1] = camera_direction_temp[1];
1905 camera_horizontal[0] = camera_horizontal_temp[0];
1906 camera_horizontal[1] = camera_horizontal_temp[1];
1909 camera_position[0] = x_min + .5 * x_dimension;
1910 camera_position[1] = y_min + .5 * y_dimension;
1912 camera_position[0] += 2. *
std::max(x_dimension, y_dimension) *
1915 camera_position[1] -= 2. *
std::max(x_dimension, y_dimension) *
1926 float min_level_min_vertex_distance = 0;
1931 (
static_cast<float>(
tria.
begin()->level()) /
1932 static_cast<float>(n_levels)) *
1933 std::max(x_dimension, y_dimension);
1936 projection_decomposition = svg_project_point(
1937 point, camera_position, camera_direction, camera_horizontal, camera_focus);
1939 x_max_perspective = projection_decomposition[0];
1940 x_min_perspective = projection_decomposition[0];
1942 y_max_perspective = projection_decomposition[1];
1943 y_min_perspective = projection_decomposition[1];
1947 point[0] = cell->vertex(0)[0];
1948 point[1] = cell->vertex(0)[1];
1955 (
static_cast<float>(cell->level()) /
static_cast<float>(n_levels)) *
1956 std::max(x_dimension, y_dimension);
1959 projection_decomposition = svg_project_point(point,
1965 if (x_max_perspective < projection_decomposition[0])
1966 x_max_perspective = projection_decomposition[0];
1967 if (x_min_perspective > projection_decomposition[0])
1968 x_min_perspective = projection_decomposition[0];
1970 if (y_max_perspective < projection_decomposition[1])
1971 y_max_perspective = projection_decomposition[1];
1972 if (y_min_perspective > projection_decomposition[1])
1973 y_min_perspective = projection_decomposition[1];
1975 point[0] = cell->vertex(1)[0];
1976 point[1] = cell->vertex(1)[1];
1978 projection_decomposition = svg_project_point(point,
1984 if (x_max_perspective < projection_decomposition[0])
1985 x_max_perspective = projection_decomposition[0];
1986 if (x_min_perspective > projection_decomposition[0])
1987 x_min_perspective = projection_decomposition[0];
1989 if (y_max_perspective < projection_decomposition[1])
1990 y_max_perspective = projection_decomposition[1];
1991 if (y_min_perspective > projection_decomposition[1])
1992 y_min_perspective = projection_decomposition[1];
1994 point[0] = cell->vertex(2)[0];
1995 point[1] = cell->vertex(2)[1];
1997 projection_decomposition = svg_project_point(point,
2003 if (x_max_perspective < projection_decomposition[0])
2004 x_max_perspective = projection_decomposition[0];
2005 if (x_min_perspective > projection_decomposition[0])
2006 x_min_perspective = projection_decomposition[0];
2008 if (y_max_perspective < projection_decomposition[1])
2009 y_max_perspective = projection_decomposition[1];
2010 if (y_min_perspective > projection_decomposition[1])
2011 y_min_perspective = projection_decomposition[1];
2013 if (cell->n_vertices() == 4)
2015 point[0] = cell->vertex(3)[0];
2016 point[1] = cell->vertex(3)[1];
2018 projection_decomposition = svg_project_point(point,
2024 if (x_max_perspective < projection_decomposition[0])
2025 x_max_perspective = projection_decomposition[0];
2026 if (x_min_perspective > projection_decomposition[0])
2027 x_min_perspective = projection_decomposition[0];
2029 if (y_max_perspective < projection_decomposition[1])
2030 y_max_perspective = projection_decomposition[1];
2031 if (y_min_perspective > projection_decomposition[1])
2032 y_min_perspective = projection_decomposition[1];
2035 if (
static_cast<unsigned int>(cell->level()) == min_level)
2036 min_level_min_vertex_distance = cell->minimum_vertex_distance();
2039 x_dimension_perspective = x_max_perspective - x_min_perspective;
2040 y_dimension_perspective = y_max_perspective - y_min_perspective;
2044 width =
static_cast<unsigned int>(
2045 .5 + height * (x_dimension_perspective / y_dimension_perspective));
2046 else if (height == 0)
2047 height =
static_cast<unsigned int>(
2048 .5 + width * (y_dimension_perspective / x_dimension_perspective));
2049 unsigned int additional_width = 0;
2051 unsigned int font_size =
2052 static_cast<unsigned int>(.5 + (height / 100.) * 1.75);
2053 cell_label_font_size =
static_cast<unsigned int>(
2055 min_level_min_vertex_distance /
std::min(x_dimension, y_dimension)));
2062 additional_width =
static_cast<unsigned int>(
2067 additional_width =
static_cast<unsigned int>(
2068 .5 + height * .175);
2080 out <<
"<svg width=\"" << width + additional_width <<
"\" height=\"" << height
2081 <<
"\" xmlns=\"http://www.w3.org/2000/svg\" version=\"1.1\">" <<
'\n'
2088 <<
" <linearGradient id=\"background_gradient\" gradientUnits=\"userSpaceOnUse\" x1=\"0\" y1=\"0\" x2=\"0\" y2=\""
2089 << height <<
"\">" <<
'\n'
2090 <<
" <stop offset=\"0\" style=\"stop-color:white\"/>" <<
'\n'
2091 <<
" <stop offset=\"1\" style=\"stop-color:lightsteelblue\"/>" <<
'\n'
2092 <<
" </linearGradient>" <<
'\n';
2098 out <<
"<!-- internal style sheet -->" <<
'\n'
2099 <<
"<style type=\"text/css\"><![CDATA[" <<
'\n';
2103 out <<
" rect.background{fill:url(#background_gradient)}" <<
'\n';
2105 out <<
" rect.background{fill:white}" <<
'\n';
2107 out <<
" rect.background{fill:none}" <<
'\n';
2110 out <<
" rect{fill:none; stroke:rgb(25,25,25); stroke-width:"
2112 <<
" text{font-family:Helvetica; text-anchor:middle; fill:rgb(25,25,25)}"
2114 <<
" line{stroke:rgb(25,25,25); stroke-width:"
2116 <<
" path{fill:none; stroke:rgb(25,25,25); stroke-width:"
2118 <<
" circle{fill:white; stroke:black; stroke-width:2}" <<
'\n'
2124 unsigned int labeling_index = 0;
2125 auto materials_it = materials.begin();
2126 auto levels_it = levels.begin();
2127 auto subdomains_it = subdomains.begin();
2128 auto level_subdomains_it = level_subdomains.begin();
2130 for (
unsigned int index = 0; index < n; ++index)
2135 h = .6 - (index / (n - 1.)) * .6;
2143 unsigned int i =
static_cast<unsigned int>(h * 6);
2145 double f = h * 6 - i;
2152 r = 255, g =
static_cast<unsigned int>(.5 + 255 * t);
2155 r =
static_cast<unsigned int>(.5 + 255 * q), g = 255;
2158 g = 255, b =
static_cast<unsigned int>(.5 + 255 * t);
2161 g =
static_cast<unsigned int>(.5 + 255 * q), b = 255;
2164 r =
static_cast<unsigned int>(.5 + 255 * t), b = 255;
2167 r = 255, b =
static_cast<unsigned int>(.5 + 255 * q);
2176 labeling_index = *materials_it++;
2179 labeling_index = *levels_it++;
2182 labeling_index = *subdomains_it++;
2185 labeling_index = *level_subdomains_it++;
2191 out <<
" path.p" << labeling_index <<
"{fill:rgb(" << r <<
',' << g
2192 <<
',' << b <<
"); "
2193 <<
"stroke:rgb(25,25,25); stroke-width:"
2196 out <<
" path.ps" << labeling_index <<
"{fill:rgb("
2197 <<
static_cast<unsigned int>(.5 + .75 * r) <<
','
2198 <<
static_cast<unsigned int>(.5 + .75 * g) <<
','
2199 <<
static_cast<unsigned int>(.5 + .75 * b) <<
"); "
2200 <<
"stroke:rgb(20,20,20); stroke-width:"
2203 out <<
" rect.r" << labeling_index <<
"{fill:rgb(" << r <<
',' << g
2204 <<
',' << b <<
"); "
2205 <<
"stroke:rgb(25,25,25); stroke-width:"
2212 out <<
"]]></style>" <<
'\n' <<
'\n';
2215 out <<
" <rect class=\"background\" width=\"" << width <<
"\" height=\""
2216 << height <<
"\"/>" <<
'\n';
2220 unsigned int x_offset = 0;
2223 x_offset =
static_cast<unsigned int>(.5 + (height / 100.) *
2224 (margin_in_percent / 2.));
2226 x_offset =
static_cast<unsigned int>(.5 + height * .025);
2229 <<
" <text x=\"" << x_offset <<
"\" y=\""
2230 <<
static_cast<unsigned int>(.5 + height * .0525) <<
'\"'
2231 <<
" style=\"font-weight:100; fill:lightsteelblue; text-anchor:start; font-family:Courier; font-size:"
2232 <<
static_cast<unsigned int>(.5 + height * .045) <<
"px\">"
2234 <<
"</text>" <<
'\n';
2253 out <<
" <!-- cells -->" <<
'\n';
2255 for (
unsigned int level_index = min_level; level_index <= max_level;
2268 out <<
" class=\"p";
2270 if (!cell->is_active() &&
2277 out << cell->material_id();
2280 out << static_cast<unsigned int>(cell->level());
2283 if (cell->is_active())
2284 out << cell->subdomain_id() + 2;
2289 out << cell->level_subdomain_id() + 2;
2300 point[0] = cell->vertex(0)[0];
2301 point[1] = cell->vertex(0)[1];
2307 (
static_cast<float>(cell->level()) /
2308 static_cast<float>(n_levels)) *
2309 std::max(x_dimension, y_dimension);
2312 projection_decomposition = svg_project_point(point,
2318 out << static_cast<unsigned int>(
2320 ((projection_decomposition[0] - x_min_perspective) /
2321 x_dimension_perspective) *
2322 (width - (width / 100.) * 2. * margin_in_percent) +
2323 ((width / 100.) * margin_in_percent))
2325 <<
static_cast<unsigned int>(
2326 .5 + height - (height / 100.) * margin_in_percent -
2327 ((projection_decomposition[1] - y_min_perspective) /
2328 y_dimension_perspective) *
2329 (height - (height / 100.) * 2. * margin_in_percent));
2333 point[0] = cell->vertex(1)[0];
2334 point[1] = cell->vertex(1)[1];
2336 projection_decomposition = svg_project_point(point,
2342 out << static_cast<unsigned int>(
2344 ((projection_decomposition[0] - x_min_perspective) /
2345 x_dimension_perspective) *
2346 (width - (width / 100.) * 2. * margin_in_percent) +
2347 ((width / 100.) * margin_in_percent))
2349 <<
static_cast<unsigned int>(
2350 .5 + height - (height / 100.) * margin_in_percent -
2351 ((projection_decomposition[1] - y_min_perspective) /
2352 y_dimension_perspective) *
2353 (height - (height / 100.) * 2. * margin_in_percent));
2357 if (cell->n_vertices() == 4)
2359 point[0] = cell->vertex(3)[0];
2360 point[1] = cell->vertex(3)[1];
2362 projection_decomposition = svg_project_point(point,
2368 out << static_cast<unsigned int>(
2370 ((projection_decomposition[0] - x_min_perspective) /
2371 x_dimension_perspective) *
2372 (width - (width / 100.) * 2. * margin_in_percent) +
2373 ((width / 100.) * margin_in_percent))
2375 <<
static_cast<unsigned int>(
2376 .5 + height - (height / 100.) * margin_in_percent -
2377 ((projection_decomposition[1] - y_min_perspective) /
2378 y_dimension_perspective) *
2379 (height - (height / 100.) * 2. * margin_in_percent));
2384 point[0] = cell->vertex(2)[0];
2385 point[1] = cell->vertex(2)[1];
2387 projection_decomposition = svg_project_point(point,
2393 out << static_cast<unsigned int>(
2395 ((projection_decomposition[0] - x_min_perspective) /
2396 x_dimension_perspective) *
2397 (width - (width / 100.) * 2. * margin_in_percent) +
2398 ((width / 100.) * margin_in_percent))
2400 <<
static_cast<unsigned int>(
2401 .5 + height - (height / 100.) * margin_in_percent -
2402 ((projection_decomposition[1] - y_min_perspective) /
2403 y_dimension_perspective) *
2404 (height - (height / 100.) * 2. * margin_in_percent));
2408 point[0] = cell->vertex(0)[0];
2409 point[1] = cell->vertex(0)[1];
2411 projection_decomposition = svg_project_point(point,
2417 out << static_cast<unsigned int>(
2419 ((projection_decomposition[0] - x_min_perspective) /
2420 x_dimension_perspective) *
2421 (width - (width / 100.) * 2. * margin_in_percent) +
2422 ((width / 100.) * margin_in_percent))
2424 <<
static_cast<unsigned int>(
2425 .5 + height - (height / 100.) * margin_in_percent -
2426 ((projection_decomposition[1] - y_min_perspective) /
2427 y_dimension_perspective) *
2428 (height - (height / 100.) * 2. * margin_in_percent));
2430 out <<
"\"/>" <<
'\n';
2437 point[0] = cell->center()[0];
2438 point[1] = cell->center()[1];
2444 (
static_cast<float>(cell->level()) /
2445 static_cast<float>(n_levels)) *
2446 std::max(x_dimension, y_dimension);
2449 const double distance_to_camera =
2451 std::pow(point[1] - camera_position[1], 2.) +
2452 std::pow(point[2] - camera_position[2], 2.));
2453 const double distance_factor =
2454 distance_to_camera / (2. *
std::max(x_dimension, y_dimension));
2456 projection_decomposition = svg_project_point(point,
2462 const unsigned int font_size_this_cell =
2463 static_cast<unsigned int>(
2465 cell_label_font_size *
2466 std::pow(.5, cell->level() - 4. + 3.5 * distance_factor));
2470 <<
static_cast<unsigned int>(
2472 ((projection_decomposition[0] - x_min_perspective) /
2473 x_dimension_perspective) *
2474 (width - (width / 100.) * 2. * margin_in_percent) +
2475 ((width / 100.) * margin_in_percent))
2477 <<
static_cast<unsigned int>(
2478 .5 + height - (height / 100.) * margin_in_percent -
2479 ((projection_decomposition[1] - y_min_perspective) /
2480 y_dimension_perspective) *
2481 (height - (height / 100.) * 2. * margin_in_percent) +
2482 0.5 * font_size_this_cell)
2483 <<
"\" style=\"font-size:" << font_size_this_cell <<
"px\">";
2487 out << cell->level();
2494 out << cell->index();
2503 <<
static_cast<std::make_signed<types::material_id>::type
>(
2504 cell->material_id());
2512 if (cell->is_active())
2514 std::make_signed<types::subdomain_id>::type
>(
2515 cell->subdomain_id());
2528 <<
static_cast<std::make_signed<types::subdomain_id>::type
>(
2529 cell->level_subdomain_id());
2532 out <<
"</text>" <<
'\n';
2539 for (
auto faceIndex : cell->face_indices())
2541 if (cell->at_boundary(faceIndex))
2543 point[0] = cell->face(faceIndex)->vertex(0)[0];
2544 point[1] = cell->face(faceIndex)->vertex(0)[1];
2550 (
static_cast<float>(cell->level()) /
2551 static_cast<float>(n_levels)) *
2552 std::max(x_dimension, y_dimension);
2555 projection_decomposition =
2556 svg_project_point(point,
2562 out <<
" <line x1=\""
2563 <<
static_cast<unsigned int>(
2565 ((projection_decomposition[0] -
2566 x_min_perspective) /
2567 x_dimension_perspective) *
2569 (width / 100.) * 2. * margin_in_percent) +
2570 ((width / 100.) * margin_in_percent))
2572 <<
static_cast<unsigned int>(
2574 (height / 100.) * margin_in_percent -
2575 ((projection_decomposition[1] -
2576 y_min_perspective) /
2577 y_dimension_perspective) *
2579 (height / 100.) * 2. * margin_in_percent));
2581 point[0] = cell->face(faceIndex)->vertex(1)[0];
2582 point[1] = cell->face(faceIndex)->vertex(1)[1];
2588 (
static_cast<float>(cell->level()) /
2589 static_cast<float>(n_levels)) *
2590 std::max(x_dimension, y_dimension);
2593 projection_decomposition =
2594 svg_project_point(point,
2601 <<
static_cast<unsigned int>(
2603 ((projection_decomposition[0] -
2604 x_min_perspective) /
2605 x_dimension_perspective) *
2607 (width / 100.) * 2. * margin_in_percent) +
2608 ((width / 100.) * margin_in_percent))
2610 <<
static_cast<unsigned int>(
2612 (height / 100.) * margin_in_percent -
2613 ((projection_decomposition[1] -
2614 y_min_perspective) /
2615 y_dimension_perspective) *
2617 (height / 100.) * 2. * margin_in_percent))
2623 const double distance_to_camera =
std::sqrt(
2624 std::pow(point[0] - camera_position[0], 2.) +
2625 std::pow(point[1] - camera_position[1], 2.) +
2626 std::pow(point[2] - camera_position[2], 2.));
2627 const double distance_factor =
2628 distance_to_camera /
2629 (2. *
std::max(x_dimension, y_dimension));
2631 const unsigned int font_size_this_edge =
2632 static_cast<unsigned int>(
2633 .5 + .5 * cell_label_font_size *
2635 cell->level() - 4. +
2636 3.5 * distance_factor));
2638 point[0] = cell->face(faceIndex)->center()[0];
2639 point[1] = cell->face(faceIndex)->center()[1];
2645 (
static_cast<float>(cell->level()) /
2646 static_cast<float>(n_levels)) *
2647 std::max(x_dimension, y_dimension);
2650 projection_decomposition =
2651 svg_project_point(point,
2657 const unsigned int xc =
static_cast<unsigned int>(
2659 ((projection_decomposition[0] - x_min_perspective) /
2660 x_dimension_perspective) *
2662 (width / 100.) * 2. * margin_in_percent) +
2663 ((width / 100.) * margin_in_percent));
2664 const unsigned int yc =
static_cast<unsigned int>(
2665 .5 + height - (height / 100.) * margin_in_percent -
2666 ((projection_decomposition[1] - y_min_perspective) /
2667 y_dimension_perspective) *
2669 (height / 100.) * 2. * margin_in_percent));
2671 out <<
" <circle cx=\"" << xc <<
"\" cy=\"" << yc
2672 <<
"\" r=\"" << font_size_this_edge <<
"\" />"
2675 out <<
" <text x=\"" << xc <<
"\" y=\"" << yc
2676 <<
"\" style=\"font-size:" << font_size_this_edge
2677 <<
"px\" dominant-baseline=\"middle\">"
2678 <<
static_cast<int>(
2679 cell->face(faceIndex)->boundary_id())
2680 <<
"</text>" <<
'\n';
2692 out <<
'\n' <<
" <!-- legend -->" <<
'\n';
2694 additional_width = 0;
2696 additional_width =
static_cast<unsigned int>(.5 + (height / 100.) * 2.5);
2704 unsigned int line_offset = 0;
2705 out <<
" <rect x=\"" << width + additional_width <<
"\" y=\""
2706 <<
static_cast<unsigned int>(.5 + (height / 100.) * margin_in_percent)
2708 <<
static_cast<unsigned int>(.5 + (height / 100.) *
2709 (40. - margin_in_percent))
2710 <<
"\" height=\"" <<
static_cast<unsigned int>(.5 + height * .215)
2713 out <<
" <text x=\""
2714 << width + additional_width +
2715 static_cast<unsigned int>(.5 + (height / 100.) * 1.25)
2717 <<
static_cast<unsigned int>(.5 +
2718 (height / 100.) * margin_in_percent +
2719 (++line_offset) * 1.5 * font_size)
2720 <<
"\" style=\"text-anchor:start; font-weight:bold; font-size:"
2721 << font_size <<
"px\">"
2723 <<
"</text>" <<
'\n';
2727 out <<
" <text x=\""
2728 << width + additional_width +
2729 static_cast<unsigned int>(.5 + (height / 100.) * 2.)
2731 <<
static_cast<unsigned int>(.5 +
2732 (height / 100.) * margin_in_percent +
2733 (++line_offset) * 1.5 * font_size)
2734 <<
"\" style=\"text-anchor:start; font-style:oblique; font-size:"
2735 << font_size <<
"px\">"
2743 out <<
"</text>" <<
'\n';
2748 out <<
" <text x=\""
2749 << width + additional_width +
2750 static_cast<unsigned int>(.5 + (height / 100.) * 2.)
2752 <<
static_cast<unsigned int>(.5 +
2753 (height / 100.) * margin_in_percent +
2754 (++line_offset) * 1.5 * font_size)
2755 <<
"\" style=\"text-anchor:start; font-style:oblique; font-size:"
2756 << font_size <<
"px\">"
2763 out <<
"</text>" <<
'\n';
2768 out <<
" <text x=\""
2769 << width + additional_width +
2770 static_cast<unsigned int>(.5 + (height / 100.) * 2.)
2772 <<
static_cast<unsigned int>(.5 +
2773 (height / 100.) * margin_in_percent +
2774 (++line_offset) * 1.5 * font_size)
2775 <<
"\" style=\"text-anchor:start; font-style:oblique; font-size:"
2776 << font_size <<
"px\">"
2783 out <<
"</text>" <<
'\n';
2788 out <<
" <text x= \""
2789 << width + additional_width +
2790 static_cast<unsigned int>(.5 + (height / 100.) * 2.)
2792 <<
static_cast<unsigned int>(.5 +
2793 (height / 100.) * margin_in_percent +
2794 (++line_offset) * 1.5 * font_size)
2795 <<
"\" style=\"text-anchor:start; font-style:oblique; font-size:"
2796 << font_size <<
"px\">"
2802 out <<
"</text>" <<
'\n';
2807 out <<
" <text x= \""
2808 << width + additional_width +
2809 static_cast<unsigned int>(.5 + (height / 100.) * 2.)
2811 <<
static_cast<unsigned int>(.5 +
2812 (height / 100.) * margin_in_percent +
2813 (++line_offset) * 1.5 * font_size)
2814 <<
"\" style=\"text-anchor:start; font-style:oblique; font-size:"
2815 << font_size <<
"px\">"
2816 <<
"level_subdomain_id"
2817 <<
"</text>" <<
'\n';
2822 out <<
" <text x=\""
2823 << width + additional_width +
2824 static_cast<unsigned int>(.5 + (height / 100.) * 1.25)
2826 <<
static_cast<unsigned int>(.5 +
2827 (height / 100.) * margin_in_percent +
2828 (++line_offset) * 1.5 * font_size)
2829 <<
"\" style=\"text-anchor:start; font-weight:bold; font-size:"
2830 << font_size <<
"px\">"
2832 <<
"</text>" <<
'\n';
2834 out <<
" <text x= \""
2835 << width + additional_width +
2836 static_cast<unsigned int>(.5 + (height / 100.) * 2.)
2838 <<
static_cast<unsigned int>(.5 +
2839 (height / 100.) * margin_in_percent +
2840 (++line_offset) * 1.5 * font_size)
2841 <<
"\" style=\"text-anchor:start; font-style:oblique; font-size:"
2842 << font_size <<
"px\">"
2844 <<
"</text>" <<
'\n';
2852 out <<
" <text x=\"" << width + additional_width <<
"\" y=\""
2853 <<
static_cast<unsigned int>(
2854 .5 + (height / 100.) * margin_in_percent + 13.75 * font_size)
2855 <<
"\" style=\"text-anchor:start; font-size:" << font_size <<
"px\">"
2864 out <<
'\n' <<
" <!-- colorbar -->" <<
'\n';
2866 out <<
" <text x=\"" << width + additional_width <<
"\" y=\""
2867 <<
static_cast<unsigned int>(
2868 .5 + (height / 100.) * (margin_in_percent + 29.) -
2870 <<
"\" style=\"text-anchor:start; font-weight:bold; font-size:"
2871 << font_size <<
"px\">";
2876 out <<
"material_id";
2879 out <<
"level_number";
2882 out <<
"subdomain_id";
2885 out <<
"level_subdomain_id";
2891 out <<
"</text>" <<
'\n';
2893 unsigned int element_height =
static_cast<unsigned int>(
2894 ((height / 100.) * (71. - 2. * margin_in_percent)) / n);
2895 unsigned int element_width =
2896 static_cast<unsigned int>(.5 + (height / 100.) * 2.5);
2898 int labeling_index = 0;
2899 auto materials_it = materials.begin();
2900 auto levels_it = levels.begin();
2901 auto subdomains_it = subdomains.begin();
2902 auto level_subdomains_it = level_subdomains.begin();
2904 for (
unsigned int index = 0; index < n; ++index)
2909 labeling_index = *materials_it++;
2912 labeling_index = *levels_it++;
2915 labeling_index = *subdomains_it++;
2918 labeling_index = *level_subdomains_it++;
2924 out <<
" <rect class=\"r" << labeling_index <<
"\" x=\""
2925 << width + additional_width <<
"\" y=\""
2926 <<
static_cast<unsigned int>(.5 + (height / 100.) *
2927 (margin_in_percent + 29)) +
2928 (n - index - 1) * element_height
2929 <<
"\" width=\"" << element_width <<
"\" height=\""
2930 << element_height <<
"\"/>" <<
'\n';
2932 out <<
" <text x=\""
2933 << width + additional_width + 1.5 * element_width <<
"\" y=\""
2934 <<
static_cast<unsigned int>(.5 + (height / 100.) *
2935 (margin_in_percent + 29)) +
2936 (n - index - 1 + .5) * element_height +
2937 static_cast<unsigned int>(.5 + font_size * .35)
2939 <<
" style=\"text-anchor:start; font-size:"
2940 <<
static_cast<unsigned int>(.5 + font_size) <<
"px";
2942 if (index == 0 || index == n - 1)
2943 out <<
"; font-weight:bold";
2945 out <<
"\">" << labeling_index;
2952 out <<
"</text>" <<
'\n';
2960 out <<
'\n' <<
"</svg>";
2976template <
int dim,
int spacedim>
2979 std::ostream & out)
const
2986 const std::time_t time1 = std::time(
nullptr);
2987 const std::tm * time = std::localtime(&time1);
2991 <<
"\n# This file was generated by the deal.II library."
2992 <<
"\n# Date = " << time->tm_year + 1900 <<
"/" << std::setfill(
'0')
2993 << std::setw(2) << time->tm_mon + 1 <<
"/" << std::setfill(
'0')
2994 << std::setw(2) << time->tm_mday <<
"\n# Time = " << std::setfill(
'0')
2995 << std::setw(2) << time->tm_hour <<
":" << std::setfill(
'0')
2996 << std::setw(2) << time->tm_min <<
":" << std::setfill(
'0')
2997 << std::setw(2) << time->tm_sec <<
"\n#"
2998 <<
"\n# For a description of the MathGL script format see the MathGL manual. "
3000 <<
"\n# Note: This file is understood by MathGL v2.1 and higher only, and can "
3001 <<
"\n# be quickly viewed in a graphical environment using \'mglview\'. "
3007 const std::string axes =
"xyz";
3022 out <<
"\nsetsize 800 800";
3023 out <<
"\nrotate 0 0";
3026 out <<
"\nsetsize 800 800";
3027 out <<
"\nrotate 60 40";
3036 <<
"\n# Vertex ordering."
3037 <<
"\n# list <vertex order> <vertex indices>"
3046 out <<
"\nlist f 0 1 2 3" <<
'\n';
3050 <<
"\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"
3059 <<
"\n# List of vertices."
3060 <<
"\n# list <id> <vertices>"
3068 for (
unsigned int i = 0; i < dim; ++i)
3075 out <<
"\nlist " << axes[i] << cell->active_cell_index() <<
" ";
3077 out << cell->vertex(j)[i] <<
" ";
3084 <<
"\n# List of cells to quadplot."
3085 <<
"\n# quadplot <vertex order> <id> <style>"
3089 out <<
"\nquadplot f ";
3090 for (
unsigned int j = 0; j < dim; ++j)
3091 out << axes[j] << i <<
" ";
3116 template <
int dim,
int spacedim,
typename ITERATOR,
typename END>
3118 generate_triangulation_patches(
3124 for (; cell != end; ++cell)
3129 patch.
data.reinit(5, cell->n_vertices());
3131 for (
const unsigned int v : cell->vertex_indices())
3133 patch.
vertices[v] = cell->vertex(v);
3134 patch.
data(0, v) = cell->level();
3136 static_cast<std::make_signed<types::manifold_id>::type
>(
3137 cell->manifold_id());
3139 static_cast<std::make_signed<types::material_id>::type
>(
3140 cell->material_id());
3141 if (cell->is_active())
3143 static_cast<std::make_signed<types::subdomain_id>::type
>(
3144 cell->subdomain_id());
3146 patch.
data(3, v) = -1;
3148 static_cast<std::make_signed<types::subdomain_id>::type
>(
3149 cell->level_subdomain_id());
3151 patches.push_back(patch);
3157 std::vector<std::string>
3158 triangulation_patch_data_names()
3160 std::vector<std::string> v(5);
3165 v[4] =
"level_subdomain";
3172 std::vector<typename Triangulation<3, 3>::active_line_iterator>
3175 std::vector<typename Triangulation<3, 3>::active_line_iterator> res;
3177 std::vector<bool> flags;
3182 for (
const auto l : face->line_indices())
3184 const auto line = face->line(l);
3185 if (line->user_flag_set() || line->has_children())
3188 line->set_user_flag();
3189 if (line->at_boundary())
3190 res.emplace_back(line);
3201 template <
int dim,
int spacedim>
3202 std::vector<typename Triangulation<dim, spacedim>::active_line_iterator>
3214 std::vector<typename Triangulation<3, 3>::active_line_iterator>
3217 std::vector<typename Triangulation<3, 3>::active_line_iterator> res;
3219 std::vector<bool> flags;
3224 for (
const auto l : face->line_indices())
3226 const auto line = face->line(l);
3227 if (line->user_flag_set() || line->has_children())
3230 line->set_user_flag();
3232 (line->boundary_id() != 0 &&
3234 res.emplace_back(line);
3244 template <
int dim,
int spacedim>
3245 std::vector<typename Triangulation<dim, spacedim>::active_line_iterator>
3256 template <
int dim,
int spacedim>
3257 std::vector<typename Triangulation<dim, spacedim>::active_face_iterator>
3260 std::vector<typename Triangulation<dim, spacedim>::active_face_iterator>
3267 res.push_back(face);
3278 template <
int dim,
int spacedim>
3279 std::vector<typename Triangulation<dim, spacedim>::active_face_iterator>
3282 std::vector<typename Triangulation<dim, spacedim>::active_face_iterator>
3289 (face->boundary_id() != 0 &&
3291 res.push_back(face);
3299template <
int dim,
int spacedim>
3302 std::ostream & out)
const
3309 const auto n_vertices =
vertices.size();
3311 out <<
"# vtk DataFile Version 3.0\n"
3312 <<
"Triangulation generated with deal.II\n"
3314 <<
"DATASET UNSTRUCTURED_GRID\n"
3315 <<
"POINTS " << n_vertices <<
" double\n";
3321 for (
unsigned int d = spacedim + 1; d <= 3; ++d)
3327 get_relevant_face_iterators(
tria) :
3328 get_boundary_face_iterators(
tria);
3330 get_relevant_edge_iterators(
tria) :
3331 get_boundary_edge_iterators(
tria);
3337 "At least one of the flags (output_cells, output_faces, output_edges) has to be enabled!"));
3354 cells_size += cell->n_vertices() + 1;
3357 for (
const auto &face : faces)
3358 cells_size += face->n_vertices() + 1;
3361 for (
const auto &edge : edges)
3362 cells_size += edge->n_vertices() + 1;
3366 out <<
"\nCELLS " << n_cells <<
' ' << cells_size <<
'\n';
3381 static const std::array<int, 8> deal_to_vtk_cell_type = {
3382 {1, 3, 5, 9, 10, 14, 13, 12}};
3383 static const std::array<unsigned int, 8> vtk_to_deal_hypercube = {
3384 {0, 1, 3, 2, 4, 5, 7, 6}};
3390 out << cell->n_vertices();
3391 for (
const unsigned int i : cell->vertex_indices())
3394 const auto reference_cell = cell->reference_cell();
3400 out << cell->vertex_index(vtk_to_deal_hypercube[i]);
3404 out << cell->vertex_index(i);
3407 static const std::array<unsigned int, 5> permutation_table{
3409 out << cell->vertex_index(permutation_table[i]);
3417 for (
const auto &face : faces)
3419 out << face->n_vertices();
3420 for (
const unsigned int i : face->vertex_indices())
3424 face->n_vertices() ?
3425 vtk_to_deal_hypercube[i] :
3431 for (
const auto &edge : edges)
3434 for (
const unsigned int i : edge->vertex_indices())
3435 out <<
' ' << edge->vertex_index(i);
3440 out <<
"\nCELL_TYPES " << n_cells <<
'\n';
3444 out << deal_to_vtk_cell_type[static_cast<int>(cell->reference_cell())]
3450 for (
const auto &face : faces)
3451 out << deal_to_vtk_cell_type[static_cast<int>(face->reference_cell())]
3457 for (
const auto &edge : edges)
3458 out << deal_to_vtk_cell_type[static_cast<int>(edge->reference_cell())]
3461 out <<
"\n\nCELL_DATA " << n_cells <<
'\n'
3462 <<
"SCALARS MaterialID int 1\n"
3463 <<
"LOOKUP_TABLE default\n";
3470 out << static_cast<std::make_signed<types::material_id>::type>(
3471 cell->material_id())
3478 for (
const auto &face : faces)
3480 out << static_cast<std::make_signed<types::boundary_id>::type>(
3481 face->boundary_id())
3488 for (
const auto &edge : edges)
3490 out << static_cast<std::make_signed<types::boundary_id>::type>(
3491 edge->boundary_id())
3496 out <<
"\n\nSCALARS ManifoldID int 1\n"
3497 <<
"LOOKUP_TABLE default\n";
3504 out << static_cast<std::make_signed<types::manifold_id>::type>(
3505 cell->manifold_id())
3512 for (
const auto &face : faces)
3514 out << static_cast<std::make_signed<types::manifold_id>::type>(
3515 face->manifold_id())
3522 for (
const auto &edge : edges)
3524 out << static_cast<std::make_signed<types::manifold_id>::type>(
3525 edge->manifold_id())
3538template <
int dim,
int spacedim>
3541 std::ostream & out)
const
3549 std::vector<DataOutBase::Patch<dim, spacedim>> patches;
3556 triangulation_patch_data_names(),
3558 std::tuple<
unsigned int,
3566 out <<
" </UnstructuredGrid>\n";
3567 out <<
"<dealiiData encoding=\"base64\">";
3568 std::stringstream outstring;
3569 boost::archive::binary_oarchive ia(outstring);
3573 out <<
"\n</dealiiData>\n";
3574 out <<
"</VTKFile>\n";
3585template <
int dim,
int spacedim>
3589 const std::string & filename_without_extension,
3590 const bool view_levels,
3591 const bool include_artificial)
const
3593 std::vector<DataOutBase::Patch<dim, spacedim>> patches;
3594 const unsigned int n_datasets = 4;
3595 std::vector<std::string> data_names;
3596 data_names.emplace_back(
"level");
3597 data_names.emplace_back(
"subdomain");
3598 data_names.emplace_back(
"level_subdomain");
3599 data_names.emplace_back(
"proc_writing");
3605 const auto &reference_cell = reference_cells[0];
3607 const unsigned int n_q_points = reference_cell.n_vertices();
3613 if (cell->has_children())
3615 if (!include_artificial &&
3619 else if (!include_artificial)
3621 if (cell->has_children() &&
3624 else if (cell->is_active() &&
3625 cell->level_subdomain_id() ==
3632 patch.
data.reinit(n_datasets, n_q_points);
3636 for (
unsigned int vertex = 0; vertex < n_q_points; ++vertex)
3638 patch.
vertices[vertex] = cell->vertex(vertex);
3639 patch.
data(0, vertex) = cell->level();
3640 if (cell->is_active())
3641 patch.
data(1, vertex) =
static_cast<double>(
3642 static_cast<std::make_signed<types::subdomain_id>::type
>(
3643 cell->subdomain_id()));
3645 patch.
data(1, vertex) = -1.0;
3646 patch.
data(2, vertex) =
static_cast<double>(
3647 static_cast<std::make_signed<types::subdomain_id>::type
>(
3648 cell->level_subdomain_id()));
3652 for (
auto f : reference_cell.face_indices())
3654 patches.push_back(patch);
3660 std::string new_file = filename_without_extension +
".vtu";
3664 new_file = filename_without_extension +
".proc" +
3669 if (tr->locally_owned_subdomain() == 0)
3671 std::vector<std::string> filenames;
3676 std::size_t pos = filename_without_extension.find_last_of(
'/');
3677 if (pos == std::string::npos)
3681 const unsigned int n_procs =
3683 for (
unsigned int i = 0; i < n_procs; ++i)
3684 filenames.push_back(filename_without_extension.substr(pos) +
3688 const std::string pvtu_filename =
3689 (filename_without_extension +
".pvtu");
3690 std::ofstream pvtu_output(pvtu_filename.c_str());
3709 std::ofstream out(new_file.c_str());
3711 std::tuple<
unsigned int,
3773template <
int dim,
int spacedim>
3778 unsigned int n_faces = 0;
3781 if ((face->at_boundary()) && (face->boundary_id() != 0))
3789template <
int dim,
int spacedim>
3796 std::vector<bool> line_flags;
3800 .clear_user_flags_line();
3802 unsigned int n_lines = 0;
3805 for (
unsigned int l = 0; l < GeometryInfo<dim>::lines_per_cell; ++l)
3806 if (cell->line(l)->at_boundary() && (cell->line(l)->boundary_id() != 0) &&
3807 (cell->line(l)->user_flag_set() ==
false))
3810 cell->line(l)->set_user_flag();
3825 const unsigned int next_element_index,
3826 std::ostream &)
const
3828 return next_element_index;
3834 const unsigned int next_element_index,
3835 std::ostream &)
const
3837 return next_element_index;
3842 const unsigned int next_element_index,
3843 std::ostream &)
const
3845 return next_element_index;
3851 const unsigned int next_element_index,
3852 std::ostream &)
const
3854 return next_element_index;
3859 const unsigned int next_element_index,
3860 std::ostream &)
const
3862 return next_element_index;
3868 const unsigned int next_element_index,
3869 std::ostream &)
const
3871 return next_element_index;
3877 const unsigned int next_element_index,
3878 std::ostream &)
const
3880 return next_element_index;
3885 const unsigned int next_element_index,
3886 std::ostream &)
const
3888 return next_element_index;
3893template <
int dim,
int spacedim>
3896 const unsigned int next_element_index,
3897 std::ostream & out)
const
3899 unsigned int current_element_index = next_element_index;
3902 if (face->at_boundary() && (face->boundary_id() != 0))
3904 out << current_element_index <<
' '
3905 << face->reference_cell().gmsh_element_type() <<
' ';
3906 out << static_cast<unsigned int>(face->boundary_id()) <<
' '
3907 <<
static_cast<unsigned int>(face->boundary_id()) <<
' '
3908 << face->n_vertices();
3910 for (
unsigned int vertex : face->vertex_indices())
3914 << face->vertex_index(
3919 out <<
' ' << face->vertex_index(vertex) + 1;
3925 ++current_element_index;
3927 return current_element_index;
3932template <
int dim,
int spacedim>
3935 const unsigned int next_element_index,
3936 std::ostream & out)
const
3938 unsigned int current_element_index = next_element_index;
3943 std::vector<bool> line_flags;
3947 .clear_user_flags_line();
3950 for (
unsigned int l = 0; l < GeometryInfo<dim>::lines_per_cell; ++l)
3951 if (cell->line(l)->at_boundary() && (cell->line(l)->boundary_id() != 0) &&
3952 (cell->line(l)->user_flag_set() ==
false))
3954 out << next_element_index <<
' '
3956 out << static_cast<unsigned int>(cell->line(l)->boundary_id()) <<
' '
3957 <<
static_cast<unsigned int>(cell->line(l)->boundary_id())
3960 for (
unsigned int vertex = 0; vertex < 2; ++vertex)
3962 << cell->line(l)->vertex_index(
3970 ++current_element_index;
3971 cell->line(l)->set_user_flag();
3979 return current_element_index;
3986 const unsigned int next_element_index,
3987 std::ostream &)
const
3989 return next_element_index;
3994 const unsigned int next_element_index,
3995 std::ostream &)
const
3997 return next_element_index;
4002 const unsigned int next_element_index,
4003 std::ostream &)
const
4005 return next_element_index;
4010 const unsigned int next_element_index,
4011 std::ostream &)
const
4013 return next_element_index;
4018 const unsigned int next_element_index,
4019 std::ostream &)
const
4021 return next_element_index;
4027 const unsigned int next_element_index,
4028 std::ostream &)
const
4030 return next_element_index;
4036 const unsigned int next_element_index,
4037 std::ostream &)
const
4039 return next_element_index;
4044 const unsigned int next_element_index,
4045 std::ostream &)
const
4047 return next_element_index;
4052template <
int dim,
int spacedim>
4055 const unsigned int next_element_index,
4056 std::ostream & out)
const
4058 unsigned int current_element_index = next_element_index;
4062 if (face->at_boundary() && (face->boundary_id() != 0))
4064 out << current_element_index <<
" "
4065 <<
static_cast<unsigned int>(face->boundary_id()) <<
" ";
4078 for (
unsigned int vertex = 0;
4079 vertex < GeometryInfo<dim>::vertices_per_face;
4081 out << face->vertex_index(
4087 ++current_element_index;
4089 return current_element_index;
4094template <
int dim,
int spacedim>
4097 const unsigned int next_element_index,
4098 std::ostream & out)
const
4100 unsigned int current_element_index = next_element_index;
4105 std::vector<bool> line_flags;
4109 .clear_user_flags_line();
4112 for (
unsigned int l = 0; l < GeometryInfo<dim>::lines_per_cell; ++l)
4113 if (cell->line(l)->at_boundary() && (cell->line(l)->boundary_id() != 0) &&
4114 (cell->line(l)->user_flag_set() ==
false))
4116 out << current_element_index <<
" "
4117 <<
static_cast<unsigned int>(cell->line(l)->boundary_id())
4120 for (
unsigned int vertex = 0; vertex < 2; ++vertex)
4121 out << cell->line(l)->vertex_index(
4130 ++current_element_index;
4131 cell->line(l)->set_user_flag();
4138 return current_element_index;
4155 template <
int spacedim>
4159 while (points.size() > 2)
4162 first_difference /= first_difference.
norm();
4164 second_difference /= second_difference.
norm();
4166 if ((first_difference - second_difference).norm() < 1e-10)
4167 points.erase(points.begin() + 1);
4175 template <
int spacedim>
4187 out <<
"# cell " << cell <<
'\n';
4189 out << cell->vertex(0) <<
' ' << cell->level() <<
' '
4190 << cell->material_id() <<
'\n'
4191 << cell->vertex(1) <<
' ' << cell->level() <<
' '
4192 << cell->material_id() <<
'\n'
4205 template <
int spacedim>
4216 const unsigned int n_additional_points =
4218 const unsigned int n_points = 2 + n_additional_points;
4223 std::vector<
Point<dim - 1>> boundary_points;
4224 if (mapping !=
nullptr)
4226 boundary_points.resize(n_points);
4227 boundary_points[0][0] = 0;
4228 boundary_points[n_points - 1][0] = 1;
4229 for (
unsigned int i = 1; i < n_points - 1; ++i)
4230 boundary_points[i](0) = 1. * i / (n_points - 1);
4232 std::vector<double> dummy_weights(n_points, 1. / n_points);
4233 Quadrature<dim - 1> quadrature(boundary_points, dummy_weights);
4242 out <<
"# cell " << cell <<
'\n';
4244 if (mapping ==
nullptr ||
4256 << cell->level() <<
' ' << cell->material_id() <<
'\n';
4257 out << cell->vertex(0) <<
' ' << cell->level() <<
' '
4258 << cell->material_id() <<
'\n'
4266 for (
const unsigned int face_no :
4269 const typename ::Triangulation<dim,
4270 spacedim>::face_iterator
4271 face = cell->face(face_no);
4272 if (dim != spacedim || face->at_boundary() ||
4278 std::vector<Point<spacedim>> line_points;
4281 const unsigned int offset = face_no * n_points;
4282 for (
unsigned int i = 0; i < n_points; ++i)
4283 line_points.push_back(
4285 cell, q_projector.
point(offset + i)));
4286 internal::remove_colinear_points(line_points);
4289 out <<
point <<
' ' << cell->level() <<
' '
4290 << cell->material_id() <<
'\n';
4292 out <<
'\n' <<
'\n';
4298 out << face->vertex(0) <<
' ' << cell->level() <<
' '
4299 << cell->material_id() <<
'\n'
4300 << face->vertex(1) <<
' ' << cell->level() <<
' '
4301 << cell->material_id() <<
'\n'
4317 template <
int spacedim>
4328 const unsigned int n_additional_points =
4330 const unsigned int n_points = 2 + n_additional_points;
4334 std::unique_ptr<Quadrature<dim>> q_projector;
4335 std::vector<Point<1>> boundary_points;
4336 if (mapping !=
nullptr)
4338 boundary_points.resize(n_points);
4339 boundary_points[0][0] = 0;
4340 boundary_points[n_points - 1][0] = 1;
4341 for (
unsigned int i = 1; i < n_points - 1; ++i)
4342 boundary_points[i](0) = 1. * i / (n_points - 1);
4344 std::vector<double> dummy_weights(n_points, 1. / n_points);
4348 QIterated<dim - 1> quadrature(quadrature1d, 1);
4349 q_projector = std::make_unique<Quadrature<dim>>(
4356 out <<
"# cell " << cell <<
'\n';
4358 if (mapping ==
nullptr || n_points == 2 ||
4359 (!cell->has_boundary_lines() &&
4365 out << cell->vertex(0) <<
' ' << cell->level() <<
' '
4366 << cell->material_id() <<
'\n'
4367 << cell->vertex(1) <<
' ' << cell->level() <<
' '
4368 << cell->material_id() <<
'\n'
4369 << cell->vertex(5) <<
' ' << cell->level() <<
' '
4370 << cell->material_id() <<
'\n'
4371 << cell->vertex(4) <<
' ' << cell->level() <<
' '
4372 << cell->material_id() <<
'\n'
4373 << cell->vertex(0) <<
' ' << cell->level() <<
' '
4374 << cell->material_id() <<
'\n'
4377 out << cell->vertex(2) <<
' ' << cell->level() <<
' '
4378 << cell->material_id() <<
'\n'
4379 << cell->vertex(3) <<
' ' << cell->level() <<
' '
4380 << cell->material_id() <<
'\n'
4381 << cell->vertex(7) <<
' ' << cell->level() <<
' '
4382 << cell->material_id() <<
'\n'
4383 << cell->vertex(6) <<
' ' << cell->level() <<
' '
4384 << cell->material_id() <<
'\n'
4385 << cell->vertex(2) <<
' ' << cell->level() <<
' '
4386 << cell->material_id() <<
'\n'
4390 out << cell->vertex(0) <<
' ' << cell->level() <<
' '
4391 << cell->material_id() <<
'\n'
4392 << cell->vertex(2) <<
' ' << cell->level() <<
' '
4393 << cell->material_id() <<
'\n'
4395 out << cell->vertex(1) <<
' ' << cell->level() <<
' '
4396 << cell->material_id() <<
'\n'
4397 << cell->vertex(3) <<
' ' << cell->level() <<
' '
4398 << cell->material_id() <<
'\n'
4400 out << cell->vertex(5) <<
' ' << cell->level() <<
' '
4401 << cell->material_id() <<
'\n'
4402 << cell->vertex(7) <<
' ' << cell->level() <<
' '
4403 << cell->material_id() <<
'\n'
4405 out << cell->vertex(4) <<
' ' << cell->level() <<
' '
4406 << cell->material_id() <<
'\n'
4407 << cell->vertex(6) <<
' ' << cell->level() <<
' '
4408 << cell->material_id() <<
'\n'
4414 for (
const unsigned int v : {0, 1, 2, 0, 3, 2})
4416 out << cell->vertex(v) <<
' ' << cell->level() <<
' '
4417 << cell->material_id() <<
'\n';
4421 for (
const unsigned int v : {3, 1})
4423 out << cell->vertex(v) <<
' ' << cell->level() <<
' '
4424 << cell->material_id() <<
'\n';
4435 for (
const unsigned int v : {0, 1, 2, 0, 3, 4, 5, 3})
4437 out << cell->vertex(v) <<
' ' << cell->level() <<
' '
4438 << cell->material_id() <<
'\n';
4442 for (
const unsigned int v : {1, 4})
4444 out << cell->vertex(v) <<
' ' << cell->level() <<
' '
4445 << cell->material_id() <<
'\n';
4449 for (
const unsigned int v : {2, 5})
4451 out << cell->vertex(v) <<
' ' << cell->level() <<
' '
4452 << cell->material_id() <<
'\n';
4459 for (
const unsigned int v : {0, 1, 3, 2, 0, 4, 1})
4461 out << cell->vertex(v) <<
' ' << cell->level() <<
' '
4462 << cell->material_id() <<
'\n';
4466 for (
const unsigned int v : {2, 4, 3})
4468 out << cell->vertex(v) <<
' ' << cell->level() <<
' '
4469 << cell->material_id() <<
'\n';
4480 for (
const unsigned int face_no :
4483 const typename ::Triangulation<dim,
4484 spacedim>::face_iterator
4485 face = cell->face(face_no);
4487 if (face->at_boundary() &&
4490 const unsigned int offset = face_no * n_points * n_points;
4491 for (
unsigned int i = 0; i < n_points - 1; ++i)
4492 for (
unsigned int j = 0; j < n_points - 1; ++j)
4497 q_projector->point(offset + i * n_points + j));
4498 out << p0 <<
' ' << cell->level() <<
' '
4499 << cell->material_id() <<
'\n';
4503 offset + (i + 1) * n_points + j)))
4504 <<
' ' << cell->level() <<
' '
4505 << cell->material_id() <<
'\n';
4509 offset + (i + 1) * n_points + j + 1)))
4510 <<
' ' << cell->level() <<
' '
4511 << cell->material_id() <<
'\n';
4514 q_projector->point(offset + i * n_points +
4516 <<
' ' << cell->level() <<
' '
4517 << cell->material_id() <<
'\n';
4519 out << p0 <<
' ' << cell->level() <<
' '
4520 << cell->material_id() <<
'\n';
4521 out <<
'\n' <<
'\n';
4526 for (
unsigned int l = 0;
4527 l < GeometryInfo<dim>::lines_per_face;
4530 const typename ::Triangulation<dim, spacedim>::
4531 line_iterator line = face->line(l);
4534 &
v1 = line->vertex(1);
4535 if (line->at_boundary() ||
4541 std::vector<Point<spacedim>> line_points;
4550 for (
unsigned int i = 0; i < n_points; ++i)
4551 line_points.push_back(
4554 (1 - boundary_points[i][0]) * u0 +
4555 boundary_points[i][0] * u1));
4556 internal::remove_colinear_points(line_points);
4558 out <<
point <<
' ' << cell->level() <<
' '
4559 <<
static_cast<unsigned int>(
4560 cell->material_id())
4564 out <<
v0 <<
' ' << cell->level() <<
' '
4565 << cell->material_id() <<
'\n'
4566 <<
v1 <<
' ' << cell->level() <<
' '
4567 << cell->material_id() <<
'\n';
4569 out <<
'\n' <<
'\n';
4586template <
int dim,
int spacedim>
4610 const unsigned int l)
4630 write_eps(const ::Triangulation<1, 2> &,
4640 write_eps(const ::Triangulation<1, 3> &,
4650 write_eps(const ::Triangulation<2, 3> &,
4661 template <
int dim,
int spacedim>
4669 using LineList = std::list<LineEntry>;
4710 for (
const unsigned int line_no : cell->line_indices())
4712 typename ::Triangulation<dim, spacedim>::line_iterator
4713 line = cell->line(line_no);
4725 if (!line->has_children() &&
4726 (mapping ==
nullptr || !line->at_boundary()))
4746 line_list.emplace_back(
4747 Point<2>(line->vertex(0)(0), line->vertex(0)(1)),
4748 Point<2>(line->vertex(1)(0), line->vertex(1)(1)),
4749 line->user_flag_set(),
4759 if (mapping !=
nullptr)
4766 std::vector<
Point<dim - 1>> boundary_points(n_points);
4768 for (
unsigned int i = 0; i < n_points; ++i)
4769 boundary_points[i](0) = 1. * (i + 1) / (n_points + 1);
4771 Quadrature<dim - 1> quadrature(boundary_points);
4780 for (
const unsigned int face_no :
4783 const typename ::Triangulation<dim, spacedim>::
4784 face_iterator face = cell->face(face_no);
4786 if (face->at_boundary())
4796 const unsigned int offset = face_no * n_points;
4797 for (
unsigned int i = 0; i < n_points; ++i)
4801 cell, q_projector.point(offset + i)));
4802 const Point<2> p1(p1_dim(0), p1_dim(1));
4804 line_list.emplace_back(p0,
4806 face->user_flag_set(),
4813 const Point<2> p1(p1_dim(0), p1_dim(1));
4814 line_list.emplace_back(p0,
4816 face->user_flag_set(),
4850 const double z_angle = eps_flags_3.azimut_angle;
4851 const double turn_angle = eps_flags_3.turn_angle;
4853 -
std::sin(z_angle * 2. * pi / 360.) *
4854 std::sin(turn_angle * 2. * pi / 360.),
4855 +
std::sin(z_angle * 2. * pi / 360.) *
4856 std::cos(turn_angle * 2. * pi / 360.),
4857 -
std::cos(z_angle * 2. * pi / 360.));
4865 ((
Point<dim>(0, 0, 1) * view_direction) * view_direction);
4874 ((
Point<dim>(1, 0, 0) * view_direction) * view_direction) -
4875 ((
Point<dim>(1, 0, 0) * unit_vector1) * unit_vector1));
4880 for (
const unsigned int line_no : cell->line_indices())
4882 typename ::Triangulation<dim, spacedim>::line_iterator
4883 line = cell->line(line_no);
4884 line_list.emplace_back(
4885 Point<2>(line->vertex(0) * unit_vector2,
4886 line->vertex(0) * unit_vector1),
4887 Point<2>(line->vertex(1) * unit_vector2,
4888 line->vertex(1) * unit_vector1),
4889 line->user_flag_set(),
4906 double x_max = x_min;
4908 double y_max = y_min;
4909 unsigned int max_level = line_list.begin()->level;
4911 for (LineList::const_iterator line = line_list.begin();
4912 line != line_list.end();
4915 x_min =
std::min(x_min, line->first(0));
4916 x_min =
std::min(x_min, line->second(0));
4918 x_max =
std::max(x_max, line->first(0));
4919 x_max =
std::max(x_max, line->second(0));
4921 y_min =
std::min(y_min, line->first(1));
4922 y_min =
std::min(y_min, line->second(1));
4924 y_max =
std::max(y_max, line->first(1));
4925 y_max =
std::max(y_max, line->second(1));
4927 max_level =
std::max(max_level, line->level);
4935 const double scale =
4936 (eps_flags_base.
size /
4947 std::time_t time1 = std::time(
nullptr);
4948 std::tm * time = std::localtime(&time1);
4949 out <<
"%!PS-Adobe-2.0 EPSF-1.2" <<
'\n'
4950 <<
"%%Title: deal.II Output" <<
'\n'
4951 <<
"%%Creator: the deal.II library" <<
'\n'
4952 <<
"%%Creation Date: " << time->tm_year + 1900 <<
"/"
4953 << time->tm_mon + 1 <<
"/" << time->tm_mday <<
" - "
4954 << time->tm_hour <<
":" << std::setw(2) << time->tm_min <<
":"
4955 << std::setw(2) << time->tm_sec <<
'\n'
4956 <<
"%%BoundingBox: "
4960 <<
static_cast<unsigned int>(
4961 std::floor(((x_max - x_min) * scale) + 1))
4963 <<
static_cast<unsigned int>(
4964 std::floor(((y_max - y_min) * scale) + 1))
4973 out <<
"/m {moveto} bind def" <<
'\n'
4974 <<
"/x {lineto stroke} bind def" <<
'\n'
4975 <<
"/b {0 0 0 setrgbcolor} def" <<
'\n'
4976 <<
"/r {1 0 0 setrgbcolor} def" <<
'\n';
4983 out <<
"/l { neg " << (max_level) <<
" add "
4984 << (0.66666 /
std::max(1U, (max_level - 1)))
4985 <<
" mul 1 0.8 sethsbcolor} def" <<
'\n';
4995 if ((dim == 2) && (eps_flags_2.write_cell_numbers ||
4996 eps_flags_2.write_vertex_numbers))
4999 << (
"/R {rmoveto} bind def\n"
5000 "/Symbol-Oblique /Symbol findfont [1 0 .167 1 0 0] makefont\n"
5001 "dup length dict begin {1 index /FID eq {pop pop} {def} ifelse} forall\n"
5002 "currentdict end definefont\n"
5003 "/MFshow {{dup dup 0 get findfont exch 1 get scalefont setfont\n"
5004 "[ currentpoint ] exch dup 2 get 0 exch rmoveto dup dup 5 get exch 4 get\n"
5005 "{show} {stringwidth pop 0 rmoveto}ifelse dup 3 get\n"
5006 "{2 get neg 0 exch rmoveto pop} {pop aload pop moveto}ifelse} forall} bind def\n"
5007 "/MFwidth {0 exch {dup 3 get{dup dup 0 get findfont exch 1 get scalefont setfont\n"
5008 "5 get stringwidth pop add}\n"
5009 "{pop} ifelse} forall} bind def\n"
5010 "/MCshow { currentpoint stroke m\n"
5011 "exch dup MFwidth -2 div 3 -1 roll R MFshow } def\n")
5015 out <<
"%%EndProlog" <<
'\n' <<
'\n';
5018 out << eps_flags_base.
line_width <<
" setlinewidth" <<
'\n';
5022 const Point<2> offset(x_min, y_min);
5024 for (LineList::const_iterator line = line_list.begin();
5025 line != line_list.end();
5028 out << line->level <<
" l " << (line->first - offset) * scale <<
" m "
5029 << (line->second - offset) *
scale <<
" x" <<
'\n';
5034 << (line->first - offset) * scale <<
" m "
5035 << (line->second - offset) *
scale <<
" x" <<
'\n';
5039 if ((dim == 2) && (eps_flags_2.write_cell_numbers ==
true))
5041 out <<
"(Helvetica) findfont 140 scalefont setfont" <<
'\n';
5045 out << (cell->center()(0) - offset(0)) * scale <<
' '
5046 << (cell->center()(1) - offset(1)) * scale <<
" m" <<
'\n'
5047 <<
"[ [(Helvetica) 12.0 0.0 true true (";
5048 if (eps_flags_2.write_cell_number_level)
5051 out << cell->index();
5054 <<
"] -6 MCshow" <<
'\n';
5059 if ((dim == 2) && (eps_flags_2.write_vertex_numbers ==
true))
5061 out <<
"(Helvetica) findfont 140 scalefont setfont" <<
'\n';
5067 std::set<unsigned int> treated_vertices;
5069 for (
const unsigned int vertex_no : cell->vertex_indices())
5070 if (treated_vertices.find(cell->vertex_index(vertex_no)) ==
5071 treated_vertices.end())
5073 treated_vertices.insert(cell->vertex_index(vertex_no));
5075 out << (cell->vertex(vertex_no)(0) - offset(0)) * scale <<
' '
5076 << (cell->vertex(vertex_no)(1) - offset(1)) * scale
5078 <<
"[ [(Helvetica) 10.0 0.0 true true ("
5079 << cell->vertex_index(vertex_no) <<
")] "
5080 <<
"] -6 MCshow" <<
'\n';
5084 out <<
"showpage" <<
'\n';
5096template <
int dim,
int spacedim>
5106template <
int dim,
int spacedim>
5113 switch (output_format)
5163template <
int dim,
int spacedim>
5174#include "grid_out.inst"
void attach_triangulation(const Triangulation< dim, spacedim > &)
void add_data_vector(const VectorType &data, const std::vector< std::string > &names, const DataVectorType type=type_automatic, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation={})
virtual void build_patches(const unsigned int n_subdivisions=0)
unsigned int n_boundary_faces(const Triangulation< dim, spacedim > &tria) const
GridOutFlags::Vtu vtu_flags
GridOutFlags::Eps< 2 > eps_flags_2
unsigned int write_ucd_lines(const Triangulation< dim, spacedim > &tria, const unsigned int next_element_index, std::ostream &out) const
void parse_parameters(ParameterHandler ¶m)
void write_svg(const Triangulation< 2, 2 > &tria, std::ostream &out) const
void set_flags(const GridOutFlags::DX &flags)
GridOutFlags::Eps< 1 > eps_flags_1
GridOutFlags::XFig xfig_flags
unsigned int write_ucd_faces(const Triangulation< dim, spacedim > &tria, const unsigned int next_element_index, std::ostream &out) const
void write_mathgl(const Triangulation< dim, spacedim > &tria, std::ostream &out) const
std::string default_suffix() const
static void declare_parameters(ParameterHandler ¶m)
static OutputFormat parse_output_format(const std::string &format_name)
void write_vtk(const Triangulation< dim, spacedim > &tria, std::ostream &out) const
GridOutFlags::Gnuplot gnuplot_flags
void write_msh(const Triangulation< dim, spacedim > &tria, std::ostream &out) const
unsigned int write_msh_lines(const Triangulation< dim, spacedim > &tria, const unsigned int next_element_index, std::ostream &out) const
void write_eps(const Triangulation< dim, spacedim > &tria, std::ostream &out, const Mapping< dim, spacedim > *mapping=nullptr) const
static std::string get_output_format_names()
GridOutFlags::Eps< 3 > eps_flags_3
void write(const Triangulation< dim, spacedim > &tria, std::ostream &out, const OutputFormat output_format, const Mapping< dim, spacedim > *mapping=nullptr) const
void write_vtu(const Triangulation< dim, spacedim > &tria, std::ostream &out) const
GridOutFlags::Ucd ucd_flags
void write_dx(const Triangulation< dim, spacedim > &tria, std::ostream &out) const
std::size_t memory_consumption() const
unsigned int n_boundary_lines(const Triangulation< dim, spacedim > &tria) const
GridOutFlags::DX dx_flags
unsigned int write_msh_faces(const Triangulation< dim, spacedim > &tria, const unsigned int next_element_index, std::ostream &out) const
void write_ucd(const Triangulation< dim, spacedim > &tria, std::ostream &out) const
GridOutFlags::Svg svg_flags
OutputFormat default_format
void write_xfig(const Triangulation< dim, spacedim > &tria, std::ostream &out, const Mapping< dim, spacedim > *mapping=nullptr) const
void write_gnuplot(const Triangulation< dim, spacedim > &tria, std::ostream &out, const Mapping< dim, spacedim > *mapping=nullptr) const
GridOutFlags::Vtk vtk_flags
GridOutFlags::Msh msh_flags
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
@ vtk
write() calls write_vtk()
@ eps
write() calls write_eps()
@ msh
write() calls write_msh()
@ xfig
write() calls write_xfig()
@ dx
write() calls write_dx()
@ ucd
write() calls write_ucd()
@ gnuplot
write() calls write_gnuplot()
@ mathgl
write() calls write_mathgl()
@ svg
write() calls write_svg()
@ none
Do nothing in write()
@ vtu
write() calls write_vtu()
GridOutFlags::MathGL mathgl_flags
Abstract base class for mapping classes.
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const =0
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const =0
long int get_integer(const std::string &entry_string) const
bool get_bool(const std::string &entry_name) const
void declare_entry(const std::string &entry, const std::string &default_value, const Patterns::PatternBase &pattern=Patterns::Anything(), const std::string &documentation="", const bool has_to_be_set=false)
std::string get(const std::string &entry_string) const
double get_double(const std::string &entry_name) const
void enter_subsection(const std::string &subsection)
static Quadrature< dim > project_to_all_faces(const Quadrature< dim - 1 > &quadrature)
const Point< dim > & point(const unsigned int i) const
unsigned int gmsh_element_type() const
numbers::NumberTraits< Number >::real_type norm() const
void save_user_flags_line(std::ostream &out) const
void save(Archive &ar, const unsigned int version) const
cell_iterator begin(const unsigned int level=0) const
virtual types::subdomain_id locally_owned_subdomain() const
unsigned int n_active_cells() const
const std::vector< Point< spacedim > > & get_vertices() const
unsigned int n_used_vertices() const
const std::vector< ReferenceCell > & get_reference_cells() const
cell_iterator end() const
const std::vector< bool > & get_used_vertices() const
active_cell_iterator begin_active(const unsigned int level=0) const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
IteratorRange< active_face_iterator > active_face_iterators() const
IteratorRange< active_cell_iterator > active_cell_iterators() const
IteratorRange< cell_iterator > cell_iterators_on_level(const unsigned int level) const
IteratorRange< cell_iterator > cell_iterators() const
void write_pvtu_record(std::ostream &out, const std::vector< std::string > &piece_names) const
static ::ExceptionBase & ExcIO()
ReferenceCell reference_cell
static void declare_parameters(ParameterHandler &prm)
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcInvalidState()
#define Assert(cond, exc)
Point< spacedim > vertices[GeometryInfo< dim >::vertices_per_cell]
#define AssertDimension(dim1, dim2)
void parse_parameters(const ParameterHandler &prm)
static ::ExceptionBase & ExcInternalError()
unsigned int n_subdivisions
bool points_are_available
std::array< unsigned int, GeometryInfo< dim >::faces_per_cell > neighbors
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
DataComponentInterpretation
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, DataComponentInterpretation::DataComponentInterpretation > > &nonscalar_data_ranges, const EpsFlags &flags, std::ostream &out)
void write_vtu_header(std::ostream &out, const VtkFlags &flags)
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, DataComponentInterpretation::DataComponentInterpretation > > &nonscalar_data_ranges, const VtkFlags &flags, std::ostream &out)
void write_vtu_main(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, DataComponentInterpretation::DataComponentInterpretation > > &nonscalar_data_ranges, const VtkFlags &flags, std::ostream &out)
void write_vtu_footer(std::ostream &out)
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, DataComponentInterpretation::DataComponentInterpretation > > &nonscalar_data_ranges, const GnuplotFlags &flags, std::ostream &out)
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
constexpr const ReferenceCell Tetrahedron
constexpr const ReferenceCell Quadrilateral
constexpr const ReferenceCell Wedge
constexpr const ReferenceCell Pyramid
constexpr const ReferenceCell Triangle
constexpr const ReferenceCell Hexahedron
constexpr const ReferenceCell Vertex
constexpr const ReferenceCell Line
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
std::string encode_base64(const std::vector< unsigned char > &binary_input)
std::string compress(const std::string &input)
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
const types::boundary_id invalid_boundary_id
static constexpr double PI
const types::boundary_id internal_face_boundary_id
const types::subdomain_id artificial_subdomain_id
static const unsigned int invalid_unsigned_int
const types::manifold_id flat_manifold_id
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > cos(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sin(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
static void declare_parameters(ParameterHandler ¶m)
void parse_parameters(ParameterHandler ¶m)
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)
static void declare_parameters(ParameterHandler ¶m)
bool color_lines_on_user_flag
unsigned int n_boundary_face_points
void parse_parameters(ParameterHandler ¶m)
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)
void parse_parameters(ParameterHandler ¶m)
void parse_parameters(ParameterHandler ¶m)
void parse_parameters(ParameterHandler ¶m)
unsigned int n_extra_curved_line_points
void parse_parameters(ParameterHandler ¶m)
Gnuplot(const bool write_cell_number=false, const unsigned int n_extra_curved_line_points=2, const bool curved_inner_cells=false, const bool write_additional_boundary_lines=true)
static void declare_parameters(ParameterHandler ¶m)
bool write_additional_boundary_lines
void parse_parameters(ParameterHandler ¶m)
static void declare_parameters(ParameterHandler ¶m)
Msh(const bool write_faces=false, const bool write_lines=false)
void parse_parameters(ParameterHandler ¶m)
static void declare_parameters(ParameterHandler ¶m)
bool label_level_subdomain_id
unsigned int line_thickness
bool convert_level_number_to_height
Svg(const unsigned int line_thickness=2, const unsigned int boundary_line_thickness=4, const 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=false, const bool label_cell_index=false, const bool label_material_id=false, const bool label_subdomain_id=false, const bool draw_colorbar=false, const bool draw_legend=false, const bool label_boundary_id=false)
@ level_subdomain_id
Convert the level subdomain id into the cell color.
@ subdomain_id
Convert the subdomain id into the cell color.
@ material_id
Convert the material id into the cell color (default)
@ level_number
Convert the level number into the cell color.
unsigned int boundary_line_thickness
float level_height_factor
Ucd(const bool write_preamble=false, const bool write_faces=false, const bool write_lines=false)
static void declare_parameters(ParameterHandler ¶m)
void parse_parameters(ParameterHandler ¶m)
bool output_only_relevant
bool serialize_triangulation
unsigned int n_boundary_face_points
@ level_number
Convert the level into the cell color.
@ material_id
Convert the material id into the cell color.
@ level_subdomain_id
Convert the level subdomain id into the cell color.
@ subdomain_id
Convert the global subdomain id into the cell color.
void parse_parameters(ParameterHandler ¶m)
static void declare_parameters(ParameterHandler ¶m)
enum GridOutFlags::XFig::Coloring color_by
constexpr Tensor< 1, dim, typename ProductType< Number1, Number2 >::type > cross_product_3d(const Tensor< 1, dim, Number1 > &src1, const Tensor< 1, dim, Number2 > &src2)
const ::Triangulation< dim, spacedim > & tria