Hi all and Wolfgang,

I did some digging in the archives and found a response from Anna Schneebeli
to almost the same question:

On Sat, 18 Oct 2003, Anna Schneebeli wrote:
>
> In fact, for the usual variational problems where  H(curl)-conforming Nedelec
> elements (those currently implemented in dealii) may be used for
> FE-discretization, it does not make sense to impose a boundary condition on 
> the
> normal component of the field - the normal trace of an H(curl)-function is not
> well defined (in a suitable space of boundary functions), e.g. it can jump
> across a boundary.
> An appropriate Dirichlet boundary condition would be  n x u = 0 , a Neumann
> condition  n x curl u = 0.
> Anna.

She continues to explain how to impose the tangential
boundary conditions

> As for imposing boundary conditions, the 
> VectorTools::interpolate_boundary_function
> cannot (yet) be used for Nedelec elements.
> The way I currenlty impose non-homogeneous Dirichlet bc in my Nedelec codes 
> is:
> loop over all boundary edges e and approximate the value of the current dof 
> (which
> is the line integral over e of some scalar function g, g being the bc for the
> tangential component of u on e) by a (line) quadrature (actually, the 
> midpoint rule
> is accurate enough).
> Have a look at the attached patch for details !
> (Also included is a possible assembling routine for the model problem  curl 
> curl u +
> u = f in a 2D domain)
> Hope it helps!
> Anna

This was helpful.  Thanks to everyone and the archive.

The nedelec element is in the polynomial space that I need,
is discontinuous in the normal direction and continuous
in the tangential directions.  These are the conditions
that I need to satisfy the inf-sup condition.  That was the
reason it was chosen.

The Nedelec dofs seem to prohibit the boundary conditions
from being enforced.  I will place
my follow-up question in a different e-mail.
Thanks for your time.

Dan


On Thu, Jun 23, 2011 at 7:53 PM, Wolfgang Bangerth
<[email protected]>wrote:

>
> Daniel,
>
> > I am trying to combine the previous e-mails into one.  I have attached
> > 3 slides from a slideshow for the MHD equations and the linear form.
> > The Nedelec element used is for the current J.
>
> If I see this right on page 3, then the only derivatives on the variable J
> is
> the divergence operator. For the divergence operator to make sense, you
> need
> continuity of the *normal* component of a finite element across element
> faces
> (something the Raviart-Thomas element provides but  the Nedelec element
> does
> not) while you don't care about continuity of the *tangential* component
> (something the Nedelec element provides but not the RT element).
>
> This is also reflected in the fact that you impose boundary conditions for
> n.J
> (which again requires continuity of the normal component), not for n \times
> J
> (which would require continuity of the tangential component). What
> motivates
> your choice of elements?
>
>
> > From reading through Anna Schneebeli's paper, I understand that
> > the shape functions don't really have support points.  Thus
> > setting boundary conditions for Nedelec elements would
> > be different from settting the boundary conditions for Lagrangian
> elements.
>
> They do have a sort of support points, but these support points are for
> vector
> components tangential to the face these points lie on. So they're no help
> for
> what you want here.
>
>
> > I was wondering how boundary conditions could be set for a variable
> > associated with a Nedelec element.  In the other e-mail you mentioned
> > Baerbel Janssen working with Hermite-type elements.  I was thinking
> > about a simple linear Hermite element on an interval, where the
> > degrees of freedom are defined using the value and
> > the derivative at the midpoint.  How would someone set the boundary
> > conditions for this element, if it has only the midpoint as a "support"
> > point.
>
> That's difficult. People use Hermite elements for the biharmonic equation
> for
> which you need boundary conditions for both the value and the derivative,
> and
> that's exactly the degrees of freedom the element defines.
>
>
> > Finally, letting x = "generalized support point", I was wondering if it
> > would make any sense to put in a constraint matrix
> > a form of "J(x) . n(x) = sum_{k=1} J_{k} phi_{k}(x) . n(x) =
> > boundary_value".
> > Where k runs through the Nedelec shape functions phi_{k}, and J_{k}
> > are the solution coefficients.  Thanks for everyone's time.
>
> Let's first figure out whether the Nedelec element is really what you want
> :-)
>
> Cheers
>  W.
>
> -------------------------------------------------------------------------
> Wolfgang Bangerth                email:            [email protected]
>                                 www: http://www.math.tamu.edu/~bangerth/
>
_______________________________________________
dealii mailing list http://poisson.dealii.org/mailman/listinfo/dealii

Reply via email to