Feel++ 0.91.0
Writing Variational Formulations

Standard formulation: the Laplacian case

Mathematical formulation

In this example, we would like to solve for the following problem in 2D find $u$ such that

\[ -\Delta u = f\ \text{in}\ \Omega = [-1;1]^2 \]

with

\[ f= 2 \pi^2 g \]

and $g$ is the exact solution

\[ g=\sin(\pi x) \cos(\pi y) \]

The following boundary conditions apply

\[ u=g_{|x=\pm 1}, \quad \frac{\partial u}{\partial n} = 0_{|y=\pm 1} \]

We propose here two possible variational formulations. The first one, handles the Dirichlet boundary conditions strongly, that is to say the condition is {incorporated} into the function space definitions. The second one handles the Dirichlet condition {weakly} and hence we have a uniform treatment for all types of boundary conditions.

Strong Dirichlet conditions

The variational formulation reads as follows, we introduce the spaces

\[ \mathcal{X} = \Big\{ v \in H_1(\Omega) \text{ such that } v=g_{|x=-1,x=1} \Big\} \]

and

\[ \mathcal{V} = \Big\{ v \in H_1(\Omega) \text{ such that } v=0_{|x=-1,x=1} \Big\} \]

We multiply (eq:1}) by $v \in \mathcal{V}$ then integrate over $\Omega$ and obtain

\begin{equation} \label{eq:13} \int_\Omega -\Delta u v = \int_\Omega f v \end{equation}

We integrate by parts and reformulate the problem as follows: we look for $u \in \mathcal{X}$ such that for all $v \in \mathcal{V}$

\[ \int_\Omega \nabla u \cdot \nabla v = \int_\Omega f v \]

In the present space setting~(eq:12}) and boundary conditions~(eq:4}), we have the boundary term from the integration by parts which is dropped being equal to 0

\[ \int_{\partial \Omega} \frac{\partial u}{\partial n} v = 0, \]

recalling that

\[ \label{eq:21} \frac{\partial u}{\partial n} \stackrel{\text{def}}{=} \nabla u \cdot n \]

where $n$ is the outward normal to $\partial \Omega$ by convention.We now discretize the problem, we create a mesh out of $\Omega$, we have

\[ \label{eq:10} \Omega = \cup_{e=1}^\nel \Omega^e \]

where $\Omega^e$ can be segments, triangles or tetrahedra depending on $d$ and we have $\nel$ of them. We introduce the finite dimensional spaces of continuous piecewise polynomial of degree $N$ functions

\[ \label{eq:17} X_h = \Big\{ v_h \in C^0(\Omega),\ {v_h}_{|\Omega^e} \in \mathbb{P}_N( \Omega^e ),\ v_h=g_{|x=-1,x=1}\Big\} \]

and

\[ \label{eq:18} V_h = \Big\{ v_h \in C^0(\Omega),\ {v_h}_{|\Omega^e} \in \mathbb{P}_N( \Omega^e ),\ v_h=0_{|x=-1,x=1}\Big\} \]

which are out trial and test function spaces respectively. We now have the problem we seek to solve which reads in our continuous Galerkin framework

we look for $u_h \in X_h \subset \mathcal{X}$ such that for all $v \in V_h \subset \mathcal{V}$

\[ \int_\Omega \nabla u_h \cdot \nabla v_h = \int_\Omega f v_h \]

Weak Dirichlet conditions

There is an alternative formulation which allows to treat weakly Dirichlet(Essential) boundary conditions similarly to Neumann(Natural) and Robin conditions. Following a similar development as in the previous section, the problem reads

We look for $u \in X_h \subset H_1(\Omega)$ such that for all $v \in X_h$

\[ \int_\Omega \nabla u \cdot \nabla v + \int_{|x=-1,x=1} -\frac{\partial u}{\partial n} v - u \frac{\partial v}{\partial n} + \frac{\mu}{h} u v = \int_\Omega f v + \int_{|x=-1,x=1} - g \frac{\partial v}{\partial n} + \frac{\mu}{h} g v \]

where

\[ \label{eq:19} X_h = \Big\{ v_h \in C^0(\Omega),\ {v_h}_{|\Omega^e} \in \mathbb{P}_N( \Omega^e ) \Big\} \]

In (eq:16}), $g$ is defined by (eq:3}). $\mu$ serves as a penalisation parameter which should be $> 0$, e.g. between 2 and 10, and $h$ is the size of the face. The inconvenient of this formulation is the introduction of the parameter $\mu$, but the advantage is the {weak} treatment of the Dirichlet condition.

Feel formulation

{sec:feel-formulation-1}

First we define the $f$ and $g$. To do that we use the new C++ auto keyword and associate to f and g their corresponding expressions

    //# marker1 #
    value_type pi = M_PI;
    auto g = sin(pi*Px())*cos(pi*Py())*cos(pi*Pz());
    gproj = vf::project( Xh, elements(mesh), g );

    auto f = pi*pi*Dim*g;
    //# endmarker1 #

where M_PI is defined in the header cmath. Using auto allows to defined f and g --- which are moderately complex object --- without having to know the actual type. auto determines automatically the type of the expression using the __typeof__ keyword internally.

Then we form the right hand side by defining a linear form whose algebraic representation will be stored in a vector_ptrtype which is provided by the chosen linear algebra backend. The linear form is equated with an integral expression defining our right hand side.

    //# marker2 #
    vector_ptrtype F( M_backend->newVector( Xh ) );
    form1( _test=Xh, _vector=F, _init=true ) =
        integrate( elements(mesh), f*id(v) )+
        integrate( markedfaces( mesh, mesh->markerName("Neumann") ), nu*gradv(gproj)*vf::N()*id(v) );
    //# endmarker2 #

form1 generates an instance of the object representing linear forms, that is to say it mimics the mathematical object $\ell$ such that

\[ \begin{array}{rccl} \ell: & X_h & \mapsto & \mathbb{R}\\ & v_h & \rightarrow &\ell(v_h)=\int_\Omega f v_h \end{array} \]

which is represented algebraically in the code by the vector F using the argument _vector. The last argument _init, if set to true, will zero-out the entries of the vector F.

Note:
_init is set to false by default.

We now turn to the left hand side and define the bilinear form using the form2 helper function which is passed (i) the trial function space using the _trial option, (ii) the test function space using the _test option, (iii) the algebraic representation using _matrix, i.e. a sparse matrix whose type is derived from one of the linear algebra backends and (iv) whether the associated matrix should initialized using _init.

    //# marker3 #
    sparse_matrix_ptrtype D( M_backend->newMatrix( Xh, Xh ) );


    form2( Xh, Xh, D, _init=true ) =
        integrate( elements(mesh), nu*gradt(u)*trans(grad(v)) );
    //# endmarker3 #

Finally, we deal with the boundary condition, we implement both formulation described in appendix~sec:vari-form-1}. For a {strong} treatment of the Dirichlet condition, we use the on() keyword of FEEL++ as follows

            //# marker5 #
            D->close();
            form2( Xh, Xh, D ) +=
                on( markedfaces(mesh, mesh->markerName("Dirichlet")), u, F, g );
            //# endmarker5 #

Notice that we add, using +=, the Dirichlet contribution for the bilinear form. The first argument is the set of boundary faces to apply the condition: in gmsh the points satisfying $x=\pm 1$ are marked using the flags $1$ and $3$ ( $x=-1$ and $x=1$ respectively.)

To implement the weak Dirichlet boundary condition, we add the following contributions to the left and right hand side:

            //# marker41 #
            form1( _test=Xh, _vector=F ) +=
                integrate( markedfaces(mesh,mesh->markerName("Dirichlet")), g*(-grad(v)*vf::N()+penaldir*id(v)/hFace()) );
            //# endmarker41 #

            //# marker10 #
            form2( Xh, Xh, D ) +=
                integrate( markedfaces(mesh,mesh->markerName("Dirichlet")),
                           -(gradt(u)*vf::N())*id(v)
                           -(grad(v)*vf::N())*idt(u)
                           +penaldir*id(v)*idt(u)/hFace());
            D->close();
            //# endmarker10 #

Note that we use the command line option --weakdir set to 1 by default to decide between weak/strong Dirichlet handling. Apart the uniform treatment of boundary conditions, the weak Dirichlet formulation has the advantage to work also in a parallel environment.

Next we solve the linear system

\[ D u = F \]

where the solve function is implemented as follows

Finally we check for the $L_2$ error in our approximation by computing

\begin{equation} \label{eq:7} \|u-u_h\|_{L_2}\ =\ \sqrt{\int_\Omega (u-u_h)^2} = \sqrt{\int_\Omega (g-u_h)^2} \end{equation}

where $u$ is the exact solution and is equal to $g$ and $u_h$ is the numerical solution of the problem and the components of $u_h$ in the $P_2$ Lagrange basis are given by solving the linear system above.

The code reads

    //# marker7 #
    double L2error2 =integrate(elements(mesh),
                               (idv(u)-g)*(idv(u)-g) ).evaluate()(0,0);
    double L2error =   math::sqrt( L2error2 );


    Log() << "||error||_L2=" << L2error << "\n";
    //# endmarker7 #

You can now verify that the $L_2$ error norm behaves like $h^{-(N+1)}$ where $h$ is the mesh size and $N$ the polynomial order. The $H_1$ error norm would be checked similarly in $h^{-N}$. The figure~fig:2} displays the results using Paraview.

laplacian.png
Solution of the Laplacian problem
laplacian_warp.png
Warped solution of the Laplacian problem

Mixed formulation: the Stokes case

We are now interested in solving the Stokes equations, we would like to solve for the following problem in 2D

find $(\mathbf{u},p)$ such that {