]> https://gitweb.dealii.org/ - dealii.git/commitdiff
deprecate VtkFalgs::ZlibCompressionLevel
authorTimo Heister <timo.heister@gmail.com>
Thu, 14 Jul 2022 01:58:24 +0000 (21:58 -0400)
committerTimo Heister <timo.heister@gmail.com>
Thu, 14 Jul 2022 01:58:24 +0000 (21:58 -0400)
examples/step-23/step-23.cc
examples/step-24/step-24.cc
examples/step-25/step-25.cc
examples/step-9/step-9.cc
include/deal.II/base/data_out_base.h
source/base/data_out_base.cc
source/base/data_out_base.inst.in
tests/simplex/step-23.cc

index 5b44c4f2217f9425a17174e567ed5222946b6f30..444a94784a7a580eb4ed957e06278e1cdc60536a 100644 (file)
@@ -412,8 +412,7 @@ namespace Step23
     // zlib compression algorithm that is optimized for speed instead of disk
     // usage since otherwise plotting the output becomes a bottleneck:
     DataOutBase::VtkFlags vtk_flags;
-    vtk_flags.compression_level =
-      DataOutBase::VtkFlags::ZlibCompressionLevel::best_speed;
+    vtk_flags.compression_level = DataOutBase::best_speed;
     data_out.set_flags(vtk_flags);
     std::ofstream output(filename);
     data_out.write_vtu(output);
index 448b44d89291c1df0d274691e70f79e812dd762a..ee7d7477a34d64f7ee53bf0a825935d09346dc4c 100644 (file)
@@ -431,8 +431,7 @@ namespace Step24
     const std::string filename =
       "solution-" + Utilities::int_to_string(timestep_number, 3) + ".vtu";
     DataOutBase::VtkFlags vtk_flags;
-    vtk_flags.compression_level =
-      DataOutBase::VtkFlags::ZlibCompressionLevel::best_speed;
+    vtk_flags.compression_level = DataOutBase::best_speed;
     std::ofstream output(filename);
     data_out.write_vtu(output);
   }
index ee2cabd13deca95a0332cf71b22519d4d7ac6713..5c0c3fa073d7590043d55d76ad077ca229a9159b 100644 (file)
@@ -559,8 +559,7 @@ namespace Step25
     const std::string filename =
       "solution-" + Utilities::int_to_string(timestep_number, 3) + ".vtu";
     DataOutBase::VtkFlags vtk_flags;
-    vtk_flags.compression_level =
-      DataOutBase::VtkFlags::ZlibCompressionLevel::best_speed;
+    vtk_flags.compression_level = DataOutBase::best_speed;
     data_out.set_flags(vtk_flags);
     std::ofstream output(filename);
     data_out.write_vtu(output);
index c00b03a26c717d9d61a838fc559f59394096ff5d..26d2db469e202570ad01bc58624d0a5f24c06b22 100644 (file)
@@ -844,8 +844,7 @@ namespace Step9
       // disk. Here we ask ZLib, a compression library, to compress the data
       // in a way that maximizes throughput.
       DataOutBase::VtkFlags vtk_flags;
-      vtk_flags.compression_level =
-        DataOutBase::VtkFlags::ZlibCompressionLevel::best_speed;
+      vtk_flags.compression_level = DataOutBase::best_speed;
       data_out.set_flags(vtk_flags);
 
       std::ofstream output("solution-" + std::to_string(cycle) + ".vtu");
index 2fa678252f15804246c5e0852f6214326798155a..75f316d514f7f0cc5d6623f37504d45a14832cbb 100644 (file)
@@ -220,6 +220,33 @@ class XDMFEntry;
  */
 namespace DataOutBase
 {
+  /**
+   * An enum for different levels of compression used in several places
+   * to determine zlib compression levels.
+   */
+  enum CompressionLevel
+  {
+    /**
+     * Do not use any compression.
+     */
+    no_compression,
+    /**
+     * Use the fastest available compression algorithm.
+     */
+    best_speed,
+    /**
+     * Use the algorithm which results in the smallest compressed
+     * files. This is the default flag.
+     */
+    best_compression,
+    /**
+     * Use the default compression algorithm. This is a compromise between
+     * speed and file size.
+     */
+    default_compression
+  };
+
+
   /**
    * Data structure describing a patch of data in <tt>dim</tt> space
    * dimensions.
@@ -1124,34 +1151,26 @@ namespace DataOutBase
     /**
      * A data type providing the different possible zlib compression
      * levels. These map directly to constants defined by zlib.
+     *
+     * @deprecated Use DataOutBase::CompressionLevel instead.
      */
-    enum ZlibCompressionLevel
-    {
-      /**
-       * Do not use any compression.
-       */
-      no_compression,
-      /**
-       * Use the fastest available compression algorithm.
-       */
-      best_speed,
-      /**
-       * Use the algorithm which results in the smallest compressed
-       * files. This is the default flag.
-       */
-      best_compression,
-      /**
-       * Use the default compression algorithm. This is a compromise between
-       * speed and file size.
-       */
-      default_compression
-    };
+    using ZlibCompressionLevel DEAL_II_DEPRECATED =
+      DataOutBase::CompressionLevel;
+
+    DEAL_II_DEPRECATED static const DataOutBase::CompressionLevel
+      no_compression = dealii::DataOutBase::no_compression;
+    DEAL_II_DEPRECATED static const DataOutBase::CompressionLevel
+      best_compression = dealii::DataOutBase::best_compression;
+    DEAL_II_DEPRECATED static const DataOutBase::CompressionLevel best_speed =
+      dealii::DataOutBase::best_speed;
+    DEAL_II_DEPRECATED static const DataOutBase::CompressionLevel
+      default_compression = dealii::DataOutBase::default_compression;
 
     /**
      * Flag determining the compression level at which zlib, if available, is
      * run. The default is <tt>best_compression</tt>.
      */
-    ZlibCompressionLevel compression_level;
+    DataOutBase::CompressionLevel compression_level;
 
     /**
      * Flag determining whether to write patches as linear cells
@@ -1200,11 +1219,11 @@ namespace DataOutBase
      * to the argument names of this function.
      */
     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_higher_order_cells = false,
+      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 CompressionLevel compression_level        = best_compression,
+      const bool             write_higher_order_cells = false,
       const std::map<std::string, std::string> &physical_units = {});
   };
 
@@ -2361,11 +2380,11 @@ namespace DataOutBase
                  unsigned int,
                  std::string,
                  DataComponentInterpretation::DataComponentInterpretation>>
-      &                                  nonscalar_data_ranges,
-    const Deal_II_IntermediateFlags &    flags,
-    const std::string &                  filename,
-    const MPI_Comm &                     comm,
-    const VtkFlags::ZlibCompressionLevel compression);
+      &                              nonscalar_data_ranges,
+    const Deal_II_IntermediateFlags &flags,
+    const std::string &              filename,
+    const MPI_Comm &                 comm,
+    const CompressionLevel           compression);
 
   /**
    * Write the data in @p data_filter to a single HDF5 file containing both the
@@ -2858,9 +2877,9 @@ public:
    */
   void
   write_deal_II_intermediate_in_parallel(
-    const std::string &                               filename,
-    const MPI_Comm &                                  comm,
-    const DataOutBase::VtkFlags::ZlibCompressionLevel compression) const;
+    const std::string &                 filename,
+    const MPI_Comm &                    comm,
+    const DataOutBase::CompressionLevel compression) const;
 
   /**
    * Create an XDMFEntry based on the data in the data_filter. This assumes
index 9caf6a6033f65226712a24660a6b60ff74600c86..701f7d7d890c7927e15e10f108151fc169165aae 100644 (file)
@@ -94,18 +94,17 @@ namespace
    * constant defined by zlib.
    */
   int
-  get_zlib_compression_level(
-    const DataOutBase::VtkFlags::ZlibCompressionLevel level)
+  get_zlib_compression_level(const DataOutBase::CompressionLevel level)
   {
     switch (level)
       {
-        case (DataOutBase::VtkFlags::no_compression):
+        case (DataOutBase::no_compression):
           return Z_NO_COMPRESSION;
-        case (DataOutBase::VtkFlags::best_speed):
+        case (DataOutBase::best_speed):
           return Z_BEST_SPEED;
-        case (DataOutBase::VtkFlags::best_compression):
+        case (DataOutBase::best_compression):
           return Z_BEST_COMPRESSION;
-        case (DataOutBase::VtkFlags::default_compression):
+        case (DataOutBase::default_compression):
           return Z_DEFAULT_COMPRESSION;
         default:
           Assert(false, ExcNotImplemented());
@@ -113,29 +112,30 @@ namespace
       }
   }
 
+#  ifdef DEAL_II_WITH_MPI
   /**
    * Convert between the enum specified inside VtkFlags and the preprocessor
    * constant defined by boost::iostreams::zlib.
    */
   int
-  get_boost_zlib_compression_level(
-    const DataOutBase::VtkFlags::ZlibCompressionLevel level)
+  get_boost_zlib_compression_level(const DataOutBase::CompressionLevel level)
   {
     switch (level)
       {
-        case (DataOutBase::VtkFlags::no_compression):
+        case (DataOutBase::no_compression):
           return boost::iostreams::zlib::no_compression;
-        case (DataOutBase::VtkFlags::best_speed):
+        case (DataOutBase::best_speed):
           return boost::iostreams::zlib::best_speed;
-        case (DataOutBase::VtkFlags::best_compression):
+        case (DataOutBase::best_compression):
           return boost::iostreams::zlib::best_compression;
-        case (DataOutBase::VtkFlags::default_compression):
+        case (DataOutBase::default_compression):
           return boost::iostreams::zlib::default_compression;
         default:
           Assert(false, ExcNotImplemented());
           return boost::iostreams::zlib::no_compression;
       }
   }
+#  endif
 
   /**
    * Do a zlib compression followed by a base64 encoding of the given data. The
@@ -2727,11 +2727,11 @@ namespace DataOutBase
 
 
 
-  VtkFlags::VtkFlags(const double                         time,
-                     const unsigned int                   cycle,
-                     const bool                           print_date_and_time,
-                     const VtkFlags::ZlibCompressionLevel compression_level,
-                     const bool write_higher_order_cells,
+  VtkFlags::VtkFlags(const double           time,
+                     const unsigned int     cycle,
+                     const bool             print_date_and_time,
+                     const CompressionLevel compression_level,
+                     const bool             write_higher_order_cells,
                      const std::map<std::string, std::string> &physical_units)
     : time(time)
     , cycle(cycle)
@@ -7602,11 +7602,11 @@ namespace DataOutBase
                  unsigned int,
                  std::string,
                  DataComponentInterpretation::DataComponentInterpretation>>
-      &                                  nonscalar_data_ranges,
-    const Deal_II_IntermediateFlags &    flags,
-    const std::string &                  filename,
-    const MPI_Comm &                     comm,
-    const VtkFlags::ZlibCompressionLevel compression)
+      &                              nonscalar_data_ranges,
+    const Deal_II_IntermediateFlags &flags,
+    const std::string &              filename,
+    const MPI_Comm &                 comm,
+    const CompressionLevel           compression)
   {
 #ifndef DEAL_II_WITH_MPI
     (void)patches;
@@ -7641,7 +7641,7 @@ namespace DataOutBase
     {
       boost::iostreams::filtering_ostream f;
 
-      if (compression != VtkFlags::no_compression)
+      if (compression != no_compression)
 #  ifdef DEAL_II_WITH_ZLIB
         f.push(boost::iostreams::zlib_compressor(
           get_boost_zlib_compression_level(compression)));
@@ -8141,9 +8141,9 @@ DataOutInterface<dim, spacedim>::write_deal_II_intermediate(
 template <int dim, int spacedim>
 void
 DataOutInterface<dim, spacedim>::write_deal_II_intermediate_in_parallel(
-  const std::string &                               filename,
-  const MPI_Comm &                                  comm,
-  const DataOutBase::VtkFlags::ZlibCompressionLevel compression) const
+  const std::string &                 filename,
+  const MPI_Comm &                    comm,
+  const DataOutBase::CompressionLevel compression) const
 {
   DataOutBase::write_deal_II_intermediate_in_parallel(
     get_patches(),
@@ -9335,8 +9335,8 @@ DataOutReader<dim, spacedim>::read_whole_parallel_file(std::istream &in)
       in.read(temp_buffer.data(), chunk_sizes[n]);
 
       boost::iostreams::filtering_istreambuf f;
-      if (static_cast<DataOutBase::VtkFlags::ZlibCompressionLevel>(
-            header.compression) != DataOutBase::VtkFlags::no_compression)
+      if (static_cast<DataOutBase::CompressionLevel>(header.compression) !=
+          DataOutBase::no_compression)
 #ifdef DEAL_II_WITH_ZLIB
         f.push(boost::iostreams::zlib_decompressor());
 #else
index 56730ced65ca32aa096a79464bdf354aa08238d5..35077900a8fb0818d429a11ad5e31a6d71fc63d3 100644 (file)
@@ -213,11 +213,11 @@ for (deal_II_dimension : OUTPUT_DIMENSIONS;
                      unsigned int,
                      std::string,
                      DataComponentInterpretation::DataComponentInterpretation>>
-          &                                  nonscalar_data_ranges,
-        const Deal_II_IntermediateFlags &    flags,
-        const std::string &                  filename,
-        const MPI_Comm &                     comm,
-        const VtkFlags::ZlibCompressionLevel compression);
+          &                              nonscalar_data_ranges,
+        const Deal_II_IntermediateFlags &flags,
+        const std::string &              filename,
+        const MPI_Comm &                 comm,
+        const CompressionLevel           compression);
 
       template void
       write_hdf5_parallel(
index 7102155424662fbdb511385d67fcf958ac31d1ed..5ad1b472982cea828c8069ec023b2c7f8a62529c 100644 (file)
@@ -312,8 +312,7 @@ namespace Step23
     const std::string filename =
       "solution-" + Utilities::int_to_string(timestep_number, 3) + ".vtu";
     DataOutBase::VtkFlags vtk_flags;
-    vtk_flags.compression_level =
-      DataOutBase::VtkFlags::ZlibCompressionLevel::best_speed;
+    vtk_flags.compression_level = DataOutBase::best_speed;
     data_out.set_flags(vtk_flags);
     std::ofstream output(filename);
     data_out.write_vtu(output);

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.