Iterative Solvers

4. 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 (i.e. not linear in the number of degrees of freedom). 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
Draw(shape);

Generate a mesh, and perform uniform mesh refinement:

ngmesh = OCCGeometry(shape).GenerateMesh(maxh=0.1)
for l in range(1):
    ngmesh.Refine()
mesh = Mesh(ngmesh)
mesh.Curve(3)
Draw (mesh);
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(v*dx)

# c = Preconditioner(a, "direct", inverse="sparsecholesky")
c = Preconditioner(a, "local")
# c = Preconditioner(a, "bddc")

gfu = GridFunction(fes)
we have 175365 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

with TaskManager():
    inv = CGSolver(mat=a.mat, pre=c.mat, printrates='\r', maxiter=400)
    gfu.vec.data = inv * f.vec
CG iteration 1, residual = 0.024168211887440998     
CG iteration 2, residual = 0.03725628120863221     
CG iteration 3, residual = 0.048754089845722055     
CG iteration 4, residual = 0.04862732578235363     
CG iteration 5, residual = 0.03549336921302025     
CG iteration 6, residual = 0.0293744614962741     
CG iteration 7, residual = 0.03066123661238225     
CG iteration 8, residual = 0.02811367911184786     
CG iteration 9, residual = 0.025913075688182834     
CG iteration 10, residual = 0.02454807040274927     
CG iteration 11, residual = 0.02171624394872283     
CG iteration 12, residual = 0.018975427803729687     
CG iteration 13, residual = 0.01653563672653103     
CG iteration 14, residual = 0.014746941865921827     
CG iteration 15, residual = 0.013725806119087544     
CG iteration 16, residual = 0.012193895458597677     
CG iteration 17, residual = 0.010740865456279053     
CG iteration 18, residual = 0.009723513116453798     
CG iteration 19, residual = 0.009066323667772342     
CG iteration 20, residual = 0.008468679817862654     
CG iteration 21, residual = 0.007732868947633278     
CG iteration 22, residual = 0.006829614491921409     
CG iteration 23, residual = 0.00587887620881805     
CG iteration 24, residual = 0.004987133044414072     
CG iteration 25, residual = 0.004178862961681182     
CG iteration 26, residual = 0.003483355314083823     
CG iteration 27, residual = 0.002936043467468231     
CG iteration 28, residual = 0.002502622295264866     
CG iteration 29, residual = 0.0021086420065424176     
CG iteration 30, residual = 0.0017546757540817886     
CG iteration 31, residual = 0.0014262183013239137     
CG iteration 32, residual = 0.0011644005319283975     
CG iteration 33, residual = 0.0009572529056334913     
CG iteration 34, residual = 0.0008013011874902202     
CG iteration 35, residual = 0.000680594270109092     
CG iteration 36, residual = 0.0005984722063953777     
CG iteration 37, residual = 0.0005311901198788222     
CG iteration 38, residual = 0.00046729461555858397     
CG iteration 39, residual = 0.000400635740190301     
CG iteration 40, residual = 0.00033647246206948563     
CG iteration 41, residual = 0.0002799054780565159     
CG iteration 42, residual = 0.00023369463891243557     
CG iteration 43, residual = 0.00019687435079165068     
CG iteration 44, residual = 0.00016898180566359445     
CG iteration 45, residual = 0.00015177909548398172     
CG iteration 46, residual = 0.00013720621413633387     
CG iteration 47, residual = 0.00012634162391539467     
CG iteration 48, residual = 0.00011682056475795088     
CG iteration 49, residual = 0.00010452962317140443     
CG iteration 50, residual = 9.020361275710905e-05     
CG iteration 51, residual = 7.665929288095819e-05     
CG iteration 52, residual = 6.497579453139789e-05     
CG iteration 53, residual = 5.481480841076131e-05     
CG iteration 54, residual = 4.706511614279396e-05     
CG iteration 55, residual = 4.144599292955043e-05     
CG iteration 56, residual = 3.6842165137481845e-05     
CG iteration 57, residual = 3.346920397065269e-05     
CG iteration 58, residual = 3.043239386675799e-05     
CG iteration 59, residual = 2.768818057802067e-05     
CG iteration 60, residual = 2.4900919295111096e-05     
CG iteration 61, residual = 2.1970492260737667e-05     
CG iteration 62, residual = 1.9138742345568966e-05     
CG iteration 63, residual = 1.64010552054693e-05     
CG iteration 64, residual = 1.3998977378845031e-05     
CG iteration 65, residual = 1.1850736605355039e-05     
CG iteration 66, residual = 1.0219370689167833e-05     
CG iteration 67, residual = 8.854120691480248e-06     
CG iteration 68, residual = 7.829275296913381e-06     
CG iteration 69, residual = 7.181276033472496e-06     
CG iteration 70, residual = 6.711152690801696e-06     
CG iteration 71, residual = 6.240562920336067e-06     
CG iteration 72, residual = 5.763534695280271e-06     
CG iteration 73, residual = 5.2165678185350175e-06     
CG iteration 74, residual = 4.5975892448129214e-06     
CG iteration 75, residual = 3.935278363510332e-06     
CG iteration 76, residual = 3.3293088457885706e-06     
CG iteration 77, residual = 2.799684728168634e-06     
CG iteration 78, residual = 2.3693553893873863e-06     
CG iteration 79, residual = 2.013462341400192e-06     
CG iteration 80, residual = 1.7682280439102759e-06     
CG iteration 81, residual = 1.5845799500684987e-06     
CG iteration 82, residual = 1.4491035588078638e-06     
CG iteration 83, residual = 1.3397657037472735e-06     
CG iteration 84, residual = 1.2271050700511731e-06     
CG iteration 85, residual = 1.0965476698712107e-06     
CG iteration 86, residual = 9.57225776774287e-07     
CG iteration 87, residual = 8.13963601621464e-07     
CG iteration 88, residual = 6.815338388700204e-07     
CG iteration 89, residual = 5.693113172908058e-07     
CG iteration 90, residual = 4.801436835407917e-07     
CG iteration 91, residual = 4.056509545021183e-07     
CG iteration 92, residual = 3.4593682569318636e-07     
CG iteration 93, residual = 2.982808260429865e-07     
CG iteration 94, residual = 2.6051732730541207e-07     
CG iteration 95, residual = 2.279699164670516e-07     
CG iteration 96, residual = 1.994991762248402e-07     
CG iteration 97, residual = 1.774442730620756e-07     
CG iteration 98, residual = 1.5913631853265494e-07     
CG iteration 99, residual = 1.422784800709499e-07     
CG iteration 100, residual = 1.2604405372872597e-07     
CG iteration 101, residual = 1.0993917465018588e-07     
CG iteration 102, residual = 9.442894034242134e-08     
CG iteration 103, residual = 7.953848210527243e-08     
CG iteration 104, residual = 6.582625885819961e-08     
CG iteration 105, residual = 5.44978876263058e-08     
CG iteration 106, residual = 4.509041234169682e-08     
CG iteration 107, residual = 3.785782839322218e-08     
CG iteration 108, residual = 3.258595439192123e-08     
CG iteration 109, residual = 2.8444819524513552e-08     
CG iteration 110, residual = 2.5638729444077072e-08     
CG iteration 111, residual = 2.3297400625283055e-08     
CG iteration 112, residual = 2.1350965683831826e-08     
CG iteration 113, residual = 1.940468112338108e-08     
CG iteration 114, residual = 1.7323281099922537e-08     
CG iteration 115, residual = 1.511824448747883e-08     
CG iteration 116, residual = 1.2915124378482981e-08     
CG iteration 117, residual = 1.0843763822364343e-08     
CG iteration 118, residual = 9.055973599084382e-09     
CG iteration 119, residual = 7.634394152668532e-09     
CG iteration 120, residual = 6.532960865961823e-09     
CG iteration 121, residual = 5.753045397720498e-09     
CG iteration 122, residual = 5.1452475895349644e-09     
CG iteration 123, residual = 4.70736007926196e-09     
CG iteration 124, residual = 4.333092435360777e-09     
CG iteration 125, residual = 3.938578616912369e-09     
CG iteration 126, residual = 3.5221050030916496e-09     
CG iteration 127, residual = 3.077844439851731e-09     
CG iteration 128, residual = 2.659388308051706e-09     
CG iteration 129, residual = 2.2988837704138722e-09     
CG iteration 130, residual = 1.9978631973803915e-09     
CG iteration 131, residual = 1.74275618741554e-09     
CG iteration 132, residual = 1.5249214181422333e-09     
CG iteration 133, residual = 1.3454679355558862e-09     
CG iteration 134, residual = 1.1957957546822729e-09     
CG iteration 135, residual = 1.0641779073964508e-09     
CG iteration 136, residual = 9.481222831612558e-10     
CG iteration 137, residual = 8.400809064759393e-10     
CG iteration 138, residual = 7.399373308334826e-10     
CG iteration 139, residual = 6.544151748904927e-10     
CG iteration 140, residual = 5.797132142445199e-10     
CG iteration 141, residual = 5.073895910509104e-10     
CG iteration 142, residual = 4.423551213484354e-10     
CG iteration 143, residual = 3.862626676178589e-10     
CG iteration 144, residual = 3.3606892773110986e-10     
CG iteration 145, residual = 2.9582327537857624e-10     
CG iteration 146, residual = 2.607466151928e-10     
CG iteration 147, residual = 2.2980930688887227e-10     
CG iteration 148, residual = 2.0268441266572106e-10     
CG iteration 149, residual = 1.784665890052008e-10     
CG iteration 150, residual = 1.5514051699692682e-10     
CG iteration 151, residual = 1.331662951719842e-10     
CG iteration 152, residual = 1.1473767751056046e-10     
CG iteration 153, residual = 9.877154338175284e-11     
CG iteration 154, residual = 8.516671442082049e-11     
CG iteration 155, residual = 7.437252562100745e-11     
CG iteration 156, residual = 6.538896555655376e-11     
CG iteration 157, residual = 5.812288887153556e-11     
CG iteration 158, residual = 5.1381917363748735e-11     
CG iteration 159, residual = 4.5260558090034034e-11     
CG iteration 160, residual = 3.96724567493656e-11     
CG iteration 161, residual = 3.4487181561734893e-11     
CG iteration 162, residual = 2.9760952117297225e-11     
CG iteration 163, residual = 2.5797982890180836e-11     
CG iteration 164, residual = 2.223330678659775e-11     
CG iteration 165, residual = 1.898605770771651e-11     
CG iteration 166, residual = 1.612199384910115e-11     
CG iteration 167, residual = 1.3682658799695273e-11     
CG iteration 168, residual = 1.1739941473261357e-11     
CG iteration 169, residual = 1.021954676866916e-11     
CG iteration 170, residual = 9.055018540591983e-12     
CG iteration 171, residual = 8.043111234513032e-12     
CG iteration 172, residual = 7.085989357340513e-12     
CG iteration 173, residual = 6.2785162675460175e-12     
CG iteration 174, residual = 5.526851634260084e-12     
CG iteration 175, residual = 4.871113969694306e-12     
CG iteration 176, residual = 4.198334513242613e-12     
CG iteration 177, residual = 3.5685173939673317e-12     
CG iteration 178, residual = 3.0175637341842797e-12     
CG iteration 179, residual = 2.5332527061819794e-12     
CG iteration 180, residual = 2.15231537542762e-12     
CG iteration 181, residual = 1.8521669262291137e-12     
CG iteration 182, residual = 1.5998414700946637e-12     
CG iteration 183, residual = 1.3994657145050026e-12     
CG iteration 184, residual = 1.2250819185450333e-12     
CG iteration 185, residual = 1.0787918759357505e-12     
CG iteration 186, residual = 9.391837969108861e-13     
CG iteration 187, residual = 8.226143669793064e-13     
CG iteration 188, residual = 7.181256780761086e-13     
CG iteration 189, residual = 6.215847353457856e-13     
CG iteration 190, residual = 5.439436707594406e-13     
CG iteration 191, residual = 4.758176488303656e-13     
CG iteration 192, residual = 4.1177727173981446e-13     
CG iteration 193, residual = 3.544282247391617e-13     
CG iteration 194, residual = 3.034094925397793e-13     
CG iteration 195, residual = 2.5757226733858625e-13     
CG iteration 196, residual = 2.1805731384748688e-13     
CG iteration 197, residual = 1.8476667970341165e-13     
CG iteration 198, residual = 1.574147380185446e-13     
CG iteration 199, residual = 1.351666003441794e-13     
CG iteration 200, residual = 1.177961698000199e-13     
CG iteration 201, residual = 1.0323025347102942e-13     
CG iteration 202, residual = 9.050060618750375e-14     
CG iteration 203, residual = 8.101672958689816e-14     
CG iteration 204, residual = 7.194844949941166e-14     
CG iteration 205, residual = 6.335979713501739e-14     
CG iteration 206, residual = 5.5437077231694794e-14     
CG iteration 207, residual = 4.7486041341990957e-14     
CG iteration 208, residual = 4.03682526757403e-14     
CG iteration 209, residual = 3.38496344165522e-14     
CG iteration 210, residual = 2.8584046454354334e-14     
CG iteration 211, residual = 2.4478783614869083e-14     
CG iteration 212, residual = 2.136066574006779e-14     
CG converged in 212 iterations to residual 2.136066574006779e-14
Draw (gfu, draw_vol=False);

Exercise:

  • Experiment with differ problem sizes and preconditioners

  • How big systems can you solve on your computer (watch out for memory usage)