]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Fix a problem where we look up elements in a matrix using the row map instead of...
authorbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Mon, 9 May 2011 22:45:32 +0000 (22:45 +0000)
committerbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Mon, 9 May 2011 22:45:32 +0000 (22:45 +0000)
git-svn-id: https://svn.dealii.org/trunk@23684 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/doc/news/changes.h
deal.II/source/lac/trilinos_sparse_matrix.cc
tests/trilinos/sparse_matrix_01.cc [new file with mode: 0644]
tests/trilinos/sparse_matrix_01/cmp/generic [new file with mode: 0644]

index 315aa47ef0390065802f5876df5b8a9695056ae6..5e7fa9daacb31cb2352e3ef437be26e763f8f5ee 100644 (file)
@@ -100,6 +100,12 @@ should be fixed now.
 <h3>Specific improvements</h3>
 
 <ol>
+<li> Fixed: The TrilinosWrappers::SparseMatrix::operator() and
+TrilinosWrappers::SparseMatrix::el() functions sometimes produced
+wrong results for rectangular matrices. This is now fixed.
+<br>
+(Habib Talavatifard, Wolfgang Bangerth 2011/05/09)
+
 <li> Changed: DoFTools is now a namespace. It has long been a class that
 had only public, static member functions, making the end result semantically
 exactly equivalent to a namespace, which is also how it was used. This is
index c4457d22a2fe4d52f9f3211f152cfd398455ef9b..7449c8accfb87614393f6409f293ee828137601b 100644 (file)
@@ -770,7 +770,7 @@ namespace TrilinosWrappers
   {
                                      // Extract local indices in
                                      // the matrix.
-    int trilinos_i = matrix->LRID(i), trilinos_j = matrix->LRID(j);
+    int trilinos_i = matrix->LRID(i), trilinos_j = matrix->LCID(j);
     TrilinosScalar value = 0.;
 
                                      // If the data is not on the
@@ -844,7 +844,7 @@ namespace TrilinosWrappers
   {
                                      // Extract local indices in
                                      // the matrix.
-    int trilinos_i = matrix->LRID(i), trilinos_j = matrix->LRID(j);
+    int trilinos_i = matrix->LRID(i), trilinos_j = matrix->LCID(j);
     TrilinosScalar value = 0.;
 
                                      // If the data is not on the
diff --git a/tests/trilinos/sparse_matrix_01.cc b/tests/trilinos/sparse_matrix_01.cc
new file mode 100644 (file)
index 0000000..971ba5f
--- /dev/null
@@ -0,0 +1,62 @@
+//-----------------  trilinos_sparse_matrix_01.cc  -------------------------
+//    $Id$
+//    Version: $Name$
+//
+//    Copyright (C) 2011 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_sparse_matrix_01.cc  -------------------------
+
+
+// TrilinosWrappers::SparseMatrix::el() and operator() had problems
+// looking up indices in rectangular matrices because they
+// accidentally used the row instead of the column map. verify that
+// this is now fixed
+//
+// this testcase is reduced from one contributed by Habib Talavatifard
+
+#include "../tests.h"
+#include <base/utilities.h>
+#include <lac/trilinos_sparse_matrix.h>
+#include <lac/trilinos_sparsity_pattern.h>
+#include <fstream>
+#include <iostream>
+
+
+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::System::MPI_InitFinalize mpi_initialization (argc, argv);
+
+  IndexSet row_partitioning (3);
+  IndexSet col_partitioning (4);
+
+  row_partitioning.add_range(0, 3);
+  col_partitioning.add_range(0, 4);
+
+                                  // Add element (2,3) to the matrix
+  TrilinosWrappers::SparsityPattern sp (row_partitioning, col_partitioning);
+  sp.add (2,3);
+  sp.compress();
+
+  TrilinosWrappers::SparseMatrix A (sp);
+  A.add (2, 3, 2.0);
+  A.compress();
+
+                                  // verify that entry (2,3) is
+                                  // indeed what we expect. verify
+                                  // that both methods of accessing
+                                  // the entry work
+  Assert (A.el(2, 3) == 2, ExcInternalError());
+  Assert (A(2, 3) == 2, ExcInternalError());
+
+  deallog << "OK" << std::endl;
+}
diff --git a/tests/trilinos/sparse_matrix_01/cmp/generic b/tests/trilinos/sparse_matrix_01/cmp/generic
new file mode 100644 (file)
index 0000000..0fd8fc1
--- /dev/null
@@ -0,0 +1,2 @@
+
+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.