]> https://gitweb.dealii.org/ - dealii.git/commitdiff
fixed bug in write_mesh_per_processor_as_vtu, added test 3952/head
authortcclevenger <tcleven@clemson.edu>
Thu, 9 Feb 2017 13:38:27 +0000 (08:38 -0500)
committertcclevenger <tcleven@clemson.edu>
Thu, 9 Feb 2017 15:54:27 +0000 (10:54 -0500)
Previous implementation did not account for non-distributed multilevel
hierarchy or active cells having different level_subdomain_id and
subdomain_id.

made requested changes

typo in changes. GridOut function, not GridTools

doc/news/changes/minor/20170209ConradClevenger [new file with mode: 0644]
source/grid/grid_out.cc
tests/grid/grid_out_per_processor_vtu_02.cc [new file with mode: 0644]
tests/grid/grid_out_per_processor_vtu_02.mpirun=2.with_zlib=on.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20170209ConradClevenger b/doc/news/changes/minor/20170209ConradClevenger
new file mode 100644 (file)
index 0000000..b92ce51
--- /dev/null
@@ -0,0 +1,3 @@
+Fixed: The GridOut::write_mesh_per_processor_as_vtu() function now works for a mesh whose multilevel hierarchy is not distributed, as well as a mesh whose level_subdomain_ids do not necessarily match its subdomain_ids for every active cell.  
+<br>
+(Conrad Clevenger, 2017/02/09)
index e8eb42d96908916aa6a38d7adb82bac66829f0fd..00d3fc68ab529781b38078a9907adecf3b0ad395 100644 (file)
@@ -2548,11 +2548,24 @@ void GridOut::write_mesh_per_processor_as_vtu (const Triangulation<dim,spacedim>
   for (cell=tria.begin(), endc=tria.end();
        cell != endc; ++cell)
     {
-      if (!include_artificial && cell->level_subdomain_id() ==
-          numbers::artificial_subdomain_id)
-        continue;
-      if (!view_levels && cell->has_children())
-        continue;
+      if (!view_levels)
+        {
+          if (cell->has_children())
+            continue;
+          if (!include_artificial &&
+              cell->subdomain_id() == numbers::artificial_subdomain_id)
+            continue;
+        }
+      else if (!include_artificial)
+        {
+          if (cell->has_children() &&
+              cell->level_subdomain_id() == numbers::artificial_subdomain_id)
+            continue;
+          else if (!cell->has_children() &&
+                   cell->level_subdomain_id() == numbers::artificial_subdomain_id &&
+                   cell->subdomain_id() == numbers::artificial_subdomain_id)
+            continue;
+        }
 
       DataOutBase::Patch<dim,spacedim> patch;
       patch.data.reinit(n_datasets, n_q_points);
diff --git a/tests/grid/grid_out_per_processor_vtu_02.cc b/tests/grid/grid_out_per_processor_vtu_02.cc
new file mode 100644 (file)
index 0000000..d425450
--- /dev/null
@@ -0,0 +1,103 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2009 - 2016 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 at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+// check GriOut::write_mesh_per_processor_as_vtu() when level_subdomain_id
+// differs from subdomain_id for a particular cell
+
+#include "../tests.h"
+#include <deal.II/base/logstream.h>
+#include <deal.II/base/tensor.h>
+#include <deal.II/base/geometry_info.h>
+#include <deal.II/grid/tria.h>
+#include <deal.II/distributed/tria.h>
+#include <deal.II/distributed/shared_tria.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_out.h>
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/dofs/dof_tools.h>
+#include <deal.II/grid/tria_accessor.h>
+#include <deal.II/grid/tria_iterator.h>
+#include <deal.II/dofs/dof_accessor.h>
+#include <deal.II/grid/grid_out.h>
+#include <deal.II/numerics/data_out.h>
+#include <deal.II/fe/fe_q.h>
+
+#include <fstream>
+
+template<int dim>
+void output(const parallel::shared::Triangulation<dim> &tr,
+            const std::string                          &filename,
+            const bool                                 view_levels,
+            const bool                                 include_artificial)
+{
+  GridOut out;
+  out.write_mesh_per_processor_as_vtu(tr, filename, view_levels, include_artificial);
+
+  // copy the .pvtu and .vtu files
+  // into the logstream
+  int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
+  if (myid==0)
+    {
+      cat_file((std::string(filename) + ".pvtu").c_str());
+      cat_file((std::string(filename) + ".proc0000.vtu").c_str());
+    }
+  else if (myid==1)
+    cat_file((std::string(filename) + ".proc0001.vtu").c_str());
+  else
+    AssertThrow(false, ExcNotImplemented());
+}
+
+template<int dim>
+void test()
+{
+  unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
+  if (myid == 0)
+    deallog << "hyper_cube" << std::endl;
+
+  parallel::shared::Triangulation<dim> tr(MPI_COMM_WORLD);
+  GridGenerator::hyper_cube(tr);
+  tr.refine_global(1);
+  typename Triangulation<dim>::active_cell_iterator
+  cell=tr.begin_active(), endc=tr.end();
+  for (; cell!=endc; ++cell)
+    {
+      if (cell->index() < 2)
+        cell->set_subdomain_id(cell->index());
+      else
+        cell->set_subdomain_id(numbers::artificial_subdomain_id);
+
+      if (cell->index() == 0 || cell->index() == 2)
+        cell->set_level_subdomain_id(numbers::artificial_subdomain_id);
+      else if (cell->index() == 1)
+        cell->set_level_subdomain_id(0);
+      else if (cell->index() == 3)
+        cell->set_level_subdomain_id(1);
+    }
+
+  output(tr, "file1", true, false);
+}
+
+
+int main(int argc, char *argv[])
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization (argc, argv, 1);
+
+  MPILogInitAll init;
+
+  deallog.push("2d");
+  test<2>();
+  deallog.pop();
+}
diff --git a/tests/grid/grid_out_per_processor_vtu_02.mpirun=2.with_zlib=on.output b/tests/grid/grid_out_per_processor_vtu_02.mpirun=2.with_zlib=on.output
new file mode 100644 (file)
index 0000000..4303e9a
--- /dev/null
@@ -0,0 +1,100 @@
+
+DEAL:0:2d::hyper_cube
+<?xml version="1.0"?>
+<!--
+#This file was generated 
+-->
+<VTKFile type="PUnstructuredGrid" version="0.1" byte_order="LittleEndian">
+  <PUnstructuredGrid GhostLevel="0">
+    <PPointData Scalars="scalars">
+    <PDataArray type="Float64" Name="subdomain" format="ascii"/>
+    <PDataArray type="Float64" Name="proc_writing" format="ascii"/>
+    </PPointData>
+    <PPoints>
+      <PDataArray type="Float64" NumberOfComponents="3"/>
+    </PPoints>
+    <Piece Source="file1.proc0000.vtu"/>
+    <Piece Source="file1.proc0001.vtu"/>
+  </PUnstructuredGrid>
+</VTKFile>
+
+<?xml version="1.0" ?> 
+<!-- 
+# vtk DataFile Version 3.0
+#This file was generated 
+-->
+<VTKFile type="UnstructuredGrid" version="0.1" compressor="vtkZLibDataCompressor" byte_order="LittleEndian">
+<UnstructuredGrid>
+<Piece NumberOfPoints="16" NumberOfCells="4" >
+  <Points>
+    <DataArray type="Float64" NumberOfComponents="3" format="binary">
+AQAAAIABAACAAQAALQAAAA==eNpjYMAHPtgzkCQP4xPSBwMP7EmTh/FxiRPrflzmfLAnzV5C6gmHDwBXYxlL
+    </DataArray>
+  </Points>
+
+  <Cells>
+    <DataArray type="Int32" Name="connectivity" format="binary">
+AQAAAEAAAABAAAAAKAAAAA==eNoNwwkSwBAQALB1VIvi/7+VzCQiIlnMVh9fm5/d6fB3edxeCvAAeQ==
+    </DataArray>
+    <DataArray type="Int32" Name="offsets" format="binary">
+AQAAABAAAAAQAAAAEwAAAA==eNpjYWBg4ABiHiAWAGIAAVAAKQ==
+    </DataArray>
+    <DataArray type="UInt8" Name="types" format="binary">
+AQAAAAQAAAAEAAAADAAAAA==eNrj5OTkBAAAXgAl
+    </DataArray>
+  </Cells>
+  <PointData Scalars="scalars">
+    <DataArray type="Float64" Name="level" format="binary">
+AQAAAIAAAACAAAAAEAAAAA==eNpjYCAGfLCnFQ0AiwIONQ==    </DataArray>
+    <DataArray type="Float64" Name="subdomain" format="binary">
+AQAAAIAAAACAAAAAGgAAAA==eNpjYACBD/sZ8NKEwAd7/DTDAVw0AOJFDnk=    </DataArray>
+    <DataArray type="Float64" Name="level_subdomain" format="binary">
+AQAAAIAAAACAAAAAFAAAAA==eNpjYCAKHCBAEwAf7HHRACjbB70=    </DataArray>
+    <DataArray type="Float64" Name="proc_writing" format="binary">
+AQAAAIAAAACAAAAADAAAAA==eNpjYBhYAAAAgAAB    </DataArray>
+  </PointData>
+ </Piece>
+ </UnstructuredGrid>
+</VTKFile>
+
+
+<?xml version="1.0" ?> 
+<!-- 
+# vtk DataFile Version 3.0
+#This file was generated 
+-->
+<VTKFile type="UnstructuredGrid" version="0.1" compressor="vtkZLibDataCompressor" byte_order="LittleEndian">
+<UnstructuredGrid>
+<Piece NumberOfPoints="16" NumberOfCells="4" >
+  <Points>
+    <DataArray type="Float64" NumberOfComponents="3" format="binary">
+AQAAAIABAACAAQAALQAAAA==eNpjYMAHPtgzkCQP4xPSBwMP7EmTh/FxiRPrflzmfLAnzV5C6gmHDwBXYxlL
+    </DataArray>
+  </Points>
+
+  <Cells>
+    <DataArray type="Int32" Name="connectivity" format="binary">
+AQAAAEAAAABAAAAAKAAAAA==eNoNwwkSwBAQALB1VIvi/7+VzCQiIlnMVh9fm5/d6fB3edxeCvAAeQ==
+    </DataArray>
+    <DataArray type="Int32" Name="offsets" format="binary">
+AQAAABAAAAAQAAAAEwAAAA==eNpjYWBg4ABiHiAWAGIAAVAAKQ==
+    </DataArray>
+    <DataArray type="UInt8" Name="types" format="binary">
+AQAAAAQAAAAEAAAADAAAAA==eNrj5OTkBAAAXgAl
+    </DataArray>
+  </Cells>
+  <PointData Scalars="scalars">
+    <DataArray type="Float64" Name="level" format="binary">
+AQAAAIAAAACAAAAAEAAAAA==eNpjYCAGfLCnFQ0AiwIONQ==    </DataArray>
+    <DataArray type="Float64" Name="subdomain" format="binary">
+AQAAAIAAAACAAAAAGgAAAA==eNpjYACBD/sZ8NKEwAd7/DTDAVw0AOJFDnk=    </DataArray>
+    <DataArray type="Float64" Name="level_subdomain" format="binary">
+AQAAAIAAAACAAAAAFAAAAA==eNpjYCAKHCBAEwAf7HHRACjbB70=    </DataArray>
+    <DataArray type="Float64" Name="proc_writing" format="binary">
+AQAAAIAAAACAAAAAEQAAAA==eNpjYACBD/YMA0QDAJLsEvE=    </DataArray>
+  </PointData>
+ </Piece>
+ </UnstructuredGrid>
+</VTKFile>
+
+

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.