This page was generated from unit-5a.1-mpi/poisson_mpi.ipynb.

5.1 Poisson Equation in Parallel

NGSolve can be executed on a cluster using the MPI message passing interface. You can download poisson_mpi.py and run it as

mpirun -np 4 python3 poisson_mpi.py

This computes and saves a solution which we can visualize it with drawsolution.py

netgen drawsolution.py

For proper parallel execution, Netgen/NGSolve must be configured with -DUSE_MPI=ON. Recent binaries for Linux and Mac are built with MPI support (?). If you are unsure, when your Netgen/NGSolve supports MPI, look for outputs like “Including MPI version 3.1” during Netgen startup.

MPI-parallel execution using ipyparallel

For the MPI jupyter-tutorials we use ipyparallel module. Please consult https://ipyparallel.readthedocs.io for installation.

On my notebook I have created the profile once by the command

ipython profile create –parallel –profile=mpi

Since I have only two physical cores, I edited the file .ipython/profile_mpi/ipcluster_config.py to allow for oversubscription:

c.MPILauncher.mpi_args = [‘–oversubscribe’]

I start the cluster via

ipcluster start –engines=MPI -n 4 –profile=mpi

In jupyter, we can then connect to the cluster via a Client object of ipyparallel.

[1]:
from ipyparallel import Client
c = Client(profile='mpi')
c.ids
[1]:
[0, 1, 2, 3]

We use mpi4py https://mpi4py.readthedocs.io/ for issuing MPI calls from Python. The %%px syntax magic causes parallel execution of that cell:

[2]:
%%px
from mpi4py import MPI
comm = MPI.COMM_WORLD
print (comm.rank, comm.size)
[stdout:0] 0 4
[stdout:1] 1 4
[stdout:2] 2 4
[stdout:3] 3 4

The master process (rank==0) generates the mesh, and distributes it within the group of processes defined by the communicator. All other ranks receive a part of the mesh. The function mesh.GetNE(VOL) returns the local number of elements:

[3]:
%%px
from ngsolve import *
from netgen.geom2d import unit_square

if comm.rank == 0:
    mesh = Mesh(unit_square.GenerateMesh(maxh=0.1).Distribute(comm))
else:
    mesh = Mesh(netgen.meshing.Mesh.Receive(comm))
print (mesh.GetNE(VOL))
[stdout:0] 0
[stdout:1] 76
[stdout:2] 81
[stdout:3] 75

We can define spaces, bilinear / linear forms, and gridfunctions in the same way as in sequential mode. But now, the degrees of freedom are distributed on the cluster following the distribution of the mesh. The finite element spaces define how the dofs match together.

[4]:
%%px
fes = H1(mesh, order=3, dirichlet=".*")
u,v = fes.TnT()

a = BilinearForm(grad(u)*grad(v)*dx)
pre = Preconditioner(a, "local")
a.Assemble()

f = LinearForm(1*v*dx).Assemble()
gfu = GridFunction(fes)

from ngsolve.krylovspace import CGSolver
inv = CGSolver(a.mat, pre.mat, printing=comm.rank==0, maxsteps=200, tol=1e-8)
gfu.vec.data = inv*f.vec
[stdout:0]
WARNING: printing is deprecated, use printrates instead!
WARNING: maxsteps is deprecated, use maxiter instead!
CG iteration 1, residual = 0.0531207003733573
CG iteration 2, residual = 0.07374408660622353
CG iteration 3, residual = 0.06613575382258267
CG iteration 4, residual = 0.050564282985203636
CG iteration 5, residual = 0.048118728884873806
CG iteration 6, residual = 0.03171133576546557
CG iteration 7, residual = 0.026339954640005568
CG iteration 8, residual = 0.018867600495310295
CG iteration 9, residual = 0.013856359824532939
CG iteration 10, residual = 0.011508956583247314
CG iteration 11, residual = 0.008869710488825622
CG iteration 12, residual = 0.004854287485627616
CG iteration 13, residual = 0.0022162679463868917
CG iteration 14, residual = 0.0015615896552425455
CG iteration 15, residual = 0.0011065170423390227
CG iteration 16, residual = 0.0007171571458209737
CG iteration 17, residual = 0.000525113124024312
CG iteration 18, residual = 0.0003498383890375573
CG iteration 19, residual = 0.00022463539798289132
CG iteration 20, residual = 0.00015923181156189808
CG iteration 21, residual = 0.00012412243193313228
CG iteration 22, residual = 8.768140471148008e-05
CG iteration 23, residual = 5.820850555540184e-05
CG iteration 24, residual = 4.05456782935216e-05
CG iteration 25, residual = 3.0221249909168776e-05
CG iteration 26, residual = 2.2281511256913116e-05
CG iteration 27, residual = 1.5466555984921416e-05
CG iteration 28, residual = 1.0276039906232303e-05
CG iteration 29, residual = 6.856538971345396e-06
CG iteration 30, residual = 4.775824481043933e-06
CG iteration 31, residual = 3.1094353975491994e-06
CG iteration 32, residual = 2.3089701477002727e-06
CG iteration 33, residual = 1.6644528884647446e-06
CG iteration 34, residual = 1.1420780176193516e-06
CG iteration 35, residual = 7.6677095140601e-07
CG iteration 36, residual = 4.84508000920816e-07
CG iteration 37, residual = 3.5489460115406744e-07
CG iteration 38, residual = 2.4269630502786904e-07
CG iteration 39, residual = 1.6023817896666353e-07
CG iteration 40, residual = 1.0291374085373982e-07
CG iteration 41, residual = 6.809889912101117e-08
CG iteration 42, residual = 4.6562019125853935e-08
CG iteration 43, residual = 3.428059952977809e-08
CG iteration 44, residual = 2.3419466928291246e-08
CG iteration 45, residual = 1.4564354581586e-08
CG iteration 46, residual = 8.733217166723476e-09
CG iteration 47, residual = 5.684444002977543e-09
CG iteration 48, residual = 4.136624311978599e-09
CG iteration 49, residual = 2.97950767293741e-09
CG iteration 50, residual = 2.0124741680358736e-09
CG iteration 51, residual = 1.223790556755452e-09
CG iteration 52, residual = 7.854711298623957e-10
CG iteration 53, residual = 6.170296483786957e-10
CG iteration 54, residual = 4.909037903452311e-10
[stdout:1] WARNING: maxsteps is deprecated, use maxiter instead!
[stdout:2] WARNING: maxsteps is deprecated, use maxiter instead!
[stdout:3] WARNING: maxsteps is deprecated, use maxiter instead!

Parallel pickling allows to serialize the distributed solution and transfer it to the client. The process with rank=0 gets the whole mesh and computed solution, all other processes get the local parts of the mesh and solution:

[5]:
gfu = c[:]["gfu"]

We can now draw the whole solution using the the master process’s gfu[0].

[6]:
from ngsolve.webgui import Draw
Draw (gfu[0])
[6]:
BaseWebGuiScene

Drawing gfu[n] will draw only part of the solution that the process with rank=n possesses.

[7]:
Draw (gfu[3])
[7]:
BaseWebGuiScene

We can also visualize the sub-domains obtained by the automatic partitioning, without using any computed solution, as follows.

[8]:
%%px
fesL2 = L2(mesh, order=0)
gfL2 = GridFunction(fesL2)
gfL2.vec.local_vec[:] = comm.rank
[9]:
gfL2 = c[:]["gfL2"]
Draw (gfL2[0])
[9]:
BaseWebGuiScene