Mercurial > octave
annotate liboctave/numeric/hess.cc @ 22322:93b3cdd36854
move most f77 function decls to separate header files
* liboctave/numeric/lo-amos-proto.h,
liboctave/numeric/lo-arpack-proto.h,
liboctave/numeric/lo-blas-proto.h,
liboctave/numeric/lo-fftpack-proto.h,
liboctave/numeric/lo-lapack-proto.h,
liboctave/numeric/lo-qrupdate-proto.h,
liboctave/numeric/lo-ranlib-proto.h,
liboctave/numeric/lo-slatec-proto.h: New files.
* liboctave/numeric/module.mk: Update.
* __pchip_deriv__.cc, dot.cc, interpreter.cc, ordschur.cc, qz.cc,
CColVector.cc, CMatrix.cc, CNDArray.cc, CRowVector.cc, CSparse.cc,
dColVector.cc, dMatrix.cc, dNDArray.cc, dRowVector.cc, dSparse.cc,
fCColVector.cc, fCMatrix.cc, fCNDArray.cc, fCRowVector.cc,
fColVector.cc, fMatrix.cc, fNDArray.cc, fRowVector.cc, EIG.cc,
aepbalance.cc, chol.cc, eigs-base.cc, fEIG.cc, gepbalance.cc, gsvd.cc,
hess.cc, lo-specfun.cc, lu.cc, oct-rand.cc, qr.cc, qrp.cc,
randpoisson.cc, schur.cc, sparse-qr.cc, svd.cc:
Use new header files.
author | John W. Eaton <jwe@octave.org> |
---|---|
date | Wed, 17 Aug 2016 00:18:08 -0400 |
parents | 6ca3acf5fad8 |
children | bac0d6f07a3e |
rev | line source |
---|---|
457 | 1 /* |
2 | |
19697
4197fc428c7d
maint: Update copyright notices for 2015.
John W. Eaton <jwe@octave.org>
parents:
17744
diff
changeset
|
3 Copyright (C) 1994-2015 John W. Eaton |
457 | 4 |
5 This file is part of Octave. | |
6 | |
7 Octave is free software; you can redistribute it and/or modify it | |
8 under the terms of the GNU General Public License as published by the | |
7016 | 9 Free Software Foundation; either version 3 of the License, or (at your |
10 option) any later version. | |
457 | 11 |
12 Octave is distributed in the hope that it will be useful, but WITHOUT | |
13 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or | |
14 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License | |
15 for more details. | |
16 | |
17 You should have received a copy of the GNU General Public License | |
7016 | 18 along with Octave; see the file COPYING. If not, see |
19 <http://www.gnu.org/licenses/>. | |
457 | 20 |
21 */ | |
22 | |
21724
aba2e6293dd8
use "#if ..." consistently instead of "#ifdef" and "#ifndef"
John W. Eaton <jwe@octave.org>
parents:
21301
diff
changeset
|
23 #if defined (HAVE_CONFIG_H) |
21301
40de9f8f23a6
Use '#include "config.h"' rather than <config.h>.
Rik <rik@octave.org>
parents:
21267
diff
changeset
|
24 # include "config.h" |
457 | 25 #endif |
26 | |
21267
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
27 #include "CMatrix.h" |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
28 #include "dMatrix.h" |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
29 #include "fCMatrix.h" |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
30 #include "fMatrix.h" |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
31 #include "hess.h" |
457 | 32 #include "lo-error.h" |
22322
93b3cdd36854
move most f77 function decls to separate header files
John W. Eaton <jwe@octave.org>
parents:
22317
diff
changeset
|
33 #include "lo-lapack-proto.h" |
457 | 34 |
22317
6ca3acf5fad8
move some new numeric classes to namespace octave::math
John W. Eaton <jwe@octave.org>
parents:
22135
diff
changeset
|
35 namespace octave |
6ca3acf5fad8
move some new numeric classes to namespace octave::math
John W. Eaton <jwe@octave.org>
parents:
22135
diff
changeset
|
36 { |
6ca3acf5fad8
move some new numeric classes to namespace octave::math
John W. Eaton <jwe@octave.org>
parents:
22135
diff
changeset
|
37 namespace math |
6ca3acf5fad8
move some new numeric classes to namespace octave::math
John W. Eaton <jwe@octave.org>
parents:
22135
diff
changeset
|
38 { |
6ca3acf5fad8
move some new numeric classes to namespace octave::math
John W. Eaton <jwe@octave.org>
parents:
22135
diff
changeset
|
39 |
21267
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
40 template <> |
5275 | 41 octave_idx_type |
21267
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
42 hess<Matrix>::init (const Matrix& a) |
457 | 43 { |
5275 | 44 octave_idx_type a_nr = a.rows (); |
45 octave_idx_type a_nc = a.cols (); | |
1932 | 46 |
457 | 47 if (a_nr != a_nc) |
21267
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
48 (*current_liboctave_error_handler) ("hess: requires square matrix"); |
457 | 49 |
1932 | 50 char job = 'N'; |
51 char side = 'R'; | |
457 | 52 |
5275 | 53 octave_idx_type n = a_nc; |
54 octave_idx_type lwork = 32 * n; | |
55 octave_idx_type info; | |
56 octave_idx_type ilo; | |
57 octave_idx_type ihi; | |
457 | 58 |
1932 | 59 hess_mat = a; |
60 double *h = hess_mat.fortran_vec (); | |
457 | 61 |
11570
57632dea2446
attempt better backward compatibility for Array constructors
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
62 Array<double> scale (dim_vector (n, 1)); |
1932 | 63 double *pscale = scale.fortran_vec (); |
64 | |
4552 | 65 F77_XFCN (dgebal, DGEBAL, (F77_CONST_CHAR_ARG2 (&job, 1), |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
66 n, h, n, ilo, ihi, pscale, info |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
67 F77_CHAR_ARG_LEN (1))); |
457 | 68 |
11570
57632dea2446
attempt better backward compatibility for Array constructors
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
69 Array<double> tau (dim_vector (n-1, 1)); |
7482
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
70 double *ptau = tau.fortran_vec (); |
1932 | 71 |
11570
57632dea2446
attempt better backward compatibility for Array constructors
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
72 Array<double> work (dim_vector (lwork, 1)); |
7482
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
73 double *pwork = work.fortran_vec (); |
457 | 74 |
7482
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
75 F77_XFCN (dgehrd, DGEHRD, (n, ilo, ihi, h, n, ptau, pwork, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
76 lwork, info)); |
457 | 77 |
7482
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
78 unitary_hess_mat = hess_mat; |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
79 double *z = unitary_hess_mat.fortran_vec (); |
457 | 80 |
7482
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
81 F77_XFCN (dorghr, DORGHR, (n, ilo, ihi, z, n, ptau, pwork, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
82 lwork, info)); |
457 | 83 |
7482
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
84 F77_XFCN (dgebak, DGEBAK, (F77_CONST_CHAR_ARG2 (&job, 1), |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
85 F77_CONST_CHAR_ARG2 (&side, 1), |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
86 n, ilo, ihi, pscale, n, z, |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
87 n, info |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
88 F77_CHAR_ARG_LEN (1) |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
89 F77_CHAR_ARG_LEN (1))); |
457 | 90 |
7482
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
91 // If someone thinks of a more graceful way of doing |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
92 // this (or faster for that matter :-)), please let |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
93 // me know! |
457 | 94 |
7482
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
95 if (n > 2) |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
96 for (octave_idx_type j = 0; j < a_nc; j++) |
29980c6b8604
don't check f77_exception_encountered
John W. Eaton <jwe@octave.org>
parents:
7017
diff
changeset
|
97 for (octave_idx_type i = j+2; i < a_nr; i++) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
98 hess_mat.elem (i, j) = 0; |
457 | 99 |
100 return info; | |
101 } | |
21267
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
102 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
103 template <> |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
104 octave_idx_type |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
105 hess<FloatMatrix>::init (const FloatMatrix& a) |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
106 { |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
107 octave_idx_type a_nr = a.rows (); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
108 octave_idx_type a_nc = a.cols (); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
109 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
110 if (a_nr != a_nc) |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
111 (*current_liboctave_error_handler) ("hess: requires square matrix"); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
112 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
113 char job = 'N'; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
114 char side = 'R'; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
115 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
116 octave_idx_type n = a_nc; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
117 octave_idx_type lwork = 32 * n; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
118 octave_idx_type info; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
119 octave_idx_type ilo; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
120 octave_idx_type ihi; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
121 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
122 hess_mat = a; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
123 float *h = hess_mat.fortran_vec (); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
124 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
125 Array<float> scale (dim_vector (n, 1)); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
126 float *pscale = scale.fortran_vec (); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
127 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
128 F77_XFCN (sgebal, SGEBAL, (F77_CONST_CHAR_ARG2 (&job, 1), |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
129 n, h, n, ilo, ihi, pscale, info |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
130 F77_CHAR_ARG_LEN (1))); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
131 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
132 Array<float> tau (dim_vector (n-1, 1)); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
133 float *ptau = tau.fortran_vec (); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
134 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
135 Array<float> work (dim_vector (lwork, 1)); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
136 float *pwork = work.fortran_vec (); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
137 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
138 F77_XFCN (sgehrd, SGEHRD, (n, ilo, ihi, h, n, ptau, pwork, |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
139 lwork, info)); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
140 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
141 unitary_hess_mat = hess_mat; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
142 float *z = unitary_hess_mat.fortran_vec (); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
143 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
144 F77_XFCN (sorghr, SORGHR, (n, ilo, ihi, z, n, ptau, pwork, |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
145 lwork, info)); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
146 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
147 F77_XFCN (sgebak, SGEBAK, (F77_CONST_CHAR_ARG2 (&job, 1), |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
148 F77_CONST_CHAR_ARG2 (&side, 1), |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
149 n, ilo, ihi, pscale, n, z, |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
150 n, info |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
151 F77_CHAR_ARG_LEN (1) |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
152 F77_CHAR_ARG_LEN (1))); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
153 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
154 // If someone thinks of a more graceful way of doing |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
155 // this (or faster for that matter :-)), please let |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
156 // me know! |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
157 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
158 if (n > 2) |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
159 for (octave_idx_type j = 0; j < a_nc; j++) |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
160 for (octave_idx_type i = j+2; i < a_nr; i++) |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
161 hess_mat.elem (i, j) = 0; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
162 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
163 return info; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
164 } |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
165 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
166 template <> |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
167 octave_idx_type |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
168 hess<ComplexMatrix>::init (const ComplexMatrix& a) |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
169 { |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
170 octave_idx_type a_nr = a.rows (); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
171 octave_idx_type a_nc = a.cols (); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
172 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
173 if (a_nr != a_nc) |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
174 (*current_liboctave_error_handler) ("hess: requires square matrix"); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
175 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
176 char job = 'N'; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
177 char side = 'R'; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
178 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
179 octave_idx_type n = a_nc; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
180 octave_idx_type lwork = 32 * n; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
181 octave_idx_type info; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
182 octave_idx_type ilo; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
183 octave_idx_type ihi; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
184 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
185 hess_mat = a; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
186 Complex *h = hess_mat.fortran_vec (); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
187 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
188 Array<double> scale (dim_vector (n, 1)); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
189 double *pscale = scale.fortran_vec (); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
190 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
191 F77_XFCN (zgebal, ZGEBAL, (F77_CONST_CHAR_ARG2 (&job, 1), |
22135
407c66ae1e20
reduce warnings from GCC's link-time optimization feature (bug #48531)
John W. Eaton <jwe@octave.org>
parents:
22133
diff
changeset
|
192 n, F77_DBLE_CMPLX_ARG (h), n, ilo, ihi, pscale, info |
21267
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
193 F77_CHAR_ARG_LEN (1))); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
194 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
195 Array<Complex> tau (dim_vector (n-1, 1)); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
196 Complex *ptau = tau.fortran_vec (); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
197 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
198 Array<Complex> work (dim_vector (lwork, 1)); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
199 Complex *pwork = work.fortran_vec (); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
200 |
22135
407c66ae1e20
reduce warnings from GCC's link-time optimization feature (bug #48531)
John W. Eaton <jwe@octave.org>
parents:
22133
diff
changeset
|
201 F77_XFCN (zgehrd, ZGEHRD, (n, ilo, ihi, F77_DBLE_CMPLX_ARG (h), n, F77_DBLE_CMPLX_ARG (ptau), F77_DBLE_CMPLX_ARG (pwork), lwork, info)); |
21267
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
202 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
203 unitary_hess_mat = hess_mat; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
204 Complex *z = unitary_hess_mat.fortran_vec (); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
205 |
22135
407c66ae1e20
reduce warnings from GCC's link-time optimization feature (bug #48531)
John W. Eaton <jwe@octave.org>
parents:
22133
diff
changeset
|
206 F77_XFCN (zunghr, ZUNGHR, (n, ilo, ihi, F77_DBLE_CMPLX_ARG (z), n, F77_DBLE_CMPLX_ARG (ptau), F77_DBLE_CMPLX_ARG (pwork), |
21267
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
207 lwork, info)); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
208 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
209 F77_XFCN (zgebak, ZGEBAK, (F77_CONST_CHAR_ARG2 (&job, 1), |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
210 F77_CONST_CHAR_ARG2 (&side, 1), |
22135
407c66ae1e20
reduce warnings from GCC's link-time optimization feature (bug #48531)
John W. Eaton <jwe@octave.org>
parents:
22133
diff
changeset
|
211 n, ilo, ihi, pscale, n, F77_DBLE_CMPLX_ARG (z), n, info |
21267
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
212 F77_CHAR_ARG_LEN (1) |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
213 F77_CHAR_ARG_LEN (1))); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
214 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
215 // If someone thinks of a more graceful way of |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
216 // doing this (or faster for that matter :-)), |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
217 // please let me know! |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
218 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
219 if (n > 2) |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
220 for (octave_idx_type j = 0; j < a_nc; j++) |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
221 for (octave_idx_type i = j+2; i < a_nr; i++) |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
222 hess_mat.elem (i, j) = 0; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
223 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
224 return info; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
225 } |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
226 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
227 template <> |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
228 octave_idx_type |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
229 hess<FloatComplexMatrix>::init (const FloatComplexMatrix& a) |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
230 { |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
231 octave_idx_type a_nr = a.rows (); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
232 octave_idx_type a_nc = a.cols (); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
233 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
234 if (a_nr != a_nc) |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
235 { |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
236 (*current_liboctave_error_handler) ("hess: requires square matrix"); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
237 return -1; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
238 } |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
239 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
240 char job = 'N'; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
241 char side = 'R'; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
242 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
243 octave_idx_type n = a_nc; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
244 octave_idx_type lwork = 32 * n; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
245 octave_idx_type info; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
246 octave_idx_type ilo; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
247 octave_idx_type ihi; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
248 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
249 hess_mat = a; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
250 FloatComplex *h = hess_mat.fortran_vec (); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
251 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
252 Array<float> scale (dim_vector (n, 1)); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
253 float *pscale = scale.fortran_vec (); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
254 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
255 F77_XFCN (cgebal, CGEBAL, (F77_CONST_CHAR_ARG2 (&job, 1), |
22135
407c66ae1e20
reduce warnings from GCC's link-time optimization feature (bug #48531)
John W. Eaton <jwe@octave.org>
parents:
22133
diff
changeset
|
256 n, F77_CMPLX_ARG (h), n, ilo, ihi, pscale, info |
21267
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
257 F77_CHAR_ARG_LEN (1))); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
258 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
259 Array<FloatComplex> tau (dim_vector (n-1, 1)); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
260 FloatComplex *ptau = tau.fortran_vec (); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
261 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
262 Array<FloatComplex> work (dim_vector (lwork, 1)); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
263 FloatComplex *pwork = work.fortran_vec (); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
264 |
22135
407c66ae1e20
reduce warnings from GCC's link-time optimization feature (bug #48531)
John W. Eaton <jwe@octave.org>
parents:
22133
diff
changeset
|
265 F77_XFCN (cgehrd, CGEHRD, (n, ilo, ihi, F77_CMPLX_ARG (h), n, F77_CMPLX_ARG (ptau), F77_CMPLX_ARG (pwork), lwork, info)); |
21267
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
266 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
267 unitary_hess_mat = hess_mat; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
268 FloatComplex *z = unitary_hess_mat.fortran_vec (); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
269 |
22135
407c66ae1e20
reduce warnings from GCC's link-time optimization feature (bug #48531)
John W. Eaton <jwe@octave.org>
parents:
22133
diff
changeset
|
270 F77_XFCN (cunghr, CUNGHR, (n, ilo, ihi, F77_CMPLX_ARG (z), n, F77_CMPLX_ARG (ptau), F77_CMPLX_ARG (pwork), |
21267
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
271 lwork, info)); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
272 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
273 F77_XFCN (cgebak, CGEBAK, (F77_CONST_CHAR_ARG2 (&job, 1), |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
274 F77_CONST_CHAR_ARG2 (&side, 1), |
22135
407c66ae1e20
reduce warnings from GCC's link-time optimization feature (bug #48531)
John W. Eaton <jwe@octave.org>
parents:
22133
diff
changeset
|
275 n, ilo, ihi, pscale, n, F77_CMPLX_ARG (z), n, info |
21267
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
276 F77_CHAR_ARG_LEN (1) |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
277 F77_CHAR_ARG_LEN (1))); |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
278 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
279 // If someone thinks of a more graceful way of |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
280 // doing this (or faster for that matter :-)), |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
281 // please let me know! |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
282 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
283 if (n > 2) |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
284 for (octave_idx_type j = 0; j < a_nc; j++) |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
285 for (octave_idx_type i = j+2; i < a_nr; i++) |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
286 hess_mat.elem (i, j) = 0; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
287 |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
288 return info; |
f5b8c3aca5f8
make better use of templates for Hessenberg decomposition
John W. Eaton <jwe@octave.org>
parents:
21202
diff
changeset
|
289 } |
22317
6ca3acf5fad8
move some new numeric classes to namespace octave::math
John W. Eaton <jwe@octave.org>
parents:
22135
diff
changeset
|
290 |
6ca3acf5fad8
move some new numeric classes to namespace octave::math
John W. Eaton <jwe@octave.org>
parents:
22135
diff
changeset
|
291 } |
6ca3acf5fad8
move some new numeric classes to namespace octave::math
John W. Eaton <jwe@octave.org>
parents:
22135
diff
changeset
|
292 } |