]> https://gitweb.dealii.org/ - dealii.git/commitdiff
preconditioner_from_diagonal_inverse
authorESeNonFossiIo <esenonfossiio@gmail.com>
Wed, 20 May 2015 16:52:01 +0000 (18:52 +0200)
committerESeNonFossiIo <esenonfossiio@gmail.com>
Sat, 23 May 2015 10:16:24 +0000 (12:16 +0200)
preconditioner

inverse check added

include/deal.II/lac/linear_operator.h
tests/lac/linear_operator_14.cc [new file with mode: 0644]
tests/lac/linear_operator_14.with_cxx11=on.output [new file with mode: 0644]

index c46a0a44e53003cbcd10af8e44fec0281ca01ca7..6cadd62cd240b057bee327f3f8a62b155ba4129e 100644 (file)
@@ -2329,6 +2329,38 @@ upper_triangular_operator(BlockMatrix &a)
   return return_op;
 }
 
+
+template <typename Range = BlockVector<double>,
+          typename Domain = Range,
+          typename BlockMatrix>
+LinearOperator<Range, Domain>
+inverse_operator( BlockMatrix &a,
+                  const LinearOperator<Range, Domain> &inverse_diagonal,
+                  bool lower = true)
+{
+
+  LinearOperator<Domain, Range> op_a;
+
+  if(lower)
+  {
+    op_a = lower_triangular_operator(a);
+  }
+  else
+  {
+    op_a = upper_triangular_operator(a);
+  }
+
+  auto id     = identity_operator(op_a.reinit_range_vector);
+  auto result = identity_operator(op_a.reinit_range_vector);
+
+
+  for(unsigned int i = 0; i < a.n_block_cols() - 1; ++i)
+    result = id - inverse_diagonal * op_a * result;
+
+  return result * inverse_diagonal;
+}
+
+
 //@}
 
 DEAL_II_NAMESPACE_CLOSE
diff --git a/tests/lac/linear_operator_14.cc b/tests/lac/linear_operator_14.cc
new file mode 100644 (file)
index 0000000..a722ad1
--- /dev/null
@@ -0,0 +1,146 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2015 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.
+//
+// ---------------------------------------------------------------------
+
+// Test non symmetric variants:
+
+#include "../tests.h"
+
+#include <deal.II/lac/block_sparse_matrix.h>
+#include <deal.II/lac/block_vector.h>
+#include <deal.II/lac/dynamic_sparsity_pattern.h>
+#include <deal.II/lac/linear_operator.h>
+#include <deal.II/lac/sparse_matrix.h>
+#include <deal.II/lac/vector.h>
+
+#define PRINTME(name, var) \
+  deallog << "Block vector: " name << ":" << std::endl; \
+  for (unsigned int i = 0; i < var.n_blocks(); ++i) \
+    deallog << "[block " << i << " ]  " << var.block(i);
+
+
+using namespace dealii;
+
+int main()
+{
+  initlog();
+  deallog << std::setprecision(12);
+
+  // BlockSparseMatrix:
+  {
+    BlockDynamicSparsityPattern dsp(3, 3);
+    for (unsigned int i = 0; i < 3; ++i)
+      for (unsigned int j = 0; j < 3; ++j)
+        dsp.block(i, j).reinit (1, 1);
+    dsp.collect_sizes ();
+
+    BlockSparsityPattern sparsity_pattern;
+    sparsity_pattern.copy_from(dsp);
+    sparsity_pattern.compress();
+
+    BlockSparseMatrix<double> a (sparsity_pattern);
+    for (unsigned int i = 0; i < 3; ++i)
+    {
+      a.block(i,i).set(0, 0, i+i +1);
+      for (unsigned int j = 0; j < i; ++j)
+        a.block(i,j).set(0, 0, 10);
+    }
+
+    BlockSparseMatrix<double> d(sparsity_pattern);
+    for (unsigned int i = 0; i < 3; ++i)
+        d.set(i,i, 1.0 / (i+i +1) );
+    // std::array<decltype(a.block(0,0)), 3> diag {{a.block(0,0), a.block(1,1),a.block(2,2)}};
+
+
+    auto op_a         = linear_operator< BlockVector<double> >(a);
+    auto diagonal_inv = linear_operator< BlockVector<double> >(d);
+    // auto diagonal_inv = block_diagonal_operator(diag);
+    auto inverse_op_a = inverse_operator< BlockVector<double>, BlockVector<double>, BlockSparseMatrix<double> >(a, diagonal_inv);
+
+    auto identity = inverse_op_a * op_a;
+
+    BlockVector<double> u;
+    BlockVector<double> v;
+
+    deallog << " -- Matrix -- " << std::endl;
+    op_a.reinit_domain_vector(u, false);
+    op_a.reinit_range_vector(v, false);
+    for(unsigned int j = 0; j<3; ++j)
+    {
+      // check 1
+      for(unsigned int i = 0; i<3; ++i)
+      {
+        u.block(i)[0] = 0;;
+        v.block(i)[0] = 0;
+      }
+      u.block(j)[0] = 1;;
+
+
+      op_a.vmult(v, u);
+      // identity.Tvmult(v2, u2);
+
+      PRINTME("v", v);
+    }
+
+    deallog << " -- Inverse -- " << std::endl;
+    inverse_op_a.reinit_domain_vector(u, false);
+    inverse_op_a.reinit_range_vector(v, false);
+    for(unsigned int j = 0; j<3; ++j)
+    {
+      // check 1
+      for(unsigned int i = 0; i<3; ++i)
+      {
+        u.block(i)[0] = 0;;
+        v.block(i)[0] = 0;
+      }
+      u.block(j)[0] = 1;;
+
+
+      inverse_op_a.vmult(v, u);
+      // identity.Tvmult(v2, u2);
+
+      PRINTME("v", v);
+    }
+
+    deallog << " -- Identity -- " << std::endl;
+    identity.reinit_domain_vector(u, false);
+    identity.reinit_range_vector(v, false);
+    for(unsigned int j = 0; j<3; ++j)
+    {
+      // check 1
+      for(unsigned int i = 0; i<3; ++i)
+      {
+        u.block(i)[0] = 0;;
+        v.block(i)[0] = 0;
+      }
+      u.block(j)[0] = 1;;
+
+
+      identity.vmult(v, u);
+      // identity.Tvmult(v2, u2);
+
+      PRINTME("v", v);
+    }
+
+    // op_a.vmult(v, u);
+    // // identity.Tvmult(v2, u2);
+    //
+    // PRINTME("v", v);
+    //
+    // inverse_op_a.vmult(v, u);
+    // // identity.Tvmult(v2, u2);
+    //
+    // PRINTME("v", v);
+  }
+}
diff --git a/tests/lac/linear_operator_14.with_cxx11=on.output b/tests/lac/linear_operator_14.with_cxx11=on.output
new file mode 100644 (file)
index 0000000..e69de29

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.