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 Cluster
c = await Cluster(engines="mpi").start_and_connect(n=4, activate=True)
Starting 4 engines with <class 'ipyparallel.cluster.launcher.MPIEngineSetLauncher'>
[2]:
%%px
from ngsolve import *
from netgen.occ import unit_square
comm = MPI.COMM_WORLD

ngmesh = unit_square.GenerateMesh(maxh=0.1, comm=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())

Create PETSc-vectors fitting to the matrix

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

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);
[9]:
gfu = c[:]["gfu"]
from ngsolve.webgui import Draw
[10]:
Draw (gfu[0]);

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, "petsc", pctype="gamg", levels=10)
pre = Preconditioner(a, "gamg")
a.Assemble();

and use it in an NGSolve - CGSolver:

[12]:
%%px
from ngsolve.krylovspace import CGSolver
inv = CGSolver(a.mat, pre, printrates=comm.rank==0)
gfu.vec.data = inv * f.vec
[stdout:0] CG iteration 1, residual = 0.15037762885386297
CG iteration 2, residual = 0.03656019842885419
CG iteration 3, residual = 0.009195248450186991
CG iteration 4, residual = 0.002500983655898158
CG iteration 5, residual = 0.0006236536664564382
CG iteration 6, residual = 0.0001460265968298943
CG iteration 7, residual = 3.094822886465786e-05
CG iteration 8, residual = 6.801434237226316e-06
CG iteration 9, residual = 1.701650763012911e-06
CG iteration 10, residual = 4.2358732809906706e-07
CG iteration 11, residual = 1.0895202048710423e-07
CG iteration 12, residual = 3.2184492103496173e-08
CG iteration 13, residual = 8.340831960378986e-09
CG iteration 14, residual = 2.1161294963907244e-09
CG iteration 15, residual = 5.805454584044729e-10
CG iteration 16, residual = 1.2849699075318396e-10
CG iteration 17, residual = 3.0674821386206156e-11
CG iteration 18, residual = 7.756932758839662e-12
CG iteration 19, residual = 1.8142260127235855e-12
CG iteration 20, residual = 4.164569543469645e-13
CG iteration 21, residual = 9.450208994620436e-14

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

[ ]: