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