]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Add new but failing test.
authorbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 5 Mar 2014 01:52:29 +0000 (01:52 +0000)
committerbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 5 Mar 2014 01:52:29 +0000 (01:52 +0000)
git-svn-id: https://svn.dealii.org/trunk@32611 0785d39b-7218-0410-832d-ea1e28bc413d

tests/mpi/trilinos_sadd_01.cc [new file with mode: 0644]
tests/mpi/trilinos_sadd_01.mpirun=2.output [new file with mode: 0644]

diff --git a/tests/mpi/trilinos_sadd_01.cc b/tests/mpi/trilinos_sadd_01.cc
new file mode 100644 (file)
index 0000000..61b4cd4
--- /dev/null
@@ -0,0 +1,106 @@
+// ---------------------------------------------------------------------
+// $Id$
+//
+// Copyright (C) 2014 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 correct behavior of sadd of Trilinos vectors
+// if they have different Epetra maps 
+
+#include "../tests.h"
+#include <deal.II/base/utilities.h>
+#include <deal.II/base/index_set.h>
+#include <deal.II/lac/trilinos_vector.h>
+#include <fstream>
+#include <iostream>
+#include <vector>
+#include <algorithm>
+
+void test ()
+{
+  const int n_proc = Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD);
+  const int my_id = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD);
+
+  //All processes should own 10 entries
+  const int entries_per_process = 10;
+
+  IndexSet locally_owned(entries_per_process*n_proc);
+  const int begin_index = my_id*entries_per_process;
+  const int end_index = (my_id+1)*entries_per_process;
+  locally_owned.add_range(begin_index, end_index);
+
+  IndexSet locally_relevant(entries_per_process*n_proc);
+  const int local_begin = std::max(0, begin_index-entries_per_process/2);
+  const int local_end = entries_per_process*n_proc;
+  locally_relevant.add_range (local_begin, local_end);
+  TrilinosWrappers::MPI::Vector ghosted, distributed;
+  distributed.reinit(locally_owned, MPI_COMM_WORLD);
+  ghosted.reinit (locally_owned, locally_relevant, MPI_COMM_WORLD);
+
+  // set the 'distributed' vector to all ones and store its results in
+  // 'ghosted'
+  distributed=1.;
+  ghosted=distributed;
+
+  // then multiply 'distributed' by two and add 'ghosted' to it again
+  distributed*=2;
+  distributed.add(ghosted, true);
+
+  // assign the result, which should contain all 3s to 'ghosted'
+  ghosted = distributed;
+
+  if (my_id==0)
+    {
+      deallog << "Distributed:" << std::endl;
+      for (unsigned int i=begin_index; i<end_index; ++i)
+       deallog << i << ": " << distributed(i) << std::endl; 
+      
+      deallog << "Ghosted:" << std::endl;
+      for (unsigned int i=local_begin; i<local_end; ++i)
+       deallog << i << ": " << ghosted(i) << std::endl;
+    }
+
+  // verify correct value
+  for (unsigned int i=begin_index; i<end_index; ++i)
+    Assert(distributed(i)==3, ExcInternalError());
+      
+  for (unsigned int i=local_begin; i<local_end; ++i)
+    Assert(ghosted(i)==3, ExcInternalError());
+}
+
+
+
+int main (int argc, char **argv)
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv);
+  
+  unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
+  deallog.push(Utilities::int_to_string(myid));
+  
+  if (myid == 0)
+  {
+    std::ofstream logfile("output");
+    deallog.attach(logfile);
+    deallog << std::setprecision(4);
+    deallog.depth_console(0);
+    deallog.threshold_double(1.e-10);
+    
+    test();
+  }
+  else
+    test();
+  
+}
diff --git a/tests/mpi/trilinos_sadd_01.mpirun=2.output b/tests/mpi/trilinos_sadd_01.mpirun=2.output
new file mode 100644 (file)
index 0000000..d954ea5
--- /dev/null
@@ -0,0 +1,33 @@
+
+DEAL:0::Distributed:
+DEAL:0::0: 3.000
+DEAL:0::1: 3.000
+DEAL:0::2: 3.000
+DEAL:0::3: 3.000
+DEAL:0::4: 3.000
+DEAL:0::5: 3.000
+DEAL:0::6: 3.000
+DEAL:0::7: 3.000
+DEAL:0::8: 3.000
+DEAL:0::9: 3.000
+DEAL:0::Ghosted:
+DEAL:0::0: 3.000
+DEAL:0::1: 3.000
+DEAL:0::2: 3.000
+DEAL:0::3: 3.000
+DEAL:0::4: 3.000
+DEAL:0::5: 3.000
+DEAL:0::6: 3.000
+DEAL:0::7: 3.000
+DEAL:0::8: 3.000
+DEAL:0::9: 3.000
+DEAL:0::10: 3.000
+DEAL:0::11: 3.000
+DEAL:0::12: 3.000
+DEAL:0::13: 3.000
+DEAL:0::14: 3.000
+DEAL:0::15: 3.000
+DEAL:0::16: 3.000
+DEAL:0::17: 3.000
+DEAL:0::18: 3.000
+DEAL:0::19: 3.000

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.