← Back to team overview

dolfin team mailing list archive

Re: [FFC-dev] Higher-order problems?

 

Leaving aside any ffc issues, this formulation is a problem mathematically
-- second derivatives of piecewise Lagrange polynomials don't make sense
because they are only C^0.  There are several workarounds:
As you've observed, Hermite or other non-Lagrange elements are needed to
resolve the higher derivatives in plate bending, biharmonic, Cahn-Hilliard,
and other 4th order problems.

- Hermite is a nonconforming C^1 (it's really C^1 in 1d, but not on
triangles or tets).
- The Morley triangle is also a suitable nonconforming element
- Argyris is a fully C^1 triangle

The bug problem with all these element is that we don't support them yet.
 They are on the to-do list, but will require some modifications to ffc
because of technicalities in how they transform from a reference element.
 Anders and I have a plan :)

- You can also use a mixed formulation, writing the 4th order problem as a
system of 2nd order equations.  In strong form,

Delta Delta u = f
-->
v = Delta u
Delta v = f

Delta is the Laplacian, so you can multiply each equation by a test function
and integrate by parts to get a weak form where you never need more than one
derivative


- You can also use the "continuous/discontinuous Galerkin method, where you
add penalty terms on the first derivatives.

Garth: you've used the C/DG method, can you contribute an example?


On Thu, Jun 25, 2009 at 5:08 PM, <ndl@xxxxxxxxxxxxxx> wrote:

> First I have to say that I'm not sure that this should work...
>
> 1- Probably I'm missing/ignoring something (mathematically speaking)
> 2- I didn't find yet, in the literature, a simple example
>    (but I've started this week :) )
>    If anyone knows such an example,
>    I would by happy to make a demo with it.
> 3- all the cases I've seen are related with plate bending and requires
>    hermite elements or some other "non-standard elements"
> 4-Concerning the code: I'm not sure on the Boundary conditions syntax.
>
> I was thinking in a simple "4th order like poisson "
> for instance:
>
> -lap lap u(x,y)=-sin(x)  in the [0,2pi]x[0,2pi]
> lap u= -sin(x)  and  grad(u).n=(cos(x),0).n  in the Boundary
>
> with solution u(x,y)=sin(x)
>
> LapLap.ufl
> ------------------------------------------------------------------------
> element = FiniteElement("Lagrange", triangle, 2)
> v = TestFunction(element)
> u = TrialFunction(element)
> f = Function(element)
> g = Function(element)
> n = VectorConstant("triangle")
> a = -div(grad(v))*div(grad(u))*dx
> L = v*f*dx-f*inner(grad(v),n)*ds+div(grad(g*n[0]))*v*ds
> ----------------------------------------------------------------------
> the main.cpp is modified from the poisson demo  and follows in attachment.
> (I'm using  0.9.2  dolfin version)
> It really fails.
>
> > This should work, but I don't know of a demo testing it.
> >
> > See if you could create a very simple test case that fails.
>
> --
> Nuno David Lopes
>
> e-mail:ndl@xxxxxxxxxxxxxx <e-mail%3Andl@xxxxxxxxxxxxxx>        (FCUL/CMAF)
>           nlopes@xxxxxxxxxxxxxxx    (ISEL)
> http://ptmat.ptmat.fc.ul.pt/%7Endl/
>
> Thu Jun 25 22:19:44 WEST 2009
>
>
> _______________________________________________
> FFC-dev mailing list
> FFC-dev@xxxxxxxxxx
> http://www.fenics.org/mailman/listinfo/ffc-dev
>
>

Follow ups

References