Dear Deal.II community,

Suppose I chose FE_Nedelec<3> of some order as finite element.

1. Is there a way to query for a given <code>cell_dof_index</code> whether 
it is a face_dof and to get the <code>face_index</code>?

2. If this dof is a <code>line_dof</code> is it by definition also a 
<code>face_dof</code>? I need to exclude line_dofs and treat them 
separately.

3. Is there a convention that (if I am not using FESystem) in the set of 
all <code>cell_dofs</code> all <code>line_dofs</code> come before all 
<code>face_dofs</code> and only then the <code>volume_dofs</code>?

I am asking since I am trying to come up with a patch for RaviartThomas and 
a pattern to implement it in the current structures without interfering 
with them (too much).

Best,
Konrad
On Wednesday, December 16, 2020 at 9:33:57 PM UTC+1 Konrad Simon wrote:

> Dear Jean-Paul, dear Deal.II community,
>
> I partially solved the problem of sign flipping and permuting the degrees 
> of freedom and would like to come up with a merge-request at least for the 
> RaviartThomas elements. I have some questions before I can go on and 
> guidance would be appreciated.
>
> I decided to go a similar way that is taken for FE_Q elements. There one 
> has only a permutation of dofs on non-standard or flipped or rotated faces 
> since a sign change in the orientation does not matter. A vector of size 
> n_faces_per_cell of which each contains a table mapping each face_dof and 
> each of the 8 combinations of the flags (non-standard | flipped | rotated) 
> to an offset of that face_dof. The system behind that relies on the 
> lexicographic ordering of dofs if I understand correctly. My questions:
>
> 1. Is there a similar numbering scheme for dofs on faces for RT elements. 
> There one has face moments as dofs and it is not fully clear to me what 
> lexicographic means there (if such a principle is usable there at all).
>
> 2. If the order (and sign) of face_dofs if permuted for each of the flags 
> (non-standard | flipped | rotated) then in what order are the permutations 
> applied? (They do not commute)
>
> I would like to avoid hard-coding that (I think up to order 2 would be 
> doable -> but already 54 face_dofs). 
>
> I believe the table similar 
> to adjust_quad_dof_index_for_face_orientation_table of FE_Q_Base makes 
> sense but I guess every FE inheriting from FE_PolyTensor must implement 
> such a table then. This is a lot of work and more complicated for Nedelec 
> elements.
>
> Best,
> Konrad
> On Saturday, December 12, 2020 at 8:07:19 PM UTC+1 Konrad Simon wrote:
>
>> Hi Jean-Paul,
>>
>> On Thursday, December 10, 2020 at 11:39:08 PM UTC+1 Jean-Paul Pelteret 
>> wrote:
>>
>>> HI Konrad,
>>>
>>> I have no familiarity with the H-div elements, so I could be wrong with 
>>> this suggestion...
>>>
>>> The Fe_Nedelec element suffered from a similar issue, where adjacent 
>>> cells had to agree on the edge sense/directionality. This requirement could 
>>> not be unconditionally guaranteed for that element type, hence the 
>>> following note in its introduction 
>>> <https://dealii.org/current/doxygen/deal.II/classFE__Nedelec.html>
>>>
>>> Even if this element is implemented for two and three space dimensions, 
>>> the definition of the node values relies on consistently oriented faces in 
>>> 3D. Therefore, care should be taken on complicated meshes.
>>>
>>
>> Yes, this issue is similar. I assume I will need to check that for these 
>> elements as well since I also need them.
>>
>>
>>> The more recent FE_NedelecSZ 
>>> <https://dealii.org/current/doxygen/deal.II/classFE__NedelecSZ.html> 
>>> element resolves the issue by querying the vertex numbers that make up each 
>>> edge/face and applying a rule that gives each face and its bounding edges a 
>>> direction. This means that the orientation of the face/edges is always 
>>> consistent no matter which of the two cells that share the face you're on. 
>>> The details are more explicitly laid out in the FE_NedelecSZ 
>>> <https://dealii.org/current/doxygen/deal.II/classFE__NedelecSZ.html> 
>>> introduction, as well as the two papers that are listed therein. 
>>>
>>
>> This is a nice alternative for the Nedelec element. Often though such 
>> elements are used in conjunction with H1 or H(div) conformal elements and 
>> need to satisfy a compatibility condition (stronger than just LBB, i.e., 
>> they must form a bounded co-chain complex that is a subcomplex to something 
>> else). I am not sure that this element does satisfy this condition together 
>> with classical RT elements etc.
>>
>>
>>> Now, I'm by no means certain that underlying issue that you're seeing 
>>> here is the same as that which plagued the canonical Nedelec element 
>>> (although your comment #3 makes me suspect that it could be). But if it is, 
>>> then perhaps what was done to fix the orientation conflict in the Zaglmayr 
>>> formulation of the Nedelec element can serve as inspiration for a fix to 
>>> the BDM and RT elements?
>>>
>>> Thanks for looking into this, and I hope that you have some success at 
>>> finding a solution to the problem. Naturally, if there's anything that we 
>>> can do to help please let us know. It would be very nice to have a robust 
>>> implementation to these elements in the library, so anything that you might 
>>> be willing to contribute would be greatly appreciated!
>>>
>>  
>> Lowest order RaviartThomas is already fixed. I will hunt down this issue 
>> as good as I can. I might need some help with local and global dof 
>> numbering conventions. I will also try to move the discussion to the issue 
>> page <https://github.com/dealii/dealii/issues/7970> on github.
>>
>> Best,
>> Konrad
>>
>

-- 
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 dealii+unsubscr...@googlegroups.com.
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/c6ceaa83-318f-438b-9982-21e2d67a39e8n%40googlegroups.com.

Reply via email to