Dear Konrad,

I do not understand when a geometry gets complicated. Is a toroid inside a 
sphere, both centred at the origin, a complicated geometry? To start with, 
I will use a relatively simple mesh. I will not use local refinement, only 
global. I can do without refinement as well, i.e. making the mesh in gmsh.

Yes, I would like to know more about orientations issues on complicated 
geometries. Could please tell me about it? I would very much prefer to use 
FE_Nedelec without hierarchical shape functions. The first order 3D 
FE_Nedelec will do nicely provided the orientation issues will be 
irrelevant to the mesh.

John


On Tuesday, April 20, 2021 at 9:43:31 PM UTC+2 Konrad Simon wrote:

> Dear John,
>
> As Jean-Paul pointed out the entries of the solution vector for a Nedelec 
> element have a different meaning. The entries are weights that result from 
> edge moments (and in case of higher order also face moments and volume 
> moments), i.e, integrals on edges. Nedelec elements have a deep geometric 
> meaning: they discretize quantities that you can integrate (a physicist 
> would say "measure") along 1D-lines in 3D.
>
> Btw, if you are using Nedelec in 2d be aware that it has some difficulties 
> on complicated meshes. In 3D the lowest order is fine. NedelecSZ (for 2D 
> lowest order and 3D all orders), as Jean-Paul pointed out, are another 
> option. They are meant to by-pass certain mesh orientation issues on 
> complicated geometries (I can tell you a bit about that since currently I 
> am chasing some problems there).
>
> Best,
> Konrad
>
> On Tuesday, April 20, 2021 at 9:23:23 PM UTC+2 Jean-Paul Pelteret wrote:
>
>> Hi John,
>>
>> Unfortunately, you’ve fallen into the trap of confusing what the entries 
>> in the solution vector mean for the different element types. For Nedelec 
>> elements, they really are coefficients of a polynomial function, so you 
>> can’t simply set each coefficient to 1 to visualise the shape function 
>> (like one might be inclined to try for Lagrange type polynomials). If you 
>> want a little piece of code that will plot the Nedelec shape functions for 
>> you, then take a look over here:
>> https://github.com/dealii/dealii/issues/8634#issuecomment-595715677
>>
>> BTW. In case you’re not aware, a newer (and probably more robust) 
>> implementation of a Nedelec element is the FE_NedelecSZ 
>> <https://dealii.org/current/doxygen/deal.II/classFE__NedelecSZ.html> that 
>> uses hierarchical shape functions to form the higher order bases. The PhD 
>> thesis of S. Zaglmayr that describes the element is mentioned in the class 
>> documentation.  
>>
>> Best,
>> Jean-Paul
>>
>> On 20. Apr 2021, at 15:47, John Smith <[email protected]> wrote:
>>
>>
>> Dear Konrad,
>>
>> Thank you for your help! I have recently modified the step-3 tutorial 
>> program to save FE_Nedelec<2> fe(0) shape function into *.csv files. It is 
>> attached to this message. The result of this program is not what I would 
>> expect it to be. I have also attached a pdf file which compares my 
>> expectations (marked as “Expectations” in the pdf file) with the results of 
>> the modified step-3 program (marked as “Reality” in the pdf file). My 
>> expectations are based on the literature I have. My question is: why the 
>> shape functions of FE_Nedelec<2> fe(0) differ from my expectations? I guess 
>> my expectations are all wrong. I would also appreciate if you will share 
>> with me a reference to a paper or a book, so I can extend my modest library.
>>
>> Thanks!
>>
>> John
>>
>> P.S. I have sent you an e-mail as you have requested. Could you please 
>> check your spam folder? 
>>
>>
>>
>> On Monday, April 19, 2021 at 8:47:29 AM UTC+2 Konrad Simon wrote:
>>
>>> Dear John,
>>>
>>> I had a look at the pdf you sent. I noticed some conceptual 
>>> inconsistencies that are important for the discretization. Let me start 
>>> like this: The curl-curl-problem that you are trying to solve is - 
>>> according to your description of the gauge - actually a curl-curl + 
>>> grad-div problem and hence a Laplace-problem that describes a 
>>> Helmholtz-decomposition. In other words in order to get a well-posed 
>>> simulation problem you need more boundary conditions. You already noticed 
>>> that since your curl-curl-matrix is singular (the kernel is quite 
>>> complicated and contains all longitudinal waves). 
>>>
>>> You have of course the option to solve (0.0.2) with a Krylov-solver but 
>>> then you need to make sure that the right-hand side is in the orthogonal 
>>> complement of this kernel before each iteration which is quite difficult. I 
>>> would not recommend that option.
>>>
>>> Another point is that if you have a solution for your field A like in 
>>> (0.0.4) you can not have a similar representation for T. This is 
>>> mathematically not possible.
>>>
>>> Since you not care about the gauge let me tell you how I would tackle 
>>> this: The reason  is  - to come back to my original point - you are missing 
>>> a boundary condition that makes you gauge unique. Since you apply natural 
>>> boundary conditions (BCs) on A you must do so as well to determine div(A). 
>>> This second BC for A is applied to a different part of the system that is 
>>> usually neglected in the literature and A is partly determined from this 
>>> part (this part then describes the missing transversal waves). The 
>>> conditions you mention on curl(A) are some what contradictory to the space 
>>> in which A lives. Nedelec elements which you must use (you can not use FE_Q 
>>> to enforce the conditions) cannot generate your desired T_0.
>>>
>>> There are some principles when discretizing these problems (which are 
>>> not obvious) that you MUST adhere to (choice of finite elements, boundary 
>>> conditions, what is the exact system etc) if you want to get a stable 
>>> solution and these are only understood recently. I am solving very similar 
>>> problems (with Deal.II) in a fluid context and will be happy to share my 
>>> experiences with you. Just email me: [email protected]
>>>
>>> Best,
>>> Konrad
>>>
>>> On Monday, April 19, 2021 at 5:51:17 AM UTC+2 Wolfgang Bangerth wrote:
>>>
>>>> On 4/18/21 12:24 PM, John Smith wrote: 
>>>> > 
>>>> > Thank you for your reply. It is a bit difficult to read formulas in 
>>>> text. So I 
>>>> > have put a few questions I have in a pdf file. Formulas are better 
>>>> there. It 
>>>> > is attached to this message. May I ask you to have a look at it? 
>>>>
>>>> John: 
>>>>
>>>> If I understand your first question right, then you are given a vector 
>>>> field J 
>>>> and you are looking for a vector field T so that 
>>>> curl T = J 
>>>> I don't really have anything to offer to this kind of question. There 
>>>> are many 
>>>> vector fields T that satisfy this because of the large null space of 
>>>> curl. You 
>>>> have a number of condition you would like to "approximate" but there 
>>>> are many 
>>>> ways to "approximate" something. In essence, you have two goals: To 
>>>> satisfy 
>>>> the equation above and to approximate something. You have to come up 
>>>> with a 
>>>> way to weigh these two goals. For example, you could look to minimize 
>>>> F(T) = ||curl T - J||^2 + alpha ||T-T_desired||^2 
>>>> where T_desired is the right hand side of (0.0.5) and you have to 
>>>> determine 
>>>> what alpha should be. 
>>>>
>>>> As for your other questions: (0.0.9) is indeed called "interpolation" 
>>>>
>>>> The issue with the Nedelec element (as opposed to a Q(k)^dim) field is 
>>>> that 
>>>> for the Nedelec element, 
>>>> \vec phi(x_j) 
>>>> is not independent of 
>>>> \vec phi(x_k) 
>>>> and so you can't choose the matrix in (0.0.8) as the identity matrix. 
>>>> You 
>>>> realize this in (0.0.13). The point I was trying to make is that 
>>>> (0.0.13) 
>>>> cannot be exactly satisfied if T(x_j) is not a vector parallel to \hat 
>>>> e_j, 
>>>> which in general it will not be. You have to come up with a different 
>>>> notion 
>>>> of what it means to solve (0.0.13) because you cannot expect the left 
>>>> and 
>>>> right hand side to be equal. 
>>>>
>>>> Best 
>>>> W> 
>>>>
>>>> -- 
>>>> ------------------------------------------------------------------------
>>>>  
>>>> Wolfgang Bangerth email:  [email protected] 
>>>> www: http://www.math.colostate.edu/~bangerth/ 
>>>>
>>>>
>> -- 
>> 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 on the web visit 
>> https://groups.google.com/d/msgid/dealii/348bd8f9-0891-42e5-b150-c3723eacb84bn%40googlegroups.com
>>  
>> <https://groups.google.com/d/msgid/dealii/348bd8f9-0891-42e5-b150-c3723eacb84bn%40googlegroups.com?utm_medium=email&utm_source=footer>
>> .
>> <FE_Nedelec.pdf><step-3.cc>
>>
>>
>>

-- 
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 on the web visit 
https://groups.google.com/d/msgid/dealii/b82d4595-8508-4201-98d7-c554019f6574n%40googlegroups.com.

Reply via email to