]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Introduce a move assignment operator for Quadrature 4692/head
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Fri, 4 Aug 2017 10:40:14 +0000 (12:40 +0200)
committerDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Fri, 4 Aug 2017 11:33:41 +0000 (13:33 +0200)
doc/news/changes/minor/20170804DanielArndt [new file with mode: 0644]
include/deal.II/base/quadrature.h
tests/base/quadrature_move_01.cc [moved from tests/base/quadrature_move.cc with 100% similarity]
tests/base/quadrature_move_01.output [moved from tests/base/quadrature_move.output with 100% similarity]
tests/base/quadrature_move_02.cc [new file with mode: 0644]
tests/base/quadrature_move_02.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20170804DanielArndt b/doc/news/changes/minor/20170804DanielArndt
new file mode 100644 (file)
index 0000000..4f536ca
--- /dev/null
@@ -0,0 +1,3 @@
+New: Quadrature has a (defaulted) move assignment operator.
+<br>
+(Daniel Arndt, 2017/08/04)
index 20d0b0b65a93353f7fcd74c76a2efbee27c2f26d..136c82d5f1de20a8f9578b75b153c1aaa46a7885 100644 (file)
@@ -173,6 +173,12 @@ public:
    */
   Quadrature &operator = (const Quadrature<dim> &);
 
+  /**
+   * Move assignment operator. Moves all data from another quadrature object
+   * to this object.
+   */
+  Quadrature &operator = (Quadrature<dim> &&) = default;
+
   /**
    * Test for equality of two quadratures.
    */
diff --git a/tests/base/quadrature_move_02.cc b/tests/base/quadrature_move_02.cc
new file mode 100644 (file)
index 0000000..0c6dd81
--- /dev/null
@@ -0,0 +1,87 @@
+// ---------------------------------------------------------------------
+//
+// 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.
+//
+// ---------------------------------------------------------------------
+
+
+#include "../tests.h"
+#include <fstream>
+#include <cmath>
+
+#include <deal.II/base/quadrature_lib.h>
+
+
+template <template <int dim> class Quad, int dim, typename... Args>
+std::string check_q_assign_move(Args &&...args)
+{
+  Quad<dim> quad1(args...);
+  const unsigned int size1 = quad1.size();
+  const std::vector<double> weights1 = quad1.get_weights();
+  const std::vector<Point<dim> > points1 = quad1.get_points();
+
+  Quadrature<dim> quad2;
+  AssertThrow (quad2.size()==0, ExcInternalError());
+
+  quad2 = std::move(quad1);
+
+  AssertThrow(quad1.size() == 0, ExcInternalError());
+  AssertThrow(quad2.size() == size1, ExcInternalError());
+
+  const std::vector<double> weights2 = quad2.get_weights();
+  const std::vector<Point<dim> > points2 = quad2.get_points();
+  for (unsigned int i=0; i<size1; ++i)
+    {
+      AssertThrow(std::abs(weights1[i]-weights2[i]) < 1.e-16, ExcInternalError());
+      AssertThrow((points1[i]-points2[i]).norm() < 1.e-16, ExcInternalError());
+    }
+
+  return "OK";
+}
+
+
+template <template <int dim> class Quad, typename... Args>
+void check_quadrature_assign_move(Args &&...args)
+{
+  deallog << check_q_assign_move<Quad,1>(std::forward<Args>(args)...) << 1 << " "
+          << check_q_assign_move<Quad,2>(std::forward<Args>(args)...) << 2 << " "
+          << check_q_assign_move<Quad,3>(std::forward<Args>(args)...) << 3
+          << std::endl;
+}
+
+
+int main()
+{
+  initlog();
+
+  check_quadrature_assign_move<QMidpoint>();
+  check_quadrature_assign_move<QTrapez>();
+  check_quadrature_assign_move<QSimpson>();
+  check_quadrature_assign_move<QMilne>();
+  check_quadrature_assign_move<QWeddle>();
+
+  for (unsigned int p = 2; p < 5; ++p)
+    {
+      check_quadrature_assign_move<QGauss>(p);
+      check_quadrature_assign_move<QGaussLobatto>(p);
+    }
+
+  const auto ep = QGaussRadauChebyshev<1>::right;
+  for (unsigned int p = 2; p < 5; ++p)
+    {
+      deallog << "Gauss Log R: " << check_q_assign_move<QGaussLogR,1>(p) << std::endl;
+      deallog << "Gauss Radau Chebyshev: "
+              << check_q_assign_move<QGaussRadauChebyshev,1>(p, ep) << std::endl;
+    }
+
+  return 0;
+}
diff --git a/tests/base/quadrature_move_02.output b/tests/base/quadrature_move_02.output
new file mode 100644 (file)
index 0000000..5e60b40
--- /dev/null
@@ -0,0 +1,18 @@
+
+DEAL::OK1 OK2 OK3
+DEAL::OK1 OK2 OK3
+DEAL::OK1 OK2 OK3
+DEAL::OK1 OK2 OK3
+DEAL::OK1 OK2 OK3
+DEAL::OK1 OK2 OK3
+DEAL::OK1 OK2 OK3
+DEAL::OK1 OK2 OK3
+DEAL::OK1 OK2 OK3
+DEAL::OK1 OK2 OK3
+DEAL::OK1 OK2 OK3
+DEAL::Gauss Log R: OK
+DEAL::Gauss Radau Chebyshev: OK
+DEAL::Gauss Log R: OK
+DEAL::Gauss Radau Chebyshev: OK
+DEAL::Gauss Log R: OK
+DEAL::Gauss Radau Chebyshev: 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.