]> https://gitweb.dealii.org/ - dealii.git/commitdiff
add test for
authorTimo Heister <timo.heister@gmail.com>
Tue, 25 Sep 2012 17:52:19 +0000 (17:52 +0000)
committerTimo Heister <timo.heister@gmail.com>
Tue, 25 Sep 2012 17:52:19 +0000 (17:52 +0000)
git-svn-id: https://svn.dealii.org/trunk@26734 0785d39b-7218-0410-832d-ea1e28bc413d

tests/mpi/distribute_sp_02.cc [new file with mode: 0644]
tests/mpi/distribute_sp_02/ncpu_10/cmp/generic [new file with mode: 0644]
tests/mpi/distribute_sp_02/ncpu_4/cmp/generic [new file with mode: 0644]

diff --git a/tests/mpi/distribute_sp_02.cc b/tests/mpi/distribute_sp_02.cc
new file mode 100644 (file)
index 0000000..cb06cd4
--- /dev/null
@@ -0,0 +1,169 @@
+//---------------------------------------------------------------------------
+//    $Id: distribute_sp_01.cc 24425 2011-09-26 13:53:32Z 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.
+//
+//---------------------------------------------------------------------------
+
+
+// check SparsityTools::distribute_sparsity_pattern for BlockCompressedSimpleSP
+
+#include "../tests.h"
+#include <deal.II/base/logstream.h>
+#include <deal.II/base/utilities.h>
+
+#include <deal.II/base/index_set.h>
+#include <deal.II/lac/sparsity_tools.h>
+#include <deal.II/lac/compressed_simple_sparsity_pattern.h>
+#include <deal.II/lac/block_sparsity_pattern.h>
+
+#include <fstream>
+
+
+void test_mpi()
+{
+  Assert( Utilities::System::job_supports_mpi(), ExcInternalError());
+
+  unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
+  const unsigned int numprocs = Utilities::MPI::n_mpi_processes (MPI_COMM_WORLD);
+
+  if (myid==0)
+    deallog << "Running on " << numprocs << " CPU(s)." << std::endl;
+
+  unsigned int num_local=10;
+  unsigned int n=numprocs*num_local;
+  std::vector<unsigned int> rows_per_cpu;
+  for (unsigned int i=0;i<numprocs;++i)
+    rows_per_cpu.push_back(num_local);
+
+  IndexSet locally_rel(n);
+  locally_rel.add_range(myid*num_local, (myid+1)*num_local);
+  if (myid>0)
+    locally_rel.add_range((myid-1)*num_local, (myid+0)*num_local);
+  if (myid<numprocs-1)
+    locally_rel.add_range((myid+1)*num_local, (myid+2)*num_local);
+
+  std::vector<IndexSet> partitioning;
+  partitioning.push_back(locally_rel);
+  
+  BlockCompressedSimpleSparsityPattern csp(partitioning);
+  
+  for (unsigned int i=0;i<n;++i)
+    csp.add(i, myid);
+  
+  if (myid==0)
+    {
+      deallog << "blocks: " << csp.n_block_rows() << "x" << csp.n_block_cols() << std::endl;
+      deallog << "size: " << csp.n_rows() << "x" << csp.n_cols() << std::endl;
+    }
+  
+  SparsityTools::distribute_sparsity_pattern<>(csp,
+                                              rows_per_cpu,
+                                              MPI_COMM_WORLD,
+                                              locally_rel);
+/*  {
+    std::ofstream f((std::string("after")+Utilities::int_to_string(myid)).c_str());
+    csp.print(f);
+    }*/
+
+                                  // checking...
+   for (unsigned int r=0;r<num_local;++r)
+    {
+      unsigned int indx=r+myid*num_local;
+      unsigned int len=csp.row_length(indx);
+
+//std::cout << "myid=" << myid << " idx=" << indx << " len=" << len <<std::endl;
+
+      if (myid>0 && myid<numprocs-1)
+       Assert(len==3, ExcInternalError());
+      if (myid==0 || myid==numprocs-1)
+       Assert(len==2, ExcInternalError());
+
+      Assert(csp.exists(indx, myid), ExcInternalError());
+      if (myid>0)
+       Assert(csp.exists(indx, myid-1), ExcInternalError());
+      if (myid<numprocs-1)
+       Assert(csp.exists(indx, myid+1), ExcInternalError());
+    }
+
+  if (myid==0)
+    deallog << "part 2" << std::endl;
+
+  IndexSet bla(1);
+  if (myid==0)
+    bla.add_index(0);
+  
+  partitioning.push_back(bla);
+
+  csp.reinit(partitioning);
+  for (unsigned int i=0;i<n;++i)
+    csp.add(i, myid);
+
+  SparsityTools::distribute_sparsity_pattern<>(csp,
+                                              rows_per_cpu,
+                                              MPI_COMM_WORLD,
+                                              locally_rel);
+
+  if (myid==0)
+    {
+      deallog << "blocks: " << csp.n_block_rows() << "x" << csp.n_block_cols() << std::endl;
+      deallog << "size: " << csp.n_rows() << "x" << csp.n_cols() << std::endl;
+    }
+  
+                                  // checking...
+  for (unsigned int r=0;r<num_local;++r)
+    {
+      unsigned int indx=r+myid*num_local;
+      unsigned int len=csp.row_length(indx);
+
+      //std::cout << "myid=" << myid << " idx=" << indx << " len=" << len <<std::endl;
+
+      if (myid>0 && myid<numprocs-1)
+       Assert(len==3, ExcInternalError());
+      if (myid==numprocs-1 || myid==0)
+       Assert(len==2, ExcInternalError());
+      
+      Assert(csp.exists(indx, myid), ExcInternalError());
+      if (myid>0)
+       Assert(csp.exists(indx, myid-1), ExcInternalError());
+      if (myid<numprocs-1)
+       Assert(csp.exists(indx, myid+1), ExcInternalError());
+    }
+   
+  if (myid==0)
+    deallog << "done" << std::endl;
+
+}
+
+
+int main(int argc, char *argv[])
+{
+#ifdef DEAL_II_COMPILER_SUPPORTS_MPI
+  Utilities::MPI::MPI_InitFinalize mpi (argc, argv);
+#else
+  (void)argc;
+  (void)argv;
+  compile_time_error;
+
+#endif
+
+  if (Utilities::MPI::this_mpi_process (MPI_COMM_WORLD) == 0)
+    {
+      std::ofstream logfile(output_file_for_mpi("distribute_sp_02").c_str());
+      deallog.attach(logfile);
+      deallog.depth_console(0);
+      deallog.threshold_double(1.e-10);
+
+      deallog.push("mpi");
+      test_mpi();
+      deallog.pop();
+    }
+  else
+    test_mpi();
+}
diff --git a/tests/mpi/distribute_sp_02/ncpu_10/cmp/generic b/tests/mpi/distribute_sp_02/ncpu_10/cmp/generic
new file mode 100644 (file)
index 0000000..21f0f8a
--- /dev/null
@@ -0,0 +1,8 @@
+
+DEAL:mpi::Running on 10 CPU(s).
+DEAL:mpi::blocks: 1x1
+DEAL:mpi::size: 100x100
+DEAL:mpi::part 2
+DEAL:mpi::blocks: 2x2
+DEAL:mpi::size: 101x101
+DEAL:mpi::done
diff --git a/tests/mpi/distribute_sp_02/ncpu_4/cmp/generic b/tests/mpi/distribute_sp_02/ncpu_4/cmp/generic
new file mode 100644 (file)
index 0000000..538e892
--- /dev/null
@@ -0,0 +1,8 @@
+
+DEAL:mpi::Running on 4 CPU(s).
+DEAL:mpi::blocks: 1x1
+DEAL:mpi::size: 40x40
+DEAL:mpi::part 2
+DEAL:mpi::blocks: 2x2
+DEAL:mpi::size: 41x41
+DEAL:mpi::done

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.