Mercurial > hg > octave-lyh
annotate liboctave/dbleQR.cc @ 7559:07522d7dcdf8
fixes to QR and Cholesky updating code
author | Jaroslav Hajek <highegg@gmail.com> |
---|---|
date | Wed, 05 Mar 2008 14:23:26 -0500 |
parents | 40574114c514 |
children | 0ef0f9802a37 |
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 | |
24 #ifdef HAVE_CONFIG_H | |
1192 | 25 #include <config.h> |
457 | 26 #endif |
27 | |
28 #include "dbleQR.h" | |
1847 | 29 #include "f77-fcn.h" |
457 | 30 #include "lo-error.h" |
7553
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
31 #include "Range.h" |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
32 #include "idx-vector.h" |
457 | 33 |
34 extern "C" | |
35 { | |
4552 | 36 F77_RET_T |
5275 | 37 F77_FUNC (dgeqrf, DGEQRF) (const octave_idx_type&, const octave_idx_type&, double*, const octave_idx_type&, |
38 double*, double*, const octave_idx_type&, octave_idx_type&); | |
457 | 39 |
4552 | 40 F77_RET_T |
5275 | 41 F77_FUNC (dorgqr, DORGQR) (const octave_idx_type&, const octave_idx_type&, const octave_idx_type&, double*, |
42 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
|
43 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
44 // these come from qrupdate |
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 F77_RET_T |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
47 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
|
48 double*, double*, const double*, const double*); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
49 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
50 F77_RET_T |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
51 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
|
52 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
|
53 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
54 F77_RET_T |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
55 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
|
56 double*, const double*, double*, const octave_idx_type&); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
57 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
58 F77_RET_T |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
59 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
|
60 const double*, double*, const double*, double*, |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
61 const octave_idx_type&, const double*); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
62 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
63 F77_RET_T |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
64 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
|
65 const double*, double*, const double*, double *, |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
66 const octave_idx_type&); |
457 | 67 } |
68 | |
539 | 69 QR::QR (const Matrix& a, QR::type qr_type) |
2763 | 70 : q (), r () |
71 { | |
72 init (a, qr_type); | |
73 } | |
74 | |
75 void | |
76 QR::init (const Matrix& a, QR::type qr_type) | |
457 | 77 { |
5275 | 78 octave_idx_type m = a.rows (); |
79 octave_idx_type n = a.cols (); | |
457 | 80 |
81 if (m == 0 || n == 0) | |
82 { | |
83 (*current_liboctave_error_handler) ("QR must have non-empty matrix"); | |
84 return; | |
85 } | |
86 | |
5275 | 87 octave_idx_type min_mn = m < n ? m : n; |
1927 | 88 Array<double> tau (min_mn); |
89 double *ptau = tau.fortran_vec (); | |
1921 | 90 |
5275 | 91 octave_idx_type lwork = 32*n; |
1927 | 92 Array<double> work (lwork); |
93 double *pwork = work.fortran_vec (); | |
1921 | 94 |
5275 | 95 octave_idx_type info = 0; |
457 | 96 |
3883 | 97 Matrix A_fact = a; |
98 if (m > n && qr_type != QR::economy) | |
99 A_fact.resize (m, m, 0.0); | |
457 | 100 |
1927 | 101 double *tmp_data = A_fact.fortran_vec (); |
457 | 102 |
1927 | 103 F77_XFCN (dgeqrf, DGEQRF, (m, n, tmp_data, m, ptau, pwork, lwork, info)); |
457 | 104 |
7482
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
105 if (qr_type == QR::raw) |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
106 { |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
107 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
|
108 { |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
109 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
|
110 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
|
111 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
|
112 } |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
113 |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
114 r = A_fact; |
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 if (m > n) |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
117 r.resize (m, n); |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
118 } |
539 | 119 else |
457 | 120 { |
7482
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
121 octave_idx_type n2 = (qr_type == QR::economy) ? min_mn : m; |
3069 | 122 |
7482
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
123 if (qr_type == QR::economy && m > n) |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
124 r.resize (n, n, 0.0); |
539 | 125 else |
7482
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
126 r.resize (m, n, 0.0); |
1921 | 127 |
7482
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
128 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
|
129 { |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
130 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
|
131 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
|
132 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
|
133 } |
457 | 134 |
7482
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
135 lwork = 32 * n2; |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
136 work.resize (lwork); |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
137 double *pwork2 = work.fortran_vec (); |
457 | 138 |
7482
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
139 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
|
140 pwork2, lwork, info)); |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
141 |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
142 q = A_fact; |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
143 q.resize (m, n2); |
1921 | 144 } |
457 | 145 } |
146 | |
7553
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
147 QR::QR (const Matrix& q, const Matrix& r) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
148 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
149 if (q.columns () != r.rows ()) |
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 (*current_liboctave_error_handler) ("QR dimensions mismatch"); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
152 return; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
153 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
154 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
155 this->q = q; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
156 this->r = r; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
157 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
158 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
159 void |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
160 QR::update (const Matrix& u, const Matrix& v) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
161 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
162 octave_idx_type m = q.rows (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
163 octave_idx_type n = r.columns (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
164 octave_idx_type k = q.columns (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
165 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
166 if (u.length () == m && v.length () == n) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
167 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
|
168 u.data (), v.data ())); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
169 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
170 (*current_liboctave_error_handler) ("QR update dimensions mismatch"); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
171 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
172 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
173 void |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
174 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
|
175 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
176 octave_idx_type m = q.rows (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
177 octave_idx_type n = r.columns (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
178 octave_idx_type k = q.columns (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
179 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
180 if (u.length () != m) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
181 (*current_liboctave_error_handler) ("QR insert dimensions mismatch"); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
182 else if (j < 1 || j > n+1) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
183 (*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
|
184 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
185 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
186 Matrix r1 (m, n+1); |
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 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
|
189 r1.fortran_vec (), j, u.data ())); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
190 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
191 r = r1; |
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 } |
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 void |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
196 QR::delete_col (octave_idx_type j) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
197 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
198 octave_idx_type m = q.rows (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
199 octave_idx_type k = r.rows (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
200 octave_idx_type n = r.columns (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
201 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
202 if (k < m && k < n) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
203 (*current_liboctave_error_handler) ("QR delete dimensions mismatch"); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
204 else if (j < 1 || j > n) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
205 (*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
|
206 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
207 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
208 Matrix r1 (k, n-1); |
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 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
|
211 r1.fortran_vec (), j)); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
212 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
213 r = r1; |
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 } |
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 void |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
218 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
|
219 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
220 octave_idx_type m = r.rows (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
221 octave_idx_type n = r.columns (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
222 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
223 if (! q.is_square () || u.length () != n) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
224 (*current_liboctave_error_handler) ("QR insert dimensions mismatch"); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
225 else if (j < 1 || j > m+1) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
226 (*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
|
227 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
228 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
229 Matrix q1 (m+1, m+1); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
230 Matrix r1 (m+1, n); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
231 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
232 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
|
233 r.data (), r1.fortran_vec (), j, u.data ())); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
234 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
235 q = q1; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
236 r = r1; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
237 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
238 } |
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 void |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
241 QR::delete_row (octave_idx_type j) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
242 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
243 octave_idx_type m = r.rows (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
244 octave_idx_type n = r.columns (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
245 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
246 if (! q.is_square ()) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
247 (*current_liboctave_error_handler) ("QR insert dimensions mismatch"); |
7559
07522d7dcdf8
fixes to QR and Cholesky updating code
Jaroslav Hajek <highegg@gmail.com>
parents:
7554
diff
changeset
|
248 else if (j < 1 || j > m) |
7553
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
249 (*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
|
250 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
251 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
252 Matrix q1 (m-1, m-1); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
253 Matrix r1 (m-1, n); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
254 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
255 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
|
256 r.data (), r1.fortran_vec (), j )); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
257 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
258 q = q1; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
259 r = r1; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
260 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
261 } |
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 void |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
264 QR::economize (void) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
265 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
266 octave_idx_type r_nc = r.columns (); |
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 if (r.rows () > r_nc) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
269 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
270 q.resize (q.rows (), r_nc); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
271 r.resize (r_nc, r_nc); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
272 } |
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 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7482
diff
changeset
|
275 |
457 | 276 /* |
277 ;;; Local Variables: *** | |
278 ;;; mode: C++ *** | |
279 ;;; End: *** | |
280 */ |