]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add mesh images
authorAlexander Grayver <agrayver@erdw.ethz.ch>
Fri, 20 Mar 2020 18:57:42 +0000 (19:57 +0100)
committerAlexander Grayver <agrayver@erdw.ethz.ch>
Sat, 21 Mar 2020 07:31:54 +0000 (08:31 +0100)
examples/step-11/doc/results.dox
examples/step-11/step-11.cc

index ebb209867be7d82cedf0249f7570c232a6a66687..ec32c393779c7e3228ef499b68fb70b7c4763500 100644 (file)
@@ -42,3 +42,24 @@ of convergence but just to reduce the constant before the convergence
 order. On the other hand, using a cubic mapping only improves the
 result further insignificantly, except for the case of very coarse
 grids.
+
+We can also visualize the underlying meshes by using, for instance
+ParaView. The image below shows initial meshes for different mapping
+degrees:
+
+<img src="https://www.dealii.org/images/steps/developer/step-11-cycle_0.png" alt="">
+
+Clearly, the effect is most pronounced when we go from the linear to
+quadratic mapping. This is also reflected in the error values given
+in the table above. The effect of going from quadratic to cubic degree
+is less dramatic, but still tangible owning to a more accurace
+description of the circular boundary.
+
+Next, let's look at the meshes after three global refinements
+
+<img src="https://www.dealii.org/images/steps/developer/step-11-cycle_3.png" alt="">
+
+Here, the differences are much less visible, especially for higher order
+mappings. Indeed, at this refinement level the error values reported
+in the table are essentially identical between mappings of degree two
+and three.
index a97816dbeb6c043edae818474bc3af4c4e402a33..10625546f47a23b4ac0d531b1b144586c7651d16 100644 (file)
@@ -41,6 +41,7 @@
 #include <deal.II/fe/mapping_q.h>
 #include <deal.II/numerics/vector_tools.h>
 #include <deal.II/numerics/matrix_tools.h>
+#include <deal.II/numerics/data_out.h>
 
 // Just this one is new: it declares a class
 // DynamicSparsityPattern, which we will use and explain
@@ -83,6 +84,7 @@ namespace Step11
     void setup_system();
     void assemble_and_solve();
     void solve();
+    void write_high_order_mesh(const unsigned cycle);
 
     Triangulation<dim> triangulation;
     FE_Q<dim>          fe;
@@ -397,6 +399,57 @@ namespace Step11
 
 
 
+  // Next, we write the solution as well as the
+  // material ids to a VTU file. This is similar to what was done in many
+  // other tutorial programs. The new ingredient presented in this tutorial
+  // program is that we want to ensure that the data written to the file
+  // used for visualization is actually a faithful representation of what
+  // is used internally by deal.II. That is because most of the visualization
+  // data formats only represent cells by their vertex coordinates, but
+  // have no way of representing the curved boundaries that are used
+  // in deal.II when using higher order mappings -- in other words, what
+  // you see in the visualization tool is not actually what you are computing
+  // on. (The same, incidentally, is true when using higher order shape
+  // functions: Most visualization tools only render bilinear/trilinear
+  // representations. This is discussed in detail in DataOut::build_patches().)
+  //
+  // So we need to ensure that a high-order representation is written
+  // to the file. We need to consider two particular topics. Firstly, we tell
+  // the DataOut object via the DataOutBase::VtkFlags that we intend to
+  // interpret the subdivisions of the elements as a high-order Lagrange
+  // polynomial rather than a collection of bilinear patches.
+  // Recent visualization programs, like ParaView version 5.5
+  // or newer, can then render a high-order solution (see a <a
+  // href="https://github.com/dealii/dealii/wiki/Notes-on-visualizing-high-order-output">wiki
+  // page</a> for more details). Secondly, we need to make sure that the mapping
+  // is passed to the DataOut::build_patches() method. Finally, the DataOut
+  // class only prints curved faces for <i>boundary</i> cells by default, so we
+  // need to ensure that also inner cells are printed in a curved representation
+  // via the mapping.
+  template <int dim>
+  void LaplaceProblem<dim>::write_high_order_mesh(const unsigned cycle)
+  {
+    DataOut<dim> data_out;
+
+    DataOutBase::VtkFlags flags;
+    flags.write_higher_order_cells = true;
+    data_out.set_flags(flags);
+
+    data_out.attach_dof_handler(dof_handler);
+    data_out.add_data_vector(solution, "solution");
+
+    data_out.build_patches(mapping,
+                           mapping.get_degree(),
+                           DataOut<dim>::curved_inner_cells);
+
+    std::ofstream file(("solution-c=" + std::to_string(cycle) +
+                        ".p=" + std::to_string(mapping.get_degree()) + ".vtu")
+                         .c_str());
+
+    data_out.write_vtu(file);
+  }
+
+
   // Finally the main function controlling the different steps to be
   // performed. Its content is rather straightforward, generating a
   // triangulation of a circle, associating a boundary to it, and then doing
@@ -419,6 +472,7 @@ namespace Step11
       {
         setup_system();
         assemble_and_solve();
+        write_high_order_mesh(cycle);
 
         triangulation.refine_global();
       }

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.