Mercurial > octave-nkf
comparison libcruft/qrupdate/dqrder.f @ 7553:56be6f31dd4e
implementation of QR factorization updating
author | Jaroslav Hajek <highegg@gmail.com> |
---|---|
date | Tue, 04 Mar 2008 21:47:11 -0500 |
parents | |
children |
comparison
equal
deleted
inserted
replaced
7552:6070c3bd69c4 | 7553:56be6f31dd4e |
---|---|
1 c Copyright (C) 2008 VZLU Prague, a.s., Czech Republic | |
2 c | |
3 c Author: Jaroslav Hajek <highegg@gmail.com> | |
4 c | |
5 c This source is free software; you can redistribute it and/or modify | |
6 c it under the terms of the GNU General Public License as published by | |
7 c the Free Software Foundation; either version 2 of the License, or | |
8 c (at your option) any later version. | |
9 c | |
10 c This program is distributed in the hope that it will be useful, | |
11 c but WITHOUT ANY WARRANTY; without even the implied warranty of | |
12 c MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the | |
13 c GNU General Public License for more details. | |
14 c | |
15 c You should have received a copy of the GNU General Public License | |
16 c along with this software; see the file COPYING. If not, see | |
17 c <http://www.gnu.org/licenses/>. | |
18 c | |
19 subroutine dqrder(m,n,Q,Q1,R,R1,j) | |
20 c purpose: updates a QR factorization after deleting a row. | |
21 c i.e., given an m-by-m orthogonal matrix Q, an m-by-n | |
22 c upper trapezoidal matrix R and index j in the range | |
23 c 1:m, this subroutine forms the (m-1)-by-(m-1) matrix | |
24 c Q1 and an (m-1)-by-n matrix R1 so that Q1 is again | |
25 c orthogonal, R1 upper trapezoidal, and | |
26 c Q1*R1 = [A(1:j-1,:); A(j+1:m,:)], where A = Q*R. | |
27 c (real version) | |
28 c | |
29 c arguments: | |
30 c m (in) number of rows of the matrix R. | |
31 c n (in) number of columns of the matrix R | |
32 c Q (in) the orthogonal matrix Q | |
33 c Q1 (out) the updated matrix Q1 | |
34 c R (in) the upper trapezoidal matrix R | |
35 c R1 (out) the updated matrix R1 | |
36 c j (in) the position of the new row in R1 | |
37 c | |
38 integer m,n,j | |
39 double precision Q(m,m),Q1(m-1,m-1),R(m,n),R1(m-1,n) | |
40 double precision c | |
41 double precision s,rr,w | |
42 external xerbla,dlacpy,dlartg,drot,dscal,daxpy | |
43 integer i | |
44 c quick return if possible | |
45 if (m == 1) return | |
46 c check arguments | |
47 info = 0 | |
48 if (m < 1) then | |
49 info = 1 | |
50 else if (j < 1 .or. j > n) then | |
51 info = 7 | |
52 end if | |
53 if (info /= 0) then | |
54 call xerbla('DQRDER',info) | |
55 end if | |
56 c setup the new matrix Q1 | |
57 c permute the columns of Q and rows of R so that the deleted row ends | |
58 c up being the topmost row. | |
59 if (j > 1) then | |
60 call dlacpy('0',j-1,m-1,Q(1,2),m,Q1(1,1),m-1) | |
61 end if | |
62 if (j < m) then | |
63 call dlacpy('0',m-j,m-1,Q(j+1,2),m,Q1(j,1),m-1) | |
64 end if | |
65 c setup the new matrix R1 | |
66 call dlacpy('0',m-1,n,R(2,1),m,R1(1,1),m-1) | |
67 c eliminate Q(j,2:m) | |
68 w = Q(j,m) | |
69 do i = m-1,2,-1 | |
70 call dlartg(Q(j,i),w,c,s,rr) | |
71 w = rr | |
72 c apply rotation to rows of R1 | |
73 if (i <= n) then | |
74 call drot(n-i+1,R1(i-1,i),m-1,R1(i,i),m-1,c,s) | |
75 end if | |
76 c apply rotation to columns of Q1 | |
77 call drot(m-1,Q1(1,i-1),1,Q1(1,i),1,c,s) | |
78 end do | |
79 c the last iteration is special, as we don't have the first row of | |
80 c R and first column of Q | |
81 call dlartg(Q(j,1),w,c,s,rr) | |
82 w = rr | |
83 call dscal(n,c,R1(1,1),m-1) | |
84 call daxpy(n,-s,R(1,1),m,R1(1,1),m-1) | |
85 c apply rotation to columns of Q1 | |
86 call dscal(m-1,c,Q1(1,1),1) | |
87 if (j > 1) then | |
88 call daxpy(j-1,-s,Q(1,1),1,Q1(1,1),1) | |
89 end if | |
90 if (j < m) then | |
91 call daxpy(m-j,-s,Q(j+1,1),1,Q1(j,1),1) | |
92 end if | |
93 end |