]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add a test where we use the ILU to compute actual inverses that we can check.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Fri, 29 Feb 2008 22:04:35 +0000 (22:04 +0000)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Fri, 29 Feb 2008 22:04:35 +0000 (22:04 +0000)
git-svn-id: https://svn.dealii.org/trunk@15818 0785d39b-7218-0410-832d-ea1e28bc413d

tests/lac/Makefile
tests/lac/sparse_ilu_inverse.cc [new file with mode: 0644]
tests/lac/sparse_ilu_inverse/cmp/generic [new file with mode: 0644]
tests/lac/sparse_ilu_inverse/cmp/mips-sgi-irix6.5+MIPSpro7.4 [new file with mode: 0644]

index 7ae78a04de27dbee79c795892ac3d437636e4642..b30389e76ffd1ece6cd9d99362d65d869fee8f85 100644 (file)
@@ -1,6 +1,6 @@
 ############################################################
 # $Id$
-# Copyright (C) 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007 by the deal.II authors
+# Copyright (C) 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008 by the deal.II authors
 ############################################################
 
 ############################################################
@@ -39,8 +39,7 @@ tests_x = vector-vector vector_memory \
        solver \
        eigen \
        filtered_matrix \
-       sparse_ilu \
-       sparse_ilu_t \
+       sparse_ilu* \
        lapack_fill \
        block_minres \
        identity_matrix* \
diff --git a/tests/lac/sparse_ilu_inverse.cc b/tests/lac/sparse_ilu_inverse.cc
new file mode 100644 (file)
index 0000000..24f02cf
--- /dev/null
@@ -0,0 +1,87 @@
+//----------------------------  sparse_ilu_inverse.cc  ---------------------------
+//    $Id$
+//    Version: $Name$ 
+//
+//    Copyright (C) 2008 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.
+//
+//----------------------------  sparse_ilu_inverse.cc  ---------------------------
+
+
+// compute the inverse of a small matrix using the SparseILU with
+// infinite fill-in
+
+#include "../tests.h"
+#include <cmath>
+#include <fstream>
+#include <iomanip>
+#include <iomanip>
+#include <cstdlib>
+#include "testmatrix.h"
+#include <base/logstream.h>
+#include <lac/sparse_matrix.h>
+#include <lac/full_matrix.h>
+#include <lac/sparse_ilu.h>
+#include <lac/vector.h>
+
+
+int main()
+{
+  std::ofstream logfile("sparse_ilu_inverse/output");
+  deallog << std::fixed;
+  deallog << std::setprecision(3);
+  deallog.attach(logfile);
+  deallog.depth_console(0);
+  deallog.threshold_double(1.e-10);
+  
+
+  for (unsigned int N=1; N<5; ++N)
+    {
+      deallog << "N=" << N << std::endl;
+      
+      SparsityPattern structure(N, N, N);
+      for (unsigned int i=0; i<N; ++i)
+       for (unsigned int j=0; j<N; ++j)
+         structure.add(i,j);
+      structure.compress();
+      SparseMatrix<double>  A(structure);
+      for (unsigned int i=0; i<N; ++i)
+       {
+         A.set(i,i,2);
+         if (i>=1)
+           A.set(i,i-1, -1);
+         if (i<N-1)
+           A.set(i,i+1, -1);
+       }
+
+      SparseILU<double> ilu;
+      ilu.initialize (A, SparseILU<double>::AdditionalData());
+
+                                      // now get an explicit
+                                      // representation of the
+                                      // inverse
+      FullMatrix<double> inverse (N,N);
+      Vector<double> tmp1(N), tmp2(N);
+      for (unsigned int i=0; i<N; ++i)
+       {
+         tmp1 = 0;
+         tmp1(i) = 1;
+         ilu.vmult (tmp2, tmp1);
+         for (unsigned int j=0; j<N; ++j)
+           inverse(i,j) = tmp2(j);
+       }
+
+      deallog << "Matrix A:" << std::endl;
+      A.print_formatted (deallog.get_file_stream());
+  
+      deallog << "Matrix A^{-1}:" << std::endl;
+      inverse.print_formatted (deallog.get_file_stream());
+
+      deallog << std::endl;
+    }
+}
+
diff --git a/tests/lac/sparse_ilu_inverse/cmp/generic b/tests/lac/sparse_ilu_inverse/cmp/generic
new file mode 100644 (file)
index 0000000..d156ab1
--- /dev/null
@@ -0,0 +1,28 @@
+
+DEAL::Size 4 Unknowns 9
+DEAL::Test 0
+DEAL::Residual with test vector 0:   left=0, right=0
+DEAL::Residual with test vector 1:   left=0, right=0
+DEAL::Residual with test vector 2:   left=0, right=0
+DEAL::Test 1
+DEAL::Residual with test vector 0:   left=0.308, right=0.350
+DEAL::Residual with test vector 1:   left=0.337, right=0.371
+DEAL::Residual with test vector 2:   left=0.290, right=0.327
+DEAL::Size 8 Unknowns 49
+DEAL::Test 0
+DEAL::Residual with test vector 0:   left=0, right=0
+DEAL::Residual with test vector 1:   left=0, right=0
+DEAL::Residual with test vector 2:   left=0, right=0
+DEAL::Test 1
+DEAL::Residual with test vector 0:   left=1.886, right=1.910
+DEAL::Residual with test vector 1:   left=2.123, right=2.138
+DEAL::Residual with test vector 2:   left=1.999, right=2.004
+DEAL::Size 16 Unknowns 225
+DEAL::Test 0
+DEAL::Residual with test vector 0:   left=0, right=0
+DEAL::Residual with test vector 1:   left=0, right=0
+DEAL::Residual with test vector 2:   left=0, right=0
+DEAL::Test 1
+DEAL::Residual with test vector 0:   left=5.944, right=5.953
+DEAL::Residual with test vector 1:   left=5.619, right=5.640
+DEAL::Residual with test vector 2:   left=6.079, right=6.094
diff --git a/tests/lac/sparse_ilu_inverse/cmp/mips-sgi-irix6.5+MIPSpro7.4 b/tests/lac/sparse_ilu_inverse/cmp/mips-sgi-irix6.5+MIPSpro7.4
new file mode 100644 (file)
index 0000000..5a918f4
--- /dev/null
@@ -0,0 +1,28 @@
+
+DEAL::Size 4 Unknowns 9
+DEAL::Test 0
+DEAL::Residual with test vector 0:   left=0.000, right=0.000
+DEAL::Residual with test vector 1:   left=0.000, right=0.000
+DEAL::Residual with test vector 2:   left=0.000, right=0.000
+DEAL::Test 1
+DEAL::Residual with test vector 0:   left=0.425, right=0.447
+DEAL::Residual with test vector 1:   left=0.266, right=0.267
+DEAL::Residual with test vector 2:   left=0.297, right=0.377
+DEAL::Size 8 Unknowns 49
+DEAL::Test 0
+DEAL::Residual with test vector 0:   left=0.000, right=0.000
+DEAL::Residual with test vector 1:   left=0.000, right=0.000
+DEAL::Residual with test vector 2:   left=0.000, right=0.000
+DEAL::Test 1
+DEAL::Residual with test vector 0:   left=1.895, right=1.900
+DEAL::Residual with test vector 1:   left=1.705, right=1.711
+DEAL::Residual with test vector 2:   left=2.061, right=2.057
+DEAL::Size 16 Unknowns 225
+DEAL::Test 0
+DEAL::Residual with test vector 0:   left=0.000, right=0.000
+DEAL::Residual with test vector 1:   left=0.000, right=0.000
+DEAL::Residual with test vector 2:   left=0.000, right=0.000
+DEAL::Test 1
+DEAL::Residual with test vector 0:   left=5.900, right=5.911
+DEAL::Residual with test vector 1:   left=5.659, right=5.662
+DEAL::Residual with test vector 2:   left=5.914, right=5.921

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.