v0.8.13
Solving the Poisson equation

Table of Contents

In this tutorial is shown how to solve Poisson equation. We introduce the basic concepts in MoFEM and exploiting MoFEM::Simple interface to set-up problem. A procedure is applied to the Poisson problem but easily can be modified to other problems from elasticity, fluid mechanics, thermo-elasticity, electromagnetics, etc. Also, methodology from this tutorial can be applied to mix-problems, discontinuous Galerkin method or coupled problems.

We introduce fundamentals of MoFEM::Interface and MoFEM::Simple interfaces and show how to link MoFEM data structures with PETSc discrete manager. Also, we show how to describe fields and how to loop/iterate finite elements. This example is restricted to linear problem and body is subjected to inhomogeneous Dirichlet boundary conditions only. We do not explain here how to implement finite element data operator; this is part of tutorial Implementing operators for the Poisson equation.

While writing this tutorial we don't know what type of person you are. Do you like to see examples first, before seeing the big picture? Or do you prefer to see big picture first and then look for the details after? If you are second kind, you would like to see this first Basic design before you start with implementation. If you like to start from examples, skip the general details and go to code.

MoFEM has a bottom-up design, have several levels of abstraction, yet enables users at any level to hack levels below. That is intended to give freedom to create new methods and solve problems in new ways.

Each level up simplifies implementation, making code short and easier to follow, from the other hand limiting flexibility. In this tutorial we focus attention on functions from level 3, we are not going to describe here how finite element is implemented, this is part of the next tutorial, Implementing operators for the Poisson equation. In this one, we focus on generic steps presented in figure below in yellow boxes.

We designed code to give you a freedom to create code which does what you want and does not force you to follow beaten path. The implementation shown here is similar other types of PDEs, even for nonlinear or time-dependent problems. You don't have to understand every detail of the code here, some times we are using advanced C++. That understanding will come later. For example, if you like to change a form of your PDE, then you need to look only how to modify differential operator to solve anisotropic the Darcy equation for flow, and that is easily done by modification of one of the operators. So, relax and look what you have below and when you are intrigued, always ask questions on our Q&A forum.

poisson_tut1_fig1.png
Figure 1. Level 3 interfaces

Mathematical formulation

We solve the most basic of all PDEs, i.e. the Poisson equation,

\[ -\nabla^2 u(\mathbf{x}) = f(x),\quad \textrm{in}\,\Omega\\ \]

where \(u=u(\mathbf{x})\) is unknown function, \(f(\mathbf{x})\) is a given source term. The Poisson equation is in so-called strong form, it has to be satisfied at every point of the body without boundary. To solve this PDE we need to describe the problem on the boundary, in this simple case for the purpose of tutorial we are self-restricted to Dirichlet boundary condition where functions values are given by

\[ u(x)=\overline{u}(\mathbf{x}),\quad \textrm{on}\,\partial\Omega \]

where \(\overline{u}\) is prescribed function on the boundary.

The Laplacian \(\nabla^2\) in 3d space with Cartesian coordinates is

\[ \nabla^2 u(x,y,z) = \nabla \cdot \nabla u = \frac{\partial^2 u}{\partial x^2}+ \frac{\partial^2 u}{\partial y^2}+ \frac{\partial^2 u}{\partial z^2} \]

The Poisson equation arises in many physical contexts, including heat conduction, electrostatics, diffusion of substances, twisting of elastic rods, and describes many other physical problems. In general PDE equation emerge from conservation law, e.g. conservation of mass, electric charge or energy, i.e. one of the fundamental laws of the universe.

Finite element variational formulation

It is well-established recipe how to change the equation in the form convenient for finite element calculations. The Poisson equation is multiplied by function \(u\) and integrated by parts. Without going into mathematical details, applying this procedure and enforcing constraints with Lagrange multiplier method we get

\[ \mathcal{L}(u,\lambda) = \frac{1}{2}\int_\Omega \left( \nabla u \right)^2 \textrm{d}\Omega - \int_\Omega uf \textrm{d}\Omega + \int_{\partial\Omega} \lambda (u-\overline{u}) \textrm{d}\partial\Omega \]

where \(\lambda=\lambda(\mathbf{x})\) is Lagrange multiplayer function on body boundary \(\partial\Omega\). Unknown functions can be expressed in finite dimension approximation space, as follows

\[ u^h = \sum_i^{n_u} \phi^i(\mathbf{x}) U_i = \boldsymbol\phi\cdot\mathbf{U} , \quad \mathbf{x} \in \Omega,\\ \lambda^h = \sum_i^{n_\lambda} \psi^i(\mathbf{x}) L_i = \boldsymbol\psi\cdot\mathbf{L} , \quad \mathbf{x} \in \partial\Omega \]

where \(U_i\) and \(L_i\) are unknown coefficients, called degrees freedom. The \(\phi^u\) and \(\psi^i\) are base functions for \(u^h\) and \(\lambda^h\), respectively. The stationary point of \(\mathcal{L}(u^h,\lambda^h)\), obtained by minimalisation of Lagrangian gives Euler-Lagrange equations for discretised problem

\[ \left[ \begin{array}{cc} \mathbf{K} & \mathbf{C}^\textrm{T}\\ \mathbf{C} & \mathbf{0} \end{array} \right] \left\{ \begin{array}{c} \mathbf{U}\\ \mathbf{L} \end{array} \right\} = \left[ \begin{array}{c} \mathbf{F}\\ \mathbf{g} \end{array} \right],\\ \mathbf{K}= \int_\Omega (\nabla \boldsymbol\phi)^\textrm{T} \nabla \boldsymbol\phi \textrm{d}\Omega,\quad \mathbf{C} = \int_{\partial\Omega} \boldsymbol\psi^\textrm{T} \boldsymbol\phi \textrm{d}\partial\Omega,\\ \mathbf{F} = \int_\Omega \boldsymbol\phi^\textrm{T} f \textrm{d}\Omega,\quad \mathbf{g} = \int_{\partial\Omega} \boldsymbol\psi^\textrm{T} \overline{u} \textrm{d}\partial\Omega. \]

Applying this procedure we transformed essential boundary conditions (in this case Dirichlet boundary conditions), into problems with natural boundary conditions only. Essential boundary conditions need to be a priory satisfied by approximation base. Natural conditions are satisfied while the problem is solved. Note that Lagrange multiplier has interpretation of flux here, which is in \(H^{-\frac{1}{2}}(\partial\Omega)\) space. However in this case approximation of Lagrange multiplier by a discontinuous piecewise result in a singular system of equations, since values are not bounded for natural (Neumann) boundary conditions when prescribed on all \(\partial\Omega\). For such case, the solution is defined up to an additive constant. To remove that difficult we will use subspace of trace space, i.e. \(H^{-1}(\partial\Omega)\), that will generate as many constraints as many degrees of freedom for \(u^h\) on body boundary, and thus enforce boundary conditions exactly. For more information about spaces and deeper insight in the Poisson equation see [19], [13] and [35].

Matrices and vectors are implemented using users data operators (UDO);

We put light on the implementation of those functions in the next tutorial Implementing operators for the Poisson equation. Here, in the following sections, we focus on the generic problem and solver setup. Complete implementation of UDO for this example is here PoissonOperators.hpp.

Exact solution

This example is designed to validate the basic implementation of finite elements operators. In this particular case, we assume exact solution. We will focus attention on particular case for which analytical solution is given by polynomial, for example

\[ \overline{u}(\mathbf{x}) = 1+x^2+y^2+z^3 \quad \textrm{in}\Omega \]

then applying Laplacian to analytical solution we can obtain source term

\[ f(\mathbf{x}) = 4+6z. \]

If finite element approximation using 3rd order polynomials or higher, a solution of the Poisson problem is exact. In MoFEM you could easily increase the order of the polynomial, as it is shown below.

Code dissection

The program with comments can be accessed here analytical_poisson.cpp. In this section, we will go step-by-step through the code.

Header files

We including three header files, first with MoFEM library files and basic finite elements, second with the implementation of operators for finite elements and third with some auxiliary functions which will share tutorials about the Poisson equation.

Exact solution functions

We choose arbitrary 3rd polynomial function; we add function for the gradient to evaluate error on element and function for Laplacian to apply it as a source term in the Poisson equation.

struct ExactFunction {
double operator()(const double x,const double y,const double z) const {
return 1+x*x+y*y+z*z*z;
}
};
FTensor::Tensor1<double,3> operator()(const double x,const double y,const double z) const {
grad(0) = 2*x;
grad(1) = 2*y;
grad(2) = 3*z*z;
return grad;
}
};
double operator()(const double x,const double y,const double z) const {
return 4+6*z;
}
};

Initialization

In the following section, we initialize PETSc and construct MoAB and MoFEM databases. Starting with PETSc and MoAB

// Initialize PETSc
MoFEM::Core::Initialize(&argc,&argv,(char *)0,help);
// Create MoAB database
moab::Core moab_core; // create database
moab::Interface& moab = moab_core; // create interface to database

and read data from a command line argument

int order = 3; // default approximation order
PetscBool flg_test = PETSC_FALSE; // true check if error is numerical error
CHKERR PetscOptionsBegin(PETSC_COMM_WORLD,"", "Poisson's problem options","none");
// Set approximation order
CHKERR PetscOptionsInt("-order","approximation order","",order,&order,PETSC_NULL);
// Set testing (used by CTest)
CHKERR PetscOptionsBool("-test","if true is ctest","",flg_test,&flg_test,PETSC_NULL);
CHKERR PetscOptionsEnd();

In this test, for simplicity we apply for approximation order uniformly. However, MoFEM in principle is designed to manage heterogeneous approximation bases.

Once we have initialized PETSc and MoAB, we create MoFEM database instance and link to it MoAB. MoFEM database is accessed by MoFEM::Interface.

MoFEM::Core mofem_core(moab); // create database
MoFEM::Interface& m_field = mofem_core; // create interface to database

Also, we need to register in PETSc implementation of MoFEM Discrete Manager

CHKERR DMRegister_MoFEM("DMMOFEM"); // register MoFEM DM in PETSc

Creating Finite element instances

We start with declaring smart pointers to finite element objects. A smart pointer is like a "raw" pointer, but you do not have to remember to delete allocated memory to which it points. More about shared pointers you find here, if you don't know what is a pointer, look here. You do not have to look at the details now; you can start working with code by making modifications of solutions which are already there.

boost::shared_ptr<ForcesAndSourcesCore> domain_lhs_fe; ///< Volume element for the matrix
boost::shared_ptr<ForcesAndSourcesCore> boundary_lhs_fe; ///< Boundary element for the matrix
boost::shared_ptr<ForcesAndSourcesCore> domain_rhs_fe; ///< Volume element to assemble vector
boost::shared_ptr<ForcesAndSourcesCore> boundary_rhs_fe; ///< Volume element to assemble vector
boost::shared_ptr<ForcesAndSourcesCore> domain_error; ///< Volume element evaluate error
boost::shared_ptr<ForcesAndSourcesCore> post_proc_volume; ///< Volume element to Post-process results
boost::shared_ptr<ForcesAndSourcesCore> null; ///< Null element do nothing

Now we have to allocate memory for the finite elements, add data operators to them. We do not cover details here, this is realised by three functions, which we will use in the following tutorials

// Add problem specific operators the generic finite elements to calculate matrices and vectors.
domain_lhs_fe,boundary_lhs_fe,domain_rhs_fe,boundary_rhs_fe
);
// Add problem specific operators the generic finite elements to calculate error on elements and global error
// in H1 norm
ExactFunction(),ExactFunctionGrad(),global_error,domain_error
);
// Post-process results

Note that those implementations are problem independent, we could use them in different problems context. To make that more visible, we create those elements before fields and problem are defined and build. The purpose of finite elements is shown in Figure 2 and Figure 3, by users data operators added to a finite element. Finite element "domain_lhs_fe" is used to calculate matrix by using PoissonExample::OpK, finite element "domain_rhs_fe" is used with PoissonExample::OpF. Similarly for a boundary, operators PoissonExample::OpC and PoissonExample::Op_g are used to integrate matrices and right vector for Lagrange multipliers, with finite elements "boundary_lhs_fe" and "boundary_rhs_fe", respectively. The operators MoFEM::OpCalculateScalarFieldValues, MoFEM::OpCalculateScalarFieldGradient and PoissonExample::OpError are used to integrate error in H1 norm in domain. Finite element "post_proc_volume" is used to post-process results for ParaView (or another visualization tool), it has set of specialised operators saving results on nodes and allows to refine mesh when ho-approximation is used. The result can be stored for example on 10-nodes tetrahedra if needed. More details about implementation of PoissonExample::CreateFiniteElements and operators is in Implementing operators for the Poisson equation. Basic of user data operators are explained here Basic design.

poisson_tut1_fig2.png
Figure 2. Operators of volume finite element
poisson_tut1_fig3.png
Figure 3. Operators of boundary finite element

Simple interface

Solution of the problem

Before we start running PETSc solver we need to get access to DM manager, it can be simply done by

DM dm;
CHKERR simple_interface->getDM(&dm);

Note that user takes responsibility for deleting DM object.

Once we have DM, we add instances of previously created finite elements with their operators to DM manager, appropriately to calculate matrix and the right-hand side

// Set operators for KSP solver
dm,simple_interface->getDomainFEName(),domain_lhs_fe,null,null
);
dm,simple_interface->getBoundaryFEName(),boundary_lhs_fe,null,null
);
// Set calculation of the right hand side vetor for KSP solver
dm,simple_interface->getDomainFEName(),domain_rhs_fe,null,null
);
dm,simple_interface->getBoundaryFEName(),boundary_rhs_fe,null,null
);

At that point we can solve problem, applying exclusively functions from native PETSc interface

// Create the right hand side vector and vector of unknowns
Vec F,D;
CHKERR DMCreateGlobalVector(dm,&F);
// Create unknown vector by creating duplicate copy of F vector. only
// structure is duplicated no values.
CHKERR VecDuplicate(F,&D);
// Create solver and link it to DM
KSP solver;
CHKERR KSPCreate(PETSC_COMM_WORLD,&solver);
CHKERR KSPSetFromOptions(solver);
CHKERR KSPSetDM(solver,dm);
// Set-up solver, is type of solver and pre-conditioners
CHKERR KSPSetUp(solver);
// At solution process, KSP solver using DM creates matrices, Calculate
// values of the left hand side and the right hand side vector. then
// solves system of equations. Results are stored in vector D.
CHKERR KSPSolve(solver,F,D);

In above code, KSP solver takes control. The KSP solver is using DM creates matrix and calls DM functions to iterate over finite elements and assembles matrix and vectors. KSP solver knows how to do it since it is using MoFEM DM which we registered at the very beginning, and get an instance of particular DM from MoFEM::Simple interface. While DM is iterating over finite elements, finite elements instance executes user data operators.

Once we have the solution, DOFs values in solution vector are scattered on the mesh

CHKERR DMoFEMMeshToGlobalVector(dm,D,INSERT_VALUES,SCATTER_REVERSE);

It is evident here the role of MoFEM, having topology (mesh) in MoAB database, MoFEM is used to manage complexities related to problem discretization and construction of algebraic systems of equations (linear in this case). While, PETSc manages complexities related to algebra and solution of the problem. For more details about the interaction of MoAB, MoFEM and PETSc see MoFEM software eccosystem (see Basic design).

Finally, we can evaluate error and verify correctness of the solution

CHKERR DMoFEMLoopFiniteElements(dm,simple_interface->getDomainFEName(),domain_error);
if(flg_test == PETSC_TRUE) {
}

Note that function DMoFEMLoopFiniteElements takes finite element instance "domain_error" and applies it to finite elements entities in the problem. For each finite element entity, finite element instance is called, and sequence of operators data operators are executed on that element.

Post-processing

Finally, we will create a mesh for post-processing. Note that in MoFEM we work with spaces that not necessarily have physical DOFs at the nodes, i.e. higher order polynomials or spaces like H-div, H-curl or L2. Discontinuities on elements edges have to be shown in post-processing accurately. To resolve this issue, a set of user data operators has been created to manage those problems, see PostProcCommonOnRefMesh for details. Here, we directly run the code

CHKERR DMoFEMLoopFiniteElements(dm,simple_interface->getDomainFEName(),post_proc_volume);
// Write results
CHKERR boost::static_pointer_cast<PostProcVolumeOnRefinedMesh>(post_proc_volume)->writeFile("out_vol.h5m");

Running the program

The easiest installation is with Docker containers; please see youtube video how to run installation link. You can also look here Installation with Docker (Linux, Mac OS X and some versions of Windows). Open terminal and run Docker container

docker run --rm=true -it --volumes-from mofem_build -v $HOME/mofem-cephas/mofem:/mofem -v $HOME:$HOME -e HOSTHOME=$HOME mofem_build /bin/bash

You are now in Docker container, go to directory where Poisson problem is located, and compile code

cd /mofem_build/um/basic_finite_elements/poisson/
make -j2

It can take some time when you compile it first time since all essential functionality is compiled

Once code is compiled, we can run test

mpirun -np 2 ./analytical_poisson -file_name cube_2part.h5m -order 3

Output dissection

Initialization

MoFEM version and number of commit

lib version 0.5.75
git commit id 4977ff148919cdd1577d07b0082b3147abc778ed

Declarations

Set-up

Solution

Visualization

Once we have solution, you can convert file to VTK and open it in ParaView

mbconvert out_vol.h5m out_vol.vtk

and copy your file from docker to the desktop

cp out_vol.vtk $HOSTHOME/Desktop/

Now your file is visible on Desktop, run ParaView and open file.

poisson_tut1_fig4.png
Figure 4. ParaView and result

Mesh formats

Since we are using a polynomial of 3rd order to approximate solution, results will not change if you would use coarser or finer mesh. You can generate mesh in gMesh, Salome or generate mesh in any format on the list

Format        Name                            Read    Write   File name description
-----------------------------------------------------------------------------------
MOAB          MOAB native (HDF5)              yes    yes     h5m mhdf
EXODUS        Exodus II                       yes    yes     exo exoII exo2 g gen
NC            Climate NC                      yes    yes     nc
UNV           IDEAS format                    yes     no     unv
MESHTAL       MCNP5 format                    yes     no     meshtal
NAS           NASTRAN format                  yes     no     nas bdf
Abaqus        mesh  ABAQUS INP mesh format    yes     no     abq
Atilla        RTT Mesh  RTT Mesh Format       yes     no     rtt
VTK           Kitware VTK                     yes    yes     vtk
OBJ mesh      OBJ mesh format                 yes     no     obj
SMS           RPI SMS                         yes     no     sms
CUBIT         Cubit                           yes     no     cub
SMF           QSlim format                    yes    yes     smf
SLAC          SLAC                             no    yes     slac
GMV           GMV                              no    yes     gmv
ANSYS         Ansys                            no    yes     ans
GMSH          Gmsh mesh file                  yes    yes     msh gmsh
STL           Stereo Lithography File (STL)   yes    yes     stl
TETGEN        TetGen output files             yes     no     node ele face edge
TEMPLATE      Template input files            yes    yes

You can obtain that list executing from the command line

mbconvert -l

Exercises

Please feel free to ask question on our Q&A forum, see link.

Notes

Note
This is similar test to one from https://fenicsproject.org/pub/tutorial/html/._ftut1004.html#ch:fundamentals.