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(4):
#    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)
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)
vecmap = n2p.VectorMapping (a.mat.row_pardofs)
[6]:
%%px
psc_f = vecmap.N2P(f.vec)
psc_u = vecmap.CreatePETScVector()

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)

ksp.solve(psc_f, psc_u)

vecmap.P2N(psc_u, gfu.vec)

netgen.meshing.SetParallelPickling(True)
pass
[7]:
gfu = c[:]["gfu"]
from ngsolve.webgui import Draw
[8]:
Draw (gfu[0])
[8]:

[ ]: