From 559ec5ff44bef632cb9d944cb8ee179ea650d2b9 Mon Sep 17 00:00:00 2001 From: Daniel Arndt Date: Mon, 11 Sep 2017 17:38:50 +0200 Subject: [PATCH] Implement move constructors for FiniteElement and FESystem --- doc/news/changes/minor/20170911DanielArndt | 3 + include/deal.II/fe/fe.h | 10 ++++ include/deal.II/fe/fe_system.h | 5 ++ tests/fe/fe_move_01.cc | 58 ++++++++++++++++++++ tests/fe/fe_move_01.output | 64 ++++++++++++++++++++++ tests/fe/fe_move_02.cc | 61 +++++++++++++++++++++ tests/fe/fe_move_02.output | 64 ++++++++++++++++++++++ 7 files changed, 265 insertions(+) create mode 100644 doc/news/changes/minor/20170911DanielArndt create mode 100644 tests/fe/fe_move_01.cc create mode 100644 tests/fe/fe_move_01.output create mode 100644 tests/fe/fe_move_02.cc create mode 100644 tests/fe/fe_move_02.output diff --git a/doc/news/changes/minor/20170911DanielArndt b/doc/news/changes/minor/20170911DanielArndt new file mode 100644 index 0000000000..508dd0dd4f --- /dev/null +++ b/doc/news/changes/minor/20170911DanielArndt @@ -0,0 +1,3 @@ +New: FiniteElement and FESystem have move constructors. +
+(Daniel Arndt, 2017/09/11) diff --git a/include/deal.II/fe/fe.h b/include/deal.II/fe/fe.h index 1639b565b3..2511949f0e 100644 --- a/include/deal.II/fe/fe.h +++ b/include/deal.II/fe/fe.h @@ -754,6 +754,16 @@ public: const std::vector &restriction_is_additive_flags, const std::vector &nonzero_components); + /** + * Move constructor. + */ + FiniteElement (FiniteElement &&) = default; + + /** + * Copy constructor. + */ + FiniteElement (const FiniteElement &) = default; + /** * Virtual destructor. Makes sure that pointers to this class are deleted * properly. diff --git a/include/deal.II/fe/fe_system.h b/include/deal.II/fe/fe_system.h index e439fb71b1..6f0efeb9d0 100644 --- a/include/deal.II/fe/fe_system.h +++ b/include/deal.II/fe/fe_system.h @@ -427,6 +427,11 @@ public: */ FESystem (const FESystem &) = delete; + /** + * Move constructor. + */ + FESystem (FESystem &&) = default; + /** * Destructor. */ diff --git a/tests/fe/fe_move_01.cc b/tests/fe/fe_move_01.cc new file mode 100644 index 0000000000..e64ff9b3ea --- /dev/null +++ b/tests/fe/fe_move_01.cc @@ -0,0 +1,58 @@ +// --------------------------------------------------------------------- +// +// 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. +// +// --------------------------------------------------------------------- + + +// check that we can move a FE_Q object in a reasonable way. + +#include "../tests.h" +#include + +#include + +template +void +check(const FiniteElement &fe) +{ + deallog << "dim: " << dim << std::endl; + deallog << "components: " << fe.n_components() << std::endl; + deallog << "blocks: " << fe.n_blocks() << std::endl; + deallog << "conforms H1: " << fe.conforms(FiniteElementData::H1) << std::endl; + deallog << "n_base_elements: " << fe.n_base_elements() << std::endl; + deallog << "memory consumption: " << fe.memory_consumption() << std::endl; + deallog << std::endl; +} + +template +void +move() +{ + FE_Q fe(1); + check(fe); + FE_Q fe2 (std::move(fe)); + check(fe); + check(fe2); +} + +int +main() +{ + initlog(); + + move<1>(); + move<2>(); + move<3>(); + + return 0; +} diff --git a/tests/fe/fe_move_01.output b/tests/fe/fe_move_01.output new file mode 100644 index 0000000000..cf83559341 --- /dev/null +++ b/tests/fe/fe_move_01.output @@ -0,0 +1,64 @@ + +DEAL::dim: 1 +DEAL::components: 1 +DEAL::blocks: 1 +DEAL::conforms H1: 1 +DEAL::n_base_elements: 1 +DEAL::memory consumption: 1368 +DEAL:: +DEAL::dim: 1 +DEAL::components: 1 +DEAL::blocks: 1 +DEAL::conforms H1: 1 +DEAL::n_base_elements: 0 +DEAL::memory consumption: 736 +DEAL:: +DEAL::dim: 1 +DEAL::components: 1 +DEAL::blocks: 1 +DEAL::conforms H1: 1 +DEAL::n_base_elements: 1 +DEAL::memory consumption: 1368 +DEAL:: +DEAL::dim: 2 +DEAL::components: 1 +DEAL::blocks: 1 +DEAL::conforms H1: 1 +DEAL::n_base_elements: 1 +DEAL::memory consumption: 3260 +DEAL:: +DEAL::dim: 2 +DEAL::components: 1 +DEAL::blocks: 1 +DEAL::conforms H1: 1 +DEAL::n_base_elements: 0 +DEAL::memory consumption: 920 +DEAL:: +DEAL::dim: 2 +DEAL::components: 1 +DEAL::blocks: 1 +DEAL::conforms H1: 1 +DEAL::n_base_elements: 1 +DEAL::memory consumption: 3260 +DEAL:: +DEAL::dim: 3 +DEAL::components: 1 +DEAL::blocks: 1 +DEAL::conforms H1: 1 +DEAL::n_base_elements: 1 +DEAL::memory consumption: 8692 +DEAL:: +DEAL::dim: 3 +DEAL::components: 1 +DEAL::blocks: 1 +DEAL::conforms H1: 1 +DEAL::n_base_elements: 0 +DEAL::memory consumption: 1288 +DEAL:: +DEAL::dim: 3 +DEAL::components: 1 +DEAL::blocks: 1 +DEAL::conforms H1: 1 +DEAL::n_base_elements: 1 +DEAL::memory consumption: 8692 +DEAL:: diff --git a/tests/fe/fe_move_02.cc b/tests/fe/fe_move_02.cc new file mode 100644 index 0000000000..1e408d5ff7 --- /dev/null +++ b/tests/fe/fe_move_02.cc @@ -0,0 +1,61 @@ +// --------------------------------------------------------------------- +// +// 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. +// +// --------------------------------------------------------------------- + +// check that we can move a FESystem object in a reasonable way. + +#include "../tests.h" +#include +#include + +#include + +template +void +check(const FiniteElement &fe) +{ + deallog << fe.get_name() << std::endl; + deallog << "components: " << fe.n_components() << std::endl; + deallog << "blocks: " << fe.n_blocks() << std::endl; + deallog << "conforms H1: " << fe.conforms(FiniteElementData::H1) << std::endl; + deallog << "n_base_elements: " << fe.n_base_elements() << std::endl; + deallog << "memory consumption: " << fe.memory_consumption() << std::endl; + deallog << std::endl; +} + +template +void +move() +{ + FESystem fe (FE_Q(2), dim, FE_Q(1), 1); + check(fe); + FESystem fe2 (std::move(fe)); + check(fe); + check(fe2); +} + +int +main() +{ + initlog(); + + move<1>(); + move<2>(); + move<3>(); + + return 0; +} + + + diff --git a/tests/fe/fe_move_02.output b/tests/fe/fe_move_02.output new file mode 100644 index 0000000000..91e1fd5c3f --- /dev/null +++ b/tests/fe/fe_move_02.output @@ -0,0 +1,64 @@ + +DEAL::FESystem<1>[FE_Q<1>(2)-FE_Q<1>(1)] +DEAL::components: 2 +DEAL::blocks: 2 +DEAL::conforms H1: 1 +DEAL::n_base_elements: 2 +DEAL::memory consumption: 4608 +DEAL:: +DEAL::FESystem<1>[] +DEAL::components: 2 +DEAL::blocks: 2 +DEAL::conforms H1: 1 +DEAL::n_base_elements: 0 +DEAL::memory consumption: 1036 +DEAL:: +DEAL::FESystem<1>[FE_Q<1>(2)-FE_Q<1>(1)] +DEAL::components: 2 +DEAL::blocks: 2 +DEAL::conforms H1: 1 +DEAL::n_base_elements: 2 +DEAL::memory consumption: 4608 +DEAL:: +DEAL::FESystem<2>[FE_Q<2>(2)^2-FE_Q<2>(1)] +DEAL::components: 3 +DEAL::blocks: 3 +DEAL::conforms H1: 1 +DEAL::n_base_elements: 2 +DEAL::memory consumption: 13032 +DEAL:: +DEAL::FESystem<2>[] +DEAL::components: 3 +DEAL::blocks: 3 +DEAL::conforms H1: 1 +DEAL::n_base_elements: 0 +DEAL::memory consumption: 2600 +DEAL:: +DEAL::FESystem<2>[FE_Q<2>(2)^2-FE_Q<2>(1)] +DEAL::components: 3 +DEAL::blocks: 3 +DEAL::conforms H1: 1 +DEAL::n_base_elements: 2 +DEAL::memory consumption: 13032 +DEAL:: +DEAL::FESystem<3>[FE_Q<3>(2)^3-FE_Q<3>(1)] +DEAL::components: 4 +DEAL::blocks: 4 +DEAL::conforms H1: 1 +DEAL::n_base_elements: 2 +DEAL::memory consumption: 56040 +DEAL:: +DEAL::FESystem<3>[] +DEAL::components: 4 +DEAL::blocks: 4 +DEAL::conforms H1: 1 +DEAL::n_base_elements: 0 +DEAL::memory consumption: 8772 +DEAL:: +DEAL::FESystem<3>[FE_Q<3>(2)^3-FE_Q<3>(1)] +DEAL::components: 4 +DEAL::blocks: 4 +DEAL::conforms H1: 1 +DEAL::n_base_elements: 2 +DEAL::memory consumption: 56040 +DEAL:: -- 2.39.5