From: Denis Davydov Date: Wed, 13 Apr 2016 07:47:58 +0000 (+0200) Subject: add a new header which combines LinearOperator and PackagedOperation X-Git-Tag: v8.5.0-rc1~1041^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=5f570efc0ba6312aeda8af9a98efa8d312d65408;p=dealii.git add a new header which combines LinearOperator and PackagedOperation --- diff --git a/doc/doxygen/headers/laoperators.h b/doc/doxygen/headers/laoperators.h index 1af479d6cf..939fe8e3f6 100644 --- a/doc/doxygen/headers/laoperators.h +++ b/doc/doxygen/headers/laoperators.h @@ -100,6 +100,13 @@ * InverseMatrixRichardson, SchurMatrix, ShiftedMatrix, * ShiftedMatrixGeneralized, TransposeMatrix * + * @note As explained below, when using LinearOperator as res = op_a*x + * 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 + * . * *

Packaged Operation

* diff --git a/include/deal.II/lac/linear_operator_tools.h b/include/deal.II/lac/linear_operator_tools.h new file mode 100644 index 0000000000..1f43fba560 --- /dev/null +++ b/include/deal.II/lac/linear_operator_tools.h @@ -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 +#include +#include +#include +#include + +#endif diff --git a/tests/lac/linear_operator_04a.cc b/tests/lac/linear_operator_04a.cc new file mode 100644 index 0000000000..8b4e9f1ea1 --- /dev/null +++ b/tests/lac/linear_operator_04a.cc @@ -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 +#include +#include +#include + +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(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 index 0000000000..0fd8fc12f0 --- /dev/null +++ b/tests/lac/linear_operator_04a.with_cxx11=on.with_trilinos=on.output @@ -0,0 +1,2 @@ + +DEAL::OK