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