]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
New test to verify that indeed Trilinos bug 5609 is and remains fixed.
authorbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Tue, 26 Jun 2012 08:46:50 +0000 (08:46 +0000)
committerbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Tue, 26 Jun 2012 08:46:50 +0000 (08:46 +0000)
git-svn-id: https://svn.dealii.org/trunk@25640 0785d39b-7218-0410-832d-ea1e28bc413d

tests/mpi/trilinos_bug_5609.cc [new file with mode: 0644]
tests/mpi/trilinos_bug_5609/ncpu_2/cmp/generic [new file with mode: 0644]

diff --git a/tests/mpi/trilinos_bug_5609.cc b/tests/mpi/trilinos_bug_5609.cc
new file mode 100644 (file)
index 0000000..dc71337
--- /dev/null
@@ -0,0 +1,87 @@
+//----------------------------  trilinos_bug_5609.cc  ---------------------------
+//    $Id$
+//    Version: $Name$
+//
+//    Copyright (C) 2011, 2012 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.
+//
+//----------------------------  trilinos_bug_5609.cc  ---------------------------
+
+
+// not really a testcase for deal.II but one that verifies that the recurring
+// Trilinos bug 5609 is indeed fixed, see
+//   http://software.sandia.gov/bugzilla/show_bug.cgi?id=5609
+
+
+#include "../tests.h"
+#include <Epetra_Map.h>
+#include <Epetra_MpiComm.h>
+#include <Epetra_FEVector.h>
+
+#include <mpi.h>
+
+#include <iostream>
+
+
+void test ()
+{
+  int NumElements = 4;
+
+  Epetra_MpiComm Comm (MPI_COMM_WORLD);
+  Epetra_Map     Map(NumElements, 0, Comm);
+  Epetra_FEVector x1(Map);
+  x1.PutScalar (0);
+
+                               // let all processors set global entry 0 to 1
+  const int GID = 0;
+  const double value = 1;
+  x1.ReplaceGlobalValues(1, &GID, &value);
+  x1.GlobalAssemble (Insert);
+  if (Comm.MyPID()==0)
+    Assert (x1[0][0] == 1, ExcInternalError());
+
+                               // copy vector
+  Epetra_FEVector x2 (x1);
+
+  x2.PutScalar(0);
+
+                               // re-apply 1 to the vector, but only on the
+                               // owning processor. should be enough to set
+                               // the value (as non-local data in x1 should
+                               // be been eliminated after calling
+                               // GlobalAssemble).
+  if (Comm.MyPID()==0)
+    x2.ReplaceGlobalValues(1, &GID, &value);
+  x2.GlobalAssemble (Insert);
+
+  if (Comm.MyPID()==0)
+    Assert (x1[0][0] == 1, 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_file_for_mpi("trilinos_bug_5609").c_str());
+      deallog.attach(logfile);
+      deallog << std::setprecision(4);
+      deallog.depth_console(0);
+      deallog.threshold_double(1.e-10);
+
+      test();
+
+      deallog << "OK" << std::endl;
+    }
+  else
+    test();
+}
diff --git a/tests/mpi/trilinos_bug_5609/ncpu_2/cmp/generic b/tests/mpi/trilinos_bug_5609/ncpu_2/cmp/generic
new file mode 100644 (file)
index 0000000..be8d055
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL:0::OK

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.