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);
};
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);
};
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>
+ 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(
+ 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(
* 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>
*
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
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
+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 ------------------------- */
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// 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>();
+}
--- /dev/null
+
+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::