]> https://gitweb.dealii.org/ - dealii.git/commitdiff
add a new header which combines LinearOperator and PackagedOperation 2505/head
authorDenis Davydov <davydden@gmail.com>
Wed, 13 Apr 2016 07:47:58 +0000 (09:47 +0200)
committerDenis Davydov <davydden@gmail.com>
Wed, 13 Apr 2016 20:12:16 +0000 (22:12 +0200)
doc/doxygen/headers/laoperators.h
include/deal.II/lac/linear_operator_tools.h [new file with mode: 0644]
tests/lac/linear_operator_04a.cc [new file with mode: 0644]
tests/lac/linear_operator_04a.with_cxx11=on.with_trilinos=on.output [new file with mode: 0644]

index 1af479d6cf0b0bf3716324f489578ab5a9cf5f1f..939fe8e3f63bcf4e43311e0e5688f31861cac260 100644 (file)
  * InverseMatrixRichardson, SchurMatrix, ShiftedMatrix,
  * ShiftedMatrixGeneralized, TransposeMatrix
  *
+ * @note As explained below, when using LinearOperator as <code>res = op_a*x</code>
+ * a PackagedOperation class instance is generated behind-the-curtains.
+ * Consequently, the user program has to include header files for both classes
+ * for compilation to be successful. In an attempt to make easier the decision of which
+ * headers to include in what circumstances and to prevent hidden templates-related
+ * compiler errors, all headers relevant to LinearOperator are grouped in
+ * <deal.ii/lac/linear_operator_tools.h>.
  *
  * <h3>Packaged Operation</h3>
  *
diff --git a/include/deal.II/lac/linear_operator_tools.h b/include/deal.II/lac/linear_operator_tools.h
new file mode 100644 (file)
index 0000000..1f43fba
--- /dev/null
@@ -0,0 +1,28 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2016 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.
+//
+// ---------------------------------------------------------------------
+
+#ifndef dealii__linear_operator_tools_h
+#define dealii__linear_operator_tools_h
+
+// many usage cases lead to a combination of LinearOperator and
+// PackagedOperation. To ease the pain of reading compilation errors, just include
+// all headers we ever need to use LO and friends in one place:
+#include <deal.II/lac/linear_operator.h>
+#include <deal.II/lac/block_linear_operator.h>
+#include <deal.II/lac/packaged_operation.h>
+#include <deal.II/lac/constrained_linear_operator.h>
+#include <deal.II/lac/schur_complement.h>
+
+#endif
diff --git a/tests/lac/linear_operator_04a.cc b/tests/lac/linear_operator_04a.cc
new file mode 100644 (file)
index 0000000..8b4e9f1
--- /dev/null
@@ -0,0 +1,49 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2016 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 that it is possible to use a operator*() of LinearOperator object for
+// Trilinos matrices and vectors
+
+#include "../tests.h"
+
+#include <deal.II/lac/linear_operator_tools.h>
+#include <deal.II/lac/trilinos_block_sparse_matrix.h>
+#include <deal.II/lac/trilinos_sparse_matrix.h>
+#include <deal.II/lac/trilinos_vector.h>
+
+using namespace dealii;
+
+int main(int argc, char *argv[])
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization (argc, argv, 1);
+
+  initlog();
+  deallog << std::setprecision(10);
+
+  typedef TrilinosWrappers::MPI::Vector vector_t;
+  typedef TrilinosWrappers::SparseMatrix matrix_t;
+
+  matrix_t a;
+
+  auto op_a  = linear_operator<vector_t>(a);
+  vector_t u,res;
+  res = op_a * u;
+  // ^^ this was not working, whereas
+  // op_a.vmult(res,u) did.
+
+  deallog << "OK" << std::endl;
+
+  return 0;
+}
diff --git a/tests/lac/linear_operator_04a.with_cxx11=on.with_trilinos=on.output b/tests/lac/linear_operator_04a.with_cxx11=on.with_trilinos=on.output
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.