This page was generated from unit-2.2-eigenvalues/pinvit.ipynb.
2.2 Programming an eigenvalue solver¶
We solve the generalized eigenvalue problem
where \(A\) comes from \(\int \nabla u \nabla v\), and \(M\) from \(\int u v\), on the space \(H_0^1\).
This tutorial shows how to implement linear algebra algorithms.
[1]:
from ngsolve import *
from ngsolve.webgui import Draw
from netgen.geom2d import unit_square
import math
import scipy.linalg
from scipy import random
mesh = Mesh(unit_square.GenerateMesh(maxh=0.1))
We setup a stiffness matrix \(A\) and a mass matrix \(M\), and declare a preconditioner for \(A\):
[2]:
fes = H1(mesh, order=4, dirichlet=".*")
u = fes.TrialFunction()
v = fes.TestFunction()
a = BilinearForm(fes)
a += grad(u)*grad(v)*dx
pre = Preconditioner(a, "multigrid")
m = BilinearForm(fes)
m += u*v*dx
a.Assemble()
m.Assemble()
u = GridFunction(fes)
The inverse iteration is
where the Rayleigh quotient
converges to the smallest eigenvalue \(\lambda_1\), with rate of convergence \(\lambda_1 / \lambda_2,\) where \(\lambda_2\) is the next smallest eigenvalue.
The preconditioned inverse iteration (PINVIT), see [Knyazef+Neymeyr], replaces \(A^{-1}\) by an approximate inverse \(C^{-1}\):
The optimal step-size \(\alpha\) is found by minimizing the Rayleigh-quotient on a two-dimensional space:
This minimization problem can be solved by a small eigenvalue problem
with matrices
Then, the new iterate is
where \(y\) is the eigenvector corresponding to the smaller eigenvalue.
Implementation in NGSolve¶
First, we create some help vectors. CreateVector
creates new vectors of the same format as the existing vector, i.e., same dimension, same real/ complex type, same entry-size, and same MPI-parallel distribution if any.
[3]:
r = u.vec.CreateVector()
w = u.vec.CreateVector()
Mu = u.vec.CreateVector()
Au = u.vec.CreateVector()
Mw = u.vec.CreateVector()
Aw = u.vec.CreateVector()
Next, we pick a random initial vector, which is zeroed on the Dirichlet boundary.
Below, the FV
method (short for FlatVector
) lets us access the abstract vector’s linear memory chunk, which in turn provides a “numpy” view of the memory. The projector clears the entries at the Dirichlet boundary:
[4]:
r.FV().NumPy()[:] = random.rand(fes.ndof)
u.vec.data = Projector(fes.FreeDofs(), True) * r
Finally, we run the PINVIT algorithm. Note that the small matrices \(a\) and \(m\) defined above are called asmall
and msmall
below. They are of type Matrix
, a class provided by NGSolve for dense matrices.
[5]:
for i in range(20):
Au.data = a.mat * u.vec
Mu.data = m.mat * u.vec
auu = InnerProduct(Au, u.vec)
muu = InnerProduct(Mu, u.vec)
# Rayleigh quotient
lam = auu/muu
print (lam / (math.pi**2))
# residual
r.data = Au - lam * Mu
w.data = pre.mat * r.data
w.data = 1/Norm(w) * w
Aw.data = a.mat * w
Mw.data = m.mat * w
# setup and solve 2x2 small eigenvalue problem
asmall = Matrix(2,2)
asmall[0,0] = auu
asmall[0,1] = asmall[1,0] = InnerProduct(Au, w)
asmall[1,1] = InnerProduct(Aw, w)
msmall = Matrix(2,2)
msmall[0,0] = muu
msmall[0,1] = msmall[1,0] = InnerProduct(Mu, w)
msmall[1,1] = InnerProduct(Mw, w)
# print ("asmall =", asmall, ", msmall = ", msmall)
eval,evec = scipy.linalg.eigh(a=asmall, b=msmall)
# print (eval, evec)
u.vec.data = float(evec[0,0]) * u.vec + float(evec[1,0]) * w
Draw (u)
22.077537355229712
2.0778907736678454
2.0015128968419047
2.000220109768116
2.0000497142129694
2.0000126970085916
2.000003292606943
2.0000008634295487
2.0000002278672877
2.000000060476191
2.000000016131165
2.0000000043407535
2.000000001195152
2.0000000003539222
2.0000000001284226
2.0000000000678706
2.000000000051582
2.0000000000471947
2.0000000000460108
2.0000000000456906
[5]:
BaseWebGuiScene
Simultaneous iteration for several eigenvalues¶
Here are the steps for extending the above to num
vectors.
Declare a GridFunction
with multiple components to store several eigenvectors:
[6]:
num = 5
u = GridFunction(fes, multidim=num)
Create list of help vectors, and a set of random initial vectors in u
, with zero boundary conditions:
[7]:
r = u.vec.CreateVector()
Av = u.vec.CreateVector()
Mv = u.vec.CreateVector()
vecs = []
for i in range(2*num):
vecs.append (u.vec.CreateVector())
for v in u.vecs:
r.FV().NumPy()[:] = random.rand(fes.ndof)
v.data = Projector(fes.FreeDofs(), True) * r
Compute num
residuals, and solve a small eigenvalue problem on a 2 \(\times\) num
dimensional space:
[8]:
asmall = Matrix(2*num, 2*num)
msmall = Matrix(2*num, 2*num)
lams = num * [1]
for i in range(20):
for j in range(num):
vecs[j].data = u.vecs[j]
r.data = a.mat * vecs[j] - lams[j] * m.mat * vecs[j]
vecs[num+j].data = pre.mat * r
for j in range(2*num):
Av.data = a.mat * vecs[j]
Mv.data = m.mat * vecs[j]
for k in range(2*num):
asmall[j,k] = InnerProduct(Av, vecs[k])
msmall[j,k] = InnerProduct(Mv, vecs[k])
ev,evec = scipy.linalg.eigh(a=asmall, b=msmall)
lams[:] = ev[0:num]
print (i, ":", [lam/math.pi**2 for lam in lams])
for j in range(num):
u.vecs[j][:] = 0.0
for k in range(2*num):
u.vecs[j].data += float(evec[k,j]) * vecs[k]
Draw (u)
0 : [12.655132001889585, 66.24401295389443, 82.1934049416161, 91.55643771072484, 98.46462849308917]
1 : [2.0497625955775343, 6.662492696108604, 7.369688733949427, 13.212703015432256, 19.92435301575199]
2 : [2.001476706094831, 5.055715880892912, 5.1688856720555805, 8.956418637931591, 12.025748537883276]
3 : [2.000191522225827, 5.0034580526944925, 5.015807771584299, 8.326658161931196, 10.507101586742221]
4 : [2.0000402017873373, 5.000319808021971, 5.00183247723273, 8.115532941755893, 10.155143232149305]
5 : [2.000009291533831, 5.000047467381901, 5.000237173651957, 8.04474007521218, 10.054457542263615]
6 : [2.000002188767801, 5.000009260819194, 5.000034664666574, 8.018496773737006, 10.020078962123122]
7 : [2.0000005223924817, 5.000002028319123, 5.000005899738493, 8.00785506366804, 10.007506430228615]
8 : [2.000000125633731, 5.000000468930171, 5.0000011534071085, 8.003385756589587, 10.002815166866338]
9 : [2.0000000304201953, 5.000000114200694, 5.000000252615285, 8.001467689841123, 10.001058188789354]
10 : [2.000000007427641, 5.000000031273024, 5.00000006222601, 8.000637995903986, 10.000398489278862]
11 : [2.000000001847122, 5.000000011519675, 5.000000019066091, 8.000277770509777, 10.000150323647752]
12 : [2.0000000004868905, 5.00000000675044, 5.000000008872362, 8.00012106374381, 10.000056825057873]
13 : [2.000000000154037, 5.0000000055793965, 5.000000006414173, 8.000052817046527, 10.000021557129523]
14 : [2.000000000072314, 5.000000005283365, 5.000000005818691, 8.00002307112953, 10.00000824125092]
15 : [2.000000000052182, 5.000000005207637, 5.0000000056745835, 8.000010098765843, 10.00000320997449]
16 : [2.000000000047211, 5.000000005188527, 5.000000005638657, 8.00000443871111, 10.000001307784554]
17 : [2.0000000000459797, 5.000000005183727, 5.000000005630309, 8.000001968190848, 10.000000588273593]
18 : [2.000000000045675, 5.000000005182127, 5.000000005627685, 8.000000889478912, 10.000000316001262]
19 : [2.0000000000455995, 5.000000005182363, 5.000000005627288, 8.000000418341893, 10.000000212935692]
[8]:
BaseWebGuiScene
The multidim-component select in the Visualization dialog box allows to switch between the components of the multidim-GridFunction.
Implementation using MultiVector¶
The simultaneous iteration can be optimized by using MultiVector
s introduced in NGSolve 6.2.2007. These are arrays of vectors of the same format. You can think of a MultiVector with m components of vector size n as an \(n \times m\) matrix.
a
MultiVector
consisting of num vectors of the same format as an existing vector vec is created via MultiVector(vec, num).we can iterate over the components of a MultiVector, and the bracket operator allows to access a subset of vectors
linear operator application is optimized for MultiVector
vector operations are optimized and called as mv * densematrix: \(x = y * mat\) results in x[i] = sum_j y[j] * mat[j,i] (where x and y are ’MultiVector’s, and mat is a dense matrix)
pair-wise inner products of two MultiVectors is available, the result is a dense matrix: InnerProduct(x,y)[i,j] =
mv.Orthogonalize() uses modified Gram-Schmidt to orthogonalize the vectors. Optionally, a matrix defining the inner product can be provided.
with mv.Append(vec) we can add another vector to the array of vectos. A new vector is created, and the values are copied
mv.AppendOrthogonalize(vec) appends a new vector, and orthogonalizes it against the existing vectors, which are assumed to be orthogonal.
[9]:
uvecs = MultiVector(u.vec, num)
vecs = MultiVector(u.vec, 2*num)
for v in vecs[0:num]:
v.SetRandom()
uvecs[:] = pre * vecs[0:num]
lams = Vector(num * [1])
[10]:
for i in range(20):
vecs[0:num] = a.mat * uvecs - (m.mat * uvecs).Scale (lams)
vecs[num:2*num] = pre * vecs[0:num]
vecs[0:num] = uvecs
vecs.Orthogonalize() # m.mat)
asmall = InnerProduct (vecs, a.mat * vecs)
msmall = InnerProduct (vecs, m.mat * vecs)
ev,evec = scipy.linalg.eigh(a=asmall, b=msmall)
lams = Vector(ev[0:num])
print (i, ":", [l/math.pi**2 for l in lams])
uvecs[:] = vecs * Matrix(evec[:,0:num])
0 : [278.7772047126754, 970.7429292345021, 990.5722620191766, 1030.4805698255402, 1071.9919975311216]
1 : [2.0095311199031047, 8.53394457323619, 10.060724007498472, 18.100339847312846, 28.57610356306631]
2 : [2.0001617818167565, 5.1106995556191235, 5.312099728755062, 9.850959164709389, 11.105479296777533]
3 : [2.0000193858069317, 5.007363829980422, 5.019913885573657, 8.754041226689143, 10.144346676245988]
4 : [2.0000041322707753, 5.0006075867135555, 5.002155682413264, 8.33353759898892, 10.025865354443946]
5 : [2.0000009859606065, 5.0000816796399175, 5.000285469643729, 8.149168665875763, 10.005569349181899]
6 : [2.000000246223299, 5.000014776973109, 5.000044381956682, 8.065066455283487, 10.001468930339382]
7 : [2.0000000630622567, 5.000002823024907, 5.000008239974368, 8.028153010528772, 10.000468799219336]
8 : [2.0000000164229195, 5.000000570275028, 5.0000017493476205, 8.012116100902565, 10.000169396595176]
9 : [2.000000004349842, 5.000000125674543, 5.00000040374791, 8.005217808521643, 10.000065008173907]
10 : [2.000000001184944, 5.000000032338049, 5.000000100068931, 8.002246540719321, 10.000025616996302]
11 : [2.000000000348883, 5.000000011629599, 5.000000028422656, 8.000968961047375, 10.000010261015335]
12 : [2.000000000126624, 5.000000006847966, 5.000000011108204, 8.000418174100572, 10.000004191534279]
13 : [2.0000000000672977, 5.000000005685576, 5.000000006889738, 8.000180760795537, 10.000001777225592]
14 : [2.0000000000514087, 5.0000000053451, 5.000000005904051, 8.0000781966409, 10.000000813577659]
15 : [2.000000000047144, 5.0000000052278235, 5.000000005691839, 8.000033884876887, 10.000000428094348]
16 : [2.0000000000459983, 5.000000005193692, 5.0000000056433125, 8.00001470700946, 10.000000273665847]
17 : [2.0000000000456875, 5.0000000051847975, 5.000000005631361, 8.000006406150307, 10.000000211725494]
18 : [2.0000000000456035, 5.000000005182605, 5.000000005628225, 8.00000280831106, 10.000000186860447]
19 : [2.0000000000455813, 5.0000000051818, 5.000000005626988, 8.00000124877095, 10.000000176870486]
The operations are implemented using late evaluation. The operations return expressions, and computation is done within the assignment operator. The advantage is to avoid dynamic allocation. An exception is InnerProduct, which allows an expression in the second argument (and then needs vector allocation in every call).
[ ]:
[ ]: