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

5.3.2 NGSolve - PETSc interface

We use the ngs2petsc interface to map vectors and matrices between NGSolve and PETSc

[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)

The Python-module ngsolve.ngs2petsc provides functionality to transfer vectors and matrices between NGSolve and Python.

Make sure that the ipyparallel server can import the module, e.g. by starting the cluster in the current directory.

[3]:
%%px
import ngsolve.ngs2petsc as n2p
import petsc4py.PETSc as psc
[4]:
%%px
fes = H1(mesh, order=1, dirichlet="left|bottom")
u,v = fes.TnT()
a = BilinearForm(grad(u)*grad(v)*dx+u*v*ds).Assemble()
f = LinearForm(x*v*dx).Assemble()
gfu = GridFunction(fes)

The function CreatePETScMatrix takes an NGSolve matrix, and creates a PETSc matrix from it. A VectorMapping object can map vectors between NGSolve and PETSc.

[5]:
%%px
psc_mat = n2p.CreatePETScMatrix(a.mat, fes.FreeDofs())
vecmap = n2p.VectorMapping (a.mat.row_pardofs, fes.FreeDofs())
[stderr:1] [1f72a6a4480f:19600] Read -1, expected 33408, errno = 1

[stderr:0] [1f72a6a4480f:19599] Read -1, expected 48512, errno = 1

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

Create PETSc-vectors fitting to the matrix

[6]:
%%px
psc_f, psc_u = psc_mat.createVecs()
pass

setting up the parallel Krylov-space solver ….

[7]:
%%px
ksp = psc.KSP()
ksp.create()
ksp.setOperators(psc_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)

moving vectors between NGSolve and PETSc, and solve:

[8]:
%%px
vecmap.N2P(f.vec, psc_f)
ksp.solve(psc_f, psc_u)
vecmap.P2N(psc_u, gfu.vec)
pass
[stderr:3] [1f72a6a4480f:19602] Read -1, expected 16964, errno = 1
[1f72a6a4480f:19602] Read -1, expected 67856, errno = 1
[1f72a6a4480f:19602] Read -1, expected 6992, errno = 1

[stderr:0] [1f72a6a4480f:19599] Read -1, expected 12552, errno = 1
[1f72a6a4480f:19599] Read -1, expected 50208, errno = 1

[stderr:2] [1f72a6a4480f:19601] Read -1, expected 8204, errno = 1
[1f72a6a4480f:19601] Read -1, expected 19168, errno = 1
[1f72a6a4480f:19601] Read -1, expected 32816, errno = 1
[1f72a6a4480f:19601] Read -1, expected 76672, errno = 1
[1f72a6a4480f:19601] Read -1, expected 7144, errno = 1
[1f72a6a4480f:19601] Read -1, expected 5720, errno = 1
[1f72a6a4480f:19601] Read -1, expected 8832, errno = 1

[stderr:1] [1f72a6a4480f:19600] Read -1, expected 9540, errno = 1
[1f72a6a4480f:19600] Read -1, expected 16204, errno = 1
[1f72a6a4480f:19600] Read -1, expected 38160, errno = 1
[1f72a6a4480f:19600] Read -1, expected 64816, errno = 1
[1f72a6a4480f:19600] Read -1, expected 5512, errno = 1
[1f72a6a4480f:19600] Read -1, expected 6704, errno = 1
[1f72a6a4480f:19600] Read -1, expected 6512, errno = 1

[9]:
gfu = c[:]["gfu"]
from ngsolve.webgui import Draw
[10]:
Draw (gfu[0])
[10]:
BaseWebGuiScene

PETSc preconditioner for NGSolve

Next we create a PETSc preconditioner, and wrap it into an NGSolve preconditioner:

[11]:
%%px
a = BilinearForm(grad(u)*grad(v)*dx+u*v*ds)
pre = Preconditioner(a, "gamg")
a.Assemble()
pass
[stderr:1] [1f72a6a4480f:19600] Read -1, expected 33408, errno = 1
[1f72a6a4480f:19600] Read -1, expected 9540, errno = 1
[1f72a6a4480f:19600] Read -1, expected 16204, errno = 1
[1f72a6a4480f:19600] Read -1, expected 38160, errno = 1
[1f72a6a4480f:19600] Read -1, expected 64816, errno = 1
[1f72a6a4480f:19600] Read -1, expected 5512, errno = 1
[1f72a6a4480f:19600] Read -1, expected 6704, errno = 1
[1f72a6a4480f:19600] Read -1, expected 6512, errno = 1

[stderr:3] [1f72a6a4480f:19602] Read -1, expected 16964, errno = 1
[1f72a6a4480f:19602] Read -1, expected 67856, errno = 1
[1f72a6a4480f:19602] Read -1, expected 6992, errno = 1

[stderr:2] [1f72a6a4480f:19601] Read -1, expected 12576, errno = 1
[1f72a6a4480f:19601] Read -1, expected 8204, errno = 1
[1f72a6a4480f:19601] Read -1, expected 19168, errno = 1
[1f72a6a4480f:19601] Read -1, expected 76672, errno = 1
[1f72a6a4480f:19601] Read -1, expected 32816, errno = 1
[1f72a6a4480f:19601] Read -1, expected 7144, errno = 1
[1f72a6a4480f:19601] Read -1, expected 5720, errno = 1
[1f72a6a4480f:19601] Read -1, expected 8832, errno = 1

[stderr:0] [1f72a6a4480f:19599] Read -1, expected 48512, errno = 1
[1f72a6a4480f:19599] Read -1, expected 12552, errno = 1
[1f72a6a4480f:19599] Read -1, expected 50208, errno = 1

and use it in an NGSolve - CGSolver:

[12]:
%%px
from ngsolve.krylovspace import CGSolver
inv = CGSolver(a.mat, pre, printing=comm.rank==0)
gfu.vec.data = inv * f.vec
[stdout:0] WARNING: printing is deprecated, use printrates instead!
CG iteration 1, residual = 0.16186373090106865
CG iteration 2, residual = 0.032256334634758486
CG iteration 3, residual = 0.0041873433365436944
CG iteration 4, residual = 0.0005968449121709221
CG iteration 5, residual = 8.461918988055599e-05
CG iteration 6, residual = 1.0980502894984023e-05
CG iteration 7, residual = 1.6480728786485835e-06
CG iteration 8, residual = 2.0479915322333332e-07
CG iteration 9, residual = 2.6152153836992176e-08
CG iteration 10, residual = 3.383295610844336e-09
CG iteration 11, residual = 4.183121875015019e-10
CG iteration 12, residual = 5.703654028140915e-11
CG iteration 13, residual = 7.797932580843257e-12
CG iteration 14, residual = 9.331733543242658e-13
CG iteration 15, residual = 1.2343363523858606e-13

[13]:
gfu = c[:]["gfu"]
Draw (gfu[0])
[13]:
BaseWebGuiScene
[ ]: