This page was generated from unit-5a.3-petsc/petsc.ipynb.

5.3.1 Using PETScΒΆ

We learn how we can interface the popular parallel toolkit PETSc for solving linear equations and much more. We use the Python interface petsc4py. The whole Python file is n2p_ex1.py

[1]:
from ipyparallel import Client
c = Client(profile='mpi')
c.ids
[1]:
[0, 1, 2, 3]
[2]:
%%px
from ngsolve import *
from netgen.geom2d import unit_square
comm = MPI.COMM_WORLD
if comm.rank == 0:
    ngmesh = unit_square.GenerateMesh(maxh=0.1).Distribute(comm)
else:
    ngmesh = netgen.meshing.Mesh.Receive(comm)

for l in range(2):
    ngmesh.Refine()
mesh = Mesh(ngmesh)
[3]:
%%px
import numpy as np
import petsc4py.PETSc as psc
[4]:
%%px
fes = H1(mesh, order=1)
u,v = fes.TnT()
a = BilinearForm(grad(u)*grad(v)*dx+u*v*ds).Assemble()
f = LinearForm(x*v*dx).Assemble()
gfu = GridFunction(fes)

We convert the local sparse matrix to a local PETSc AIJ matrix:

[5]:
%%px
locmat = a.mat.local_mat
val,col,ind = locmat.CSR()
ind = np.array(ind, dtype='int32')
apsc_loc = psc.Mat().createAIJ(size=(locmat.height, locmat.width), csr=(ind,col,val), comm=MPI.COMM_SELF)

The NGSolve ParallelDofs object corresponds to a PETSc IndexSet object, which we create next. In PETSc, dofs are globally enumerated, what is not the case in NGSolve. For this purpose, a ParallelDofs class can generate a globally consistent enumeration of dofs. The generated globnums array contains the global dof-numbers of the local dofs.

[6]:
%%px
pardofs = fes.ParallelDofs()
globnums, nglob = pardofs.EnumerateGlobally()
# print (list(globnums))

iset = psc.IS().createGeneral (indices=globnums, comm=comm)
lgmap = psc.LGMap().createIS(iset)

We can now create the global matrix using our local2global map:

[7]:
%%px
mat = psc.Mat().createPython(size=nglob, comm=comm)
mat.setType(psc.Mat.Type.IS)
mat.setLGMap(lgmap)
mat.setISLocalMat(apsc_loc)
mat.assemble()
mat.convert("mpiaij")
[stderr:0] [1f72a6a4480f:19599] Read -1, expected 50640, errno = 1

[stderr:1] [1f72a6a4480f:19600] Read -1, expected 32064, errno = 1

[stderr:2] [1f72a6a4480f:19601] Read -1, expected 13008, errno = 1

Out[0:26]: <petsc4py.PETSc.Mat at 0x7f8318e5f950>
Out[1:26]: <petsc4py.PETSc.Mat at 0x7fb51a0cd180>
Out[2:26]: <petsc4py.PETSc.Mat at 0x7fca1a08de00>
Out[3:26]: <petsc4py.PETSc.Mat at 0x7f45b9fda7c0>
[8]:
%%px

f.vec.Cumulate()

v1, v2 = mat.createVecs()

v2loc = v2.getSubVector(iset)
v2loc.getArray()[:] = f.vec.FV()
v2.restoreSubVector(iset, v2loc)
[9]:
%%px
ksp = psc.KSP()
ksp.create()
ksp.setOperators(mat)
ksp.setType(psc.KSP.Type.CG)
ksp.setNormType(psc.KSP.NormType.NORM_NATURAL)
ksp.getPC().setType("gamg")
ksp.setTolerances(rtol=1e-6, atol=0, divtol=1e16, max_it=400)
ksp.solve(v2,v1)

print ("petsc-its =", ksp.its)
[stdout:1] petsc-its = 7

[stdout:3] petsc-its = 7

[stdout:0] petsc-its = 7

[stdout:2] petsc-its = 7

[stderr:1] [1f72a6a4480f:19600] Read -1, expected 11536, errno = 1
[1f72a6a4480f:19600] Read -1, expected 15788, errno = 1
[1f72a6a4480f:19600] Read -1, expected 46144, errno = 1
[1f72a6a4480f:19600] Read -1, expected 63152, errno = 1
[1f72a6a4480f:19600] Read -1, expected 4096, errno = 1
[1f72a6a4480f:19600] Read -1, expected 6112, errno = 1
[1f72a6a4480f:19600] Read -1, expected 6968, errno = 1
[1f72a6a4480f:19600] Read -1, expected 6880, errno = 1

[stderr:3] [1f72a6a4480f:19602] Read -1, expected 16972, errno = 1
[1f72a6a4480f:19602] Read -1, expected 67888, errno = 1
[1f72a6a4480f:19602] Read -1, expected 7696, errno = 1

[stderr:2] [1f72a6a4480f:19601] Read -1, expected 8580, errno = 1
[1f72a6a4480f:19601] Read -1, expected 19656, errno = 1
[1f72a6a4480f:19601] Read -1, expected 34320, errno = 1
[1f72a6a4480f:19601] Read -1, expected 78624, errno = 1
[1f72a6a4480f:19601] Read -1, expected 7400, errno = 1
[1f72a6a4480f:19601] Read -1, expected 5392, errno = 1
[1f72a6a4480f:19601] Read -1, expected 9360, errno = 1

[stderr:0] [1f72a6a4480f:19599] Read -1, expected 15344, errno = 1
[1f72a6a4480f:19599] Read -1, expected 61376, errno = 1
[1f72a6a4480f:19599] Read -1, expected 4736, errno = 1

[10]:
%%px
v1loc = v1.getSubVector(iset)
for i in range(len(gfu.vec)):
    gfu.vec.FV()[i] = v1loc.getArray()[i]
[11]:
from ngsolve.webgui import Draw
gfu = c[:]["gfu"]
Draw (gfu[0])
[11]:
BaseWebGuiScene
[ ]:

[ ]: