]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Delete the disallowed functions instead 7727/head
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Thu, 14 Feb 2019 21:50:50 +0000 (22:50 +0100)
committerDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Thu, 14 Feb 2019 21:53:01 +0000 (22:53 +0100)
include/deal.II/lac/linear_operator.h
tests/lac/linear_operator_18.cc [deleted file]
tests/lac/linear_operator_18.output [deleted file]

index 64dff9a0996016d55595ec6e44bbfc41c0ad3d55..b8e118fe112cf4a69843e8e76c217308c63ef64a 100644 (file)
@@ -1477,9 +1477,7 @@ linear_operator(const LinearOperator<Range, Domain, Payload> &operator_exemplar,
 
 //
 // 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.
+// to avoid "stack use after free".
 //
 
 template <
@@ -1491,13 +1489,7 @@ template <
   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."));
-}
+linear_operator(const OperatorExemplar &, Matrix &&) = delete;
 
 template <
   typename Range   = Vector<double>,
@@ -1506,15 +1498,12 @@ template <
   typename OperatorExemplar,
   typename Matrix,
   typename = typename std::enable_if<
-    !std::is_lvalue_reference<OperatorExemplar>::value>::type>
+    !std::is_lvalue_reference<OperatorExemplar>::value>::type,
+  typename = typename std::enable_if<
+    !std::is_same<OperatorExemplar,
+                  LinearOperator<Range, Domain, Payload>>::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."));
-}
+linear_operator(OperatorExemplar &&, const Matrix &) = delete;
 
 template <
   typename Range   = Vector<double>,
@@ -1525,16 +1514,12 @@ template <
   typename =
     typename std::enable_if<!std::is_lvalue_reference<Matrix>::value>::type,
   typename = typename std::enable_if<
-    !std::is_lvalue_reference<OperatorExemplar>::value>::type>
+    !std::is_lvalue_reference<OperatorExemplar>::value>::type,
+  typename = typename std::enable_if<
+    !std::is_same<OperatorExemplar,
+                  LinearOperator<Range, Domain, Payload>>::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."));
-}
-
+linear_operator(OperatorExemplar &&, Matrix &&) = delete;
 
 template <
   typename Range   = Vector<double>,
@@ -1544,13 +1529,8 @@ template <
   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."));
-}
+linear_operator(const LinearOperator<Range, Domain, Payload> &,
+                Matrix &&) = delete;
 
 template <
   typename Range   = Vector<double>,
@@ -1560,13 +1540,7 @@ template <
   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."));
-}
+linear_operator(Matrix &&) = delete;
 
 template <typename Payload,
           typename Solver,
@@ -1583,13 +1557,7 @@ template <typename Payload,
 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."));
-}
+                 Preconditioner &&) = delete;
 
 #endif // DOXYGEN
 
diff --git a/tests/lac/linear_operator_18.cc b/tests/lac/linear_operator_18.cc
deleted file mode 100644 (file)
index 08fbfb1..0000000
+++ /dev/null
@@ -1,55 +0,0 @@
-// ---------------------------------------------------------------------
-//
-// 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
deleted file mode 100644 (file)
index 8aad575..0000000
+++ /dev/null
@@ -1,6 +0,0 @@
-
-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.