]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Added get_names to Extractors, and moved implementation to cc file.
authorLuca Heltai <luca.heltai@sissa.it>
Tue, 19 Mar 2019 18:32:08 +0000 (19:32 +0100)
committerLuca Heltai <luca.heltai@sissa.it>
Thu, 21 Mar 2019 23:56:02 +0000 (00:56 +0100)
doc/news/changes/minor/20190319LucaHeltai [new file with mode: 0644]
include/deal.II/fe/fe_values_extractors.h
source/fe/CMakeLists.txt
source/fe/fe_values_extractors.cc [new file with mode: 0644]
tests/fe/fe_values_extractor_02.cc [new file with mode: 0644]
tests/fe/fe_values_extractor_02.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20190319LucaHeltai b/doc/news/changes/minor/20190319LucaHeltai
new file mode 100644 (file)
index 0000000..5328898
--- /dev/null
@@ -0,0 +1,4 @@
+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)
index db3729129af05020990efbafed054e9637bd3920..a7e1944bd3b1bb0440466bb3af5dc07ad8d3d4e1 100644 (file)
@@ -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
 
index 79a5bb8466f360fe9dea11ccc410077de89fd0eb..fbc00aeafc050323f0d0762daee8a7a51a2c30a3 100644 (file)
@@ -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 (file)
index 0000000..8d0c417
--- /dev/null
@@ -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 <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
diff --git a/tests/fe/fe_values_extractor_02.cc b/tests/fe/fe_values_extractor_02.cc
new file mode 100644 (file)
index 0000000..c48a224
--- /dev/null
@@ -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 <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;
+}
diff --git a/tests/fe/fe_values_extractor_02.output b/tests/fe/fe_values_extractor_02.output
new file mode 100644 (file)
index 0000000..1b5a914
--- /dev/null
@@ -0,0 +1,5 @@
+
+DEAL::Scalar(0)
+DEAL::Vector(1)
+DEAL::Tensor<2>(3)
+DEAL::SymmetricTensor<2>(4)

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.