3427
|
1 ## Copyright (C) 1993, 1998, 1999 Auburn University. All rights reserved. |
|
2 ## |
|
3 ## This file is part of Octave. |
|
4 ## |
|
5 ## Octave is free software; you can redistribute it and/or modify it |
|
6 ## under the terms of the GNU General Public License as published by the |
|
7 ## Free Software Foundation; either version 2, or (at your option) any |
|
8 ## later version. |
|
9 ## |
|
10 ## Octave is distributed in the hope that it will be useful, but WITHOUT |
|
11 ## ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or |
|
12 ## FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License |
|
13 ## for more details. |
|
14 ## |
|
15 ## You should have received a copy of the GNU General Public License |
|
16 ## along with Octave; see the file COPYING. If not, write to the Free |
5307
|
17 ## Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA |
|
18 ## 02110-1301 USA. |
3427
|
19 |
3439
|
20 ## -*- texinfo -*- |
5555
|
21 ## @deftypefn {Function File} {[@var{u}, @var{h}, @var{nu}] =} krylov (@var{a}, @var{v}, @var{k}, @var{eps1}, @var{pflg}) |
|
22 ## Construct an orthogonal basis @var{u} of block Krylov subspace |
|
23 ## |
|
24 ## @example |
|
25 ## [v a*v a^2*v ... a^(k+1)*v] |
|
26 ## @end example |
|
27 ## |
|
28 ## @noindent |
|
29 ## Using Householder reflections to guard against loss of orthogonality. |
3427
|
30 ## |
5555
|
31 ## If @var{v} is a vector, then @var{h} contains the Hessenberg matrix |
|
32 ## such that @code{a*u == u*h}. Otherwise, @var{h} is meaningless. |
|
33 ## |
|
34 ## The value of @var{nu} is the dimension of the span of the krylov |
|
35 ## subspace (based on @var{eps1}). |
|
36 ## |
|
37 ## If @var{b} is a vector and @var{k} is greater than @var{m-1}, then |
|
38 ## @var{h} contains the Hessenberg decompostion of @var{a}. |
|
39 ## |
|
40 ## The optional parameter @var{eps1} is the threshold for zero. The |
|
41 ## default value is 1e-12. |
|
42 ## |
|
43 ## If the optional parameter @var{pflg} is nonzero, row pivoting is used |
|
44 ## to improve numerical behavior. The default value is 0. |
3427
|
45 ## |
|
46 ## Reference: Hodel and Misra, "Partial Pivoting in the Computation of |
5555
|
47 ## Krylov Subspaces", to be submitted to Linear Algebra and its |
|
48 ## Applications |
3439
|
49 ## @end deftypefn |
|
50 |
|
51 ## Author: A. Scottedward Hodel <a.s.hodel@eng.auburn.edu> |
3240
|
52 |
4609
|
53 function [Uret, H, nu] = krylov (A, V, k, eps1, pflg); |
3240
|
54 |
|
55 defeps = 1e-12; |
3457
|
56 |
|
57 if (nargin < 3 || nargin > 5) |
|
58 usage ("[U, nu] = krylov (A, V, k, eps1, pflg)") |
|
59 elseif (nargin < 5) |
3240
|
60 pflg = 0; # default permutation flag |
|
61 endif |
3457
|
62 |
3240
|
63 if(nargin < 4) |
|
64 eps1 = defeps; # default tolerance parameter |
|
65 endif |
3273
|
66 |
3457
|
67 if (isempty (eps1)) |
|
68 eps1 = defeps; |
|
69 endif |
|
70 |
4030
|
71 na = issquare (A); |
3457
|
72 if (! na) |
|
73 error ("A(%d x %d) must be square", rows (A), columns (A)); |
3240
|
74 endif |
3225
|
75 |
3457
|
76 [m, kb] = size(V); |
|
77 if (m != na) |
|
78 error("A(%d x %d), V(%d x %d): argument dimensions do not match", |
|
79 na, na, m, kb) |
3240
|
80 endif |
3233
|
81 |
4030
|
82 if (! isscalar (k)) |
3457
|
83 error ("krylov: third argument must be a scalar integer"); |
|
84 endif |
|
85 |
|
86 Vnrm = norm (V, Inf); |
3273
|
87 |
3457
|
88 ## check for trivial solution |
|
89 if (Vnrm == 0) |
|
90 Uret = []; |
4609
|
91 H = []; |
3457
|
92 nu = 0; |
|
93 return; |
3211
|
94 endif |
|
95 |
3240
|
96 # identify trivial null space |
3457
|
97 abm = max (abs ([A, V]')); |
|
98 nzidx = find (abm != 0); |
|
99 zidx = find (abm == 0); |
3240
|
100 |
|
101 # set up vector of pivot points |
|
102 pivot_vec = 1:na; |
3225
|
103 |
3240
|
104 iter = 0; |
3273
|
105 alpha = []; |
3240
|
106 nh = 0; |
3457
|
107 while (length(alpha) < na) && (columns(V) > 0) && (iter < k) |
3240
|
108 iter++; |
3211
|
109 |
3457
|
110 ## get orthogonal basis of V |
3240
|
111 jj = 1; |
3457
|
112 while (jj <= columns (V) && length (alpha) < na) |
|
113 ## index of next Householder reflection |
|
114 nu = length(alpha)+1; |
3273
|
115 |
3240
|
116 short_pv = pivot_vec(nu:na); |
|
117 q = V(:,jj); |
|
118 short_q = q(short_pv); |
3211
|
119 |
3457
|
120 if (norm (short_q) < eps1) |
|
121 ## insignificant column; delete |
|
122 nv = columns (V); |
|
123 if (jj != nv) |
|
124 [V(:,jj), V(:,nv)] = swap (V(:,jj), V(:,nv)); |
|
125 ## XXX FIXME XXX -- H columns should be swapped too. Not done |
|
126 ## since Block Hessenberg structure is lost anyway. |
3240
|
127 endif |
|
128 V = V(:,1:(nv-1)); |
3457
|
129 ## one less reflection |
|
130 nu--; |
3240
|
131 else |
3457
|
132 ## new householder reflection |
|
133 if (pflg) |
|
134 ## locate max magnitude element in short_q |
|
135 asq = abs (short_q); |
|
136 maxv = max (asq); |
|
137 maxidx = find (asq == maxv); |
3240
|
138 pivot_idx = short_pv(maxidx(1)); |
3273
|
139 |
3457
|
140 ## see if need to change the pivot list |
|
141 if (pivot_idx != pivot_vec(nu)) |
3240
|
142 swapidx = maxidx(1) + (nu-1); |
3457
|
143 [pivot_vec(nu), pivot_vec(swapidx)] = ... |
|
144 swap (pivot_vec(nu), pivot_vec(swapidx)); |
3240
|
145 endif |
|
146 endif |
3273
|
147 |
3457
|
148 ## isolate portion of vector for reflection |
3240
|
149 idx = pivot_vec(nu:na); |
|
150 jdx = pivot_vec(1:nu); |
|
151 |
3457
|
152 [hv, av, z] = housh (q(idx), 1, 0); |
3240
|
153 alpha(nu) = av; |
|
154 U(idx,nu) = hv; |
3211
|
155 |
3240
|
156 # reduce V per the reflection |
|
157 V(idx,:) = V(idx,:) - av*hv*(hv' * V(idx,:)); |
|
158 if(iter > 1) |
3457
|
159 ## XXX FIXME XXX -- not done correctly for block case |
3240
|
160 H(nu,nu-1) = V(pivot_vec(nu),jj); |
|
161 endif |
|
162 |
3457
|
163 ## advance to next column of V |
|
164 jj++; |
3240
|
165 endif |
|
166 endwhile |
3211
|
167 |
3457
|
168 ## check for oversize V (due to full rank) |
|
169 if ((columns (V) > na) && (length (alpha) == na)) |
|
170 ## trim to size |
3273
|
171 V = V(:,1:na); |
3457
|
172 elseif (columns(V) > na) |
3273
|
173 krylov_V = V |
|
174 krylov_na = na |
3457
|
175 krylov_length_alpha = length (alpha) |
|
176 error ("This case should never happen; submit a bug report"); |
3273
|
177 endif |
|
178 |
3457
|
179 if (columns (V) > 0) |
|
180 ## construct next Q and multiply |
|
181 Q = zeros (size (V)); |
|
182 for kk = 1:columns (Q) |
3240
|
183 Q(pivot_vec(nu-columns(Q)+kk),kk) = 1; |
|
184 endfor |
3273
|
185 |
3457
|
186 ## apply Householder reflections |
3240
|
187 for ii = nu:-1:1 |
|
188 idx = pivot_vec(ii:na); |
|
189 hv = U(idx,ii); |
|
190 av = alpha(ii); |
|
191 Q(idx,:) = Q(idx,:) - av*hv*(hv'*Q(idx,:)); |
|
192 endfor |
3211
|
193 endif |
|
194 |
3457
|
195 ## multiply to get new vector; |
3240
|
196 V = A*Q; |
3457
|
197 ## project off of previous vectors |
|
198 nu = length (alpha); |
|
199 for i = 1:nu |
|
200 hv = U(:,i); |
|
201 av = alpha(i); |
3240
|
202 V = V - av*hv*(hv'*V); |
|
203 H(i,nu-columns(V)+(1:columns(V))) = V(pivot_vec(i),:); |
|
204 end |
|
205 |
3211
|
206 endwhile |
|
207 |
3457
|
208 ## Back out complete U matrix |
|
209 ## back out U matrix ; |
|
210 j1 = columns (U); |
|
211 for i = j1:-1:1; |
3240
|
212 idx = pivot_vec(i:na); |
|
213 hv = U(idx,i); |
|
214 av = alpha(i); |
3457
|
215 U(:,i) = zeros (na, 1); |
3240
|
216 U(idx(1),i) = 1; |
|
217 U(idx,i:j1) = U(idx,i:j1)-av*hv*(hv'*U(idx,i:j1)); |
3211
|
218 endfor |
3273
|
219 |
3457
|
220 nu = length (alpha); |
3240
|
221 Uret = U; |
3457
|
222 if (max (max (abs (Uret(zidx,:)))) > 0) |
|
223 warning ("krylov: trivial null space corrupted; set pflg = 1 or eps1 > %e", |
|
224 eps1); |
3225
|
225 endif |
|
226 |
3211
|
227 endfunction |