Getting Started

Requirements

FMG requires a functional FEniCS 1.0 installation (does currently not work with newer versions) with PETSc >=3.3 (might also work with older 3.x versions).

Installation

For TU-Chemnitz members the following installation procedure was tested successfully. It requires VPN or direct access to the MRZ:

ssh -X -l username 134.109.41.1  # connect to newton (only necessary when connecting from outside)
ssh -X riemann
source /usr/local/FEniCS/1.0.0/setup_fenics_with_common
wget www-user.tu-chemnitz.de/~feo/fmg/fmg.tar.gz
tar -xf fmg.tar.gz
cd fmg-src
sh install.sh

For non TU-Chemnitz members use the following guide:

  1. Setup your Dolfin environment (if not already):

    source /path_to_your_dolfin_installation/dolfin.conf
    
  2. Download the source code:

    wget www-user.tu-chemnitz.de/~feo/fmg/fmg.tar.gz
  3. Extract the files:

    tar -xf fmg.tar.gz
  4. Run install script:

    cd fmg-src
    sh install.sh
  5. If you get errors use the following command to create the makefile with a user specific installation path for fmg:

    set FMG_PATH=~/fmg
    cmake -DCMAKE_INSTALL_PREFIX:PATH=$FMG_PATH .

    If you get errors you probably need to:

    locate dolfin-config.cmake
    export dolfin_PATH=<DIRECTORY TO dolfin-config.cmake>

    After that run:

    make && make install
  1. Setup your environment variables:

    source $FMG_PATH/fmg.conf

    For future reference, we recommend that you add this command to your configuration (.bashrc, .profile or similar).

Linking with FMG

If you use cmake to build your projects you need to add fmg to the list of linked libraries. Therefore you need to edit your CMakeLists.txt and add the line:

target_link_libraries(${PROJECT_NAME} fmg)

For other build systems you find all required libraries and include files within $FMG_PATH.

Using FMG

Solving Poisson equation

The following example shows the modified demo for solving the Poisson equation (included in the Dolfin demos) using multigrid:

// This demo program solves Poisson's equation
//
//     - div grad u(x, y) = f(x, y)
//
// on the unit square with source f given by
//
//     f(x, y) = 10*exp(-((x - 0.5)^2 + (y - 0.5)^2) / 0.02)
//
// and boundary conditions given by
//
//     u(x, y) = 0        for x = 0 or x = 1
// du/dn(x, y) = sin(5*x) for y = 0 or y = 1

#include <dolfin.h>
#include <fmg.h>
#include "Poisson.h"

using namespace dolfin;

// Source term (right-hand side)
class Source : public Expression
{
        void eval(Array<double>& values, const Array<double>& x) const
        {
                double dx = x[0] - 0.5;
                double dy = x[1] - 0.5;
                values[0] = 10*exp(-(dx*dx + dy*dy) / 0.02);
        }
};

// Normal derivative (Neumann boundary condition)
class dUdN : public Expression
{
        void eval(Array<double>& values, const Array<double>& x) const
        {
                values[0] = sin(5*x[0]);
        }
};

// Sub domain for Dirichlet boundary condition
class DirichletBoundary : public SubDomain
{
        bool inside(const Array<double>& x, bool on_boundary) const
        {
                return x[0] < DOLFIN_EPS or x[0] > 1.0 - DOLFIN_EPS;
        }
};

int main()
{
        // Create mesh and function space
        UnitSquare mesh(4, 4);
        Poisson::FunctionSpace V(mesh);

        // Define boundary condition
        Constant u0(0.0);
        DirichletBoundary boundary;
        DirichletBC bc(V, u0, boundary);

        // Define variational forms
        Poisson::BilinearForm a(V, V);
        Poisson::LinearForm L(V);
        Source f;
        dUdN g;
        L.f = f;
        L.g = g;

        // Compute solution using 3 uniform mesh refinements
        Function u(V);
        LinearVariationalProblem problem(a == L, u, bc);
        fmg::MultigridPreconditionedKrylovSolver fmg_solver(problem, 3);
        fmg_solver.solve();

        // Get solution on the fine mesh
        Function& u_fine = *fmg_solver.solution();

        // Save solution in VTK format
        File file("poisson.pvd");
        file << u_fine;

        // Plot solution
        plot(u_fine);
        interactive();

        return 0;
}

Configuring the solver

The FMG solver can be configured with various settings. The following example uses a more detailed instantiation for the solver and also changes the smoother:

LinearVariationalProblem problem(a, L, u, bcs);

// create multigrid solver
fmg::MultigridSolver mg_solver(problem, 3);

// set some parameters (optional)
mg_solver.parameters["pre_smoother"] = "jacobi";
mg_solver.parameters["pre_smoother_relax"] = 0.6;
mg_solver.parameters["post_smoother"] = "jacobi";
mg_solver.parameters["post_smoother_relax"] = 0.6;
mg_solver.parameters["direct_solver_type"] = "lu";

// create multigrid pcg solver
fmg::MultigridPreconditionedKrylovSolver mgpcg_solver(mg_solver, "cg");

// solve problem
mgpcg_solver.solve();

// plot solution
plot(mgpcg_solver.solution());

For more information refer to the MultigridSolver class reference or have a look at the demo directory.

Performance Testing

In order to see if FMG is right for your problem. You can run a quick performance analysis on your problem:

LinearVariationalProblem problem(a, L, u, bcs);

// create tests instance
fmg::Tests tests(problem);

// set some parameters (optional)
tests.parameters("fmg_solver")["pre_smoother"] = "jacobi";
tests.parameters("fmg_solver")["pre_smoother_relax"] = 0.6;
tests.parameters("fmg_solver")["post_smoother"] = "jacobi";
tests.parameters("fmg_solver")["post_smoother_relax"] = 0.6;
tests.parameters["test_solver"] = "cg,cg+fmg";
tests.parameters["num_refinements"] = 6;
tests.parameters.parse(argc, argv);

// run the tests on problem
tests.run();
which will print out something like:

The following table computes the additional time that is needed for the fmg solver with 7 refinements to initialize in contrast to using a conventional solver (with no additional initialization time) on the fine mesh directly.

name description time loss %
t_mesh direct problem mesh creation 0.0905  
t_mesh_mg multigrid adaptive mesh creation 0.247  
t_loss_mesh loss by mesh creation 0.157 9.24
t_init direct problem initialization (without t_mesh) 0.812  
t_init_mg multigrid adaptive problem initialization (without t_mesh_mg) 1.09  
t_loss_init loss by problem initialization 0.276 16.3
t_asm_bc direct problem initialization of boundary values 1.71  
t_asm_bc_mg multigrid problem initialization of boundary values 2.26  
t_loss_asm_bc loss by initialization of boundary values 0.558 32.9
t_asm_lhs direct problem assembly of lhs 0.450  
t_asm_lhs_mg multigrid problem assembly of lhs 0.602  
t_loss_asm_lhs loss by assembly of lhs 0.152 8.97
t_asm_rhs direct problem assembly of rhs 0.327  
t_asm_rhs_mg multigrid problem assembly of rhs 0.439  
t_loss_asm_rhs loss by assembly of rhs 0.112 6.61
t_loss_pro loss by prolongation operator creation 0.321 18.9
t_loss_sr loss by system reduction 0.111 6.57
t_loss_other other losses 0.00863 0.509
t_direct total direct problem initialization (“ready to solve”) 3.39  
t_mg total multigrid problem initialization (“ready to solve”) 5.08  
t_loss total losses 1.70

The following table shows t_loss for each refinement level.

refinements dofs t_direct t_mg t_loss t_loss %
0 48 0.00183 0.00267 0.000840 45.8
1 195 0.00370 0.00812 0.00442
2 783 0.0115 0.0220 0.0105 90.8
3 3135 0.0488 0.0780 0.0291 59.7
4 12543 0.208 0.311 0.104 49.8
5 50175 0.854 1.27 0.417 48.9
6 200703 3.39 5.08 1.70 50.1

The following table shows the number of iterations and solve time versus the number of mesh refinements for all methods listed in the “test_solver” parameter. All solver timings don’t include the time for initialization of the solver. You have to add the t_loss from the previous table to the timings for all solvers using fmg to have the total time that is required when doing only a single solve.

refinements dofs cg iter cg time fmg iter fmg time cg+hypre_amg iter cg+hypre_amg time cg+fmg iter cg+fmg time
0 48 22 0.000110 1 0.000216 3 0.000444 1 0.000351
1 195 44 0.000347 9 0.00155 4 0.00135 5 0.00135
2 783 87 0.00234 10 0.00331 4 0.00451 5 0.00244
3 3135 172 0.0181 10 0.00843 4 0.0187 5 0.00617
4 12543 343 0.242 10 0.0310 4 0.0848 5 0.0202
5 50175 685 2.11 10 0.114 4 0.366 5 0.0751
6 200703 1368 16.8 10 0.425 4 1.54 5 0.330

The following table shows the number of iterations and solve time for different smoother methods (as listed in the “test_smoothers” parameter up to “max_smoothing_reps” repetitions) for the cg+fmg solver on the finest level.

smoother reps fmg iter fmg time cg+fmg iter cg+fmg time
jacobi@0.66 1 18 0.842 7 0.346
jacobi@0.66 2 10 0.709 6 0.519
fsor+bsor 1 10 0.450 5 0.342
fsor+bsor 2 5 0.590 4 0.694
ssor 1 7 0.461 4 0.375
ssor 2 4 0.816 3 0.935

The following table shows the number of iterations and solve time for different solvers for the coarse level (as listed in the “test_coarse_solver” parameter) for the cg+fmg solver on the finest level.

coarse solver fmg iter fmg time cg+fmg iter cg+fmg time
petsc 10 0.428 5 0.323
mumps 10 0.434 5 0.319
cg 10 0.433 5 0.318

The following table shows the number of iterations and solve time for different cycle patterns (as listed in the “test_cycle_patterns” parameter) for the cg+fmg solver on the finest level.

pattern fmg iter fmg time cg+fmg iter cg+fmg time
V 10 0.427 5 0.316
2/V 5 0.422 4 0.496
1/2/V 8 0.430 5 0.382
1/1/2/V 9 0.405 5 0.332

The following table shows the number of iterations and solve time for different cycle depths (the number of coarsenings relative to the finest level) for the cg+fmg solver on the finest level.

depth pattern coarse dofs fmg iter fmg time cg+fmg iter cg+fmg time
0 1 200703 1 9.31 1 9.31
1 1/1 50175 8 1.51 4 1.42
2 1/1/1 12543 9 0.509 4 0.400
3 1/1/1/1 3135 9 0.396 5 0.333
4 1/1/1/1/1 783 10 0.422 5 0.318
5 1/1/1/1/1/1 195 10 0.422 5 0.316
6 1/1/1/1/1/1/1 48 10 0.422 5 0.316

The following table shows the execution time of different operations for 7 cycles (0.28s) of the multigrid algorithm on the finest level.

operation reps min time max time total time/rep time/rep % time/iter time/iter %
calc defect 42 6.91e-06 0.00964 0.0705 0.00168 25.2 0.0101 25.3
pre smooth 42 8.82e-06 0.00957 0.0718 0.00171 25.6 0.0103 25.8
post smooth 42 7.87e-06 0.0102 0.0777 0.00185 27.7 0.0111 27.9
restrict 42 6.91e-06 0.00313 0.0246 0.000586 8.79 0.00351 8.85
prolongate 42 5.01e-06 0.00329 0.0241 0.000573 8.60 0.00344 8.66
correct error 42 2.86e-06 0.00141 0.00909 0.000217 3.25 0.00130 3.27
coarse solve 7 4.48e-05 5.29e-05 0.000332 4.74e-05 0.711 4.74e-05 0.119
test break 7 5.96e-06 7.15e-06 4.48e-05 6.40e-06 0.0961 6.40e-06 0.0161

For more information refer to the Tests class reference or have a look at the demo directory.