Mercurial > octave
annotate liboctave/dbleQR.cc @ 7553:56be6f31dd4e
implementation of QR factorization updating
author | Jaroslav Hajek <highegg@gmail.com> |
---|---|
date | Tue, 04 Mar 2008 21:47:11 -0500 |
parents | 29980c6b8604 |
children | 40574114c514 |
rev | line source |
---|---|
457 | 1 /* |
2 | |
7017 | 3 Copyright (C) 1994, 1995, 1996, 1997, 2002, 2003, 2004, 2005, 2007 |
4 John W. Eaton | |
457 | 5 |
6 This file is part of Octave. | |
7 | |
8 Octave is free software; you can redistribute it and/or modify it | |
9 under the terms of the GNU General Public License as published by the | |
7016 | 10 Free Software Foundation; either version 3 of the License, or (at your |
11 option) any later version. | |
457 | 12 |
13 Octave is distributed in the hope that it will be useful, but WITHOUT | |
14 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or | |
15 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License | |
16 for more details. | |
17 | |
18 You should have received a copy of the GNU General Public License | |
7016 | 19 along with Octave; see the file COPYING. If not, see |
20 <http://www.gnu.org/licenses/>. | |
457 | 21 |
22 */ | |
23 | |
7553
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
24 // updating/downdating by Jaroslav Hajek 2008 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
25 |
457 | 26 #ifdef HAVE_CONFIG_H |
1192 | 27 #include <config.h> |
457 | 28 #endif |
29 | |
30 #include "dbleQR.h" | |
1847 | 31 #include "f77-fcn.h" |
457 | 32 #include "lo-error.h" |
7553
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
33 #include "Range.h" |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
34 #include "idx-vector.h" |
457 | 35 |
36 extern "C" | |
37 { | |
4552 | 38 F77_RET_T |
5275 | 39 F77_FUNC (dgeqrf, DGEQRF) (const octave_idx_type&, const octave_idx_type&, double*, const octave_idx_type&, |
40 double*, double*, const octave_idx_type&, octave_idx_type&); | |
457 | 41 |
4552 | 42 F77_RET_T |
5275 | 43 F77_FUNC (dorgqr, DORGQR) (const octave_idx_type&, const octave_idx_type&, const octave_idx_type&, double*, |
44 const octave_idx_type&, double*, double*, const octave_idx_type&, octave_idx_type&); | |
7553
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
45 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
46 // these come from qrupdate |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
47 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
48 F77_RET_T |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
49 F77_FUNC (dqr1up, DQR1UP) (const octave_idx_type&, const octave_idx_type&, const octave_idx_type&, |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
50 double*, double*, const double*, const double*); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
51 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
52 F77_RET_T |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
53 F77_FUNC (dqrinc, DQRINC) (const octave_idx_type&, const octave_idx_type&, const octave_idx_type&, |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
54 double*, const double*, double*, const octave_idx_type&, const double*); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
55 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
56 F77_RET_T |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
57 F77_FUNC (dqrdec, DQRDEC) (const octave_idx_type&, const octave_idx_type&, const octave_idx_type&, |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
58 double*, const double*, double*, const octave_idx_type&); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
59 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
60 F77_RET_T |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
61 F77_FUNC (dqrinr, DQRINR) (const octave_idx_type&, const octave_idx_type&, |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
62 const double*, double*, const double*, double*, |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
63 const octave_idx_type&, const double*); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
64 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
65 F77_RET_T |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
66 F77_FUNC (dqrder, DQRDER) (const octave_idx_type&, const octave_idx_type&, |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
67 const double*, double*, const double*, double *, |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
68 const octave_idx_type&); |
457 | 69 } |
70 | |
539 | 71 QR::QR (const Matrix& a, QR::type qr_type) |
2763 | 72 : q (), r () |
73 { | |
74 init (a, qr_type); | |
75 } | |
76 | |
77 void | |
78 QR::init (const Matrix& a, QR::type qr_type) | |
457 | 79 { |
5275 | 80 octave_idx_type m = a.rows (); |
81 octave_idx_type n = a.cols (); | |
457 | 82 |
83 if (m == 0 || n == 0) | |
84 { | |
85 (*current_liboctave_error_handler) ("QR must have non-empty matrix"); | |
86 return; | |
87 } | |
88 | |
5275 | 89 octave_idx_type min_mn = m < n ? m : n; |
1927 | 90 Array<double> tau (min_mn); |
91 double *ptau = tau.fortran_vec (); | |
1921 | 92 |
5275 | 93 octave_idx_type lwork = 32*n; |
1927 | 94 Array<double> work (lwork); |
95 double *pwork = work.fortran_vec (); | |
1921 | 96 |
5275 | 97 octave_idx_type info = 0; |
457 | 98 |
3883 | 99 Matrix A_fact = a; |
100 if (m > n && qr_type != QR::economy) | |
101 A_fact.resize (m, m, 0.0); | |
457 | 102 |
1927 | 103 double *tmp_data = A_fact.fortran_vec (); |
457 | 104 |
1927 | 105 F77_XFCN (dgeqrf, DGEQRF, (m, n, tmp_data, m, ptau, pwork, lwork, info)); |
457 | 106 |
7482
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
107 if (qr_type == QR::raw) |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
108 { |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
109 for (octave_idx_type j = 0; j < min_mn; j++) |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
110 { |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
111 octave_idx_type limit = j < min_mn - 1 ? j : min_mn - 1; |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
112 for (octave_idx_type i = limit + 1; i < m; i++) |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
113 A_fact.elem (i, j) *= tau.elem (j); |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
114 } |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
115 |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
116 r = A_fact; |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
117 |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
118 if (m > n) |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
119 r.resize (m, n); |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
120 } |
539 | 121 else |
457 | 122 { |
7482
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
123 octave_idx_type n2 = (qr_type == QR::economy) ? min_mn : m; |
3069 | 124 |
7482
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
125 if (qr_type == QR::economy && m > n) |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
126 r.resize (n, n, 0.0); |
539 | 127 else |
7482
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
128 r.resize (m, n, 0.0); |
1921 | 129 |
7482
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
130 for (octave_idx_type j = 0; j < n; j++) |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
131 { |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
132 octave_idx_type limit = j < min_mn-1 ? j : min_mn-1; |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
133 for (octave_idx_type i = 0; i <= limit; i++) |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
134 r.elem (i, j) = tmp_data[m*j+i]; |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
135 } |
457 | 136 |
7482
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
137 lwork = 32 * n2; |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
138 work.resize (lwork); |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
139 double *pwork2 = work.fortran_vec (); |
457 | 140 |
7482
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
141 F77_XFCN (dorgqr, DORGQR, (m, n2, min_mn, tmp_data, m, ptau, |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
142 pwork2, lwork, info)); |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
143 |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
144 q = A_fact; |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
145 q.resize (m, n2); |
1921 | 146 } |
457 | 147 } |
148 | |
7553
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
149 QR::QR (const Matrix& q, const Matrix& r) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
150 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
151 if (q.columns () != r.rows ()) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
152 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
153 (*current_liboctave_error_handler) ("QR dimensions mismatch"); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
154 return; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
155 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
156 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
157 this->q = q; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
158 this->r = r; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
159 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
160 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
161 void |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
162 QR::update (const Matrix& u, const Matrix& v) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
163 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
164 octave_idx_type m = q.rows (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
165 octave_idx_type n = r.columns (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
166 octave_idx_type k = q.columns (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
167 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
168 if (u.length () == m && v.length () == n) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
169 F77_XFCN (dqr1up, DQR1UP, (m, n, k, q.fortran_vec (), r.fortran_vec (), |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
170 u.data (), v.data ())); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
171 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
172 (*current_liboctave_error_handler) ("QR update dimensions mismatch"); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
173 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
174 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
175 void |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
176 QR::insert_col (const Matrix& u, octave_idx_type j) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
177 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
178 octave_idx_type m = q.rows (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
179 octave_idx_type n = r.columns (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
180 octave_idx_type k = q.columns (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
181 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
182 if (u.length () != m) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
183 (*current_liboctave_error_handler) ("QR insert dimensions mismatch"); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
184 else if (j < 1 || j > n+1) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
185 (*current_liboctave_error_handler) ("QR insert index out of range"); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
186 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
187 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
188 Matrix r1 (m, n+1); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
189 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
190 F77_XFCN (dqrinc, DQRINC, (m, n, k, q.fortran_vec (), r.data (), |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
191 r1.fortran_vec (), j, u.data ())); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
192 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
193 r = r1; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
194 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
195 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
196 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
197 void |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
198 QR::delete_col (octave_idx_type j) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
199 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
200 octave_idx_type m = q.rows (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
201 octave_idx_type k = r.rows (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
202 octave_idx_type n = r.columns (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
203 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
204 if (k < m && k < n) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
205 (*current_liboctave_error_handler) ("QR delete dimensions mismatch"); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
206 else if (j < 1 || j > n) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
207 (*current_liboctave_error_handler) ("QR delete index out of range"); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
208 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
209 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
210 Matrix r1 (k, n-1); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
211 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
212 F77_XFCN (dqrdec, DQRDEC, (m, n, k, q.fortran_vec (), r.data (), |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
213 r1.fortran_vec (), j)); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
214 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
215 r = r1; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
216 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
217 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
218 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
219 void |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
220 QR::insert_row (const Matrix& u, octave_idx_type j) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
221 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
222 octave_idx_type m = r.rows (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
223 octave_idx_type n = r.columns (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
224 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
225 if (! q.is_square () || u.length () != n) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
226 (*current_liboctave_error_handler) ("QR insert dimensions mismatch"); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
227 else if (j < 1 || j > m+1) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
228 (*current_liboctave_error_handler) ("QR insert index out of range"); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
229 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
230 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
231 Matrix q1 (m+1, m+1); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
232 Matrix r1 (m+1, n); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
233 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
234 F77_XFCN (dqrinr, DQRINR, (m, n, q.data (), q1.fortran_vec (), |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
235 r.data (), r1.fortran_vec (), j, u.data ())); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
236 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
237 q = q1; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
238 r = r1; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
239 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
240 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
241 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
242 void |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
243 QR::delete_row (octave_idx_type j) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
244 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
245 octave_idx_type m = r.rows (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
246 octave_idx_type n = r.columns (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
247 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
248 if (! q.is_square ()) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
249 (*current_liboctave_error_handler) ("QR insert dimensions mismatch"); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
250 else if (j < 1 || j > n) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
251 (*current_liboctave_error_handler) ("QR delete index out of range"); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
252 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
253 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
254 Matrix q1 (m-1, m-1); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
255 Matrix r1 (m-1, n); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
256 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
257 F77_XFCN (dqrder, DQRDER, (m, n, q.data (), q1.fortran_vec (), |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
258 r.data (), r1.fortran_vec (), j )); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
259 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
260 q = q1; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
261 r = r1; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
262 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
263 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
264 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
265 void |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
266 QR::economize (void) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
267 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
268 idx_vector c (':'), i (Range (1, r.columns ())); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
269 q = Matrix (q.index (c, i)); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
270 r = Matrix (r.index (i, c)); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
271 #if 0 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
272 octave_idx_type r_nc = r.columns (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
273 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
274 if (r.rows () > r_nc) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
275 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
276 q.resize (q.rows (), r_nc); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
277 r.resize (r_nc, r_nc); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
278 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
279 #endif |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
280 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
281 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
282 |
457 | 283 /* |
284 ;;; Local Variables: *** | |
285 ;;; mode: C++ *** | |
286 ;;; End: *** | |
287 */ |