Quick Starts with singularity

1. Installation Quick Start

Using Feel++ inside Docker is the recommended and fastest way to use Feel++. The Docker chapter is dedicated to Docker and using Feel++ in Docker.

We strongly encourage you to follow these steps if you begin with Feel++ in particular as an end-user.

People who would like to develop with and in Feel++ should read through the remaining sections of this chapter.

2. Usage Start

Start the Docker container feelpp/feelpp-base or feelpp/feelpp-toolboxes as follows

> docker run -it -v $HOME/feel:/feel feelpp/feelpp-toolboxes
these steps are explained in the chapter on Feel++ containers.

Then run e.g. the Quickstart Laplacian that solves the Laplacian problem in Quickstart Laplacian sequential or in Quickstart Laplacian on 4 cores in parallel.

Quickstart Laplacian sequential
> feelpp_qs_laplacian_2d --config-file Testcases/quickstart/laplacian/feelpp2d/feelpp2d.cfg

The results are stored in Docker in

/feel/qs_laplacian/feelpp2d/np_1/exports/ensightgold/qs_laplacian/

and on your computer

$HOME/feel/qs_laplacian/feelpp2d/np_1/exports/ensightgold/qs_laplacian/

The mesh and solutions can be visualized using e.g. Parariew or Visit.

ParaView (recommended)

is an open-source, multi-platform data analysis and visualization application.

Visit

is a distributed, parallel visualization and graphical analysis tool for data defined on two- and three-dimensional (2D and 3D) meshes

Quickstart Laplacian on 4 cores in parallel
> mpirun -np 4 feelpp_qs_laplacian_2d --config-file Testcases/quickstart/laplacian/feelpp2d/feelpp2d.cfg

The results are stored in a simular place as above: just replace np_1 by np_4 in the paths above. The results should look like

ufeelpp2d

meshfeelpp2d

Solution

Mesh

3. Syntax Start

Here are some excerpts from Quickstart Laplacian that solves the Laplacian problem. It shows some of the features of Feel++ and in particular the domain specific language for Galerkin methods.

First we load the mesh, define the function space define some expressions

Laplacian problem in an arbitrary geometry, loading mesh and defining spaces
    tic();
    auto mesh = loadMesh(_mesh=new Mesh<Simplex<FEELPP_DIM,1>>);
    toc("loadMesh");

    tic();
    auto Vh = Pch<2>( mesh ); (1)
    auto u = Vh->element("u"); (2)
    auto mu = expr(soption(_name="functions.mu")); // diffusion term (3)
    auto f = expr( soption(_name="functions.f"), "f" ); (4)
    auto r_1 = expr( soption(_name="functions.a"), "a" ); // Robin left hand side expression (5)
    auto r_2 = expr( soption(_name="functions.b"), "b" ); // Robin right hand side expression (6)
    auto n = expr( soption(_name="functions.c"), "c" ); // Neumann expression (7)
    auto solution = expr( checker().solution(), "solution" ); (8)
    auto g = checker().check()?solution:expr( soption(_name="functions.g"), "g" ); (9)
    auto v = Vh->element( g, "g" ); (3)
    toc("Vh");

Second we define the linear and bilinear forms to solve the problem

Laplacian problem in an arbitrary geometry, defining forms and solving
    tic();
    auto l = form1( _test=Vh );
    l = integrate(_range=elements(mesh),
                  _expr=f*id(v));
    l+=integrate(_range=markedfaces(mesh,"Robin"), _expr=r_2*id(v));
    l+=integrate(_range=markedfaces(mesh,"Neumann"), _expr=n*id(v));
    toc("l");

    tic();
    auto a = form2( _trial=Vh, _test=Vh);
    tic();
    a = integrate(_range=elements(mesh),
                  _expr=mu*inner(gradt(u),grad(v)) );
    toc("a");
    a+=integrate(_range=markedfaces(mesh,"Robin"), _expr=r_1*idt(u)*id(v));
    a+=on(_range=markedfaces(mesh,"Dirichlet"), _rhs=l, _element=u, _expr=g );
    //! if no markers Robin Neumann or Dirichlet are present in the mesh then
    //! impose Dirichlet boundary conditions over the entire boundary
    if ( !mesh->hasAnyMarker({"Robin", "Neumann","Dirichlet"}) )
        a+=on(_range=boundaryfaces(mesh), _rhs=l, _element=u, _expr=g );
    toc("a");

More explanations are available in the Laplacian example.