Power Iteration and its Variants¶
Copyright (C) 2020 Andreas Kloeckner
MIT License
Permission is hereby granted, free of charge, to any person obtaining a copy of this software and associated documentation files (the "Software"), to deal in the Software without restriction, including without limitation the rights to use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies of the Software, and to permit persons to whom the Software is furnished to do so, subject to the following conditions:The above copyright notice and this permission notice shall be included in all copies or substantial portions of the Software.
THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.
import numpy as np
import numpy.linalg as la
np.set_printoptions(precision=3, linewidth=120)
Let's prepare a matrix with some random or deliberately chosen eigenvalues:
n = 6
if 1:
np.random.seed(70)
eigvecs = np.random.randn(n, n)
eigvals = np.sort(np.random.randn(n))
# Uncomment for near-duplicate largest-magnitude eigenvalue
# eigvals[1] = eigvals[0] + 1e-3
A = eigvecs.dot(np.diag(eigvals)).dot(la.inv(eigvecs))
print(eigvals)
else:
# Complex eigenvalues
np.random.seed(40)
A = np.random.randn(n, n)
print(la.eig(A)[0])
[-2.668 -0.958 -0.33 -0.292 -0.186 -0.144]
Let's also pick an initial vector:
x0 = np.random.randn(n)
x0
array([ 2.269, 0.664, 0.899, -0.366, 0.463, 0.08 ])
Power iteration¶
x = x0
Now implement plain power iteration.
for i in range(20):
x = A @ x
print(x)
[-7.705 22.151 2.865 -4.648 4.043 11.651] [ 19.254 -66.76 -10.275 13.745 -12.354 -31.851] [-50.872 180.931 28.982 -38.003 32.878 84.804] [ 135.725 -484.052 -78.458 103.027 -87.053 -226.183] [-362.275 1292.214 210.269 -276.615 231.404 603.437] [ 966.69 -3447.945 -561.808 739.696 -616.452 -1609.861] [-2579.07 9198.61 1499.54 -1974.989 1643.638 4294.662] [ 6880.334 -24539.336 -4001.041 5270.259 -4383.84 -11456.774] [-18354.595 65463.009 10674.137 -14060.82 11693.775 30562.786] [ 48963.909 -174633.058 -28475.596 37510.904 -31194.161 -81530.956] [-130618.864 465860.62 75963.646 -100067.48 83214.419 217496.239] [ 348445.777 -1242754.094 -202645.16 266946.529 -221986.341 -580204.159] [-929531.949 3315234.721 540587.199 -712121.531 592181.424 1547782.299] [ 2479667.045 -8843889.704 -1442098.59 1899693.011 -1579732.736 -4128943.082] [-6614886.466 23592411.656 3847016.323 -5067719.196 4214175.002 11014579.211] [ 17646208.664 -62936320.902 -10262497.484 13518907.353 -11241947.561 -29383053.282] [-4.707e+07 1.679e+08 2.738e+07 -3.606e+07 2.999e+07 7.838e+07] [ 1.256e+08 -4.479e+08 -7.303e+07 9.621e+07 -8.000e+07 -2.091e+08] [-3.350e+08 1.195e+09 1.948e+08 -2.566e+08 2.134e+08 5.578e+08] [ 8.936e+08 -3.187e+09 -5.197e+08 6.846e+08 -5.693e+08 -1.488e+09]
- What's the problem with this method?
- Does anything useful come of this?
- How do we fix it?
Normalized power iteration¶
Back to the beginning: Reset to the initial vector.
x0 = np.random.randn(n)
x = x0
Implement normalized power iteration.
for i in range(10):
x = A @ x
nrm = la.norm(x)
x = x/nrm
print(la.solve(eigvecs, x))
#print(x)
print(nrm)
[ 0.242 0.297 -0.081 0.058 -0.045 -0.038] [-0.259 -0.114 0.011 -0.007 0.003 0.002] [ 2.709e-01 4.282e-02 -1.380e-03 7.734e-04 -2.481e-04 -1.241e-04] [-2.757e-01 -1.565e-02 1.738e-04 -8.599e-05 1.763e-05 6.822e-06] [ 2.774e-01 5.653e-03 -2.165e-05 9.455e-06 -1.239e-06 -3.710e-07] [-2.780e-01 -2.035e-03 2.686e-06 -1.036e-06 8.678e-08 2.010e-08] [ 2.782e-01 7.312e-04 -3.327e-07 1.133e-07 -6.067e-09 -1.087e-09] [-2.783e-01 -2.627e-04 4.119e-08 -1.238e-08 4.239e-10 5.877e-11] [ 2.783e-01 9.433e-05 -5.099e-09 1.353e-09 -2.962e-11 -3.178e-12] [-2.783e-01 -3.388e-05 6.312e-10 -1.479e-10 2.070e-12 1.724e-13] 2.6675529342865976
Checking convergence¶
x = x0
errors = []
coeffs = la.solve(eigvecs, x0)
for i in range(10):
x = A @ x
errors.append(
la.norm(x/eigvals[0]**(i+1) - coeffs[0]*eigvecs[:,0], 2))
print("coefficients:", la.solve(eigvecs, x/la.norm(x,2)))
conv_factor = eigvals[1]/eigvals[0]
errors = np.array(errors)
for i in range(len(errors)-1):
print(f"{i=}: {errors[i]=:.6e}, {errors[i+1]/errors[i]=:.6g}, {conv_factor=:.6g}")
coefficients: [ 0.242 0.297 -0.081 0.058 -0.045 -0.038] coefficients: [-0.259 -0.114 0.011 -0.007 0.003 0.002] coefficients: [ 2.709e-01 4.282e-02 -1.380e-03 7.734e-04 -2.481e-04 -1.241e-04] coefficients: [-2.757e-01 -1.565e-02 1.738e-04 -8.599e-05 1.763e-05 6.822e-06] coefficients: [ 2.774e-01 5.653e-03 -2.165e-05 9.455e-06 -1.239e-06 -3.710e-07] coefficients: [-2.780e-01 -2.035e-03 2.686e-06 -1.036e-06 8.678e-08 2.010e-08] coefficients: [ 2.782e-01 7.312e-04 -3.327e-07 1.133e-07 -6.067e-09 -1.087e-09] coefficients: [-2.783e-01 -2.627e-04 4.119e-08 -1.238e-08 4.239e-10 5.877e-11] coefficients: [ 2.783e-01 9.433e-05 -5.099e-09 1.353e-09 -2.962e-11 -3.178e-12] coefficients: [-2.783e-01 -3.388e-05 6.312e-10 -1.479e-10 2.071e-12 1.736e-13] i=0: errors[i]=3.370824e+00, errors[i+1]/errors[i]=0.446913, conv_factor=0.359107 i=1: errors[i]=1.506465e+00, errors[i+1]/errors[i]=0.374136, conv_factor=0.359107 i=2: errors[i]=5.636223e-01, errors[i+1]/errors[i]=0.361846, conv_factor=0.359107 i=3: errors[i]=2.039444e-01, errors[i+1]/errors[i]=0.359503, conv_factor=0.359107 i=4: errors[i]=7.331854e-02, errors[i+1]/errors[i]=0.35911, conv_factor=0.359107 i=5: errors[i]=2.632940e-02, errors[i+1]/errors[i]=0.359074, conv_factor=0.359107 i=6: errors[i]=9.454210e-03, errors[i+1]/errors[i]=0.359087, conv_factor=0.359107 i=7: errors[i]=3.394882e-03, errors[i+1]/errors[i]=0.359097, conv_factor=0.359107 i=8: errors[i]=1.219094e-03, errors[i+1]/errors[i]=0.359103, conv_factor=0.359107
- Now try the matrix variants above.
Inverse iteration¶
What if we want the eigenvalue closest to a give value $\sigma$?
Once again, reset to the beginning.
x = x0/la.norm(x0)
sigma = 1
A_sigma = A-sigma*np.eye(A.shape[0])
for i in range(30):
x = la.solve(A_sigma, x)
nrm = la.norm(x)
x = x/nrm
print(x)
[ 0.011 -0.809 -0.191 0.223 -0.186 -0.474] [-0.127 0.775 0.149 -0.232 0.16 0.531] [ 0.192 -0.728 -0.123 0.246 -0.135 -0.582] [-0.246 0.674 0.103 -0.261 0.108 0.628] [ 0.295 -0.616 -0.085 0.271 -0.08 -0.669] [-0.338 0.556 0.07 -0.277 0.052 0.702] [ 0.376 -0.497 -0.057 0.277 -0.026 -0.728] [-0.41 0.441 0.047 -0.273 0.002 0.749] [ 0.44 -0.387 -0.038 0.264 0.019 -0.765] [-0.465 0.338 0.031 -0.252 -0.039 0.777] [ 0.488 -0.291 -0.026 0.237 0.056 -0.786] [-0.508 0.249 0.022 -0.219 -0.071 0.792] [ 0.525 -0.209 -0.019 0.199 0.084 -0.796] [-0.54 0.172 0.017 -0.178 -0.095 0.798] [ 0.553 -0.138 -0.016 0.155 0.105 -0.799] [-0.565 0.107 0.016 -0.132 -0.114 0.799] [ 0.575 -0.078 -0.016 0.108 0.122 -0.798] [-0.584 0.05 0.016 -0.084 -0.129 0.796] [ 0.591 -0.025 -0.017 0.06 0.135 -0.792] [-0.597 0.001 0.018 -0.035 -0.14 0.789] [ 0.603 0.021 -0.02 0.011 0.144 -0.784] [-0.607 -0.041 0.022 0.013 -0.148 0.779] [ 0.611 0.061 -0.024 -0.036 0.151 -0.774] [-0.614 -0.079 0.025 0.059 -0.154 0.768] [ 0.616 0.096 -0.027 -0.082 0.156 -0.761] [-0.618 -0.112 0.03 0.103 -0.158 0.755] [ 0.619 0.126 -0.032 -0.124 0.16 -0.748] [-0.619 -0.14 0.034 0.145 -0.161 0.741] [ 0.619 0.154 -0.036 -0.164 0.162 -0.734] [-0.619 -0.166 0.038 0.183 -0.163 0.726]
Rayleigh quotient iteration¶
Can we feed an estimate of the current approximate eigenvalue back into the calculation? (Hint: Rayleigh quotient)
Reset once more.
x = x0/la.norm(x0)
Run this cell in-place (Ctrl-Enter) many times.
for i in range(10):
sigma = np.dot(x, np.dot(A, x))/np.dot(x, x)
x = la.solve(A-sigma*np.eye(n), x)
x = x/la.norm(x)
print(x)
[ 0.063 -0.792 -0.173 0.23 -0.173 -0.505] [-0.191 0.726 0.126 -0.248 0.133 0.585] [ 0.521 -0.254 -0.051 0.121 0.058 -0.802] [-0.544 -0.35 0.081 0.533 -0.15 0.519] [ 0.52 0.378 -0.086 -0.58 0.146 -0.472] [-0.518 -0.379 0.086 0.583 -0.145 0.469] [ 0.518 0.379 -0.086 -0.583 0.145 -0.469] [-0.518 -0.379 0.086 0.583 -0.145 0.469] [-0.518 -0.379 0.086 0.583 -0.145 0.469] [ 0.518 0.379 -0.086 -0.583 0.145 -0.469]
- What's a reasonable stopping criterion?
- Computational downside of this iteration?