Quoting Ronan Perrussel <[EMAIL PROTECTED]>:
> Dear getfem++ users,
>
> I would like to use Nedelec finite elements for some computations. I
> made some small two dimensional experiments for evaluating the use of
> these finite element family but I do not manage to make all what I wish.
> For assembling a term like "\int_\Omega \curl(E) \curl(E') dx", it
> seems to be fine with the following piece of code :
> ""
> getfem::generic_assembly assem;
> assem.push_mi(mi);
> assem.push_mf(mf);
> assem.push_mat(K);
> assem.set("M(#1,#1)+=comp(vGrad(#1).vGrad(#1))(:,2,1,:,2,1)"
> "+comp(vGrad(#1).vGrad(#1))(:,1,2,:,1,2)"
> "-comp(vGrad(#1).vGrad(#1))(:,2,1,:,1,2)"
> "-comp(vGrad(#1).vGrad(#1))(:,1,2,:,2,1)");
> assem.assembly();
> ""
> but for assembling a term like "\int_\Omega E . E' dx", the simple
> piece of code for the mass matrix does not give the good answer :
> ""
> getfem::asm_mass_matrix(M, mi, mf);
> ""
> and I do not understand exactly why.
>
> Does somebody have a suggestion for assembling the term " \int_\Omega E
> . E' dx"?
>
> Thank you in advance for the answer,
> Best regards,
> Ronan Perrussel
Ok, after checking the implementation of Nedelec elements in Getfem++,
I see that there may be a difference with the common definition of
these elements.
In order to have a general implementation (for curved elements for
instance) and a maximum of computation done only once, they are
defined via a reference element.
What seems to be litigious is that that the base functions are
transported via the geometric transformation like vectors. That is,
they are multiplied to the left by the gradient of the transformation.
I tried to change this and transport them like gradients (i.e.
multiplied by the transposed inverse of the gradient of the
transformation) and I have the same results as your python program.
It is a little bit annoying because this implies that there is two
kind of intrinsic vectorial elements, those which are transported as
vector and thos which are transported as gradients.
My question is : Is there a big difference in the capability of
approximation between the two versions ? The one which is presently
implemented in getfem is nevertheless a a valid finite element (not
really the Nedelec one). Should I give two versions ?
I am not a so much specialist in Nedelec elements, so I would
appreciate your advice.
Yves.
_______________________________________________
Getfem-users mailing list
[email protected]
https://mail.gna.org/listinfo/getfem-users