Dear Lukas,

would it not help to talk to Sebastian Kinnewig who implemented recently 


a lit of this stuff for quads and hexes? 

Of course for your question 

Both FE_PolyTensor and FE_Simplex_Polyseem like reasonable starting points.
However, from my initial review, FE_Simplex_Poly only takes unit_support_points in its constructor, and I am unsure how it supports generalized support points.

For the FE_PolyTensor approach, I could refer to FE_Nedelec, but I am not certain whether this approach is compatible with a local reorientation strategy.


I would also assume that the core deal.II developers have a reasonable suggestion.

Best regards, Thomas

+++++
Prof. Dr. Thomas Wick
www.thomaswick.org
Sent from my iPhone

On 30. Oct 2025, at 14:21, Siarhei Uzunbajakau <[email protected]> wrote:

Dear Lukas,

I am also interested in the Nedelec finite elements. Unfortunately, I have absolutely no experience in implementing new finite elements, so I cannot help you with it.  I strongly suspect  that I will have to learn how to do that in the future as I would like to make a few finite elements that will be handy in combination with FE_NedelecSZ.

I have never tried the finite elements on triangular and tetrahedral cells in deal.II so I cannot tell you what is going on with cell orientation there. I know that the deal.II library compensates for edge/face orientation in the case of the FE_Nedelec finite elements on the quadrilateral and hexahedral cells in 2D and 3D if the order of the finite elements is less than five. It does the compensation by remembering the orientation of each edge and face. The mismatch in edge orientation is compensated by flipping the sign of the shape function. The mismatch in face orientation is compensated by a table lookup. There are eight possible face orientations and eight entries in the table. Each entry encodes which shape functions to swap and which shape functions needs to flip sign. The orientation compensation  of the 2D FE_Nedelec finite elements is done in source/fe/fe_poly_tensor.cc. The 3D FE_Nedelec finite elements  are compensated in source/fe/fe_nedelec.cc. See comments in these two files for more details.

I prefer not to think in terms of global cell orientation. Local orientation is good enough for me. Cells that can map themselves onto the reference cell correctly is all we need. In my opinion the cells are mapped correctly onto the reference cell in 2D and 3D if element degree <5, see test file tests/fe/fe_nedelec_face_edge_orientation.h . I do not experience any problems with the global mesh orientation and sign conflicts under these conditions.

The classical FE_Nedelec finite elements have a disadvantage: they do not work quite well (upset De Rham complex) if the degree of the finite elements varies from cell to cell (p- adaptation), see for example introduction to this article,  DOI:10.1016/S0898-1221(00)00062-6, or first pages in the thesis of Sabine Zaglmayr (listed in the deal.II bibliography page). I guess the main reason why people create new Nedelec finite elements, such as FE_NedelecSZ, is the desire to overcome this disadvantage.  Although, I must admit that “incorporating the global orientation into the finite element itself and resolving sign conflicts” is a cool feature of FE_NedelecSZ. The FE_Nedelec implements this cool feature by compensation of the mismatch in the edge/face orientation as described  above.  From the user’s stand point, the real significant difference between FE_Nedelec and FE_NedelecSZ is that FE_NedelecSZ allows a nonuniform distribution of the degree of the finite elements over the mesh. The degree of the FE_Nedelec finite elements must be uniform all over the mesh.

Do not be discouraged by the above mentioned disadvantage of the FE_Nedelec finite elements. In my humble opinion, one can have a lot of fun even if the degree of the finite elements is uniform all over the mesh.

Please, let me know if you have made any progress in creating new finite elements. I am interested in this topic and will be more than happy to learn from your experience.

Kind regards,
Siarhei
On Monday, October 27, 2025 at 11:44:52 AM UTC+1 Lukas Dreyer wrote:

Dear deal.II community,

I am currently exploring the implementation of a Nédélec finite element for triangles and tetrahedra.

My focus is on linear elements, where all DoFs are located on the edges of the element. For triangles, the basis functions are defined as:

p_1 = (1 - y,  x),  p_2 = ( -y ,  x),  p_3 = (  y , 1 - x)

For the construction of the linear functionals N_i, I define the tangential vectors t_i and edge midpoints q_i as follows:

t_1 = (1, 0),  t_2 = (-1, 1),  t_3 = (0, 1)

q_1 = (0.5, 0),  q_2 = (0.5, 0.5),  q_3 = (0, 0.5)

Then, with N_i(p) = p(q_i) · t_i, we obtain N_i(p_j) = δ_ij.

These functionals require generalized support points, as the functionals are not point evaluations.

I have found two existing implementations for Nédélec elements on hypercubes: FE_Nedelec and FE_NedelecSZ.
>From my understanding, FE_Nedelec requires globally oriented meshes, where each local edge orientation matches the global orientation.
In contrast, FE_NedelecSZ incorporates the global orientation into the finite element itself and resolves sign conflicts by multiplying incorrectly oriented DoFs by -1.

I am currently deciding which base class would be most suitable for my implementation.
Is there a better option than deriving directly from FiniteElement?

Both FE_PolyTensor and FE_Simplex_Poly seem like reasonable starting points.
However, from my initial review, FE_Simplex_Poly only takes unit_support_points in its constructor, and I am unsure how it supports generalized support points.

For the FE_PolyTensor approach, I could refer to FE_Nedelec, but I am not certain whether this approach is compatible with a local reorientation strategy.

I would greatly appreciate any insights or recommendations on how to implement the simplicial linear nedelec element with minimal code duplication.

Best regards,
Lukas Dreyer



--
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see https://groups.google.com/d/forum/dealii?hl=en
---
You received this message because you are subscribed to the Google Groups "deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email to [email protected].
To view this discussion visit https://groups.google.com/d/msgid/dealii/daed4470-8e65-49d8-8f3a-284f80a038c6n%40googlegroups.com.

--
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see https://groups.google.com/d/forum/dealii?hl=en
---
You received this message because you are subscribed to the Google Groups "deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email to [email protected].
To view this discussion visit https://groups.google.com/d/msgid/dealii/DAB46113-8F14-47C8-BA6C-7C47780CD77B%40gmail.com.

Reply via email to