From 5fdaaaf6c1ee8ff5152dff86bb7fec3a7ea6ce72 Mon Sep 17 00:00:00 2001 From: Luca Heltai Date: Tue, 19 Mar 2019 19:32:08 +0100 Subject: [PATCH] Added get_names to Extractors, and moved implementation to cc file. --- doc/news/changes/minor/20190319LucaHeltai | 4 ++ include/deal.II/fe/fe_values_extractors.h | 24 ++++++++ source/fe/CMakeLists.txt | 1 + source/fe/fe_values_extractors.cc | 67 +++++++++++++++++++++++ tests/fe/fe_values_extractor_02.cc | 40 ++++++++++++++ tests/fe/fe_values_extractor_02.output | 5 ++ 6 files changed, 141 insertions(+) create mode 100644 doc/news/changes/minor/20190319LucaHeltai create mode 100644 source/fe/fe_values_extractors.cc create mode 100644 tests/fe/fe_values_extractor_02.cc create mode 100644 tests/fe/fe_values_extractor_02.output diff --git a/doc/news/changes/minor/20190319LucaHeltai b/doc/news/changes/minor/20190319LucaHeltai new file mode 100644 index 0000000000..5328898f64 --- /dev/null +++ b/doc/news/changes/minor/20190319LucaHeltai @@ -0,0 +1,4 @@ +New: FEValuesExtractor classes now have a new method get_name() that returns a unique string identifier +for each extractor. +
+(Luca Heltai, 2019/03/19) diff --git a/include/deal.II/fe/fe_values_extractors.h b/include/deal.II/fe/fe_values_extractors.h index db3729129a..a7e1944bd3 100644 --- a/include/deal.II/fe/fe_values_extractors.h +++ b/include/deal.II/fe/fe_values_extractors.h @@ -112,6 +112,12 @@ namespace FEValuesExtractors * Constructor. Take the selected vector component as argument. */ Scalar(const unsigned int component); + + /** + * Return a string that uniquely identifies this finite element extractor. + */ + std::string + get_name() const; }; @@ -162,6 +168,12 @@ namespace FEValuesExtractors * FEValues object as argument. */ Vector(const unsigned int first_vector_component); + + /** + * Return a string that uniquely identifies this finite element extractor. + */ + std::string + get_name() const; }; @@ -205,6 +217,12 @@ namespace FEValuesExtractors * FEValues object as argument. */ SymmetricTensor(const unsigned int first_tensor_component); + + /** + * Return a string that uniquely identifies this finite element extractor. + */ + std::string + get_name() const; }; @@ -248,6 +266,12 @@ namespace FEValuesExtractors * FEValues object as argument. */ Tensor(const unsigned int first_tensor_component); + + /** + * Return a string that uniquely identifies this finite element extractor. + */ + std::string + get_name() const; }; } // namespace FEValuesExtractors diff --git a/source/fe/CMakeLists.txt b/source/fe/CMakeLists.txt index 79a5bb8466..fbc00aeafc 100644 --- a/source/fe/CMakeLists.txt +++ b/source/fe/CMakeLists.txt @@ -63,6 +63,7 @@ SET(_unity_include_src SET(_separate_src fe_values.cc + fe_values_extractors.cc fe_values_inst2.cc fe_values_inst3.cc fe_values_inst4.cc diff --git a/source/fe/fe_values_extractors.cc b/source/fe/fe_values_extractors.cc new file mode 100644 index 0000000000..8d0c417f6e --- /dev/null +++ b/source/fe/fe_values_extractors.cc @@ -0,0 +1,67 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2019 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. +// +// --------------------------------------------------------------------- + +#include + +#include + +DEAL_II_NAMESPACE_OPEN + +namespace FEValuesExtractors +{ + std::string + Scalar::get_name() const + { + return "Scalar(" + Utilities::int_to_string(component) + ")"; + } + + + std::string + Vector::get_name() const + { + return "Vector(" + Utilities::int_to_string(first_vector_component) + ")"; + } + + + template + std::string + Tensor::get_name() const + { + return "Tensor<" + Utilities::int_to_string(rank) + ">(" + + Utilities::int_to_string(first_tensor_component) + ")"; + } + + + template + std::string + SymmetricTensor::get_name() const + { + return "SymmetricTensor<" + Utilities::int_to_string(rank) + ">(" + + Utilities::int_to_string(first_tensor_component) + ")"; + } + + // Explicit instantiations + template struct Tensor<1>; + template struct Tensor<2>; + template struct Tensor<3>; + template struct Tensor<4>; + template struct SymmetricTensor<2>; + template struct SymmetricTensor<4>; + +} // namespace FEValuesExtractors + + + +DEAL_II_NAMESPACE_CLOSE \ No newline at end of file diff --git a/tests/fe/fe_values_extractor_02.cc b/tests/fe/fe_values_extractor_02.cc new file mode 100644 index 0000000000..c48a224ced --- /dev/null +++ b/tests/fe/fe_values_extractor_02.cc @@ -0,0 +1,40 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2019 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 the FEValuesExtractors::get_name() work as expected + +#include +#include + +#include "../tests.h" + +int +main() +{ + initlog(); + + FEValuesExtractors::Scalar scalar(0); + FEValuesExtractors::Vector vector(1); + FEValuesExtractors::Tensor<2> tensor(3); + FEValuesExtractors::SymmetricTensor<2> symmetric_tensor(4); + + + deallog << scalar.get_name() << std::endl + << vector.get_name() << std::endl + << tensor.get_name() << std::endl + << symmetric_tensor.get_name() << std::endl; +} diff --git a/tests/fe/fe_values_extractor_02.output b/tests/fe/fe_values_extractor_02.output new file mode 100644 index 0000000000..1b5a914406 --- /dev/null +++ b/tests/fe/fe_values_extractor_02.output @@ -0,0 +1,5 @@ + +DEAL::Scalar(0) +DEAL::Vector(1) +DEAL::Tensor<2>(3) +DEAL::SymmetricTensor<2>(4) -- 2.39.5