]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
DataOut was not prepared to deal with processors that have no locally owned cells...
authorbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 16 Dec 2010 18:56:32 +0000 (18:56 +0000)
committerbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 16 Dec 2010 18:56:32 +0000 (18:56 +0000)
git-svn-id: https://svn.dealii.org/trunk@22976 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/source/base/data_out_base.cc
deal.II/source/numerics/data_out.cc
tests/mpi/p4est_data_out_01.cc [new file with mode: 0644]
tests/mpi/p4est_data_out_01/ncpu_10/cmp/generic [new file with mode: 0644]
tests/mpi/p4est_data_out_01/ncpu_4/cmp/generic [new file with mode: 0644]

index 2ae5a1074ad3ad5e0bca5431b8be14cddaac05fc..aefae526a6b3861dc0269d6582a9bcb062916404 100644 (file)
@@ -1909,8 +1909,26 @@ void DataOutBase::write_ucd (const std::vector<Patch<dim,spacedim> > &patches,
 {
   AssertThrow (out, ExcIO());
 
+#ifndef DEAL_II_COMPILER_SUPPORTS_MPI
+                                  // verify that there are indeed
+                                  // patches to be written out. most
+                                  // of the times, people just forget
+                                  // to call build_patches when there
+                                  // are no patches, so a warning is
+                                  // in order. that said, the
+                                  // assertion is disabled if we
+                                  // support MPI since then it can
+                                  // happen that on the coarsest
+                                  // mesh, a processor simply has no
+                                  // cells it actually owns, and in
+                                  // that case it is legit if there
+                                  // are no patches
   Assert (patches.size() > 0, ExcNoPatches());
-
+#else
+  if (patches.size() == 0)
+    return;
+#endif
+  
   const unsigned int n_data_sets = data_names.size();
 
   UcdStream ucd_out(out, flags);
@@ -1990,7 +2008,25 @@ void DataOutBase::write_dx (const std::vector<Patch<dim,spacedim> > &patches,
 {
   AssertThrow (out, ExcIO());
 
+#ifndef DEAL_II_COMPILER_SUPPORTS_MPI
+                                  // verify that there are indeed
+                                  // patches to be written out. most
+                                  // of the times, people just forget
+                                  // to call build_patches when there
+                                  // are no patches, so a warning is
+                                  // in order. that said, the
+                                  // assertion is disabled if we
+                                  // support MPI since then it can
+                                  // happen that on the coarsest
+                                  // mesh, a processor simply has no
+                                  // cells it actually owns, and in
+                                  // that case it is legit if there
+                                  // are no patches
   Assert (patches.size() > 0, ExcNoPatches());
+#else
+  if (patches.size() == 0)
+    return;
+#endif
                                   // Stream with special features for dx output
   DXStream dx_out(out, flags);
 
@@ -2256,8 +2292,26 @@ void DataOutBase::write_gnuplot (const std::vector<Patch<dim,spacedim> > &patche
 {
   AssertThrow (out, ExcIO());
 
+#ifndef DEAL_II_COMPILER_SUPPORTS_MPI
+                                  // verify that there are indeed
+                                  // patches to be written out. most
+                                  // of the times, people just forget
+                                  // to call build_patches when there
+                                  // are no patches, so a warning is
+                                  // in order. that said, the
+                                  // assertion is disabled if we
+                                  // support MPI since then it can
+                                  // happen that on the coarsest
+                                  // mesh, a processor simply has no
+                                  // cells it actually owns, and in
+                                  // that case it is legit if there
+                                  // are no patches
   Assert (patches.size() > 0, ExcNoPatches());
-
+#else
+  if (patches.size() == 0)
+    return;
+#endif
+  
   const unsigned int n_data_sets = data_names.size();
 
                                   // write preamble
@@ -2473,7 +2527,25 @@ void DataOutBase::write_povray (const std::vector<Patch<dim,spacedim> > &patches
 {
   AssertThrow (out, ExcIO());
 
+#ifndef DEAL_II_COMPILER_SUPPORTS_MPI
+                                  // verify that there are indeed
+                                  // patches to be written out. most
+                                  // of the times, people just forget
+                                  // to call build_patches when there
+                                  // are no patches, so a warning is
+                                  // in order. that said, the
+                                  // assertion is disabled if we
+                                  // support MPI since then it can
+                                  // happen that on the coarsest
+                                  // mesh, a processor simply has no
+                                  // cells it actually owns, and in
+                                  // that case it is legit if there
+                                  // are no patches
   Assert (patches.size() > 0, ExcNoPatches());
+#else
+  if (patches.size() == 0)
+    return;
+#endif
   Assert (dim==2, ExcNotImplemented());        // only for 2-D surfaces on a 2-D plane
   Assert (spacedim==2, ExcNotImplemented());
 
@@ -2819,8 +2891,26 @@ void DataOutBase::write_eps (const std::vector<Patch<dim,spacedim> > &patches,
 {
   Assert (out, ExcIO());
 
+#ifndef DEAL_II_COMPILER_SUPPORTS_MPI
+                                  // verify that there are indeed
+                                  // patches to be written out. most
+                                  // of the times, people just forget
+                                  // to call build_patches when there
+                                  // are no patches, so a warning is
+                                  // in order. that said, the
+                                  // assertion is disabled if we
+                                  // support MPI since then it can
+                                  // happen that on the coarsest
+                                  // mesh, a processor simply has no
+                                  // cells it actually owns, and in
+                                  // that case it is legit if there
+                                  // are no patches
   Assert (patches.size() > 0, ExcNoPatches());
-
+#else
+  if (patches.size() == 0)
+    return;
+#endif
+  
                                   // Do not allow volume rendering
   AssertThrow (dim==2, ExcNotImplemented());
 
@@ -3181,8 +3271,27 @@ void DataOutBase::write_gmv (const std::vector<Patch<dim,spacedim> > &patches,
 {
   Assert(dim<=3, ExcNotImplemented());
   AssertThrow (out, ExcIO());
+  
+#ifndef DEAL_II_COMPILER_SUPPORTS_MPI
+                                  // verify that there are indeed
+                                  // patches to be written out. most
+                                  // of the times, people just forget
+                                  // to call build_patches when there
+                                  // are no patches, so a warning is
+                                  // in order. that said, the
+                                  // assertion is disabled if we
+                                  // support MPI since then it can
+                                  // happen that on the coarsest
+                                  // mesh, a processor simply has no
+                                  // cells it actually owns, and in
+                                  // that case it is legit if there
+                                  // are no patches
   Assert (patches.size() > 0, ExcNoPatches());
-
+#else
+  if (patches.size() == 0)
+    return;
+#endif
+  
   GmvStream gmv_out(out, flags);
   const unsigned int n_data_sets = data_names.size();
                                   // check against # of data sets in
@@ -3317,8 +3426,26 @@ void DataOutBase::write_tecplot (const std::vector<Patch<dim,spacedim> > &patche
 {
   AssertThrow (out, ExcIO());
 
+#ifndef DEAL_II_COMPILER_SUPPORTS_MPI
+                                  // verify that there are indeed
+                                  // patches to be written out. most
+                                  // of the times, people just forget
+                                  // to call build_patches when there
+                                  // are no patches, so a warning is
+                                  // in order. that said, the
+                                  // assertion is disabled if we
+                                  // support MPI since then it can
+                                  // happen that on the coarsest
+                                  // mesh, a processor simply has no
+                                  // cells it actually owns, and in
+                                  // that case it is legit if there
+                                  // are no patches
   Assert (patches.size() > 0, ExcNoPatches());
-
+#else
+  if (patches.size() == 0)
+    return;
+#endif
+  
   TecplotStream tecplot_out(out, flags);
 
   const unsigned int n_data_sets = data_names.size();
@@ -3590,8 +3717,26 @@ void DataOutBase::write_tecplot_binary (const std::vector<Patch<dim,spacedim> >
 
   AssertThrow (out, ExcIO());
 
+#ifndef DEAL_II_COMPILER_SUPPORTS_MPI
+                                  // verify that there are indeed
+                                  // patches to be written out. most
+                                  // of the times, people just forget
+                                  // to call build_patches when there
+                                  // are no patches, so a warning is
+                                  // in order. that said, the
+                                  // assertion is disabled if we
+                                  // support MPI since then it can
+                                  // happen that on the coarsest
+                                  // mesh, a processor simply has no
+                                  // cells it actually owns, and in
+                                  // that case it is legit if there
+                                  // are no patches
   Assert (patches.size() > 0, ExcNoPatches());
-
+#else
+  if (patches.size() == 0)
+    return;
+#endif
+  
   const unsigned int n_data_sets = data_names.size();
                                   // check against # of data sets in
                                   // first patch. checks against all
@@ -3883,8 +4028,26 @@ DataOutBase::write_vtk (const std::vector<Patch<dim,spacedim> > &patches,
 {
   AssertThrow (out, ExcIO());
 
+#ifndef DEAL_II_COMPILER_SUPPORTS_MPI
+                                  // verify that there are indeed
+                                  // patches to be written out. most
+                                  // of the times, people just forget
+                                  // to call build_patches when there
+                                  // are no patches, so a warning is
+                                  // in order. that said, the
+                                  // assertion is disabled if we
+                                  // support MPI since then it can
+                                  // happen that on the coarsest
+                                  // mesh, a processor simply has no
+                                  // cells it actually owns, and in
+                                  // that case it is legit if there
+                                  // are no patches
   Assert (patches.size() > 0, ExcNoPatches());
-
+#else
+  if (patches.size() == 0)
+    return;
+#endif
+  
   VtkStream vtk_out(out, flags);
 
   const unsigned int n_data_sets = data_names.size();
@@ -4116,8 +4279,26 @@ DataOutBase::write_vtu (const std::vector<Patch<dim,spacedim> > &patches,
 {
   AssertThrow (out, ExcIO());
 
+#ifndef DEAL_II_COMPILER_SUPPORTS_MPI
+                                  // verify that there are indeed
+                                  // patches to be written out. most
+                                  // of the times, people just forget
+                                  // to call build_patches when there
+                                  // are no patches, so a warning is
+                                  // in order. that said, the
+                                  // assertion is disabled if we
+                                  // support MPI since then it can
+                                  // happen that on the coarsest
+                                  // mesh, a processor simply has no
+                                  // cells it actually owns, and in
+                                  // that case it is legit if there
+                                  // are no patches
   Assert (patches.size() > 0, ExcNoPatches());
-
+#else
+  if (patches.size() == 0)
+    return;
+#endif
+  
   VtuStream vtu_out(out, flags);
 
   const unsigned int n_data_sets = data_names.size();
@@ -5125,6 +5306,7 @@ merge (const DataOutReader<dim,spacedim> &source)
 {
   typedef typename dealii::DataOutBase::Patch<dim,spacedim> Patch;
 
+
   const std::vector<Patch> source_patches = source.get_patches ();
   Assert (patches.size () != 0,        ExcNoPatches ());
   Assert (source_patches.size () != 0, ExcNoPatches ());
index d72ff30abb18038e1a5acfd12abf792ddb6ffeae..b1467ebd1c219bb66bd1cd86405e34dafe9dd8c6 100644 (file)
@@ -983,13 +983,17 @@ void DataOut<dim,DH>::build_patches (const Mapping<DH::dimension,DH::space_dimen
   std::vector<std::pair<cell_iterator, unsigned int> > all_cells;
   {
                                     // set the index of the first
-                                    // cell. if first_locally_owned_cell/next_locally_owned_cell
+                                    // cell. if
+                                    // first_locally_owned_cell /
+                                    // next_locally_owned_cell
                                     // returns non-active cells, then
                                     // the index is not usable
                                     // anyway, but otherwise we
                                     // should keep track where we are
     unsigned int index;
-    if (first_locally_owned_cell()->has_children())
+    if ((first_locally_owned_cell() == this->dofs->end())
+       ||
+       (first_locally_owned_cell()->has_children()))
       index = 0;
     else
       index = std::distance (this->dofs->begin_active(),
diff --git a/tests/mpi/p4est_data_out_01.cc b/tests/mpi/p4est_data_out_01.cc
new file mode 100644 (file)
index 0000000..23c15fb
--- /dev/null
@@ -0,0 +1,107 @@
+//---------------------------------------------------------------------------
+//    $Id: 2d_coarse_grid_01.cc 17444 2008-10-31 19:35:14Z bangerth $
+//    Version: $Name$
+//
+//    Copyright (C) 2009, 2010 by the deal.II authors
+//
+//    This file is subject to QPL and may not be  distributed
+//    without copyright and license information. Please refer
+//    to the file deal.II/doc/license.html for the  text  and
+//    further information on this license.
+//
+//---------------------------------------------------------------------------
+
+
+// create a parallel DoFHandler and output data on a single
+// cell. DataOut was not prepared to handle situations where a
+// processor has no active cells at all.
+
+#include "../tests.h"
+#include "coarse_grid_common.h"
+#include <base/logstream.h>
+#include <base/tensor.h>
+#include <grid/tria.h>
+#include <distributed/tria.h>
+#include <grid/grid_generator.h>
+#include <grid/grid_out.h>
+#include <dofs/dof_handler.h>
+#include <dofs/dof_handler.h>
+#include <dofs/dof_tools.h>
+#include <grid/tria_accessor.h>
+#include <grid/tria_iterator.h>
+#include <dofs/dof_accessor.h>
+#include <numerics/data_out.h>
+
+#include <fe/fe_q.h>
+#include <lac/trilinos_vector.h>
+
+#include <fstream>
+
+
+template<int dim>
+void test()
+{
+  unsigned int myid = Utilities::System::get_this_mpi_process (MPI_COMM_WORLD);
+
+  if (Utilities::System::get_this_mpi_process (MPI_COMM_WORLD) == 0)
+    deallog << "hyper_cube" << std::endl;
+
+  parallel::distributed::Triangulation<dim> tr(MPI_COMM_WORLD);
+
+  GridGenerator::hyper_cube(tr);
+  DoFHandler<dim> dofh(tr);
+
+  static const FE_Q<dim> fe(2);
+  dofh.distribute_dofs (fe);
+
+  DataOut<dim> data_out;
+  data_out.attach_dof_handler (dofh);
+
+  TrilinosWrappers::MPI::Vector x;
+  x.reinit(dofh.locally_owned_dofs(), MPI_COMM_WORLD);
+  x=2.0;
+
+  data_out.add_data_vector (x, "x");
+  data_out.build_patches ();
+    
+  if (myid==0)
+    {
+      for (unsigned int i=0; i<dofh.n_locally_owned_dofs_per_processor().size(); ++i)
+       deallog << dofh.n_locally_owned_dofs_per_processor()[i] << std::endl;
+      data_out.write_vtu (deallog.get_file_stream());
+    }
+}
+
+
+int main(int argc, char *argv[])
+{
+#ifdef DEAL_II_COMPILER_SUPPORTS_MPI
+  MPI_Init (&argc,&argv);
+#else
+  (void)argc;
+  (void)argv;
+#endif
+
+  unsigned int myid = Utilities::System::get_this_mpi_process (MPI_COMM_WORLD);
+
+
+  deallog.push(Utilities::int_to_string(myid));
+
+  if (myid == 0)
+    {
+      std::ofstream logfile(output_file_for_mpi("p4est_data_out_01").c_str());
+      deallog.attach(logfile);
+      deallog.depth_console(0);
+      deallog.threshold_double(1.e-10);
+
+      deallog.push("2d");
+      test<2>();
+      deallog.pop();
+    }
+  else
+    test<2>();
+
+#ifdef DEAL_II_COMPILER_SUPPORTS_MPI
+  MPI_Finalize();
+#endif
+}
diff --git a/tests/mpi/p4est_data_out_01/ncpu_10/cmp/generic b/tests/mpi/p4est_data_out_01/ncpu_10/cmp/generic
new file mode 100644 (file)
index 0000000..3cdf9a9
--- /dev/null
@@ -0,0 +1,12 @@
+
+DEAL:0:2d::hyper_cube
+DEAL:0:2d::0
+DEAL:0:2d::0
+DEAL:0:2d::0
+DEAL:0:2d::0
+DEAL:0:2d::0
+DEAL:0:2d::0
+DEAL:0:2d::0
+DEAL:0:2d::0
+DEAL:0:2d::0
+DEAL:0:2d::9
diff --git a/tests/mpi/p4est_data_out_01/ncpu_4/cmp/generic b/tests/mpi/p4est_data_out_01/ncpu_4/cmp/generic
new file mode 100644 (file)
index 0000000..f884b80
--- /dev/null
@@ -0,0 +1,6 @@
+
+DEAL:0:2d::hyper_cube
+DEAL:0:2d::0
+DEAL:0:2d::0
+DEAL:0:2d::0
+DEAL:0:2d::9

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.