--- /dev/null
+New: FEValuesExtractor classes now have a new method get_name() that returns a unique string identifier
+for each extractor.
+<br>
+(Luca Heltai, 2019/03/19)
* 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;
};
* 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;
};
* 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;
};
* 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
SET(_separate_src
fe_values.cc
+ fe_values_extractors.cc
fe_values_inst2.cc
fe_values_inst3.cc
fe_values_inst4.cc
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// 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 <deal.II/base/utilities.h>
+
+#include <deal.II/fe/fe_values_extractors.h>
+
+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 <int rank>
+ std::string
+ Tensor<rank>::get_name() const
+ {
+ return "Tensor<" + Utilities::int_to_string(rank) + ">(" +
+ Utilities::int_to_string(first_tensor_component) + ")";
+ }
+
+
+ template <int rank>
+ std::string
+ SymmetricTensor<rank>::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
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// 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 <deal.II/fe/fe_values.h>
+#include <deal.II/fe/fe_values_extractors.h>
+
+#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;
+}
--- /dev/null
+
+DEAL::Scalar(0)
+DEAL::Vector(1)
+DEAL::Tensor<2>(3)
+DEAL::SymmetricTensor<2>(4)