]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Put in three tests that illustrate a dirty bug exactly.
authoryoung <young@0785d39b-7218-0410-832d-ea1e28bc413d>
Sun, 20 Oct 2013 15:58:53 +0000 (15:58 +0000)
committeryoung <young@0785d39b-7218-0410-832d-ea1e28bc413d>
Sun, 20 Oct 2013 15:58:53 +0000 (15:58 +0000)
git-svn-id: https://svn.dealii.org/branches/branch_petscscalar_complex@31340 0785d39b-7218-0410-832d-ea1e28bc413d

tests/petsc_scalar_complex/sparse_matrix_01.cc [new file with mode: 0644]
tests/petsc_scalar_complex/sparse_matrix_01/cmp/generic [new file with mode: 0644]
tests/petsc_scalar_complex/sparse_matrix_02.cc [new file with mode: 0644]
tests/petsc_scalar_complex/sparse_matrix_02/cmp/generic [new file with mode: 0644]
tests/petsc_scalar_complex/vector_02.cc [new file with mode: 0644]
tests/petsc_scalar_complex/vector_02/cmp/generic [new file with mode: 0644]

diff --git a/tests/petsc_scalar_complex/sparse_matrix_01.cc b/tests/petsc_scalar_complex/sparse_matrix_01.cc
new file mode 100644 (file)
index 0000000..ed35aa4
--- /dev/null
@@ -0,0 +1,83 @@
+// ---------------------------------------------------------------------
+// $Id sparse_matrix_01.cc
+//
+// Copyright (C) 2013 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 SparseMatrix::add(other, factor)
+
+#include "../tests.h"
+#include <deal.II/lac/petsc_vector.h>
+#include <deal.II/lac/petsc_sparse_matrix.h>
+#include <fstream>
+#include <iostream>
+#include <vector>
+
+
+void test ()
+{
+  const unsigned int s = 10;
+  PETScWrappers::SparseMatrix m(s,s,s);
+  for (unsigned int k=0; k<m.m(); ++k)
+    for (unsigned int l=0; l<=k; ++l)
+      m.set (k,l, k+2*l);
+
+  m.compress (VectorOperation::add);
+
+  PETScWrappers::SparseMatrix m2(s,s,s);
+  m2.set(0,1,5.0);
+  for (unsigned int k=0; k<m2.m(); ++k)
+    m2.set(k,k,PetscScalar(1.0+k,-1.0-k));
+  m2.compress(VectorOperation::add);
+
+  // we now print the matrix m, it is all real. Then print after
+  // adding m2 which is complex, and then subtract m2 again to get the
+  // original matrix back.
+
+  deallog << "before: " << m(0,1) << std::endl;
+  for (unsigned int k=0; k<s; ++k)
+    deallog << m(k,k) << " ";
+  deallog << std::endl;
+
+  m.add(m2,1.0);
+
+  deallog << "after: " << m(0,1) << std::endl;
+  for (unsigned int k=0; k<s; ++k)
+    deallog << m(k,k) << " ";
+  deallog << std::endl;
+
+  m.add(m2,-1.0);
+
+  deallog << "back to original: " << m(0,1) << std::endl;
+  for (unsigned int k=0; k<s; ++k)
+    deallog << m(k,k) << " ";
+  deallog << std::endl;
+
+  deallog << "OK" << std::endl;
+}
+
+
+
+int main (int argc, char **argv)
+{
+  std::ofstream logfile("sparse_matrix_01/output");
+  deallog.attach(logfile);
+  deallog.depth_console(0);
+  deallog.threshold_double(1.e-10);
+
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+  test ();
+
+}
diff --git a/tests/petsc_scalar_complex/sparse_matrix_01/cmp/generic b/tests/petsc_scalar_complex/sparse_matrix_01/cmp/generic
new file mode 100644 (file)
index 0000000..0fac50f
--- /dev/null
@@ -0,0 +1,8 @@
+
+DEAL::before: (0.00000,0.00000)
+DEAL::(0.00000,0.00000) (3.00000,0.00000) (6.00000,0.00000) (9.00000,0.00000) (12.0000,0.00000) (15.0000,0.00000) (18.0000,0.00000) (21.0000,0.00000) (24.0000,0.00000) (27.0000,0.00000) 
+DEAL::after: (5.00000,0.00000)
+DEAL::(1.00000,-1.00000) (5.00000,-2.00000) (9.00000,-3.00000) (13.0000,-4.00000) (17.0000,-5.00000) (21.0000,-6.00000) (25.0000,-7.00000) (29.0000,-8.00000) (33.0000,-9.00000) (37.0000,-10.0000) 
+DEAL::back to original: (0.00000,0.00000)
+DEAL::(0.00000,0.00000) (3.00000,0.00000) (6.00000,0.00000) (9.00000,0.00000) (12.0000,0.00000) (15.0000,0.00000) (18.0000,0.00000) (21.0000,0.00000) (24.0000,0.00000) (27.0000,0.00000) 
+DEAL::OK
diff --git a/tests/petsc_scalar_complex/sparse_matrix_02.cc b/tests/petsc_scalar_complex/sparse_matrix_02.cc
new file mode 100644 (file)
index 0000000..9b7c83d
--- /dev/null
@@ -0,0 +1,94 @@
+// ---------------------------------------------------------------------
+// $Id sparse_matrix_01.cc
+//
+// Copyright (C) 2013 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 SparseMatrix::add(other, factor)
+
+// This is the same as sparse_matrix_02.cc, but uses integer input to
+// complex numbers. AT the time of writing, this documents the element
+// access problem for matrices.
+//
+// TODO: Fix this before going further with petsc-scalar=complex.
+
+#include "../tests.h"
+#include <deal.II/lac/petsc_vector.h>
+#include <deal.II/lac/petsc_sparse_matrix.h>
+#include <fstream>
+#include <iostream>
+#include <vector>
+
+
+void test ()
+{
+  const unsigned int s = 10;
+  PETScWrappers::SparseMatrix m(s,s,s);
+  for (unsigned int k=0; k<m.m(); ++k)
+    for (unsigned int l=0; l<=k; ++l)
+      m.set (k,l, k+2*l);
+
+  m.compress (VectorOperation::add);
+
+  PETScWrappers::SparseMatrix m2(s,s,s);
+  m2.set(0,1,5.0);
+  for (unsigned int k=0; k<m2.m(); ++k)
+    {
+      // This fails
+      m2.set(k,k,PetscScalar(k,-k));
+      // This is ok
+      // m2.set(k,k,PetscScalar(k,-1.0*k));
+    }
+  m2.compress(VectorOperation::add);
+      
+  // we now print the matrix m, it is all real. Then print after
+  // adding m2 which is complex, and then subtract m2 again to get the
+  // original matrix back.
+
+  deallog << "before: " << m(0,1) << std::endl;
+  for (unsigned int k=0; k<s; ++k)
+    deallog << m(k,k) << " ";
+  deallog << std::endl;
+
+  m.add(m2,1.0);
+
+  deallog << "after: " << m(0,1) << std::endl;
+  for (unsigned int k=0; k<s; ++k)
+    deallog << m(k,k) << " ";
+  deallog << std::endl;
+
+  m.add(m2,-1.0);
+
+  deallog << "back to original: " << m(0,1) << std::endl;
+  for (unsigned int k=0; k<s; ++k)
+    deallog << m(k,k) << " ";
+  deallog << std::endl;
+
+  deallog << "OK" << std::endl;
+}
+
+
+
+int main (int argc, char **argv)
+{
+  std::ofstream logfile("sparse_matrix_02/output");
+  deallog.attach(logfile);
+  deallog.depth_console(0);
+  deallog.threshold_double(1.e-10);
+
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+  test ();
+
+}
diff --git a/tests/petsc_scalar_complex/sparse_matrix_02/cmp/generic b/tests/petsc_scalar_complex/sparse_matrix_02/cmp/generic
new file mode 100644 (file)
index 0000000..5383659
--- /dev/null
@@ -0,0 +1,8 @@
+
+DEAL::before: (0.00000,0.00000)
+DEAL::(0.00000,0.00000) (3.00000,0.00000) (6.00000,0.00000) (9.00000,0.00000) (12.0000,0.00000) (15.0000,0.00000) (18.0000,0.00000) (21.0000,0.00000) (24.0000,0.00000) (27.0000,0.00000) 
+DEAL::after: (5.00000,0.00000)
+DEAL::(0.00000,0.00000) (4.00000,-1.00000) (8.00000,-2.00000) (12.0000,-3.00000) (16.0000,-4.00000) (20.0000,-5.00000) (24.0000,-6.00000) (28.0000,-7.00000) (32.0000,-8.00000) (36.0000,-9.00000) 
+DEAL::back to original: (0.00000,0.00000)
+DEAL::(0.00000,0.00000) (3.00000,0.00000) (6.00000,0.00000) (9.00000,0.00000) (12.0000,0.00000) (15.0000,0.00000) (18.0000,0.00000) (21.0000,0.00000) (24.0000,0.00000) (27.0000,0.00000) 
+DEAL::OK
diff --git a/tests/petsc_scalar_complex/vector_02.cc b/tests/petsc_scalar_complex/vector_02.cc
new file mode 100644 (file)
index 0000000..c602a65
--- /dev/null
@@ -0,0 +1,90 @@
+// ---------------------------------------------------------------------
+// $Id sparse_matrix_01.cc
+//
+// Copyright (C) 2013 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 Vector
+
+// This is the similar to sparse_matrix_02.cc as it uses integer input
+// to complex numbers. This purposefully illustrates that this
+// assignment DOES NOT works for PETScWrappers::Vector.
+
+#include "../tests.h"
+#include <deal.II/lac/petsc_vector.h>
+#include <fstream>
+#include <iostream>
+#include <vector>
+
+
+void test ()
+{
+  const unsigned int s = 10;
+  PETScWrappers::Vector v(s);
+  for (unsigned int k=0; k<v.size(); ++k)
+    v(k) = k;
+
+  v.compress (VectorOperation::add);
+
+  PETScWrappers::Vector v2(s);
+  for (unsigned int k=0; k<v2.size(); ++k)
+    {
+      // This fails
+      v2(k) = PetscScalar (k,-k);
+      // This is ok
+      // v2(k) = PetscScalar (k,-1.0*k);
+    }
+
+  v2.compress(VectorOperation::add);
+
+  // we now print the vector v, it is all real. Then print after
+  // adding v2 which is complex, and then subtract v2 again to get the
+  // original vector back.
+
+   deallog << "before: " << std::endl;
+   for (unsigned int k=0; k<s; ++k)
+     deallog << "(" << v(k).real () << "," << v(k).imag () << "i) ";
+   deallog << std::endl;
+
+   v.add(1.0,v2);
+
+   deallog << "after: " << std::endl;
+   for (unsigned int k=0; k<s; ++k)
+     deallog << "(" << v(k).real () << "," << v(k).imag () << "i) ";
+   deallog << std::endl;
+
+  v.add(-1.0,v2);
+
+  deallog << "back to original: " << std::endl;
+  for (unsigned int k=0; k<s; ++k)
+    deallog << "(" << v(k).real () << "," << v(k).imag () << "i) ";
+  deallog << std::endl;
+  
+  deallog << "OK" << std::endl;
+}
+
+
+
+int main (int argc, char **argv)
+{
+  std::ofstream logfile("vector_02/output");
+  deallog.attach(logfile);
+  deallog.depth_console(0);
+  deallog.threshold_double(1.e-10);
+
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+  test ();
+
+}
diff --git a/tests/petsc_scalar_complex/vector_02/cmp/generic b/tests/petsc_scalar_complex/vector_02/cmp/generic
new file mode 100644 (file)
index 0000000..d20636e
--- /dev/null
@@ -0,0 +1,8 @@
+
+DEAL::before: 
+DEAL::(0,0i) (1.00000,0i) (2.00000,0i) (3.00000,0i) (4.00000,0i) (5.00000,0i) (6.00000,0i) (7.00000,0i) (8.00000,0i) (9.00000,0i) 
+DEAL::after: 
+DEAL::(0,0i) (2.00000,-1.00000i) (4.00000,-2.00000i) (6.00000,-3.00000i) (8.00000,-4.00000i) (10.0000,-5.00000i) (12.0000,-6.00000i) (14.0000,-7.00000i) (16.0000,-8.00000i) (18.0000,-9.00000i) 
+DEAL::back to original: 
+DEAL::(0,0i) (1.00000,0i) (2.00000,0i) (3.00000,0i) (4.00000,0i) (5.00000,0i) (6.00000,0i) (7.00000,0i) (8.00000,0i) (9.00000,0i) 
+DEAL::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.