]> https://gitweb.dealii.org/ - dealii.git/commitdiff
make DiagonalMatrix usable as a linear operator
authorDenis Davydov <davydden@gmail.com>
Tue, 1 Aug 2017 20:28:05 +0000 (22:28 +0200)
committerDenis Davydov <davydden@gmail.com>
Tue, 1 Aug 2017 20:28:05 +0000 (22:28 +0200)
doc/news/changes/minor/20170801DenisDavydov-b [new file with mode: 0644]
include/deal.II/lac/diagonal_matrix.h
tests/lac/diagonal_matrix_03.cc [new file with mode: 0644]
tests/lac/diagonal_matrix_03.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20170801DenisDavydov-b b/doc/news/changes/minor/20170801DenisDavydov-b
new file mode 100644 (file)
index 0000000..31ba64f
--- /dev/null
@@ -0,0 +1,3 @@
+Fixed: DiagonalMatrix<VectorType> can now be used as a linear operator.
+<br>
+(Denis Davydov, 2017/08/01)
index 54966bb2937af3f0563f4329850f4810ba9415a5..621ece07099be4e6410e2fc6741a99187381856f 100644 (file)
@@ -178,6 +178,12 @@ public:
   void Tvmult_add (VectorType       &dst,
                    const VectorType &src) const;
 
+  /**
+   * Initialize vector @p dst. This is a part of the interface required
+   * by linear operator.
+   */
+  void initialize_dof_vector(VectorType &dst) const;
+
   /**
    * Return the memory consumption of this object.
    */
@@ -221,6 +227,15 @@ DiagonalMatrix<VectorType>::reinit(const VectorType &vec)
 
 
 
+template <typename VectorType>
+void
+DiagonalMatrix<VectorType>::initialize_dof_vector(VectorType &dst) const
+{
+  dst.reinit(diagonal);
+}
+
+
+
 template <typename VectorType>
 void
 DiagonalMatrix<VectorType>::compress (VectorOperation::values operation)
diff --git a/tests/lac/diagonal_matrix_03.cc b/tests/lac/diagonal_matrix_03.cc
new file mode 100644 (file)
index 0000000..18696b3
--- /dev/null
@@ -0,0 +1,108 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2017 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.
+//
+// ---------------------------------------------------------------------
+
+
+// Make sure DiagonalMatrix can be used with linear_operator
+// Otherwise this test is exact copy of diagonal_matrix.cc
+
+
+#include "../tests.h"
+#include <deal.II/base/logstream.h>
+#include <deal.II/lac/diagonal_matrix.h>
+#include <deal.II/lac/vector.h>
+#include <deal.II/lac/linear_operator.h>
+
+
+
+void
+check()
+{
+  const unsigned int size = 10;
+  Vector<double> vec(size);
+  for (unsigned int i=0; i<size; ++i)
+    vec(i) = 2;
+  DiagonalMatrix<Vector<double> > mat;
+  mat.reinit(vec);
+  const auto op = linear_operator<Vector<double>>(mat);
+
+  Vector<double> in (size), out(size), exact(size);
+  for (unsigned int i=0; i<size; ++i)
+    in(i) = 0.5;
+  exact = 1;
+
+  op.vmult(out, in);
+  out -= exact;
+  deallog << "Error vmult: " << out.linfty_norm() << std::endl;
+
+  op.Tvmult(out, in);
+  out -= exact;
+  deallog << "Error Tvmult: " << out.linfty_norm() << std::endl;
+
+  out = exact;
+  op.vmult_add(out, in);
+  out -= exact;
+  out -= exact;
+  deallog << "Error vmult_add: " << out.linfty_norm() << std::endl;
+
+  out = exact;
+  op.Tvmult_add(out, in);
+  out -= exact;
+  out -= exact;
+  deallog << "Error Tvmult_add: " << out.linfty_norm() << std::endl;
+
+  for (unsigned int i=0; i<size; ++i)
+    {
+      mat(i,i) = i+1;
+      in(i) = 1./mat(i,i);
+    }
+
+  op.vmult(out, in);
+  out -= exact;
+  deallog << "Error vmult set 1: " << out.linfty_norm() << std::endl;
+
+  for (unsigned int i=0; i<size; ++i)
+    {
+      mat.get_vector()(i) = 1./(i+1);
+      in(i) = i+1;
+    }
+  op.vmult(out, in);
+  out -= exact;
+  deallog << "Error vmult set 2: " << out.linfty_norm() << std::endl;
+
+  for (unsigned int i=0; i<size; ++i)
+    {
+      const double mat_entry =
+        const_cast<const DiagonalMatrix<Vector<double> > &>(mat)(i,i);
+      mat(i,i) = in(i);
+      in(i) = mat_entry;
+    }
+  op.vmult(out, in);
+  out -= exact;
+  deallog << "Error vmult set 3: " << out.linfty_norm() << std::endl;
+}
+
+
+int main()
+{
+  std::ofstream logfile("output");
+  deallog << std::fixed;
+  deallog << std::setprecision(2);
+  deallog.attach(logfile);
+  deallog.threshold_double(1.e-10);
+
+  check();
+
+  return 0;
+}
diff --git a/tests/lac/diagonal_matrix_03.output b/tests/lac/diagonal_matrix_03.output
new file mode 100644 (file)
index 0000000..97ce615
--- /dev/null
@@ -0,0 +1,8 @@
+
+DEAL::Error vmult: 0
+DEAL::Error Tvmult: 0
+DEAL::Error vmult_add: 0
+DEAL::Error Tvmult_add: 0
+DEAL::Error vmult set 1: 0
+DEAL::Error vmult set 2: 0
+DEAL::Error vmult set 3: 0

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.