Iterative Solvers

20.5. Iterative Solvers#

So far we have used direct solvers to solve the linear system of equations. Although a direct solver can profit from the sparse matrix, it’s arithmetic complexity is sub-optimal. For large-scale problems iterative solvers are a must.

The conjugate gradient (cg) method is the standard method for symmetric and positive definite matrices. It’s convergence rate depends on a preconditioner, what is a cheap approximative inverse to the matrix.

from ngsolve import *
from ngsolve.webgui import Draw

We generate a 3D geometry and mesh using the OCC constructive solid geometry (CSG) modeler:

from netgen.occ import *
cube = Box((0,0,0),(1,1,1))
cyl = Cylinder((0,0.5,0.5),X, r=0.2, h=1)
cube.faces.name = "outer"
cyl.faces.name = "cyl"
shape = cube-cyl
ea = { "euler_angles" : [-180, -80, 165] }
Draw(shape, **ea);
ngmesh = OCCGeometry(shape).GenerateMesh(maxh=0.1)

mesh = Mesh(ngmesh).Curve(3)

for l in range(0):
    mesh.Refine()

Draw (mesh, **ea);
fes = H1(mesh, order=3, dirichlet="outer", wb_withedges=False)
print ("we have", fes.ndof, "unknowns")

u,v = fes.TnT()

a = BilinearForm(grad(u)*grad(v)*dx)
f = LinearForm(1*v*dx)

# pre = preconditioners.MultiGrid(a)
pre = preconditioners.Local(a)
# pre = preconditioners.BDDC(a)

gfu = GridFunction(fes)
we have 23787 unknowns

assemble system and setup preconditioner in parallel:

with TaskManager():
    a.Assemble()
    f.Assemble()

solve the system using the preconditioned conjugate gradient method:

from ngsolve.krylovspace import CGSolver
from time import time

ts = time()
with TaskManager():
    # inv = a.mat.Inverse(inverse="sparsecholesky", freedofs=fes.FreeDofs())
    inv = CGSolver(mat=a.mat, pre=pre, printrates='\r', maxiter=1000)
    gfu.vec.data = inv * f.vec
te = time()
print ("needed", te-ts, "seconds for", fes.ndof, "unknowns")
CG iteration 1, residual = 0.047739949075105856     
CG iteration 2, residual = 0.05294727938543792     
CG iteration 3, residual = 0.05500960698375211     
CG iteration 4, residual = 0.04606877270694542     
CG iteration 5, residual = 0.029876107300354112     
CG iteration 6, residual = 0.02338270673926799     
CG iteration 7, residual = 0.018047146455209776     
CG iteration 8, residual = 0.012650995831141891     
CG iteration 9, residual = 0.008947539939651404     
CG iteration 10, residual = 0.0071342097470362235     
CG iteration 11, residual = 0.006486603723518196     
CG iteration 12, residual = 0.0064453373377096005     
CG iteration 13, residual = 0.00667923520810839     
CG iteration 14, residual = 0.00693806713138742     
CG iteration 15, residual = 0.006705024574677978     
CG iteration 16, residual = 0.005583967413619985     
CG iteration 17, residual = 0.004207964726985902     
CG iteration 18, residual = 0.003077708743241626     
CG iteration 19, residual = 0.002292553710133818     
CG iteration 20, residual = 0.0017210286974316956     
CG iteration 21, residual = 0.00139188918625965     
CG iteration 22, residual = 0.0012867256293030039     
CG iteration 23, residual = 0.0012207824425805747     
CG iteration 24, residual = 0.0011508422442719902     
CG iteration 25, residual = 0.0010326259360892084     
CG iteration 26, residual = 0.0008650872416098801     
CG iteration 27, residual = 0.000683412255381723     
CG iteration 28, residual = 0.0005203227711120677     
CG iteration 29, residual = 0.0003856893581649341     
CG iteration 30, residual = 0.00030289016916732186     
CG iteration 31, residual = 0.0002468748185900819     
CG iteration 32, residual = 0.0002005777049766653     
CG iteration 33, residual = 0.0001639970190888681     
CG iteration 34, residual = 0.00013686401056622847     
CG iteration 35, residual = 0.00011413361303568408     
CG iteration 36, residual = 9.451369431613458e-05     
CG iteration 37, residual = 8.065214654944001e-05     
CG iteration 38, residual = 7.150818336039979e-05     
CG iteration 39, residual = 6.345362194690667e-05     
CG iteration 40, residual = 5.598287673765576e-05     
CG iteration 41, residual = 4.817075802670546e-05     
CG iteration 42, residual = 3.927699677793061e-05     
CG iteration 43, residual = 3.084502899259499e-05     
CG iteration 44, residual = 2.3441691859581678e-05     
CG iteration 45, residual = 1.806809166896807e-05     
CG iteration 46, residual = 1.4349833606446255e-05     
CG iteration 47, residual = 1.1766739852763842e-05     
CG iteration 48, residual = 9.964987268310969e-06     
CG iteration 49, residual = 9.056719106388718e-06     
CG iteration 50, residual = 8.561745619842749e-06     
CG iteration 51, residual = 8.026554471089201e-06     
CG iteration 52, residual = 7.09753603319757e-06     
CG iteration 53, residual = 5.944747573670306e-06     
CG iteration 54, residual = 4.555910099875375e-06     
CG iteration 55, residual = 3.3819070629162298e-06     
CG iteration 56, residual = 2.4711107258406207e-06     
CG iteration 57, residual = 1.8477611726611034e-06     
CG iteration 58, residual = 1.4857145554662001e-06     
CG iteration 59, residual = 1.2825964453132414e-06     
CG iteration 60, residual = 1.1972381968942671e-06     
CG iteration 61, residual = 1.1728671909923245e-06     
CG iteration 62, residual = 1.126179049197469e-06     
CG iteration 63, residual = 9.894511688700755e-07     
CG iteration 64, residual = 7.958772966268733e-07     
CG iteration 65, residual = 5.901602467726177e-07     
CG iteration 66, residual = 4.3534259787806927e-07     
CG iteration 67, residual = 3.3237064477404487e-07     
CG iteration 68, residual = 2.698792505812133e-07     
CG iteration 69, residual = 2.3028612227108508e-07     
CG iteration 70, residual = 2.083853711175276e-07     
CG iteration 71, residual = 1.9228490301377265e-07     
CG iteration 72, residual = 1.7336077681162059e-07     
CG iteration 73, residual = 1.4925977570562424e-07     
CG iteration 74, residual = 1.2342270369943355e-07     
CG iteration 75, residual = 1.0181258758281511e-07     
CG iteration 76, residual = 8.450738566286675e-08     
CG iteration 77, residual = 7.132426095176737e-08     
CG iteration 78, residual = 6.050788355810563e-08     
CG iteration 79, residual = 4.996349766921651e-08     
CG iteration 80, residual = 4.009749246323375e-08     
CG iteration 81, residual = 3.144267392739601e-08     
CG iteration 82, residual = 2.5066840926766224e-08     
CG iteration 83, residual = 2.0327720064246213e-08     
CG iteration 84, residual = 1.7509128806904694e-08     
CG iteration 85, residual = 1.608237196625368e-08     
CG iteration 86, residual = 1.507860069173036e-08     
CG iteration 87, residual = 1.372484911187188e-08     
CG iteration 88, residual = 1.2121807998454798e-08     
CG iteration 89, residual = 1.00888339100624e-08     
CG iteration 90, residual = 7.890534948262154e-09     
CG iteration 91, residual = 5.87579891950059e-09     
CG iteration 92, residual = 4.400928225527029e-09     
CG iteration 93, residual = 3.428025435380982e-09     
CG iteration 94, residual = 2.907898412691934e-09     
CG iteration 95, residual = 2.6245354827297815e-09     
CG iteration 96, residual = 2.4637850210760234e-09     
CG iteration 97, residual = 2.2859495940564325e-09     
CG iteration 98, residual = 2.131693730347244e-09     
CG iteration 99, residual = 1.8863068454439993e-09     
CG iteration 100, residual = 1.5392597513245563e-09     
CG iteration 101, residual = 1.1752103715671209e-09     
CG iteration 102, residual = 8.699962262971321e-10     
CG iteration 103, residual = 6.432581480113498e-10     
CG iteration 104, residual = 4.872851375917496e-10     
CG iteration 105, residual = 3.907432860671821e-10     
CG iteration 106, residual = 3.394412281261527e-10     
CG iteration 107, residual = 3.1514292983610024e-10     
CG iteration 108, residual = 3.0195942984448364e-10     
CG iteration 109, residual = 2.8036642753790477e-10     
CG iteration 110, residual = 2.438539949321038e-10     
CG iteration 111, residual = 1.9483936951231735e-10     
CG iteration 112, residual = 1.4896146572263817e-10     
CG iteration 113, residual = 1.1169992483961141e-10     
CG iteration 114, residual = 8.262811653519696e-11     
CG iteration 115, residual = 6.514626468032496e-11     
CG iteration 116, residual = 5.35320933937949e-11     
CG iteration 117, residual = 4.736050413544766e-11     
CG iteration 118, residual = 4.348247742525626e-11     
CG iteration 119, residual = 3.98365173150979e-11     
CG iteration 120, residual = 3.6271237423315454e-11     
CG iteration 121, residual = 3.2001632650508726e-11     
CG iteration 122, residual = 2.684155206797013e-11     
CG iteration 123, residual = 2.1632016046157746e-11     
CG iteration 124, residual = 1.6963671373212918e-11     
CG iteration 125, residual = 1.3089006529540631e-11     
CG iteration 126, residual = 1.0272985452754067e-11     
CG iteration 127, residual = 8.389378504760537e-12     
CG iteration 128, residual = 6.730311501819364e-12     
CG iteration 129, residual = 5.718089083371147e-12     
CG iteration 130, residual = 5.215433426170504e-12     
CG iteration 131, residual = 4.996278179654494e-12     
CG iteration 132, residual = 4.68654879210119e-12     
CG iteration 133, residual = 4.188276580230564e-12     
CG iteration 134, residual = 3.4449302708412135e-12     
CG iteration 135, residual = 2.6179376365126524e-12     
CG iteration 136, residual = 1.9001535179908234e-12     
CG iteration 137, residual = 1.367394511656965e-12     
CG iteration 138, residual = 9.948833458483642e-13     
CG iteration 139, residual = 7.754672454287577e-13     
CG iteration 140, residual = 6.747101554413715e-13     
CG iteration 141, residual = 6.420861407288143e-13     
CG iteration 142, residual = 6.333470268391048e-13     
CG iteration 143, residual = 5.859147516185604e-13     
CG iteration 144, residual = 5.100310500151748e-13     
CG iteration 145, residual = 4.076744753954105e-13     
CG iteration 146, residual = 3.043110403520709e-13     
CG iteration 147, residual = 2.1878948873045242e-13     
CG iteration 148, residual = 1.6124246728221954e-13     
CG iteration 149, residual = 1.2393746680490095e-13     
CG iteration 150, residual = 1.0400757617816995e-13     
CG iteration 151, residual = 9.535998649686239e-14     
CG iteration 152, residual = 9.072494303508726e-14     
CG iteration 153, residual = 8.75599210331743e-14     
CG iteration 154, residual = 8.170994888761618e-14     
CG iteration 155, residual = 7.021016375172058e-14     
CG iteration 156, residual = 5.5940620103756005e-14     
CG iteration 157, residual = 4.200651678096323e-14     
CG converged in 157 iterations to residual 4.200651678096323e-14
needed 0.060534000396728516 seconds for 23787 unknowns
Draw (gfu, **ea);