]> https://gitweb.dealii.org/ - dealii.git/commitdiff
lac/linear_operator.h: Prevent users from creating LinearOperators with temoporaries
authorMatthias Maier <tamiko@43-1.org>
Thu, 14 Feb 2019 19:12:06 +0000 (13:12 -0600)
committerMatthias Maier <tamiko@43-1.org>
Thu, 14 Feb 2019 19:51:54 +0000 (13:51 -0600)
include/deal.II/lac/linear_operator.h
tests/lac/linear_operator_18.cc [new file with mode: 0644]
tests/lac/linear_operator_18.output [new file with mode: 0644]

index f0115ff2578dc3e8ff8713458cd1c76aea8daaf0..64dff9a0996016d55595ec6e44bbfc41c0ad3d55 100644 (file)
@@ -1473,6 +1473,125 @@ linear_operator(const LinearOperator<Range, Domain, Payload> &operator_exemplar,
 
 //@}
 
+#ifndef DOXYGEN
+
+//
+// Ensure that we never capture a reference to a temporary by accident.
+// This ensures that instead of silently allowing a "stack use after free",
+// we at least bail out with a runtime error message that is slightly more
+// understandable.
+//
+
+template <
+  typename Range   = Vector<double>,
+  typename Domain  = Range,
+  typename Payload = internal::LinearOperatorImplementation::EmptyPayload,
+  typename OperatorExemplar,
+  typename Matrix,
+  typename =
+    typename std::enable_if<!std::is_lvalue_reference<Matrix>::value>::type>
+LinearOperator<Range, Domain, Payload>
+linear_operator(const OperatorExemplar &, Matrix &&)
+{
+  Assert(false,
+         ExcMessage(
+           "You are trying to construct a linear operator from a temporary "
+           "matrix object."));
+}
+
+template <
+  typename Range   = Vector<double>,
+  typename Domain  = Range,
+  typename Payload = internal::LinearOperatorImplementation::EmptyPayload,
+  typename OperatorExemplar,
+  typename Matrix,
+  typename = typename std::enable_if<
+    !std::is_lvalue_reference<OperatorExemplar>::value>::type>
+LinearOperator<Range, Domain, Payload>
+linear_operator(OperatorExemplar &&, const Matrix &)
+{
+  Assert(false,
+         ExcMessage(
+           "You are trying to construct a linear operator with a temporary "
+           "operator_exemplar object."));
+}
+
+template <
+  typename Range   = Vector<double>,
+  typename Domain  = Range,
+  typename Payload = internal::LinearOperatorImplementation::EmptyPayload,
+  typename OperatorExemplar,
+  typename Matrix,
+  typename =
+    typename std::enable_if<!std::is_lvalue_reference<Matrix>::value>::type,
+  typename = typename std::enable_if<
+    !std::is_lvalue_reference<OperatorExemplar>::value>::type>
+LinearOperator<Range, Domain, Payload>
+linear_operator(OperatorExemplar &&, Matrix &&)
+{
+  Assert(false,
+         ExcMessage(
+           "You are trying to construct a linear operator from a temporary "
+           "matrix object."));
+}
+
+
+template <
+  typename Range   = Vector<double>,
+  typename Domain  = Range,
+  typename Payload = internal::LinearOperatorImplementation::EmptyPayload,
+  typename Matrix,
+  typename =
+    typename std::enable_if<!std::is_lvalue_reference<Matrix>::value>::type>
+LinearOperator<Range, Domain, Payload>
+linear_operator(const LinearOperator<Range, Domain, Payload> &, Matrix &&)
+{
+  Assert(false,
+         ExcMessage(
+           "You are trying to construct a linear operator from a temporary "
+           "matrix object."));
+}
+
+template <
+  typename Range   = Vector<double>,
+  typename Domain  = Range,
+  typename Payload = internal::LinearOperatorImplementation::EmptyPayload,
+  typename Matrix,
+  typename =
+    typename std::enable_if<!std::is_lvalue_reference<Matrix>::value>::type>
+LinearOperator<Range, Domain, Payload>
+linear_operator(Matrix &&)
+{
+  Assert(false,
+         ExcMessage(
+           "You are trying to construct a linear operator from a temporary "
+           "matrix object."));
+}
+
+template <typename Payload,
+          typename Solver,
+          typename Preconditioner,
+          typename Range  = typename Solver::vector_type,
+          typename Domain = Range,
+          typename        = typename std::enable_if<
+            !std::is_lvalue_reference<Preconditioner>::value>::type,
+          typename = typename std::enable_if<
+            !std::is_same<Preconditioner, PreconditionIdentity>::value>::type,
+          typename = typename std::enable_if<
+            !std::is_same<Preconditioner,
+                          LinearOperator<Range, Domain, Payload>>::value>::type>
+LinearOperator<Domain, Range, Payload>
+inverse_operator(const LinearOperator<Range, Domain, Payload> &,
+                 Solver &,
+                 Preconditioner &&)
+{
+  Assert(false,
+         ExcMessage(
+           "You are trying to construct an inverse operator with a temporary "
+           "preconditioner object."));
+}
+
+#endif // DOXYGEN
 
 DEAL_II_NAMESPACE_CLOSE
 
diff --git a/tests/lac/linear_operator_18.cc b/tests/lac/linear_operator_18.cc
new file mode 100644 (file)
index 0000000..08fbfb1
--- /dev/null
@@ -0,0 +1,55 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2015 - 2018 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.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+// Test that we cannot accidentally create a linear operator from a
+// temporary object
+
+#include <deal.II/lac/linear_operator.h>
+#include <deal.II/lac/solver_cg.h>
+#include <deal.II/lac/sparse_matrix.h>
+
+#include "../tests.h"
+
+#define CATCH(call)                                           \
+  try                                                         \
+    {                                                         \
+      call;                                                   \
+      deallog << "Error: Something went wrong!" << std::endl; \
+    }                                                         \
+  catch (ExceptionBase & e)                                   \
+    {                                                         \
+      deallog << e.get_exc_name() << std::endl;               \
+    }
+
+using namespace dealii;
+
+int
+main()
+{
+  initlog();
+  deal_II_exceptions::disable_abort_on_exception();
+
+  SparseMatrix<double> A;
+  const auto           op_A = linear_operator(A);
+
+  CATCH(linear_operator(SparseMatrix<double>()))
+  CATCH(linear_operator(SparseMatrix<double>(), A))
+  CATCH(linear_operator(A, SparseMatrix<double>()))
+  CATCH(linear_operator(SparseMatrix<double>(), SparseMatrix<double>()))
+
+  SolverControl            solver_control;
+  SolverCG<Vector<double>> solver(solver_control);
+  CATCH(inverse_operator(op_A, solver, SparseMatrix<double>()))
+}
diff --git a/tests/lac/linear_operator_18.output b/tests/lac/linear_operator_18.output
new file mode 100644 (file)
index 0000000..8aad575
--- /dev/null
@@ -0,0 +1,6 @@
+
+DEAL::ExcMessage( "You are trying to construct a linear operator from a temporary " "matrix object.")
+DEAL::ExcMessage( "You are trying to construct a linear operator with a temporary " "operator_exemplar object.")
+DEAL::ExcMessage( "You are trying to construct a linear operator from a temporary " "matrix object.")
+DEAL::ExcMessage( "You are trying to construct a linear operator from a temporary " "matrix object.")
+DEAL::ExcMessage( "You are trying to construct an inverse operator with a temporary " "preconditioner object.")

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.