462
|
1 /* |
|
2 |
2847
|
3 Copyright (C) 1996, 1997 John W. Eaton |
462
|
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 |
|
9 Free Software Foundation; either version 2, or (at your option) any |
|
10 later version. |
|
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 |
|
18 along with Octave; see the file COPYING. If not, write to the Free |
1315
|
19 Software Foundation, 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. |
462
|
20 |
|
21 */ |
|
22 |
4192
|
23 #if defined (__GNUG__) && defined (USE_PRAGMA_INTERFACE_IMPLEMENTATION) |
1296
|
24 #pragma implementation |
|
25 #endif |
|
26 |
462
|
27 #ifdef HAVE_CONFIG_H |
1192
|
28 #include <config.h> |
462
|
29 #endif |
|
30 |
|
31 #include "EIG.h" |
2815
|
32 #include "dColVector.h" |
1847
|
33 #include "f77-fcn.h" |
462
|
34 #include "lo-error.h" |
|
35 |
|
36 extern "C" |
|
37 { |
4552
|
38 F77_RET_T |
|
39 F77_FUNC (dgeev, DGEEV) (F77_CONST_CHAR_ARG_DECL, |
|
40 F77_CONST_CHAR_ARG_DECL, |
|
41 const int&, double*, const int&, double*, |
|
42 double*, double*, const int&, double*, |
|
43 const int&, double*, const int&, int& |
|
44 F77_CHAR_ARG_LEN_DECL |
|
45 F77_CHAR_ARG_LEN_DECL); |
462
|
46 |
4552
|
47 F77_RET_T |
|
48 F77_FUNC (zgeev, ZGEEV) (F77_CONST_CHAR_ARG_DECL, |
|
49 F77_CONST_CHAR_ARG_DECL, |
|
50 const int&, Complex*, const int&, Complex*, |
|
51 Complex*, const int&, Complex*, const int&, |
|
52 Complex*, const int&, double*, int& |
|
53 F77_CHAR_ARG_LEN_DECL |
|
54 F77_CHAR_ARG_LEN_DECL); |
2815
|
55 |
4552
|
56 F77_RET_T |
|
57 F77_FUNC (dsyev, DSYEV) (F77_CONST_CHAR_ARG_DECL, |
|
58 F77_CONST_CHAR_ARG_DECL, |
|
59 const int&, double*, const int&, double*, |
|
60 double*, const int&, int& |
|
61 F77_CHAR_ARG_LEN_DECL |
|
62 F77_CHAR_ARG_LEN_DECL); |
2815
|
63 |
4552
|
64 F77_RET_T |
|
65 F77_FUNC (zheev, ZHEEV) (F77_CONST_CHAR_ARG_DECL, |
|
66 F77_CONST_CHAR_ARG_DECL, |
|
67 const int&, Complex*, const int&, double*, |
|
68 Complex*, const int&, double*, int& |
|
69 F77_CHAR_ARG_LEN_DECL |
|
70 F77_CHAR_ARG_LEN_DECL); |
462
|
71 } |
|
72 |
|
73 int |
4725
|
74 EIG::init (const Matrix& a, bool calc_ev) |
462
|
75 { |
2815
|
76 if (a.is_symmetric ()) |
4725
|
77 return symmetric_init (a, calc_ev); |
2815
|
78 |
1934
|
79 int n = a.rows (); |
|
80 |
|
81 if (n != a.cols ()) |
462
|
82 { |
|
83 (*current_liboctave_error_handler) ("EIG requires square matrix"); |
|
84 return -1; |
|
85 } |
|
86 |
2815
|
87 int info = 0; |
462
|
88 |
1934
|
89 Matrix atmp = a; |
|
90 double *tmp_data = atmp.fortran_vec (); |
|
91 |
|
92 Array<double> wr (n); |
|
93 double *pwr = wr.fortran_vec (); |
|
94 |
|
95 Array<double> wi (n); |
|
96 double *pwi = wi.fortran_vec (); |
|
97 |
4725
|
98 int nvr = calc_ev ? n : 0; |
|
99 Matrix vr (nvr, nvr); |
462
|
100 double *pvr = vr.fortran_vec (); |
1934
|
101 |
2815
|
102 // XXX FIXME XXX -- it might be possible to choose a better value of |
|
103 // lwork that would result in more efficient computations. |
|
104 |
462
|
105 int lwork = 8*n; |
1934
|
106 Array<double> work (lwork); |
|
107 double *pwork = work.fortran_vec (); |
462
|
108 |
1365
|
109 double *dummy = 0; |
462
|
110 int idummy = 1; |
|
111 |
4552
|
112 F77_XFCN (dgeev, DGEEV, (F77_CONST_CHAR_ARG2 ("N", 1), |
4725
|
113 F77_CONST_CHAR_ARG2 (calc_ev ? "V" : "N", 1), |
4552
|
114 n, tmp_data, n, pwr, pwi, dummy, |
|
115 idummy, pvr, n, pwork, lwork, info |
|
116 F77_CHAR_ARG_LEN (1) |
|
117 F77_CHAR_ARG_LEN (1))); |
462
|
118 |
2815
|
119 if (f77_exception_encountered || info < 0) |
1934
|
120 (*current_liboctave_error_handler) ("unrecoverable error in dgeev"); |
|
121 else |
462
|
122 { |
2815
|
123 if (info > 0) |
|
124 (*current_liboctave_error_handler) ("dgeev failed to converge"); |
|
125 else |
|
126 { |
|
127 lambda.resize (n); |
4725
|
128 v.resize (nvr, nvr); |
1934
|
129 |
2815
|
130 for (int j = 0; j < n; j++) |
462
|
131 { |
2815
|
132 if (wi.elem (j) == 0.0) |
|
133 { |
|
134 lambda.elem (j) = Complex (wr.elem (j)); |
4725
|
135 for (int i = 0; i < nvr; i++) |
2815
|
136 v.elem (i, j) = vr.elem (i, j); |
|
137 } |
|
138 else |
1934
|
139 { |
2815
|
140 if (j+1 >= n) |
|
141 { |
|
142 (*current_liboctave_error_handler) |
|
143 ("EIG: internal error"); |
|
144 return -1; |
|
145 } |
|
146 |
3178
|
147 lambda.elem(j) = Complex (wr.elem(j), wi.elem(j)); |
|
148 lambda.elem(j+1) = Complex (wr.elem(j+1), wi.elem(j+1)); |
|
149 |
4725
|
150 for (int i = 0; i < nvr; i++) |
2815
|
151 { |
|
152 double real_part = vr.elem (i, j); |
|
153 double imag_part = vr.elem (i, j+1); |
|
154 v.elem (i, j) = Complex (real_part, imag_part); |
|
155 v.elem (i, j+1) = Complex (real_part, -imag_part); |
|
156 } |
|
157 j++; |
1934
|
158 } |
2815
|
159 } |
|
160 } |
|
161 } |
|
162 |
|
163 return info; |
|
164 } |
|
165 |
|
166 int |
4725
|
167 EIG::symmetric_init (const Matrix& a, bool calc_ev) |
2815
|
168 { |
|
169 int n = a.rows (); |
1934
|
170 |
2815
|
171 if (n != a.cols ()) |
|
172 { |
|
173 (*current_liboctave_error_handler) ("EIG requires square matrix"); |
|
174 return -1; |
|
175 } |
|
176 |
|
177 int info = 0; |
|
178 |
|
179 Matrix atmp = a; |
|
180 double *tmp_data = atmp.fortran_vec (); |
|
181 |
3585
|
182 ColumnVector wr (n); |
2815
|
183 double *pwr = wr.fortran_vec (); |
|
184 |
|
185 // XXX FIXME XXX -- it might be possible to choose a better value of |
|
186 // lwork that would result in more efficient computations. |
|
187 |
|
188 int lwork = 8*n; |
|
189 Array<double> work (lwork); |
|
190 double *pwork = work.fortran_vec (); |
|
191 |
4725
|
192 F77_XFCN (dsyev, DSYEV, (F77_CONST_CHAR_ARG2 (calc_ev ? "V" : "N", 1), |
4552
|
193 F77_CONST_CHAR_ARG2 ("U", 1), |
|
194 n, tmp_data, n, pwr, pwork, lwork, info |
|
195 F77_CHAR_ARG_LEN (1) |
|
196 F77_CHAR_ARG_LEN (1))); |
2815
|
197 |
|
198 if (f77_exception_encountered || info < 0) |
|
199 (*current_liboctave_error_handler) ("unrecoverable error in dsyev"); |
3585
|
200 else if (info > 0) |
|
201 (*current_liboctave_error_handler) ("dsyev failed to converge"); |
2815
|
202 else |
|
203 { |
3585
|
204 lambda = ComplexColumnVector (wr); |
4725
|
205 v = calc_ev ? ComplexMatrix (atmp) : ComplexMatrix (); |
462
|
206 } |
|
207 |
|
208 return info; |
|
209 } |
|
210 |
|
211 int |
4725
|
212 EIG::init (const ComplexMatrix& a, bool calc_ev) |
462
|
213 { |
2815
|
214 if (a.is_hermitian ()) |
4725
|
215 return hermitian_init (a, calc_ev); |
2815
|
216 |
1934
|
217 int n = a.rows (); |
|
218 |
|
219 if (n != a.cols ()) |
462
|
220 { |
|
221 (*current_liboctave_error_handler) ("EIG requires square matrix"); |
|
222 return -1; |
|
223 } |
|
224 |
2815
|
225 int info = 0; |
462
|
226 |
1934
|
227 ComplexMatrix atmp = a; |
|
228 Complex *tmp_data = atmp.fortran_vec (); |
462
|
229 |
2815
|
230 ComplexColumnVector w (n); |
|
231 Complex *pw = w.fortran_vec (); |
|
232 |
4725
|
233 int nvr = calc_ev ? n : 0; |
|
234 ComplexMatrix vtmp (nvr, nvr); |
2815
|
235 Complex *pv = vtmp.fortran_vec (); |
|
236 |
|
237 // XXX FIXME XXX -- it might be possible to choose a better value of |
|
238 // lwork that would result in more efficient computations. |
|
239 |
462
|
240 int lwork = 8*n; |
1934
|
241 Array<Complex> work (lwork); |
|
242 Complex *pwork = work.fortran_vec (); |
|
243 |
2815
|
244 int lrwork = 2*n; |
|
245 Array<double> rwork (lrwork); |
1934
|
246 double *prwork = rwork.fortran_vec (); |
462
|
247 |
1365
|
248 Complex *dummy = 0; |
462
|
249 int idummy = 1; |
|
250 |
4552
|
251 F77_XFCN (zgeev, ZGEEV, (F77_CONST_CHAR_ARG2 ("N", 1), |
4725
|
252 F77_CONST_CHAR_ARG2 (calc_ev ? "V" : "N", 1), |
4552
|
253 n, tmp_data, n, pw, dummy, idummy, |
|
254 pv, n, pwork, lwork, prwork, info |
|
255 F77_CHAR_ARG_LEN (1) |
|
256 F77_CHAR_ARG_LEN (1))); |
2815
|
257 |
|
258 if (f77_exception_encountered || info < 0) |
|
259 (*current_liboctave_error_handler) ("unrecoverable error in zgeev"); |
|
260 else if (info > 0) |
|
261 (*current_liboctave_error_handler) ("zgeev failed to converge"); |
|
262 else |
|
263 { |
|
264 lambda = w; |
|
265 v = vtmp; |
|
266 } |
|
267 |
|
268 return info; |
|
269 } |
|
270 |
|
271 int |
4725
|
272 EIG::hermitian_init (const ComplexMatrix& a, bool calc_ev) |
2815
|
273 { |
|
274 int n = a.rows (); |
|
275 |
|
276 if (n != a.cols ()) |
|
277 { |
|
278 (*current_liboctave_error_handler) ("EIG requires square matrix"); |
|
279 return -1; |
|
280 } |
|
281 |
|
282 int info = 0; |
462
|
283 |
2815
|
284 ComplexMatrix atmp = a; |
|
285 Complex *tmp_data = atmp.fortran_vec (); |
|
286 |
3585
|
287 ColumnVector wr (n); |
|
288 double *pwr = wr.fortran_vec (); |
2815
|
289 |
|
290 // XXX FIXME XXX -- it might be possible to choose a better value of |
|
291 // lwork that would result in more efficient computations. |
|
292 |
|
293 int lwork = 8*n; |
|
294 Array<Complex> work (lwork); |
|
295 Complex *pwork = work.fortran_vec (); |
|
296 |
|
297 int lrwork = 3*n; |
|
298 Array<double> rwork (lrwork); |
|
299 double *prwork = rwork.fortran_vec (); |
|
300 |
4725
|
301 F77_XFCN (zheev, ZHEEV, (F77_CONST_CHAR_ARG2 (calc_ev ? "V" : "N", 1), |
4552
|
302 F77_CONST_CHAR_ARG2 ("U", 1), |
|
303 n, tmp_data, n, pwr, pwork, lwork, prwork, info |
|
304 F77_CHAR_ARG_LEN (1) |
|
305 F77_CHAR_ARG_LEN (1))); |
2815
|
306 |
|
307 if (f77_exception_encountered || info < 0) |
|
308 (*current_liboctave_error_handler) ("unrecoverable error in zheev"); |
|
309 else if (info > 0) |
|
310 (*current_liboctave_error_handler) ("zheev failed to converge"); |
|
311 else |
|
312 { |
3585
|
313 lambda = ComplexColumnVector (wr); |
4725
|
314 v = calc_ev ? ComplexMatrix (atmp) : ComplexMatrix (); |
2815
|
315 } |
462
|
316 |
|
317 return info; |
|
318 } |
|
319 |
|
320 /* |
|
321 ;;; Local Variables: *** |
|
322 ;;; mode: C++ *** |
|
323 ;;; End: *** |
|
324 */ |