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