]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Make all template specialization of FPE compile for float 14155/head
authorPeter Munch <peterrmuench@gmail.com>
Wed, 20 Jul 2022 18:52:03 +0000 (20:52 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Thu, 21 Jul 2022 19:06:16 +0000 (21:06 +0200)
include/deal.II/base/derivative_form.h
include/deal.II/matrix_free/fe_point_evaluation.h
tests/matrix_free/point_evaluation_18.cc [new file with mode: 0644]
tests/matrix_free/point_evaluation_18.output [new file with mode: 0644]

index 10d6545ed1c3bd8a3a9bdcda9933dbccaacf04ce..bb4fa64e3c0ea2eb177ff2b7a03d5797add866d8 100644 (file)
@@ -108,12 +108,14 @@ public:
    * Number>. In particular, if order == 1 and the derivative is the Jacobian of
    * $\mathbf F(\mathbf x)$, then Tensor[i] = $\nabla F_i(\mathbf x)$.
    */
-  operator Tensor<order + 1, dim, Number>() const;
+  template <typename Number2>
+  operator Tensor<order + 1, dim, Number2>() const;
 
   /**
    * Converts a DerivativeForm<1, dim, 1, Number> to Tensor<1, dim, Number>.
    */
-  operator Tensor<1, dim, Number>() const;
+  template <typename Number2>
+  operator Tensor<1, dim, Number2>() const;
 
   /**
    * Return the transpose of a rectangular DerivativeForm,
@@ -273,8 +275,9 @@ DerivativeForm<order, dim, spacedim, Number>::operator[](
 
 
 template <int order, int dim, int spacedim, typename Number>
+template <typename Number2>
 inline DerivativeForm<order, dim, spacedim, Number>::
-operator Tensor<1, dim, Number>() const
+operator Tensor<1, dim, Number2>() const
 {
   Assert((1 == spacedim) && (order == 1),
          ExcMessage("Only allowed for spacedim==1."));
@@ -285,12 +288,13 @@ operator Tensor<1, dim, Number>() const
 
 
 template <int order, int dim, int spacedim, typename Number>
+template <typename Number2>
 inline DerivativeForm<order, dim, spacedim, Number>::
-operator Tensor<order + 1, dim, Number>() const
+operator Tensor<order + 1, dim, Number2>() const
 {
   Assert((dim == spacedim), ExcMessage("Only allowed when dim==spacedim."));
 
-  Tensor<order + 1, dim, Number> t;
+  Tensor<order + 1, dim, Number2> t;
 
   if (dim == spacedim)
     for (unsigned int j = 0; j < dim; ++j)
index f8e3269b24b1fd1b78b8c795fbd97ad014944353..4e64cfc2df80ec449ffec76b5f43e08dabea68fb 100644 (file)
@@ -947,15 +947,11 @@ FEPointEvaluation<n_components, dim, spacedim, Number>::evaluate(
                     Number>::set_gradient(val_and_grad.second,
                                           j,
                                           unit_gradients[i + j]);
-                  gradients[i + j] =
-                    static_cast<typename internal::FEPointEvaluation::
-                                  EvaluatorTypeTraits<dim,
-                                                      n_components,
-                                                      Number>::gradient_type>(
-                      apply_transformation(mapping_info->get_mapping_data()
-                                             .inverse_jacobians[i + j]
-                                             .transpose(),
-                                           unit_gradients[i + j]));
+                  gradients[i + j] = static_cast<gradient_type>(
+                    apply_transformation(mapping_info->get_mapping_data()
+                                           .inverse_jacobians[i + j]
+                                           .transpose(),
+                                         unit_gradients[i + j]));
                 }
             }
         }
@@ -1092,9 +1088,7 @@ FEPointEvaluation<n_components, dim, spacedim, Number>::integrate(
             for (unsigned int j = 0; j < n_lanes && i + j < n_points; ++j)
               {
                 gradients[i + j] =
-                  static_cast<typename internal::FEPointEvaluation::
-                                EvaluatorTypeTraits<dim, n_components, Number>::
-                                  gradient_type>(apply_transformation(
+                  static_cast<gradient_type>(apply_transformation(
                     mapping_info->get_mapping_data().inverse_jacobians[i + j],
                     gradients[i + j]));
                 internal::FEPointEvaluation::
diff --git a/tests/matrix_free/point_evaluation_18.cc b/tests/matrix_free/point_evaluation_18.cc
new file mode 100644 (file)
index 0000000..ab96667
--- /dev/null
@@ -0,0 +1,76 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2022 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.
+//
+// ---------------------------------------------------------------------
+
+
+// Check that all template specializations of FEPointEvaluation are
+// compiling.
+
+#include <deal.II/matrix_free/fe_point_evaluation.h>
+
+#include "../tests.h"
+
+template <int n_components, int dim, int spacedim, typename Number>
+void
+test()
+{
+  return; // nothing to do, since we are only interested if the code
+          // compiles
+
+  std::unique_ptr<Mapping<dim, spacedim>>       mapping;
+  std::unique_ptr<FiniteElement<dim, spacedim>> fe;
+
+  FEPointEvaluation<n_components, dim, spacedim, Number> fpe(
+    *mapping, *fe, UpdateFlags::update_default);
+
+  Triangulation<dim, spacedim> tria;
+
+  fpe.reinit(tria.begin(), ArrayView<const Point<dim>>());
+
+  fpe.evaluate(ArrayView<const Number>(), EvaluationFlags::values);
+
+  fpe.integrate(ArrayView<Number>(), EvaluationFlags::values);
+}
+
+int
+main()
+{
+  initlog();
+
+  test<1, 1, 1, double>();
+  test<2, 1, 1, double>();
+
+  test<1, 2, 2, double>();
+  test<2, 2, 2, double>();
+  test<3, 2, 2, double>();
+
+  test<1, 3, 3, double>();
+  test<2, 3, 3, double>();
+  test<3, 3, 3, double>();
+  test<4, 3, 3, double>();
+
+  test<1, 1, 1, float>();
+  test<2, 1, 1, float>();
+
+  test<1, 2, 2, float>();
+  test<2, 2, 2, float>();
+  test<3, 2, 2, float>();
+
+  test<1, 3, 3, float>();
+  test<2, 3, 3, float>();
+  test<3, 3, 3, float>();
+  test<4, 3, 3, float>();
+
+  deallog << "OK!" << std::endl;
+}
diff --git a/tests/matrix_free/point_evaluation_18.output b/tests/matrix_free/point_evaluation_18.output
new file mode 100644 (file)
index 0000000..5cfb783
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL::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.