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 |
|
23 #ifdef HAVE_CONFIG_H |
1192
|
24 #include <config.h> |
462
|
25 #endif |
|
26 |
|
27 #include "EIG.h" |
2815
|
28 #include "dColVector.h" |
1847
|
29 #include "f77-fcn.h" |
462
|
30 #include "lo-error.h" |
|
31 |
|
32 extern "C" |
|
33 { |
4552
|
34 F77_RET_T |
|
35 F77_FUNC (dgeev, DGEEV) (F77_CONST_CHAR_ARG_DECL, |
|
36 F77_CONST_CHAR_ARG_DECL, |
|
37 const int&, double*, const int&, double*, |
|
38 double*, double*, const int&, double*, |
|
39 const int&, double*, const int&, int& |
|
40 F77_CHAR_ARG_LEN_DECL |
|
41 F77_CHAR_ARG_LEN_DECL); |
462
|
42 |
4552
|
43 F77_RET_T |
|
44 F77_FUNC (zgeev, ZGEEV) (F77_CONST_CHAR_ARG_DECL, |
|
45 F77_CONST_CHAR_ARG_DECL, |
|
46 const int&, Complex*, const int&, Complex*, |
|
47 Complex*, const int&, Complex*, const int&, |
|
48 Complex*, const int&, double*, int& |
|
49 F77_CHAR_ARG_LEN_DECL |
|
50 F77_CHAR_ARG_LEN_DECL); |
2815
|
51 |
4552
|
52 F77_RET_T |
|
53 F77_FUNC (dsyev, DSYEV) (F77_CONST_CHAR_ARG_DECL, |
|
54 F77_CONST_CHAR_ARG_DECL, |
|
55 const int&, double*, const int&, double*, |
|
56 double*, const int&, int& |
|
57 F77_CHAR_ARG_LEN_DECL |
|
58 F77_CHAR_ARG_LEN_DECL); |
2815
|
59 |
4552
|
60 F77_RET_T |
|
61 F77_FUNC (zheev, ZHEEV) (F77_CONST_CHAR_ARG_DECL, |
|
62 F77_CONST_CHAR_ARG_DECL, |
|
63 const int&, Complex*, const int&, double*, |
|
64 Complex*, const int&, double*, int& |
|
65 F77_CHAR_ARG_LEN_DECL |
|
66 F77_CHAR_ARG_LEN_DECL); |
462
|
67 } |
|
68 |
|
69 int |
4725
|
70 EIG::init (const Matrix& a, bool calc_ev) |
462
|
71 { |
2815
|
72 if (a.is_symmetric ()) |
4725
|
73 return symmetric_init (a, calc_ev); |
2815
|
74 |
1934
|
75 int n = a.rows (); |
|
76 |
|
77 if (n != a.cols ()) |
462
|
78 { |
|
79 (*current_liboctave_error_handler) ("EIG requires square matrix"); |
|
80 return -1; |
|
81 } |
|
82 |
2815
|
83 int info = 0; |
462
|
84 |
1934
|
85 Matrix atmp = a; |
|
86 double *tmp_data = atmp.fortran_vec (); |
|
87 |
|
88 Array<double> wr (n); |
|
89 double *pwr = wr.fortran_vec (); |
|
90 |
|
91 Array<double> wi (n); |
|
92 double *pwi = wi.fortran_vec (); |
|
93 |
4933
|
94 volatile int nvr = calc_ev ? n : 0; |
4725
|
95 Matrix vr (nvr, nvr); |
462
|
96 double *pvr = vr.fortran_vec (); |
1934
|
97 |
4800
|
98 int lwork = -1; |
|
99 double dummy_work; |
462
|
100 |
1365
|
101 double *dummy = 0; |
462
|
102 int idummy = 1; |
|
103 |
4552
|
104 F77_XFCN (dgeev, DGEEV, (F77_CONST_CHAR_ARG2 ("N", 1), |
4725
|
105 F77_CONST_CHAR_ARG2 (calc_ev ? "V" : "N", 1), |
4552
|
106 n, tmp_data, n, pwr, pwi, dummy, |
4800
|
107 idummy, pvr, n, &dummy_work, lwork, info |
4552
|
108 F77_CHAR_ARG_LEN (1) |
|
109 F77_CHAR_ARG_LEN (1))); |
462
|
110 |
4800
|
111 if (! f77_exception_encountered && info == 0) |
462
|
112 { |
4800
|
113 lwork = static_cast<int> (dummy_work); |
|
114 Array<double> work (lwork); |
|
115 double *pwork = work.fortran_vec (); |
|
116 |
|
117 F77_XFCN (dgeev, DGEEV, (F77_CONST_CHAR_ARG2 ("N", 1), |
|
118 F77_CONST_CHAR_ARG2 (calc_ev ? "V" : "N", 1), |
|
119 n, tmp_data, n, pwr, pwi, dummy, |
|
120 idummy, pvr, n, pwork, lwork, info |
|
121 F77_CHAR_ARG_LEN (1) |
|
122 F77_CHAR_ARG_LEN (1))); |
|
123 |
|
124 if (f77_exception_encountered || info < 0) |
|
125 { |
|
126 (*current_liboctave_error_handler) ("unrecoverable error in dgeev"); |
|
127 return info; |
|
128 } |
|
129 |
2815
|
130 if (info > 0) |
|
131 { |
4800
|
132 (*current_liboctave_error_handler) ("dgeev failed to converge"); |
|
133 return info; |
|
134 } |
1934
|
135 |
4800
|
136 lambda.resize (n); |
|
137 v.resize (nvr, nvr); |
|
138 |
|
139 for (int j = 0; j < n; j++) |
|
140 { |
|
141 if (wi.elem (j) == 0.0) |
462
|
142 { |
4800
|
143 lambda.elem (j) = Complex (wr.elem (j)); |
|
144 for (int i = 0; i < nvr; i++) |
|
145 v.elem (i, j) = vr.elem (i, j); |
|
146 } |
|
147 else |
|
148 { |
|
149 if (j+1 >= n) |
2815
|
150 { |
4800
|
151 (*current_liboctave_error_handler) ("EIG: internal error"); |
|
152 return -1; |
2815
|
153 } |
4800
|
154 |
|
155 lambda.elem(j) = Complex (wr.elem(j), wi.elem(j)); |
|
156 lambda.elem(j+1) = Complex (wr.elem(j+1), wi.elem(j+1)); |
2815
|
157 |
4800
|
158 for (int i = 0; i < nvr; i++) |
|
159 { |
|
160 double real_part = vr.elem (i, j); |
|
161 double imag_part = vr.elem (i, j+1); |
|
162 v.elem (i, j) = Complex (real_part, imag_part); |
|
163 v.elem (i, j+1) = Complex (real_part, -imag_part); |
1934
|
164 } |
4800
|
165 j++; |
2815
|
166 } |
|
167 } |
|
168 } |
4800
|
169 else |
|
170 (*current_liboctave_error_handler) ("dgeev workspace query failed"); |
2815
|
171 |
|
172 return info; |
|
173 } |
|
174 |
|
175 int |
4725
|
176 EIG::symmetric_init (const Matrix& a, bool calc_ev) |
2815
|
177 { |
|
178 int n = a.rows (); |
1934
|
179 |
2815
|
180 if (n != a.cols ()) |
|
181 { |
|
182 (*current_liboctave_error_handler) ("EIG requires square matrix"); |
|
183 return -1; |
|
184 } |
|
185 |
|
186 int info = 0; |
|
187 |
|
188 Matrix atmp = a; |
|
189 double *tmp_data = atmp.fortran_vec (); |
|
190 |
3585
|
191 ColumnVector wr (n); |
2815
|
192 double *pwr = wr.fortran_vec (); |
|
193 |
4800
|
194 int lwork = -1; |
|
195 double dummy_work; |
2815
|
196 |
4725
|
197 F77_XFCN (dsyev, DSYEV, (F77_CONST_CHAR_ARG2 (calc_ev ? "V" : "N", 1), |
4552
|
198 F77_CONST_CHAR_ARG2 ("U", 1), |
4800
|
199 n, tmp_data, n, pwr, &dummy_work, lwork, info |
4552
|
200 F77_CHAR_ARG_LEN (1) |
|
201 F77_CHAR_ARG_LEN (1))); |
2815
|
202 |
4800
|
203 if (! f77_exception_encountered && info == 0) |
2815
|
204 { |
4800
|
205 lwork = static_cast<int> (dummy_work); |
|
206 Array<double> work (lwork); |
|
207 double *pwork = work.fortran_vec (); |
|
208 |
|
209 F77_XFCN (dsyev, DSYEV, (F77_CONST_CHAR_ARG2 (calc_ev ? "V" : "N", 1), |
|
210 F77_CONST_CHAR_ARG2 ("U", 1), |
|
211 n, tmp_data, n, pwr, pwork, lwork, info |
|
212 F77_CHAR_ARG_LEN (1) |
|
213 F77_CHAR_ARG_LEN (1))); |
|
214 |
|
215 if (f77_exception_encountered || info < 0) |
|
216 { |
|
217 (*current_liboctave_error_handler) ("unrecoverable error in dsyev"); |
|
218 return info; |
|
219 } |
|
220 |
|
221 if (info > 0) |
|
222 { |
|
223 (*current_liboctave_error_handler) ("dsyev failed to converge"); |
|
224 return info; |
|
225 } |
|
226 |
3585
|
227 lambda = ComplexColumnVector (wr); |
4725
|
228 v = calc_ev ? ComplexMatrix (atmp) : ComplexMatrix (); |
462
|
229 } |
4800
|
230 else |
|
231 (*current_liboctave_error_handler) ("dsyev workspace query failed"); |
462
|
232 |
|
233 return info; |
|
234 } |
|
235 |
|
236 int |
4725
|
237 EIG::init (const ComplexMatrix& a, bool calc_ev) |
462
|
238 { |
2815
|
239 if (a.is_hermitian ()) |
4725
|
240 return hermitian_init (a, calc_ev); |
2815
|
241 |
1934
|
242 int n = a.rows (); |
|
243 |
|
244 if (n != a.cols ()) |
462
|
245 { |
|
246 (*current_liboctave_error_handler) ("EIG requires square matrix"); |
|
247 return -1; |
|
248 } |
|
249 |
2815
|
250 int info = 0; |
462
|
251 |
1934
|
252 ComplexMatrix atmp = a; |
|
253 Complex *tmp_data = atmp.fortran_vec (); |
462
|
254 |
2815
|
255 ComplexColumnVector w (n); |
|
256 Complex *pw = w.fortran_vec (); |
|
257 |
4725
|
258 int nvr = calc_ev ? n : 0; |
|
259 ComplexMatrix vtmp (nvr, nvr); |
2815
|
260 Complex *pv = vtmp.fortran_vec (); |
|
261 |
4800
|
262 int lwork = -1; |
|
263 Complex dummy_work; |
1934
|
264 |
2815
|
265 int lrwork = 2*n; |
|
266 Array<double> rwork (lrwork); |
1934
|
267 double *prwork = rwork.fortran_vec (); |
462
|
268 |
1365
|
269 Complex *dummy = 0; |
462
|
270 int idummy = 1; |
|
271 |
4552
|
272 F77_XFCN (zgeev, ZGEEV, (F77_CONST_CHAR_ARG2 ("N", 1), |
4725
|
273 F77_CONST_CHAR_ARG2 (calc_ev ? "V" : "N", 1), |
4552
|
274 n, tmp_data, n, pw, dummy, idummy, |
4800
|
275 pv, n, &dummy_work, lwork, prwork, info |
4552
|
276 F77_CHAR_ARG_LEN (1) |
|
277 F77_CHAR_ARG_LEN (1))); |
2815
|
278 |
4800
|
279 if (! f77_exception_encountered && info == 0) |
2815
|
280 { |
4800
|
281 lwork = static_cast<int> (dummy_work.real ()); |
|
282 Array<Complex> work (lwork); |
|
283 Complex *pwork = work.fortran_vec (); |
|
284 |
|
285 F77_XFCN (zgeev, ZGEEV, (F77_CONST_CHAR_ARG2 ("N", 1), |
|
286 F77_CONST_CHAR_ARG2 (calc_ev ? "V" : "N", 1), |
|
287 n, tmp_data, n, pw, dummy, idummy, |
|
288 pv, n, pwork, lwork, prwork, info |
|
289 F77_CHAR_ARG_LEN (1) |
|
290 F77_CHAR_ARG_LEN (1))); |
|
291 |
|
292 if (f77_exception_encountered || info < 0) |
|
293 { |
|
294 (*current_liboctave_error_handler) ("unrecoverable error in zgeev"); |
|
295 return info; |
|
296 } |
|
297 |
|
298 if (info > 0) |
|
299 { |
|
300 (*current_liboctave_error_handler) ("zgeev failed to converge"); |
|
301 return info; |
|
302 } |
|
303 |
2815
|
304 lambda = w; |
|
305 v = vtmp; |
|
306 } |
4800
|
307 else |
|
308 (*current_liboctave_error_handler) ("zgeev workspace query failed"); |
2815
|
309 |
|
310 return info; |
|
311 } |
|
312 |
|
313 int |
4725
|
314 EIG::hermitian_init (const ComplexMatrix& a, bool calc_ev) |
2815
|
315 { |
|
316 int n = a.rows (); |
|
317 |
|
318 if (n != a.cols ()) |
|
319 { |
|
320 (*current_liboctave_error_handler) ("EIG requires square matrix"); |
|
321 return -1; |
|
322 } |
|
323 |
|
324 int info = 0; |
462
|
325 |
2815
|
326 ComplexMatrix atmp = a; |
|
327 Complex *tmp_data = atmp.fortran_vec (); |
|
328 |
3585
|
329 ColumnVector wr (n); |
|
330 double *pwr = wr.fortran_vec (); |
2815
|
331 |
4800
|
332 int lwork = -1; |
|
333 Complex dummy_work; |
2815
|
334 |
|
335 int lrwork = 3*n; |
|
336 Array<double> rwork (lrwork); |
|
337 double *prwork = rwork.fortran_vec (); |
|
338 |
4725
|
339 F77_XFCN (zheev, ZHEEV, (F77_CONST_CHAR_ARG2 (calc_ev ? "V" : "N", 1), |
4552
|
340 F77_CONST_CHAR_ARG2 ("U", 1), |
4800
|
341 n, tmp_data, n, pwr, &dummy_work, lwork, |
|
342 prwork, info |
4552
|
343 F77_CHAR_ARG_LEN (1) |
|
344 F77_CHAR_ARG_LEN (1))); |
2815
|
345 |
4800
|
346 if (! f77_exception_encountered && info == 0) |
2815
|
347 { |
4800
|
348 lwork = static_cast<int> (dummy_work.real ()); |
|
349 Array<Complex> work (lwork); |
|
350 Complex *pwork = work.fortran_vec (); |
|
351 |
|
352 F77_XFCN (zheev, ZHEEV, (F77_CONST_CHAR_ARG2 (calc_ev ? "V" : "N", 1), |
|
353 F77_CONST_CHAR_ARG2 ("U", 1), |
|
354 n, tmp_data, n, pwr, pwork, lwork, prwork, info |
|
355 F77_CHAR_ARG_LEN (1) |
|
356 F77_CHAR_ARG_LEN (1))); |
|
357 |
|
358 if (f77_exception_encountered || info < 0) |
|
359 { |
|
360 (*current_liboctave_error_handler) ("unrecoverable error in zheev"); |
|
361 return info; |
|
362 } |
|
363 |
|
364 if (info > 0) |
|
365 { |
|
366 (*current_liboctave_error_handler) ("zheev failed to converge"); |
|
367 return info; |
|
368 } |
|
369 |
3585
|
370 lambda = ComplexColumnVector (wr); |
4725
|
371 v = calc_ev ? ComplexMatrix (atmp) : ComplexMatrix (); |
2815
|
372 } |
4800
|
373 else |
|
374 (*current_liboctave_error_handler) ("zheev workspace query failed"); |
462
|
375 |
|
376 return info; |
|
377 } |
|
378 |
|
379 /* |
|
380 ;;; Local Variables: *** |
|
381 ;;; mode: C++ *** |
|
382 ;;; End: *** |
|
383 */ |