]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Change default compression_level of vtu output 14983/head
authorChristoph Schmidt <christoph.schmidt@tum.de>
Tue, 28 Mar 2023 15:34:27 +0000 (17:34 +0200)
committerChristoph Schmidt <christoph.schmidt@tum.de>
Sat, 1 Apr 2023 09:29:36 +0000 (11:29 +0200)
from best_compression to best_speed, as discussed in pull request #14958.
The compression flags of the tests is set to best_compression to conserve
the current output behavior of the tests.

20 files changed:
doc/news/changes/minor/20230328Schmidt [new file with mode: 0644]
include/deal.II/base/data_out_base.h
source/grid/grid_out.cc
tests/data_out/data_out_12.cc
tests/data_out/data_out_13.cc
tests/data_out/data_out_base_vtu_02.cc
tests/data_out/data_out_base_vtu_04.cc
tests/data_out/data_out_base_vtu_05.cc
tests/data_out/data_out_postprocessor_tensor_02.cc
tests/grid/grid_out_07.cc
tests/grid/grid_out_08.cc
tests/grid/grid_out_per_processor_vtu_01.cc
tests/grid/grid_out_per_processor_vtu_02.cc
tests/grid/grid_out_per_processor_vtu_03.cc
tests/hp/solution_transfer_03.cc
tests/hp/solution_transfer_05.cc
tests/mpi/codim_01.cc
tests/mpi/parallel_vtu_01.cc
tests/simplex/data_out_write_vtu_01.cc
tests/simplex/write_mesh_per_processor_as_vtu_01.cc

diff --git a/doc/news/changes/minor/20230328Schmidt b/doc/news/changes/minor/20230328Schmidt
new file mode 100644 (file)
index 0000000..f071317
--- /dev/null
@@ -0,0 +1,3 @@
+Changed: The default `compression_level` for vtu output is changed from `best_compression` to `best_speed`.
+<br>
+(Christoph Schmidt, 2023/03/28)
index 5479903734369211009241107b15e8406ce1bfcc..a28ee2e52708934b1a55a09dd5627fcb32bb932a 100644 (file)
@@ -1188,7 +1188,7 @@ namespace DataOutBase
 
     /**
      * Flag determining the compression level at which zlib, if available, is
-     * run. The default is <tt>best_compression</tt>.
+     * run. The default is <tt>best_speed</tt>.
      */
     DataOutBase::CompressionLevel compression_level;
 
@@ -1242,9 +1242,8 @@ namespace DataOutBase
       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 =
-        CompressionLevel::best_compression,
-      const bool write_higher_order_cells                      = false,
+      const CompressionLevel compression_level   = CompressionLevel::best_speed,
+      const bool             write_higher_order_cells          = false,
       const std::map<std::string, std::string> &physical_units = {});
   };
 
index 7dcfcb56581067b31e8bacb2f82bf578be872a18..0525c8ddafe976c8ed533328f105e33aa281ef9d 100644 (file)
@@ -3718,9 +3718,9 @@ GridOut::write_mesh_per_processor_as_vtu(
                unsigned int,
                std::string,
                DataComponentInterpretation::DataComponentInterpretation>>
-                        vector_data_ranges;
-  DataOutBase::VtkFlags flags;
-  DataOutBase::write_vtu(patches, data_names, vector_data_ranges, flags, out);
+    vector_data_ranges;
+  DataOutBase::write_vtu(
+    patches, data_names, vector_data_ranges, vtu_flags, out);
 }
 
 
index e9968b9180cadcd550dce2d0c7270763d614fb2f..179956783f2710d07e4697a81f366478e592517a 100644 (file)
@@ -166,6 +166,9 @@ check()
   std::vector<std::string> component_name(dim + dim * dim, "tensor");
   std::fill(component_name.begin(), component_name.begin() + dim, "vector");
 
+  DataOutBase::VtkFlags vtk_flags;
+  vtk_flags.compression_level = DataOutBase::CompressionLevel::best_compression;
+
   DataOut<dim> data_out;
   data_out.attach_dof_handler(dof_handler);
   data_out.add_data_vector(v,
@@ -173,6 +176,7 @@ check()
                            DataOut<dim>::type_dof_data,
                            data_component_interpretation);
   data_out.build_patches();
+  data_out.set_flags(vtk_flags);
 
   /*
   std::ofstream out("output_" + Utilities::int_to_string(dim) + "d.vtu");
index 4ca79d3ead7150ad588dd7e4707f6aef177c02a9..174b79b348111f850cad47650a4457e0fea69c58 100644 (file)
@@ -166,6 +166,9 @@ check()
   std::vector<std::string> component_name(dim + dim * dim, "tensor");
   std::fill(component_name.begin(), component_name.begin() + dim, "vector");
 
+  DataOutBase::VtkFlags vtk_flags;
+  vtk_flags.compression_level = DataOutBase::CompressionLevel::best_compression;
+
   DataOut<dim> data_out;
   data_out.attach_dof_handler(dof_handler);
   data_out.add_data_vector(v,
@@ -173,6 +176,7 @@ check()
                            DataOut<dim>::type_dof_data,
                            data_component_interpretation);
   data_out.build_patches();
+  data_out.set_flags(vtk_flags);
 
   std::vector<std::string> filenames;
   filenames.push_back("output_" + Utilities::int_to_string(dim) + "d.vtu");
index 7ea1bcbbd7dec1ea9d8257fd652a79f56e8736b8..2c8b29990b21b7cfdab83cc41fbbb46fef146278 100644 (file)
@@ -67,6 +67,7 @@ check_all(std::ostream &log)
 
   char                  name[100];
   DataOutBase::VtkFlags flags;
+  flags.compression_level = DataOutBase::CompressionLevel::best_compression;
   if (true)
     {
       sprintf(name, "%d%d.vtu", dim, spacedim);
index 2ead86a3d9c1853c528797c8b8836029a1e6aef0..f54b320690d69235b964e2f50ef8990421f1db29 100644 (file)
@@ -61,6 +61,7 @@ check(std::ostream &log, unsigned cell_order)
 
   DataOutBase::VtkFlags flags;
   flags.write_higher_order_cells = true;
+  flags.compression_level = DataOutBase::CompressionLevel::best_compression;
 
   DataOut<dim> data_out;
   data_out.set_flags(flags);
index 7d15e30d4c96fd8149a2b4c649c1c85801c2cd3f..ab387484b5b29ccf19d17bd3f374a2ed5c5c511d 100644 (file)
@@ -53,6 +53,7 @@ check(std::ostream &log, unsigned cell_order)
 
   DataOutBase::VtkFlags flags;
   flags.write_higher_order_cells = true;
+  flags.compression_level = DataOutBase::CompressionLevel::best_compression;
 
   DataOut<dim> data_out;
   data_out.set_flags(flags);
index 163157ad38b81b5386adb3759f0461e14683300a..1b2ea58cd14cdb1daed9cefae9a4786ede47acfb 100644 (file)
@@ -362,6 +362,10 @@ namespace Step8
 
     StrainPostprocessor<dim> grad_u;
 
+    DataOutBase::VtkFlags vtk_flags;
+    vtk_flags.compression_level =
+      DataOutBase::CompressionLevel::best_compression;
+
     DataOut<dim> data_out;
     data_out.attach_dof_handler(dof_handler);
 
@@ -376,6 +380,7 @@ namespace Step8
                              data_component_interpretation);
     data_out.add_data_vector(solution, grad_u);
     data_out.build_patches();
+    data_out.set_flags(vtk_flags);
     data_out.write_vtu(deallog.get_file_stream());
   }
 
index b4131967c1de342ac41293a14cd443797c5951c1..2b666eb3cb4fb91b6ac278c6b0b9602d683f63a1 100644 (file)
@@ -44,7 +44,10 @@ test(std::ostream &logfile)
   tria.begin_active()->set_refine_flag();
   tria.execute_coarsening_and_refinement();
 
-  GridOut grid_out;
+  GridOut           grid_out;
+  GridOutFlags::Vtu vtu_flags;
+  vtu_flags.compression_level = DataOutBase::CompressionLevel::best_compression;
+  grid_out.set_flags(vtu_flags);
   grid_out.write_vtu(tria, logfile);
 }
 
index f98eaa571d66d2855fc0a5a2269a88d2f6e71d11..813978b9b27c5ca83b63bdfbf1745b5378803656 100644 (file)
@@ -71,7 +71,10 @@ test()
 
   tria.refine_global(1);
 
-  GridOut grid_out;
+  GridOut           grid_out;
+  GridOutFlags::Vtu vtu_flags;
+  vtu_flags.compression_level = DataOutBase::CompressionLevel::best_compression;
+  grid_out.set_flags(vtu_flags);
   grid_out.write_vtu(tria, deallog.get_file_stream());
 }
 
index fd68938a5e134218c6677564baadca1c1458143b..52622a5ddad1e716f4e0ac1d3f66bd4a504f832d 100644 (file)
@@ -46,7 +46,10 @@ output(const parallel::distributed::Triangulation<dim> &tr,
        const bool                                       view_levels,
        const bool                                       include_artificial)
 {
-  GridOut out;
+  GridOut           out;
+  GridOutFlags::Vtu vtu_flags;
+  vtu_flags.compression_level = DataOutBase::CompressionLevel::best_compression;
+  out.set_flags(vtu_flags);
   out.write_mesh_per_processor_as_vtu(tr,
                                       filename,
                                       view_levels,
index 699db930bd928c250460d2594f3190dd96f0f2e7..020bce3e788a9be60f4bcfe8fbbd921f83282426 100644 (file)
@@ -47,7 +47,10 @@ output(const parallel::shared::Triangulation<dim> &tr,
        const bool                                  view_levels,
        const bool                                  include_artificial)
 {
-  GridOut out;
+  GridOut           out;
+  GridOutFlags::Vtu vtu_flags;
+  vtu_flags.compression_level = DataOutBase::CompressionLevel::best_compression;
+  out.set_flags(vtu_flags);
   out.write_mesh_per_processor_as_vtu(tr,
                                       filename,
                                       view_levels,
index aa3f9d85126abaeb60136cc86ba9964441d2719a..e82926272024f4e11174d394b83f4e741f4ac7fd 100644 (file)
@@ -49,8 +49,11 @@ test()
   GridGenerator::hyper_cube(tr);
   tr.refine_global(3);
 
-  std::string filename = "file" + Utilities::int_to_string(dim);
-  GridOut     grid_out;
+  std::string       filename = "file" + Utilities::int_to_string(dim);
+  GridOut           grid_out;
+  GridOutFlags::Vtu vtu_flags;
+  vtu_flags.compression_level = DataOutBase::CompressionLevel::best_compression;
+  grid_out.set_flags(vtu_flags);
   grid_out.write_mesh_per_processor_as_vtu(tr, filename, true);
 
   cat_file((std::string(filename) + ".vtu").c_str());
index acbc107c5428aa69aaeebbe8c96f46f05a66737e..fc0a09615e7ecee1620d49adcc2563e584a87822 100644 (file)
@@ -78,12 +78,16 @@ main()
   Vector<double> solution(dof_handler.n_dofs());
   solution = 1.0;
 
+  // Define compression level for output data
+  DataOutBase::VtkFlags vtk_flags;
+  vtk_flags.compression_level = DataOutBase::CompressionLevel::best_compression;
 
   // Save output
   DataOut<2> data_out;
   data_out.attach_dof_handler(dof_handler);
   data_out.add_data_vector(solution, "Solution");
   data_out.build_patches();
+  data_out.set_flags(vtk_flags);
   data_out.write_vtu(deallog.get_file_stream());
 
 
@@ -109,5 +113,6 @@ main()
   data_out2.attach_dof_handler(dof_handler);
   data_out2.add_data_vector(new_solution, "Solution");
   data_out2.build_patches();
+  data_out2.set_flags(vtk_flags);
   data_out2.write_vtu(deallog.get_file_stream());
 }
index 48d1acedb114ad0b9bbb4b5a45dc0b512dcdd7e9..b7c49e8b77ee93403569fbd33b5b748b9130109f 100644 (file)
@@ -88,12 +88,17 @@ main()
   Vector<double> new_solution(dof_handler.n_dofs());
   solultion_trans.interpolate(solution, new_solution);
 
+  // Define compression level for output data
+  DataOutBase::VtkFlags vtk_flags;
+  vtk_flags.compression_level = DataOutBase::CompressionLevel::best_compression;
+
   // a follow-up error to the one fixed with _04 was that DataOut also got
   // itself confused
   DataOut<2> data_out2;
   data_out2.attach_dof_handler(dof_handler);
   data_out2.add_data_vector(new_solution, "Solution");
   data_out2.build_patches();
+  data_out2.set_flags(vtk_flags);
   data_out2.write_vtu(deallog.get_file_stream());
 
   // we are good if we made it to here
index 10dd0a8653cf2645ca7491994e346973f11f60b3..22cb4964fade3def6ae3234b61f42748f6b486f0 100644 (file)
@@ -89,12 +89,16 @@ test(std::ostream & /*out*/)
   dh.distribute_dofs(fe);
   deallog << "dofs " << dh.n_dofs() << std::endl;
 
+  DataOutBase::VtkFlags vtk_flags;
+  vtk_flags.compression_level = DataOutBase::CompressionLevel::best_compression;
+
   DataOut<dim, spacedim> data_out;
   data_out.attach_triangulation(tr);
   Vector<float> subdomain(tr.n_active_cells());
   for (unsigned int i = 0; i < subdomain.size(); ++i)
     subdomain(i) = tr.locally_owned_subdomain();
   data_out.add_data_vector(subdomain, "subdomain");
+  data_out.set_flags(vtk_flags);
 
   std::string name = "f0.vtu";
   name[1] += tr.locally_owned_subdomain();
index 01ebc5030848bfc8b127ac339083739349f23b56..77f08c7f2d012e00c775016aa5ff2ada2a30b8fb 100644 (file)
@@ -64,10 +64,14 @@ test()
   x.reinit(dofh.locally_owned_dofs(), MPI_COMM_WORLD);
   x = 2.0;
 
+  DataOutBase::VtkFlags vtk_flags;
+  vtk_flags.compression_level = DataOutBase::CompressionLevel::best_compression;
+
   DataOut<dim> data_out;
   data_out.attach_dof_handler(dofh);
   data_out.add_data_vector(x, "x");
   data_out.build_patches();
+  data_out.set_flags(vtk_flags);
 
   data_out.write_vtu_in_parallel("output.vtu", MPI_COMM_WORLD);
   MPI_Barrier(MPI_COMM_WORLD);
index b1c67611960b48e3461a2ac5509aaa14045e9801..813518fb15c238720fdab855a029f366b44cd41f 100644 (file)
@@ -73,12 +73,16 @@ test(const FiniteElement<dim, spacedim> &fe, const unsigned int n_components)
                            RightHandSideFunction<dim>(n_components),
                            solution);
 
+  DataOutBase::VtkFlags vtk_flags;
+  vtk_flags.compression_level = DataOutBase::CompressionLevel::best_compression;
+
   for (unsigned int n_subdivisions = 1; n_subdivisions <= 2; ++n_subdivisions)
     {
       DataOut<dim> data_out;
 
       data_out.attach_dof_handler(dof_handler);
       data_out.add_data_vector(solution, "solution");
+      data_out.set_flags(vtk_flags);
 
 
       data_out.build_patches(mapping, n_subdivisions);
index d8babf853904409ad4e3e8c3cc572540f36dec52..7e04942a3f66f1413a424632fed129835279e950 100644 (file)
@@ -42,7 +42,10 @@ output(const Triangulation<dim> &tr,
        const bool                view_levels,
        const bool                include_artificial)
 {
-  GridOut out;
+  GridOut           out;
+  GridOutFlags::Vtu vtu_flags;
+  vtu_flags.compression_level = DataOutBase::CompressionLevel::best_compression;
+  out.set_flags(vtu_flags);
   out.write_mesh_per_processor_as_vtu(tr,
                                       filename,
                                       view_levels,

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.