From: Wolfgang Bangerth Date: Sun, 12 May 2024 18:11:12 +0000 (-0600) Subject: Link FE_Q, FE_SimplexP, and friends in the documentation. X-Git-Tag: v9.6.0-rc1~272^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=refs%2Fpull%2F17011%2Fhead;p=dealii.git Link FE_Q, FE_SimplexP, and friends in the documentation. --- diff --git a/include/deal.II/fe/fe_pyramid_p.h b/include/deal.II/fe/fe_pyramid_p.h index af8cdea728..9593376113 100644 --- a/include/deal.II/fe/fe_pyramid_p.h +++ b/include/deal.II/fe/fe_pyramid_p.h @@ -54,7 +54,9 @@ public: /** * Implementation of a scalar Lagrange finite element on a pyramid that yields * the finite element space of continuous, piecewise polynomials of - * degree $k$. + * degree $k$. The corresponding element on simplex (triangular or tetrahedral) + * cells is FE_SimplexP, on hypercube cells it is FE_Q, and + * on wedges it is FE_WedgeP. * * @note Currently, only linear polynomials (degree=1) are implemented. See * also the documentation of ScalarLagrangePolynomialPyramid. diff --git a/include/deal.II/fe/fe_q.h b/include/deal.II/fe/fe_q.h index aa4c2efa83..0e9e948e63 100644 --- a/include/deal.II/fe/fe_q.h +++ b/include/deal.II/fe/fe_q.h @@ -30,11 +30,15 @@ DEAL_II_NAMESPACE_OPEN */ /** - * Implementation of a scalar Lagrange finite element @p Qp that yields the + * This class provides an implementation of the scalar Lagrange $Q_p$ finite + * element on hypercube (line segments, quadrilaterals, or hexahedra) cells that + * yields the * finite element space of continuous, piecewise polynomials of degree @p p in * each coordinate direction. This class is realized using tensor product * polynomials based on 1d Lagrange polynomials with equidistant (degree up to - * 2), Gauss-Lobatto (starting from degree 3), or given support points. + * 2), Gauss-Lobatto (starting from degree 3), or given support points. The + * corresponding element on simplex (triangular or tetrahedral) cells is + * FE_SimplexP, whereas on other cell kinds it is FE_WedgeP and FE_PyramidP. * * The standard constructor of this class takes the degree @p p of this finite * element. Alternatively, it can take a quadrature formula @p points defining diff --git a/include/deal.II/fe/fe_simplex_p.h b/include/deal.II/fe/fe_simplex_p.h index de1df4cf8e..05f55e4569 100644 --- a/include/deal.II/fe/fe_simplex_p.h +++ b/include/deal.II/fe/fe_simplex_p.h @@ -122,7 +122,8 @@ protected: /** * Implementation of a scalar Lagrange finite element $P_k$ that yields * the finite element space of continuous, piecewise polynomials of - * degree $k$. + * degree $k$. The corresponding element on hypercube cells is FE_Q, on + * wegdes it is FE_WedgeP, and on pyramids it is FE_PyramidP. * * Also see * @ref simplex "Simplex support". diff --git a/include/deal.II/fe/fe_wedge_p.h b/include/deal.II/fe/fe_wedge_p.h index 9006d3ecbb..75d4f28b1f 100644 --- a/include/deal.II/fe/fe_wedge_p.h +++ b/include/deal.II/fe/fe_wedge_p.h @@ -54,7 +54,9 @@ public: /** * Implementation of a scalar Lagrange finite element on a wedge that yields * the finite element space of continuous, piecewise polynomials of - * degree $k$. + * degree $k$. The corresponding element on simplex (triangular or tetrahedral) + * cells is FE_SimplexP, on hypercube cells it is FE_Q, and + * on pyramids it is FE_PyramidP. * * @note Currently, only linear (degree=1) and quadratic polynomials * (degree=2) are implemented. See also the documentation of