s_grad_grad(3) grad_grad-like operator for the Stokes stream function computation

SYNOPSIS


form(const space V, const space& V, "s_grad_grad");


DESCRIPTION

Assembly the form associated to the -div(grad) variant operator on a finite element space V. The V space may be a either P1 or P2 finite element space. See also form(2) and space(2). On cartesian coordinate systems, the form coincide with the "grad_grad" one (see grad_grad(3)):

               /
               |
      a(u,v) = |  grad(u).grad(v) dx
               |
               / Omega

The stream function on tri-dimensionnal cartesian coordinate systems is such that

       u = curl psi
       div psi = 0

where u is the velocity field. Taking the curl of the first relation, using the identity:

        curl(curl(psi)) = -div(grad(psi)) + grad(div(psi))

and using the div(psi)=0 relation leads to:

        -div(grad(psi)) = curl(u)

This relation leads to a variational formulation involving the the "grad_grad" and the "curl" forms (see grad_grad(3), curl(3)).

In the axisymmetric case, the stream function psi is scalar ans is defined from the velocity field u=(ur,uz) by (see Batchelor, 6th ed., 1967, p 543):

                 d psi                       d psi
      uz = (1/r) -----    and   ur = - (1/r) -----
                  d r                         d r

See also http://en.wikipedia.org/wiki/Stokes_stream_function . Multiplying by rot(xi)=(d xi/dr, -d xi/dz), and integrating with r dr dz, we get a well-posed variationnal problem:

        a(psi,xi) = b(xi,u)

with

                  /
                  | (d psi d xi   d psi d xi)
      a(psi,xi) = | (----- ---- + ----- ----) dr dz
                  | ( d r  d r     d z  d z )
                  / Omega

and

                /
                | (d xi      d xi   )
      b(xi,u) = | (---- ur - ---- uz) r dr dz
                | (d z       d r    )
                / Omega

Notice that a is symmetric definite positive, but without the 'r' weight as is is usual for axisymmetric standard forms. The b form is named "s_curl", for the Stokes curl variant of the "curl" operator (see s_curl(3)) as it is closely related to the "curl" operator, but differs by the r and 1/r factors, as:

                   (       d (r xi)     d xi )
        curl(xi) = ( (1/r) -------- ; - -----)
                   (         d r        d z  )

while

                   ( d xi       d xi )
      s_curl(xi) = ( ----  ;  - ---- )
                   ( d r        d z  )

EXAMPLE

The following piece of code build the form associated to the P1 approximation:

        geo g("square");
        space V(g, "P1");
        form a(V, V, "s_grad_grad");