From 735493a1aaec08914b8788a658975a91fa360262 Mon Sep 17 00:00:00 2001 From: Wolfgang Bangerth Date: Wed, 15 Sep 2010 11:24:58 +0000 Subject: [PATCH] Add documentation on the order. git-svn-id: https://svn.dealii.org/trunk@21976 0785d39b-7218-0410-832d-ea1e28bc413d --- deal.II/deal.II/include/fe/fe_nedelec.h | 23 ++++++++++++++++++----- 1 file changed, 18 insertions(+), 5 deletions(-) diff --git a/deal.II/deal.II/include/fe/fe_nedelec.h b/deal.II/deal.II/include/fe/fe_nedelec.h index 2e36dbf7e6..f368f280d0 100644 --- a/deal.II/deal.II/include/fe/fe_nedelec.h +++ b/deal.II/deal.II/include/fe/fe_nedelec.h @@ -37,11 +37,24 @@ template class MappingQ; * space Hcurl. These elements generate vector fields with * tangential components continuous between mesh cells. * - * We follow the usual definition of the degree of Nédélec elements, - * which denotes the polynomial degree of the lowest complete polynomial - * subspace contained in the Nédélec space. Then, approximation order of - * the function itself is degree. In this scheme, the lowest - * order element would be created by the call FE_Nedelec(0). + * We follow the convention that the degree of Nédélec elements + * denotes the polynomial degree of the largest complete polynomial + * subspace contained in the Nédélec space. This leads to the consistently + * numbered sequence of spaces + * @f[ + * Q_{k+1} + * \stackrel{\text{grad}}{\rightarrow} + * \text{Nedelec}_k + * \stackrel{\text{curl}}{\rightarrow} + * \text{RaviartThomas}_k + * \stackrel{\text{div}}{\rightarrow} + * DGQ_{k} + * @f] + * Consequently, approximation order of + * the Nedelec space equals the value degree given to the constructor. + * In this scheme, the lowest order element would be created by the call + * FE_Nedelec(0). Note that this follows the convention of Brezzi and + * Raviart, though not the one used in the original paper by Nedelec. * * This class is not implemented for the codimension one case * (spacedim != dim). -- 2.39.5