]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Extend GridOut::write_svg() for simplices 11090/head
authorElias Dejene <elias.dejene@tum.de>
Fri, 23 Oct 2020 08:23:26 +0000 (10:23 +0200)
committerElias Dejene <elias.dejene@tum.de>
Wed, 28 Oct 2020 13:26:32 +0000 (14:26 +0100)
Add output for ctest. Replace forgotten write_eps by write_svg.

Add comment for quad-case. Avoid writing extra point in case of simplices.

Update ctest source-file.

Replace post-increment by pre-increment.

source/grid/grid_out.cc
tests/simplex/grid_out_svg.cc [new file with mode: 0644]
tests/simplex/grid_out_svg.with_simplex_support=on.output [new file with mode: 0644]

index 16e1a7f5e8f73b3449ce7c7e093feade75e515ee..73e2ea34b79d1b052e7e4687eb6043b435cdf0a1 100644 (file)
@@ -1584,7 +1584,8 @@ GridOut::write_svg(const Triangulation<2, 2> &tria, std::ostream &out) const
   // (, and level subdomain id).
   for (const auto &cell : tria.cell_iterators())
     {
-      for (unsigned int vertex_index = 0; vertex_index < 4; vertex_index++)
+      for (unsigned int vertex_index = 0; vertex_index < cell->n_vertices();
+           ++vertex_index)
         {
           if (cell->vertex(vertex_index)[0] < x_min)
             x_min = cell->vertex(vertex_index)[0];
@@ -1834,24 +1835,27 @@ GridOut::write_svg(const Triangulation<2, 2> &tria, std::ostream &out) const
       if (y_min_perspective > projection_decomposition[1])
         y_min_perspective = projection_decomposition[1];
 
-      point[0] = cell->vertex(3)[0];
-      point[1] = cell->vertex(3)[1];
+      if (cell->n_vertices() == 4) // in case of quadrilateral
+        {
+          point[0] = cell->vertex(3)[0];
+          point[1] = cell->vertex(3)[1];
 
-      projection_decomposition = svg_project_point(point,
-                                                   camera_position,
-                                                   camera_direction,
-                                                   camera_horizontal,
-                                                   camera_focus);
+          projection_decomposition = svg_project_point(point,
+                                                       camera_position,
+                                                       camera_direction,
+                                                       camera_horizontal,
+                                                       camera_focus);
 
-      if (x_max_perspective < projection_decomposition[0])
-        x_max_perspective = projection_decomposition[0];
-      if (x_min_perspective > projection_decomposition[0])
-        x_min_perspective = projection_decomposition[0];
+          if (x_max_perspective < projection_decomposition[0])
+            x_max_perspective = projection_decomposition[0];
+          if (x_min_perspective > projection_decomposition[0])
+            x_min_perspective = projection_decomposition[0];
 
-      if (y_max_perspective < projection_decomposition[1])
-        y_max_perspective = projection_decomposition[1];
-      if (y_min_perspective > projection_decomposition[1])
-        y_min_perspective = projection_decomposition[1];
+          if (y_max_perspective < projection_decomposition[1])
+            y_max_perspective = projection_decomposition[1];
+          if (y_min_perspective > projection_decomposition[1])
+            y_min_perspective = projection_decomposition[1];
+        }
 
       if (static_cast<unsigned int>(cell->level()) == min_level)
         min_level_min_vertex_distance = cell->minimum_vertex_distance();
@@ -2175,29 +2179,32 @@ GridOut::write_svg(const Triangulation<2, 2> &tria, std::ostream &out) const
 
           out << " L ";
 
-          point[0] = cell->vertex(3)[0];
-          point[1] = cell->vertex(3)[1];
+          if (cell->n_vertices() == 4) // in case of quadrilateral
+            {
+              point[0] = cell->vertex(3)[0];
+              point[1] = cell->vertex(3)[1];
 
-          projection_decomposition = svg_project_point(point,
-                                                       camera_position,
-                                                       camera_direction,
-                                                       camera_horizontal,
-                                                       camera_focus);
+              projection_decomposition = svg_project_point(point,
+                                                           camera_position,
+                                                           camera_direction,
+                                                           camera_horizontal,
+                                                           camera_focus);
 
-          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))
-              << ' '
-              << 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));
+              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))
+                  << ' '
+                  << 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));
 
-          out << " L ";
+              out << " L ";
+            }
 
           point[0] = cell->vertex(2)[0];
           point[1] = cell->vertex(2)[1];
@@ -2354,8 +2361,7 @@ GridOut::write_svg(const Triangulation<2, 2> &tria, std::ostream &out) const
           // the additional boundary line
           if (svg_flags.boundary_line_thickness)
             {
-              for (const unsigned int faceIndex :
-                   GeometryInfo<2>::face_indices())
+              for (auto faceIndex : cell->face_indices())
                 {
                   if (cell->at_boundary(faceIndex))
                     {
diff --git a/tests/simplex/grid_out_svg.cc b/tests/simplex/grid_out_svg.cc
new file mode 100644 (file)
index 0000000..73c795f
--- /dev/null
@@ -0,0 +1,60 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+// Write a file in the SVG format with (linear) triangular and tetrahedral
+// elements created by the GMSH program.
+
+#include <deal.II/grid/grid_in.h>
+#include <deal.II/grid/grid_out.h>
+#include <deal.II/grid/tria.h>
+
+#include <fstream>
+
+#include "../tests.h"
+
+template <int dim, int spacedim>
+void
+check_file(const std::string &file_name)
+{
+  Triangulation<dim, spacedim> tria;
+
+  GridIn<dim, spacedim> grid_in;
+  grid_in.attach_triangulation(tria);
+  std::ifstream input_file(file_name);
+  grid_in.read_vtk(input_file);
+
+  GridOut grid_out;
+#if false
+  std::ofstream out("mesh-" + std::to_string(spacedim) + ".svg");
+  grid_out.write_svg(tria, out);
+#else
+  grid_out.write_svg(tria, deallog.get_file_stream());
+#endif
+
+  deallog << "OK!" << std::endl;
+}
+
+
+int
+main()
+{
+  initlog();
+  // TRIANGULAR ELEMENTS
+  // dim = spacedim = 2
+  deallog.push("triangluar_elements_dim2_spacedim2: ");
+  check_file<2, 2>(std::string(SOURCE_DIR "/grid_in_vtk/tri.vtk"));
+  deallog.pop();
+}
diff --git a/tests/simplex/grid_out_svg.with_simplex_support=on.output b/tests/simplex/grid_out_svg.with_simplex_support=on.output
new file mode 100644 (file)
index 0000000..92f3b86
--- /dev/null
@@ -0,0 +1,30 @@
+
+<svg width="1000" height="1000" xmlns="http://www.w3.org/2000/svg" version="1.1">
+
+
+<!-- internal style sheet -->
+<style type="text/css"><![CDATA[
+ rect.background{fill:white}
+ rect{fill:none; stroke:rgb(25,25,25); stroke-width:2}
+ text{font-family:Helvetica; text-anchor:middle; fill:rgb(25,25,25)}
+ line{stroke:rgb(25,25,25); stroke-width:4}
+ path{fill:none; stroke:rgb(25,25,25); stroke-width:2}
+ circle{fill:white; stroke:black; stroke-width:2}
+
+ path.p0{fill:rgb(0,102,255); stroke:rgb(25,25,25); stroke-width:2}
+ path.ps0{fill:rgb(0,77,191); stroke:rgb(20,20,20); stroke-width:2}
+ rect.r0{fill:rgb(0,102,255); stroke:rgb(25,25,25); stroke-width:2}
+]]></style>
+
+ <rect class="background" width="1000" height="1000"/>
+  <!-- cells -->
+  <path class="p0" d="M 80 920 L 920 920 L 500 500 L 80 920"/>
+  <line x1="80" y1="920" x2="920" y2="920"/>
+  <path class="p0" d="M 80 80 L 80 920 L 500 500 L 80 80"/>
+  <line x1="80" y1="80" x2="80" y2="920"/>
+  <path class="p0" d="M 920 920 L 920 80 L 500 500 L 920 920"/>
+  <line x1="920" y1="920" x2="920" y2="80"/>
+  <path class="p0" d="M 920 80 L 80 80 L 500 500 L 920 80"/>
+  <line x1="920" y1="80" x2="80" y2="80"/>
+
+</svg>DEAL:triangluar_elements_dim2_spacedim2: ::OK!

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.