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.
We start an MPI server:
[1]:
from ipyparallel import Cluster
c = await Cluster(engines="mpi").start_and_connect(n=4, activate=True)
Starting 4 engines with <class 'ipyparallel.cluster.launcher.MPIEngineSetLauncher'>
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:2] 2 4
[stdout:1] 1 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:1] 73
[stdout:3] 78
[stdout:0] 0
[stdout:2] 79
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:3] WARNING: maxsteps is deprecated, use maxiter instead!
[stdout:1] WARNING: maxsteps is deprecated, use maxiter instead!
[stdout:2] WARNING: maxsteps is deprecated, use maxiter instead!
[stdout:0] WARNING: printing is deprecated, use printrates instead!
WARNING: maxsteps is deprecated, use maxiter instead!
CG iteration 1, residual = 0.053126041648229295
CG iteration 2, residual = 0.07344287617289662
CG iteration 3, residual = 0.06110786204339784
CG iteration 4, residual = 0.049294272666440334
CG iteration 5, residual = 0.04641237052420435
CG iteration 6, residual = 0.030446095751156126
CG iteration 7, residual = 0.0257639814575267
CG iteration 8, residual = 0.016975392921489894
CG iteration 9, residual = 0.012989664894467531
CG iteration 10, residual = 0.011308228515895663
CG iteration 11, residual = 0.007800282946410079
CG iteration 12, residual = 0.003524698587117302
CG iteration 13, residual = 0.0017401033251837994
CG iteration 14, residual = 0.0011366587144063154
CG iteration 15, residual = 0.0008035602552022153
CG iteration 16, residual = 0.0005243481976976642
CG iteration 17, residual = 0.0003741983982428156
CG iteration 18, residual = 0.0002749727554288431
CG iteration 19, residual = 0.00017229471964983172
CG iteration 20, residual = 0.00012857735127167917
CG iteration 21, residual = 8.997393301231609e-05
CG iteration 22, residual = 6.195695703998971e-05
CG iteration 23, residual = 4.865877020420745e-05
CG iteration 24, residual = 3.387607816378412e-05
CG iteration 25, residual = 2.2199678719709072e-05
CG iteration 26, residual = 1.5023055764552838e-05
CG iteration 27, residual = 1.058100460987856e-05
CG iteration 28, residual = 7.355751857291855e-06
CG iteration 29, residual = 5.352395923296746e-06
CG iteration 30, residual = 3.633926663895672e-06
CG iteration 31, residual = 2.04443377406657e-06
CG iteration 32, residual = 1.575126704942929e-06
CG iteration 33, residual = 1.3159605053736425e-06
CG iteration 34, residual = 8.346360401328236e-07
CG iteration 35, residual = 4.604383800956418e-07
CG iteration 36, residual = 2.9450046172586994e-07
CG iteration 37, residual = 1.9185646830432075e-07
CG iteration 38, residual = 1.3661842726540912e-07
CG iteration 39, residual = 9.49786272702141e-08
CG iteration 40, residual = 6.613160428743896e-08
CG iteration 41, residual = 4.2282308657154504e-08
CG iteration 42, residual = 2.7713291458243708e-08
CG iteration 43, residual = 1.8141562264855626e-08
CG iteration 44, residual = 1.2727903304187594e-08
CG iteration 45, residual = 8.504452240352474e-09
CG iteration 46, residual = 5.56435823744953e-09
CG iteration 47, residual = 3.6807848441792823e-09
CG iteration 48, residual = 2.3751444017405335e-09
CG iteration 49, residual = 1.6993731519085797e-09
CG iteration 50, residual = 1.0657829476235226e-09
CG iteration 51, residual = 7.154229522411755e-10
CG iteration 52, residual = 5.014272480699724e-10
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