const unsigned int base_element = 0);
/**
- * Return which kinds of elements are supported by MatrixFree.
+ * Return whether an element is supported by MatrixFree.
+ *
+ * The following scalar elements are supported:
+ * - FENothing, FE_DGP, and FE_Q_DG0
+ * - polynomial tensor-product elements based on
+ * Polynomials::Polynomial (FE_Q, FE_DG_Q, FE_DGQArbitraryNodes,
+ * FE_DGQHermite, FE_DGQLegendre) or Polynomials::PiecewisePolynomial
+ * (FE_Q_iso_Q1).
+ * - elements for simplex, pyramids, and wedges (FE_SimplexP,
+ * FE_SimplexDGP, FE_PyramidP, FE_PyramidDGP, FE_WedgeP, FE_WedgeDGP)
+ *
+ * In the case of vectorial elements, FE_RaviartThomasNodal
+ * and FESystem with base elements from the scalar elements
+ * listed above are supported.
*/
template <int dim, int spacedim>
static bool
#include <deal.II/fe/fe.h>
#include <deal.II/fe/fe_dgp.h>
#include <deal.II/fe/fe_dgq.h>
+#include <deal.II/fe/fe_nothing.h>
#include <deal.II/fe/fe_poly.h>
#include <deal.II/fe/fe_pyramid_p.h>
#include <deal.II/fe/fe_q.h>
return;
}
else if (quad_in.is_tensor_product() == false ||
- dynamic_cast<const FE_SimplexP<dim> *>(
- &fe_in.base_element(base_element_number)) ||
- dynamic_cast<const FE_SimplexDGP<dim> *>(
- &fe_in.base_element(base_element_number)) ||
- dynamic_cast<const FE_WedgeP<dim> *>(
- &fe_in.base_element(base_element_number)) ||
- dynamic_cast<const FE_PyramidP<dim> *>(
- &fe_in.base_element(base_element_number)))
+ dynamic_cast<const FE_SimplexPoly<dim, dim> *>(
+ &fe_in.base_element(base_element_number)) != nullptr ||
+ dynamic_cast<const FE_WedgePoly<dim, dim> *>(
+ &fe_in.base_element(base_element_number)) != nullptr ||
+ dynamic_cast<const FE_PyramidPoly<dim, dim> *>(
+ &fe_in.base_element(base_element_number)) != nullptr)
{
// specialization for arbitrary finite elements and quadrature rules
// as needed in the context, e.g., of simplices
if (dim != spacedim)
return false;
+ if (dynamic_cast<const FE_RaviartThomasNodal<dim> *>(&fe))
+ return true;
+
for (unsigned int base = 0; base < fe.n_base_elements(); ++base)
{
const FiniteElement<dim, spacedim> *fe_ptr = &(fe.base_element(base));
dynamic_cast<const FE_Poly<dim, spacedim> *>(fe_ptr);
// Simplices are a special case since the polynomial family is not
// indicative of their support
- if (dynamic_cast<const FE_SimplexP<dim> *>(fe_poly_ptr) ||
- dynamic_cast<const FE_SimplexDGP<dim> *>(fe_poly_ptr) ||
- dynamic_cast<const FE_WedgeP<dim> *>(fe_poly_ptr) ||
- dynamic_cast<const FE_PyramidP<dim> *>(fe_poly_ptr))
+ if (dynamic_cast<const FE_SimplexPoly<dim, dim> *>(fe_poly_ptr) !=
+ nullptr ||
+ dynamic_cast<const FE_WedgePoly<dim, dim> *>(fe_poly_ptr) !=
+ nullptr ||
+ dynamic_cast<const FE_PyramidPoly<dim, dim> *>(fe_poly_ptr) !=
+ nullptr)
return true;
- if (dynamic_cast<const TensorProductPolynomials<dim> *>(
+ if (dynamic_cast<const TensorProductPolynomials<
+ dim,
+ Polynomials::Polynomial<double>> *>(
&fe_poly_ptr->get_poly_space()) == nullptr &&
dynamic_cast<const TensorProductPolynomials<
dim,
nullptr)
return false;
}
+ else if (dynamic_cast<const FE_Nothing<dim, spacedim> *>(fe_ptr) !=
+ nullptr)
+ return true;
else
return false;
}
DEAL::FE_FaceQ<2>(1) supported by MatrixFree: false
DEAL::FE_P1NC supported by MatrixFree: false
DEAL::FE_Nedelec<2>(1) supported by MatrixFree: false
-DEAL::FE_Nothing<2>() supported by MatrixFree: false
+DEAL::FE_Nothing<2>() supported by MatrixFree: true
DEAL::FE_Q<1,2>(1) supported by MatrixFree: false
DEAL::FE_Q<2>(1) supported by MatrixFree: true
DEAL::FE_Q_Bubbles<2>(1) supported by MatrixFree: false
DEAL::FE_Q_iso_Q1<2>(1) supported by MatrixFree: true
DEAL::FE_RannacherTurek<2>(0, 2) supported by MatrixFree: false
DEAL::FE_RaviartThomas<2>(1) supported by MatrixFree: false
-DEAL::FE_RaviartThomasNodal<2>(1) supported by MatrixFree: false
+DEAL::FE_RaviartThomasNodal<2>(1) supported by MatrixFree: true
DEAL::FE_TraceQ<2>(1) supported by MatrixFree: false