]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Address comments
authorAlexander Grayver <agrayver@erdw.ethz.ch>
Thu, 26 Jul 2018 20:31:52 +0000 (22:31 +0200)
committerAlexander Grayver <agrayver@erdw.ethz.ch>
Fri, 27 Jul 2018 07:45:35 +0000 (09:45 +0200)
include/deal.II/base/data_out_base.h
source/base/data_out_base.cc
tests/data_out/data_out_base_vtk_03.cc [new file with mode: 0644]
tests/data_out/data_out_base_vtk_03.output [new file with mode: 0644]

index c448e9185e33fa0a63d30c93f32444d45c921bb6..8ecc24f01fa0be8cd7eb3485a94ac4f6be3d4a17 100644 (file)
@@ -1154,7 +1154,7 @@ namespace DataOutBase
      *
      * Default is <tt>false</tt>.
      */
-    bool write_lagrange_cells;
+    bool write_higher_order_cells;
 
     /**
      * Constructor.
@@ -1162,9 +1162,9 @@ namespace DataOutBase
     VtkFlags(
       const double       time  = std::numeric_limits<double>::min(),
       const unsigned int cycle = std::numeric_limits<unsigned int>::min(),
-      const bool         print_date_and_time          = true,
-      const ZlibCompressionLevel compression_level    = best_compression,
-      const bool                 write_lagrange_cells = false);
+      const bool         print_date_and_time              = true,
+      const ZlibCompressionLevel compression_level        = best_compression,
+      const bool                 write_higher_order_cells = false);
   };
 
 
index daa7f29e5c4d7dccafb2beb1c00c9b58976d20a4..08a8d5c4126ff423d16e129cb6818dbe3f948665 100644 (file)
@@ -792,16 +792,16 @@ namespace
   }
 
   /**
-   * Given (i,j,k) coordinates within the Lagrange cell, return an
-   * offset into the local connectivity array
+   * Given (i,j,k) coordinates within the Lagrange quadrilateral, return an
+   * offset into the local connectivity array.
    *
    * Modified from
-   * https://github.com/Kitware/VTK/blob/265ca48a79a36538c95622c237da11133608bbe5/Common/DataModel/vtkLagrangeQuadrilateral.cxx#L558
+   * https://github.com/Kitware/VTK/blob/265ca48a/Common/DataModel/vtkLagrangeQuadrilateral.cxx#L558
    */
   int
-  vtk_point_index_from_ijk(int i,
-                           int j,
-                           int,
+  vtk_point_index_from_ijk(unsigned i,
+                           unsigned j,
+                           unsigned,
                            const std::array<unsigned, 2> &order)
   {
     bool ibdy = (i == 0 || i == order[0]);
@@ -810,7 +810,7 @@ namespace
     int nbdy = (ibdy ? 1 : 0) + (jbdy ? 1 : 0);
 
     if (nbdy == 2) // Vertex DOF
-      { // ijk is a corner node. Return the proper index (somewhere in [0,7]):
+      { // ijk is a corner node. Return the proper index (somewhere in [0,3]):
         return (i ? (j ? 2 : 1) : (j ? 3 : 0));
       }
 
@@ -836,15 +836,16 @@ namespace
   }
 
   /**
-   * Same as above, but for hexes.
+   * Given (i,j,k) coordinates within the Lagrange hexahedron, return an
+   * offset into the local connectivity array.
    *
    * Modified from
-   * https://github.com/Kitware/VTK/blob/265ca48a79a36538c95622c237da11133608bbe5/Common/DataModel/vtkLagrangeHexahedron.cxx#L734
+   * https://github.com/Kitware/VTK/blob/265ca48a/Common/DataModel/vtkLagrangeHexahedron.cxx#L734
    */
   int
-  vtk_point_index_from_ijk(int                            i,
-                           int                            j,
-                           int                            k,
+  vtk_point_index_from_ijk(unsigned                       i,
+                           unsigned                       j,
+                           unsigned                       k,
                            const std::array<unsigned, 3> &order)
   {
     bool ibdy = (i == 0 || i == order[0]);
@@ -907,16 +908,22 @@ namespace
   }
 
   int
-  vtk_point_index_from_ijk(int, int, int, const std::array<unsigned, 0> &)
+  vtk_point_index_from_ijk(unsigned,
+                           unsigned,
+                           unsigned,
+                           const std::array<unsigned, 0> &)
   {
-    ExcNotImplemented();
+    Assert(false, ExcNotImplemented());
     return 0;
   }
 
   int
-  vtk_point_index_from_ijk(int, int, int, const std::array<unsigned, 1> &)
+  vtk_point_index_from_ijk(unsigned,
+                           unsigned,
+                           unsigned,
+                           const std::array<unsigned, 1> &)
   {
-    ExcNotImplemented();
+    Assert(false, ExcNotImplemented());
     return 0;
   }
 
@@ -2383,12 +2390,12 @@ namespace DataOutBase
                      const unsigned int                   cycle,
                      const bool                           print_date_and_time,
                      const VtkFlags::ZlibCompressionLevel compression_level,
-                     const bool                           write_lagrange_cells)
+                     const bool write_higher_order_cells)
     : time(time)
     , cycle(cycle)
     , print_date_and_time(print_date_and_time)
     , compression_level(compression_level)
-    , write_lagrange_cells(write_lagrange_cells)
+    , write_higher_order_cells(write_higher_order_cells)
   {}
 
 
@@ -2547,9 +2554,9 @@ namespace DataOutBase
         const unsigned int n2 = (dim > 1) ? n_subdivisions : 1;
         const unsigned int n3 = (dim > 2) ? n_subdivisions : 1;
         // Offsets of outer loops
-        unsigned int d1 = 1;
-        unsigned int d2 = n;
-        unsigned int d3 = n * n;
+        const unsigned int d1 = 1;
+        const unsigned int d2 = n;
+        const unsigned int d3 = n * n;
         for (unsigned int i3 = 0; i3 < n3; ++i3)
           for (unsigned int i2 = 0; i2 < n2; ++i2)
             for (unsigned int i1 = 0; i1 < n1; ++i1)
@@ -2596,17 +2603,17 @@ namespace DataOutBase
         const unsigned int n2 = (dim > 1) ? n_subdivisions : 0;
         const unsigned int n3 = (dim > 2) ? n_subdivisions : 0;
         // Offsets of outer loops
-        unsigned int d1 = 1;
-        unsigned int d2 = n;
-        unsigned int d3 = n * n;
+        const unsigned int d1 = 1;
+        const unsigned int d2 = n;
+        const unsigned int d3 = n * n;
         for (unsigned int i3 = 0; i3 <= n3; ++i3)
           for (unsigned int i2 = 0; i2 <= n2; ++i2)
             for (unsigned int i1 = 0; i1 <= n1; ++i1)
               {
-                const unsigned int local_offset = i3 * d3 + i2 * d2 + i1 * d1;
-                const unsigned int connectivity_idx =
+                const unsigned int local_index = i3 * d3 + i2 * d2 + i1 * d1;
+                const unsigned int connectivity_index =
                   vtk_point_index_from_ijk(i1, i2, i3, cell_order);
-                connectivity[connectivity_idx] = local_offset;
+                connectivity[connectivity_index] = local_index;
               }
 
         out.template write_lagrange_cell<dim>(count++,
@@ -5274,11 +5281,11 @@ namespace DataOutBase
     unsigned int n_cells;
     compute_sizes<dim, spacedim>(patches, n_nodes, n_cells);
 
-    // If a user set to out Lagrange cells, we treat n_subdivions
-    // as cell order and adjust variables accordingly, otherwise
+    // If a user set to output Lagrange cells, we treat n_subdivisions
+    // as cell order and adjust variables accordingly, otherwise
     // each patch is written as a linear cell.
     unsigned int n_points_per_cell = GeometryInfo<dim>::vertices_per_cell;
-    if (flags.write_lagrange_cells)
+    if (flags.write_higher_order_cells)
       {
         n_cells           = patches.size();
         n_points_per_cell = n_nodes / n_cells;
@@ -5315,7 +5322,7 @@ namespace DataOutBase
     // now for the cells
     out << "CELLS " << n_cells << ' ' << n_cells * (n_points_per_cell + 1)
         << '\n';
-    if (flags.write_lagrange_cells)
+    if (flags.write_higher_order_cells)
       write_lagrange_cells(patches, vtk_out);
     else
       write_cells(patches, vtk_out);
@@ -5325,7 +5332,7 @@ namespace DataOutBase
     out << "CELL_TYPES " << n_cells << '\n';
 
     // need to distinguish between linear and Lagrange cells
-    const unsigned int vtk_cell_id = flags.write_lagrange_cells ?
+    const unsigned int vtk_cell_id = flags.write_higher_order_cells ?
                                        vtk_lagrange_cell_type[dim] :
                                        vtk_cell_type[dim];
     for (unsigned int i = 0; i < n_cells; ++i)
diff --git a/tests/data_out/data_out_base_vtk_03.cc b/tests/data_out/data_out_base_vtk_03.cc
new file mode 100644 (file)
index 0000000..98c111f
--- /dev/null
@@ -0,0 +1,53 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2006 - 2017 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.
+//
+// ---------------------------------------------------------------------
+
+// Test high-order Lagrange VTK output for on a unit cube.
+// Added by Alexander Grayver
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria.h>
+
+#include <deal.II/numerics/data_out.h>
+
+#include <string>
+
+#include "../tests.h"
+
+template <int dim>
+void
+check(std::ostream &log, unsigned cell_order)
+{
+  Triangulation<dim> triangulation;
+  GridGenerator::hyper_cube(triangulation);
+
+  DataOutBase::VtkFlags flags;
+  flags.write_higher_order_cells = true;
+
+  DataOut<dim> data_out;
+  data_out.set_flags(flags);
+  data_out.attach_triangulation(triangulation);
+  data_out.build_patches(cell_order);
+  data_out.write_vtk(log);
+}
+
+int
+main()
+{
+  std::ofstream logfile("output");
+
+  unsigned cell_order = 4;
+  check<2>(logfile, cell_order);
+  check<3>(logfile, cell_order);
+}
diff --git a/tests/data_out/data_out_base_vtk_03.output b/tests/data_out/data_out_base_vtk_03.output
new file mode 100644 (file)
index 0000000..0af34c3
--- /dev/null
@@ -0,0 +1,176 @@
+# vtk DataFile Version 3.0
+#This file was generated 
+ASCII
+DATASET UNSTRUCTURED_GRID
+
+POINTS 25 double
+0 0 0
+0.25 0 0
+0.5 0 0
+0.75 0 0
+1 0 0
+0 0.25 0
+0.25 0.25 0
+0.5 0.25 0
+0.75 0.25 0
+1 0.25 0
+0 0.5 0
+0.25 0.5 0
+0.5 0.5 0
+0.75 0.5 0
+1 0.5 0
+0 0.75 0
+0.25 0.75 0
+0.5 0.75 0
+0.75 0.75 0
+1 0.75 0
+0 1 0
+0.25 1 0
+0.5 1 0
+0.75 1 0
+1 1 0
+
+CELLS 1 26
+25     0       4       24      20      1       2       3       9       14      19      21      22      23      5       10      15      6       7       8       11      12      13      16      17      18
+
+CELL_TYPES 1
+ 70
+POINT_DATA 25
+# vtk DataFile Version 3.0
+#This file was generated 
+ASCII
+DATASET UNSTRUCTURED_GRID
+
+POINTS 125 double
+0 0 0
+0.25 0 0
+0.5 0 0
+0.75 0 0
+1 0 0
+0 0.25 0
+0.25 0.25 0
+0.5 0.25 0
+0.75 0.25 0
+1 0.25 0
+0 0.5 0
+0.25 0.5 0
+0.5 0.5 0
+0.75 0.5 0
+1 0.5 0
+0 0.75 0
+0.25 0.75 0
+0.5 0.75 0
+0.75 0.75 0
+1 0.75 0
+0 1 0
+0.25 1 0
+0.5 1 0
+0.75 1 0
+1 1 0
+0 0 0.25
+0.25 0 0.25
+0.5 0 0.25
+0.75 0 0.25
+1 0 0.25
+0 0.25 0.25
+0.25 0.25 0.25
+0.5 0.25 0.25
+0.75 0.25 0.25
+1 0.25 0.25
+0 0.5 0.25
+0.25 0.5 0.25
+0.5 0.5 0.25
+0.75 0.5 0.25
+1 0.5 0.25
+0 0.75 0.25
+0.25 0.75 0.25
+0.5 0.75 0.25
+0.75 0.75 0.25
+1 0.75 0.25
+0 1 0.25
+0.25 1 0.25
+0.5 1 0.25
+0.75 1 0.25
+1 1 0.25
+0 0 0.5
+0.25 0 0.5
+0.5 0 0.5
+0.75 0 0.5
+1 0 0.5
+0 0.25 0.5
+0.25 0.25 0.5
+0.5 0.25 0.5
+0.75 0.25 0.5
+1 0.25 0.5
+0 0.5 0.5
+0.25 0.5 0.5
+0.5 0.5 0.5
+0.75 0.5 0.5
+1 0.5 0.5
+0 0.75 0.5
+0.25 0.75 0.5
+0.5 0.75 0.5
+0.75 0.75 0.5
+1 0.75 0.5
+0 1 0.5
+0.25 1 0.5
+0.5 1 0.5
+0.75 1 0.5
+1 1 0.5
+0 0 0.75
+0.25 0 0.75
+0.5 0 0.75
+0.75 0 0.75
+1 0 0.75
+0 0.25 0.75
+0.25 0.25 0.75
+0.5 0.25 0.75
+0.75 0.25 0.75
+1 0.25 0.75
+0 0.5 0.75
+0.25 0.5 0.75
+0.5 0.5 0.75
+0.75 0.5 0.75
+1 0.5 0.75
+0 0.75 0.75
+0.25 0.75 0.75
+0.5 0.75 0.75
+0.75 0.75 0.75
+1 0.75 0.75
+0 1 0.75
+0.25 1 0.75
+0.5 1 0.75
+0.75 1 0.75
+1 1 0.75
+0 0 1
+0.25 0 1
+0.5 0 1
+0.75 0 1
+1 0 1
+0 0.25 1
+0.25 0.25 1
+0.5 0.25 1
+0.75 0.25 1
+1 0.25 1
+0 0.5 1
+0.25 0.5 1
+0.5 0.5 1
+0.75 0.5 1
+1 0.5 1
+0 0.75 1
+0.25 0.75 1
+0.5 0.75 1
+0.75 0.75 1
+1 0.75 1
+0 1 1
+0.25 1 1
+0.5 1 1
+0.75 1 1
+1 1 1
+
+CELLS 1 126
+125    0       4       24      20      100     104     124     120     1       2       3       9       14      19      21      22      23      5       10      15      101     102     103     109     114     119     121     122     123     105     110     115     25      50      75      29      54      79      45      70      95      49      74      99      30      35      40      55      60      65      80      85      90      34      39      44      59      64      69      84      89      94      26      27      28      51      52      53      76      77      78      46      47      48      71      72      73      96      97      98      6       7       8       11      12      13      16      17      18      106     107     108     111     112     113     116     117     118     31      32      33      36      37      38      41      42      43      56      57      58      61      62      63      66      67      68      81      82      83      86      87      88      91      92      93
+
+CELL_TYPES 1
+ 72
+POINT_DATA 125

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.