]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add FEEvaluation::fast_evaluation_supported() 12128/head
authorPeter Munch <peterrmuench@gmail.com>
Mon, 3 May 2021 19:17:24 +0000 (21:17 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Tue, 4 May 2021 20:10:19 +0000 (22:10 +0200)
include/deal.II/matrix_free/evaluation_template_factory.h
include/deal.II/matrix_free/evaluation_template_factory.templates.h
include/deal.II/matrix_free/fe_evaluation.h
tests/matrix_free/fe_evaluation_print_01.cc [new file with mode: 0644]
tests/matrix_free/fe_evaluation_print_01.output [new file with mode: 0644]

index 3052af2dffd40d673753dbcc538558fdc699a795..35dc8ccb1dcb5f0185a89c0f46827d17b337c819 100644 (file)
@@ -60,6 +60,10 @@ namespace internal
       VectorizedArrayType *gradients_quad,
       VectorizedArrayType *scratch_data,
       const bool           sum_into_values_array);
+
+    static bool
+    fast_evaluation_supported(const unsigned int given_degree,
+                              const unsigned int n_q_points_1d);
   };
 
 
@@ -143,6 +147,10 @@ namespace internal
       const std::array<unsigned int, VectorizedArrayType::size()>
                                     face_orientations,
       const Table<2, unsigned int> &orientation_map);
+
+    static bool
+    fast_evaluation_supported(const unsigned int given_degree,
+                              const unsigned int n_q_points_1d);
   };
 
 
index 1a5c2548999eaa2c584e043631eda2a0faf6aafe..67a3deac9bc4e675451338c77a2c9809072fd2e4 100644 (file)
@@ -63,6 +63,16 @@ namespace internal
       return EvaluatorType::template run<-1, 0>(args...);
   }
 
+  struct FastEvaluationSupported
+  {
+    template <int fe_degree, int n_q_points_1d>
+    static bool
+    run()
+    {
+      return fe_degree != -1;
+    }
+  };
+
 
 
   template <int dim, typename Number, typename VectorizedArrayType>
@@ -123,6 +133,18 @@ namespace internal
 
 
 
+  template <int dim, typename Number, typename VectorizedArrayType>
+  bool
+  FEEvaluationFactory<dim, Number, VectorizedArrayType>::
+    fast_evaluation_supported(const unsigned int given_degree,
+                              const unsigned int n_q_points_1d)
+  {
+    return instantiation_helper_run<1, FastEvaluationSupported>(given_degree,
+                                                                n_q_points_1d);
+  }
+
+
+
   template <int dim, typename Number, typename VectorizedArrayType>
   void
   FEFaceEvaluationFactory<dim, Number, VectorizedArrayType>::evaluate(
@@ -307,6 +329,18 @@ namespace internal
 
 
 
+  template <int dim, typename Number, typename VectorizedArrayType>
+  bool
+  FEFaceEvaluationFactory<dim, Number, VectorizedArrayType>::
+    fast_evaluation_supported(const unsigned int given_degree,
+                              const unsigned int n_q_points_1d)
+  {
+    return instantiation_helper_run<1, FastEvaluationSupported>(given_degree,
+                                                                n_q_points_1d);
+  }
+
+
+
   template <int dim, typename Number, typename VectorizedArrayType>
   void
   CellwiseInverseMassFactory<dim, Number, VectorizedArrayType>::apply(
index 5f4be8d93b40474053538285e29243c5953891c9..024ae5734ae5fac324fad5ba9d47254c1692c4b9 100644 (file)
@@ -2320,7 +2320,10 @@ protected:
  * setting the macro `FE_EVAL_FACTORY_DEGREE_MAX` to the desired integer and
  * instantiating the classes FEEvaluationFactory and FEFaceEvaluationFactory
  * (the latter for FEFaceEvaluation) creates paths to templated functions for
- * a possibly larger set of degrees.
+ * a possibly larger set of degrees. You can check if fast
+ * evaluation/integration for a given degree/n_quadrature_points pair by calling
+ * FEEvaluation::fast_evaluation_supported() or
+ * FEFaceEvaluation::fast_evaluation_supported().
  *
  * <h3>Handling multi-component systems</h3>
  *
@@ -2763,6 +2766,13 @@ public:
   void
   reinit(const typename Triangulation<dim>::cell_iterator &cell);
 
+  /**
+   * Check if face evaluation/integration is supported.
+   */
+  static bool
+  fast_evaluation_supported(const unsigned int given_degree,
+                            const unsigned int give_n_q_points_1d);
+
   /**
    * Evaluate the function values, the gradients, and the Hessians of the
    * polynomial interpolation from the DoF values in the input vector to the
@@ -3175,6 +3185,13 @@ public:
   void
   reinit(const unsigned int cell_batch_number, const unsigned int face_number);
 
+  /**
+   * Check if face evaluation/integration is supported.
+   */
+  static bool
+  fast_evaluation_supported(const unsigned int given_degree,
+                            const unsigned int give_n_q_points_1d);
+
   /**
    * Evaluates the function values, the gradients, and the Laplacians of the
    * FE function given at the DoF values stored in the internal data field
@@ -9538,6 +9555,54 @@ FEFaceEvaluation<dim,
 
 
 
+template <int dim,
+          int fe_degree,
+          int n_q_points_1d,
+          int n_components_,
+          typename Number,
+          typename VectorizedArrayType>
+bool
+FEEvaluation<dim,
+             fe_degree,
+             n_q_points_1d,
+             n_components_,
+             Number,
+             VectorizedArrayType>::
+  fast_evaluation_supported(const unsigned int given_degree,
+                            const unsigned int give_n_q_points_1d)
+{
+  return fe_degree == -1 ?
+           internal::FEEvaluationFactory<dim, Number, VectorizedArrayType>::
+             fast_evaluation_supported(given_degree, give_n_q_points_1d) :
+           true;
+}
+
+
+
+template <int dim,
+          int fe_degree,
+          int n_q_points_1d,
+          int n_components_,
+          typename Number,
+          typename VectorizedArrayType>
+bool
+FEFaceEvaluation<dim,
+                 fe_degree,
+                 n_q_points_1d,
+                 n_components_,
+                 Number,
+                 VectorizedArrayType>::
+  fast_evaluation_supported(const unsigned int given_degree,
+                            const unsigned int give_n_q_points_1d)
+{
+  return fe_degree == -1 ?
+           internal::FEFaceEvaluationFactory<dim, Number, VectorizedArrayType>::
+             fast_evaluation_supported(given_degree, give_n_q_points_1d) :
+           true;
+}
+
+
+
 /*------------------------- end FEFaceEvaluation ------------------------- */
 
 
diff --git a/tests/matrix_free/fe_evaluation_print_01.cc b/tests/matrix_free/fe_evaluation_print_01.cc
new file mode 100644 (file)
index 0000000..e6f1239
--- /dev/null
@@ -0,0 +1,66 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2021 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 FEEvaluation::fast_evaluation_supported()
+
+#include <deal.II/dofs/dof_handler.h>
+
+#include <deal.II/fe/fe_q.h>
+
+#include <deal.II/grid/grid_generator.h>
+
+#include <deal.II/lac/affine_constraints.h>
+
+#include <deal.II/matrix_free/fe_evaluation.h>
+#include <deal.II/matrix_free/matrix_free.h>
+
+#include "../tests.h"
+
+
+
+template <int dim>
+void
+test()
+{
+  for (unsigned int i = 1; i < 25; i++)
+    {
+      for (unsigned int j = 1; j < 25; j++)
+        if (FEEvaluation<dim, -1, 0, 1>::fast_evaluation_supported(i, j))
+          deallog << 1 << " ";
+        else
+          deallog << "  ";
+      deallog << std::endl;
+    }
+  for (unsigned int i = 1; i < 25; i++)
+    {
+      for (unsigned int j = 1; j < 25; j++)
+        if (FEFaceEvaluation<dim, -1, 0, 1>::fast_evaluation_supported(i, j))
+          deallog << 1 << " ";
+        else
+          deallog << "  ";
+      deallog << std::endl;
+    }
+}
+
+
+
+int
+main()
+{
+  initlog();
+  test<1>();
+}
diff --git a/tests/matrix_free/fe_evaluation_print_01.output b/tests/matrix_free/fe_evaluation_print_01.output
new file mode 100644 (file)
index 0000000..742ac8d
--- /dev/null
@@ -0,0 +1,49 @@
+
+DEAL::1 1 1                                           
+DEAL::  1 1 1                                         
+DEAL::    1 1 1                                       
+DEAL::      1 1 1 1                                   
+DEAL::        1 1 1 1                                 
+DEAL::          1 1 1   1                             
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::1 1 1                                           
+DEAL::  1 1 1                                         
+DEAL::    1 1 1                                       
+DEAL::      1 1 1 1                                   
+DEAL::        1 1 1 1                                 
+DEAL::          1 1 1   1                             
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                
+DEAL::                                                

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.