]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Make FE_EVAL_FACTORY_DEGREE_MAX more robust
authorMartin Kronbichler <martin.kronbichler@uni-a.de>
Mon, 19 Jun 2023 16:58:34 +0000 (18:58 +0200)
committerMartin Kronbichler <martin.kronbichler@uni-a.de>
Mon, 19 Jun 2023 16:58:34 +0000 (18:58 +0200)
include/deal.II/matrix_free/evaluation_kernels.h
include/deal.II/matrix_free/evaluation_template_face_factory.templates.h
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

index 0f083d655c485a47d3c83b73b43fcce3d3bd8d3d..422c597c0a45c569eae5464e111c8b9894e5347a 100644 (file)
@@ -5856,6 +5856,20 @@ namespace internal
     }
   };
 
+  /**
+   * This struct is used to implement
+   * FEEvaluation::fast_evaluation_supported() and
+   * FEFaceEvaluation::fast_evaluation_supported().
+   */
+  struct FastEvaluationSupported
+  {
+    template <int fe_degree, int n_q_points_1d>
+    static bool
+    run()
+    {
+      return fe_degree != -1;
+    }
+  };
 } // end of namespace internal
 
 
index fa068963bad7209d3178d3a82a3510cbf6128e6f..e9c541136b4ae673ba25c13e144df64406549cca 100644 (file)
@@ -118,6 +118,21 @@ namespace internal
       fe_eval);
   }
 
+
+
+  // It is important that this file sits in the same compilation unit as the
+  // evaluate() and integrate() calls of this class, to ensure that all
+  // options choose the same code path when compiling FEFaceEvaluationFactory
+  // outside of deal.II.
+  template <int dim, typename Number>
+  bool
+  FEFaceEvaluationFactory<dim, Number>::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);
+  }
 } // end of namespace internal
 
 DEAL_II_NAMESPACE_CLOSE
index a179402e54712d5e908fbb446fd04714be2c008f..0d8e036c33145f506de67dbae635b0972ad6b088 100644 (file)
@@ -71,6 +71,10 @@ namespace internal
               const EvaluationFlags::EvaluationFlags integration_flag,
               Number *                               values_dofs,
               FEEvaluationData<dim, Number, true> &  fe_eval);
+
+    static bool
+    fast_evaluation_supported(const unsigned int given_degree,
+                              const unsigned int n_q_points_1d);
   };
 
 
index bd0f37fda627d0f766911f27b4a63308c6c0d519..51b98ee4d0fe2569b6ee47ad4e9fefac573b58ab 100644 (file)
@@ -29,18 +29,6 @@ DEAL_II_NAMESPACE_OPEN
 
 namespace internal
 {
-  struct FastEvaluationSupported
-  {
-    template <int fe_degree, int n_q_points_1d>
-    static bool
-    run()
-    {
-      return fe_degree != -1;
-    }
-  };
-
-
-
   template <int dim, typename Number>
   void
   FEEvaluationFactory<dim, Number>::evaluate(
@@ -81,6 +69,10 @@ namespace internal
 
 
 
+  // It is important that this file sits in the same compilation unit as the
+  // evaluate() and integrate() calls of this class, to ensure that all
+  // options choose the same code path when compiling FEFaceEvaluationFactory
+  // outside of deal.II.
   template <int dim, typename Number>
   bool
   FEEvaluationFactory<dim, Number>::fast_evaluation_supported(
index a171bcda110ea2627a6b2a34ee57a2715accd092..59b8f75c6b21f35cf9e2eab03549d2659610db34 100644 (file)
@@ -1709,8 +1709,16 @@ 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. You can check if fast
- * evaluation/integration for a given degree/n_quadrature_points pair by calling
+ * a possibly larger set of degrees. This can both be set when configuring
+ * deal.II by passing the flag `-D FE_EVAL_FACTORY_DEGREE_MAX=8` (in case you
+ * want to compile all degrees up to eight; recommended setting) or by
+ * compiling `evaluation_template_factory.templates.h` and
+ * `evaluation_template_face_factory.templates.h` with the
+ * `FE_EVAL_FACTORY_DEGREE_MAX` overriden to the desired value. In the second
+ * option, symbols will be available twice, and it depends on your linker and
+ * dynamic library loader whether the user-specified setting takes
+ * precendence. 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().
  *
@@ -9249,7 +9257,7 @@ FEFaceEvaluation<dim,
                             const unsigned int given_n_q_points_1d)
 {
   return fe_degree == -1 ?
-           internal::FEEvaluationFactory<dim, VectorizedArrayType>::
+           internal::FEFaceEvaluationFactory<dim, VectorizedArrayType>::
              fast_evaluation_supported(given_degree, given_n_q_points_1d) :
            true;
 }

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.