3
|
1 // Extra Matrix manipulations. -*- C++ -*- |
|
2 /* |
|
3 |
|
4 Copyright (C) 1992 John W. Eaton |
|
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 |
|
10 Free Software Foundation; either version 2, or (at your option) any |
|
11 later version. |
|
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 |
|
19 along with Octave; see the file COPYING. If not, write to the Free |
|
20 Software Foundation, 675 Mass Ave, Cambridge, MA 02139, USA. |
|
21 |
|
22 */ |
|
23 |
|
24 #ifdef __GNUG__ |
|
25 #pragma implementation |
|
26 #endif |
|
27 |
|
28 #include "Matrix.h" |
|
29 #include "mx-inlines.cc" |
227
|
30 #include "lo-error.h" |
233
|
31 #include "f77-uscore.h" |
|
32 |
|
33 // Fortran functions we call. |
|
34 |
|
35 extern "C" |
|
36 { |
|
37 int F77_FCN (dgesv) (const int*, const int*, double*, const int*, |
|
38 int*, double*, const int*, int*); |
|
39 |
|
40 int F77_FCN (dgeqrf) (const int*, const int*, double*, const int*, |
|
41 double*, double*, const int*, int*); |
|
42 |
|
43 int F77_FCN (dorgqr) (const int*, const int*, const int*, double*, |
|
44 const int*, double*, double*, const int*, int*); |
|
45 |
|
46 int F77_FCN (dgeev) (const char*, const char*, const int*, double*, |
|
47 const int*, double*, double*, double*, |
|
48 const int*, double*, const int*, double*, |
|
49 const int*, int*, long, long); |
|
50 |
|
51 int F77_FCN (dgeesx) (const char*, const char*, int (*)(), const char*, |
|
52 const int*, double*, const int*, int*, double*, |
|
53 double*, double*, const int*, double*, double*, |
|
54 double*, const int*, int*, const int*, int*, |
|
55 int*, long, long); |
|
56 |
|
57 int F77_FCN (dgebal) (const char*, const int*, double*, |
|
58 const int*, int*, int*, double*, |
|
59 int*, long, long); |
|
60 |
|
61 int F77_FCN (dgebak) (const char*, const char*, const int*, const int*, |
|
62 const int*, double*, const int*, double*, const int*, |
|
63 int*, long, long); |
|
64 |
|
65 int F77_FCN (dgehrd) (const int*, const int*, const int*, |
|
66 double*, const int*, double*, double*, |
|
67 const int*, int*, long, long); |
|
68 |
|
69 int F77_FCN (dorghr) (const int*, const int*, const int*, |
|
70 double*, const int*, double*, double*, |
|
71 const int*, int*, long, long); |
|
72 |
|
73 int F77_FCN (dgesvd) (const char*, const char*, const int*, |
|
74 const int*, double*, const int*, double*, |
|
75 double*, const int*, double*, const int*, |
|
76 double*, const int*, int*, long, long); |
|
77 |
|
78 int F77_FCN (dpotrf) (const char*, const int*, double*, const int*, |
|
79 int*, long); |
|
80 |
|
81 // |
|
82 // fortran functions for generalized eigenvalue problems |
|
83 // |
|
84 int F77_FCN (reduce) (const int*, const int*, double*, |
|
85 const int*, double*, |
|
86 int*, int*, double*, double*); |
|
87 |
|
88 int F77_FCN (scaleg) (const int*, const int*, double*, |
|
89 const int*, double*, |
|
90 const int*, const int*, double*, double*, double*); |
|
91 |
|
92 int F77_FCN (gradeq) (const int*, const int*, double*, |
|
93 const int*, double*, |
|
94 int*, int*, double*, double*); |
|
95 |
|
96 /* |
|
97 * f2c translates complex*16 as |
|
98 * |
|
99 * typedef struct { doublereal re, im; } doublecomplex; |
|
100 * |
|
101 * and Complex.h from libg++ uses |
|
102 * |
|
103 * protected: |
|
104 * double re; |
|
105 * double im; |
|
106 * |
|
107 * as the only data members, so this should work (fingers crossed that |
|
108 * things don't change). |
|
109 */ |
|
110 |
|
111 int F77_FCN (zgesv) (const int*, const int*, Complex*, const int*, |
|
112 int*, Complex*, const int*, int*); |
|
113 |
|
114 int F77_FCN (zgeqrf) (const int*, const int*, Complex*, const int*, |
|
115 Complex*, Complex*, const int*, int*); |
|
116 |
|
117 int F77_FCN (zgeesx) (const char*, const char*, int (*)(), const char*, |
|
118 const int*, Complex*, const int*, int*, |
|
119 Complex*, Complex*, const int*, double*, double*, |
|
120 Complex*, const int*, double*, int*, int*, |
|
121 long, long); |
|
122 |
|
123 int F77_FCN (zgebal) (const char*, const int*, Complex*, const int*, |
|
124 int*, int*, double*, int*, long, long); |
|
125 |
|
126 int F77_FCN (zgebak) (const char*, const char*, const int*, const int*, |
|
127 const int*, double*, const int*, Complex*, |
|
128 const int*, int*, long, long); |
|
129 |
|
130 int F77_FCN (zgehrd) (const int*, const int*, const int*, Complex*, |
|
131 const int*, Complex*, Complex*, const int*, |
|
132 int*, long, long); |
|
133 |
|
134 int F77_FCN (zunghr) (const int*, const int*, const int*, Complex*, |
|
135 const int*, Complex*, Complex*, const int*, |
|
136 int*, long, long); |
|
137 |
|
138 int F77_FCN (zungqr) (const int*, const int*, const int*, Complex*, |
|
139 const int*, Complex*, Complex*, const int*, int*); |
|
140 |
|
141 int F77_FCN (zgeev) (const char*, const char*, const int*, Complex*, |
|
142 const int*, Complex*, Complex*, const int*, |
|
143 Complex*, const int*, Complex*, const int*, |
|
144 double*, int*, long, long); |
|
145 |
|
146 int F77_FCN (zgesvd) (const char*, const char*, const int*, |
|
147 const int*, Complex*, const int*, double*, |
|
148 Complex*, const int*, Complex*, const int*, |
|
149 Complex*, const int*, double*, int*, long, long); |
|
150 |
|
151 int F77_FCN (zpotrf) (const char*, const int*, Complex*, const int*, |
|
152 int*, long); |
|
153 } |
3
|
154 |
|
155 /* |
22
|
156 * AEPBALANCE operations |
|
157 */ |
|
158 |
|
159 int |
|
160 AEPBALANCE::init (const Matrix& a, const char *balance_job) |
|
161 { |
|
162 if (a.nr != a.nc) |
227
|
163 { |
|
164 (*current_liboctave_error_handler) ("AEPBALANCE requires square matrix"); |
|
165 return -1; |
|
166 } |
22
|
167 |
|
168 int n = a.nc; |
|
169 |
|
170 // Parameters for balance call. |
|
171 |
|
172 int info; |
|
173 int ilo; |
|
174 int ihi; |
|
175 double *scale = new double [n]; |
|
176 |
|
177 // Copy matrix into local structure. |
|
178 |
|
179 balanced_mat = a; |
|
180 |
|
181 F77_FCN (dgebal) (balance_job, &n, balanced_mat.fortran_vec (), |
|
182 &n, &ilo, &ihi, scale, &info, 1L, 1L); |
|
183 |
|
184 // Initialize balancing matrix to identity. |
|
185 |
|
186 balancing_mat = Matrix (n, n, 0.0); |
|
187 for (int i = 0; i < n; i++) |
|
188 balancing_mat.elem (i ,i) = 1.0; |
|
189 |
|
190 F77_FCN (dgebak) (balance_job, "R", &n, &ilo, &ihi, scale, &n, |
|
191 balancing_mat.fortran_vec (), &n, &info, 1L, 1L); |
|
192 |
|
193 delete [] scale; |
|
194 |
|
195 return info; |
|
196 } |
|
197 |
|
198 int |
|
199 ComplexAEPBALANCE::init (const ComplexMatrix& a, const char *balance_job) |
|
200 { |
|
201 |
|
202 int n = a.nc; |
|
203 |
|
204 // Parameters for balance call. |
|
205 |
|
206 int info; |
|
207 int ilo; |
|
208 int ihi; |
|
209 double *scale = new double [n]; |
|
210 |
|
211 // Copy matrix into local structure. |
|
212 |
|
213 balanced_mat = a; |
|
214 |
|
215 F77_FCN (zgebal) (balance_job, &n, balanced_mat.fortran_vec (), |
|
216 &n, &ilo, &ihi, scale, &info, 1L, 1L); |
|
217 |
|
218 // Initialize balancing matrix to identity. |
|
219 |
|
220 balancing_mat = Matrix (n, n, 0.0); |
|
221 for (int i = 0; i < n; i++) |
|
222 balancing_mat (i, i) = 1.0; |
|
223 |
|
224 F77_FCN (zgebak) (balance_job, "R", &n, &ilo, &ihi, scale, &n, |
|
225 balancing_mat.fortran_vec(), &n, &info, 1L, 1L); |
|
226 |
|
227 delete [] scale; |
|
228 |
|
229 return info; |
|
230 } |
|
231 |
|
232 /* |
|
233 * GEPBALANCE operations |
|
234 */ |
|
235 |
|
236 int |
|
237 GEPBALANCE::init (const Matrix& a, const Matrix& b, const char *balance_job) |
|
238 { |
|
239 if (a.nr != a.nc || a.nr != b.nr || b.nr != b.nc) |
227
|
240 { |
|
241 (*current_liboctave_error_handler) |
|
242 ("GEPBALANCE requires square matrices of the same size"); |
|
243 return -1; |
|
244 } |
22
|
245 |
|
246 int n = a.nc; |
|
247 |
|
248 // Parameters for balance call. |
|
249 |
|
250 int info; |
|
251 int ilo; |
|
252 int ihi; |
|
253 double *cscale = new double [n]; |
|
254 double *cperm = new double [n]; |
|
255 Matrix wk (n, 6, 0.0); |
|
256 |
|
257 // Back out the permutations: |
|
258 // |
|
259 // cscale contains the exponents of the column scaling factors in its |
|
260 // ilo through ihi locations and the reducing column permutations in |
|
261 // its first ilo-1 and its ihi+1 through n locations. |
|
262 // |
|
263 // cperm contains the column permutations applied in grading the a and b |
|
264 // submatrices in its ilo through ihi locations. |
|
265 // |
|
266 // wk contains the exponents of the row scaling factors in its ilo |
|
267 // through ihi locations, the reducing row permutations in its first |
|
268 // ilo-1 and its ihi+1 through n locations, and the row permutations |
|
269 // applied in grading the a and b submatrices in its n+ilo through |
|
270 // n+ihi locations. |
|
271 |
|
272 // Copy matrices into local structure. |
|
273 |
|
274 balanced_a_mat = a; |
|
275 balanced_b_mat = b; |
|
276 |
|
277 // Initialize balancing matrices to identity. |
|
278 |
|
279 left_balancing_mat = Matrix(n,n,0.0); |
|
280 for (int i = 0; i < n; i++) |
|
281 left_balancing_mat (i, i) = 1.0; |
|
282 |
|
283 right_balancing_mat = left_balancing_mat; |
|
284 |
|
285 // Check for permutation option. |
|
286 |
|
287 if (*balance_job == 'P' || *balance_job == 'B') |
|
288 { |
233
|
289 F77_FCN (reduce) (&n, &n, balanced_a_mat.fortran_vec (), |
|
290 &n, balanced_b_mat.fortran_vec (), &ilo, &ihi, |
|
291 cscale, wk.fortran_vec ()); |
22
|
292 } |
|
293 else |
|
294 { |
|
295 |
|
296 // Set up for scaling later. |
|
297 |
|
298 ilo = 1; |
|
299 ihi = n; |
|
300 } |
|
301 |
|
302 // Check for scaling option. |
|
303 |
|
304 if ((*balance_job == 'S' || *balance_job == 'B') && ilo != ihi) |
|
305 { |
233
|
306 F77_FCN (scaleg) (&n, &n, balanced_a_mat.fortran_vec (), |
|
307 &n, balanced_b_mat.fortran_vec (), &ilo, &ihi, |
|
308 cscale, cperm, wk.fortran_vec ()); |
22
|
309 } |
|
310 else |
|
311 { |
|
312 |
|
313 // Set scaling data to 0's. |
|
314 |
|
315 for (int tmp = ilo-1; tmp < ihi; tmp++) |
|
316 { |
|
317 cscale[tmp] = 0.0; |
|
318 wk.elem(tmp,0) = 0.0; |
|
319 } |
|
320 } |
|
321 |
|
322 // Scaleg returns exponents, not values, so... |
|
323 |
|
324 for (int tmp = ilo-1; tmp < ihi; tmp++) |
|
325 { |
|
326 cscale[tmp] = pow(2.0,cscale[tmp]); |
|
327 wk.elem(tmp,0) = pow(2.0,-wk.elem(tmp,0)); |
|
328 } |
|
329 |
|
330 // Column permutations/scaling. |
|
331 |
|
332 F77_FCN (dgebak) (balance_job, "R", &n, &ilo, &ihi, cscale, &n, |
|
333 right_balancing_mat.fortran_vec (), &n, &info, 1L, |
|
334 1L); |
|
335 |
|
336 // Row permutations/scaling. |
|
337 |
|
338 F77_FCN (dgebak) (balance_job, "L", &n, &ilo, &ihi, &wk.elem (0, 0), &n, |
|
339 left_balancing_mat.fortran_vec (), &n, &info, 1L, 1L); |
|
340 |
|
341 // XXX FIXME XXX --- these four lines need to be added and debugged. |
|
342 // GEPBALANCE::init will work without them, though, so here they are. |
|
343 |
|
344 #if 0 |
|
345 if ((*balance_job == 'P' || *balance_job == 'B') && ilo != ihi) |
|
346 { |
|
347 F77_FCN (gradeq) (&n, &n, balanced_a_mat.fortran_vec (), |
|
348 &n, balanced_b_mat.fortran_vec (), &ilo, &ihi, |
|
349 cperm, &wk.elem (0, 1)); |
|
350 } |
|
351 #endif |
|
352 |
|
353 // Transpose for aa = cc*a*dd convention... |
|
354 left_balancing_mat = left_balancing_mat.transpose (); |
|
355 |
|
356 delete [] cscale; |
|
357 delete [] cperm; |
|
358 |
|
359 return info; |
|
360 } |
|
361 |
|
362 /* |
182
|
363 * CHOL stuff |
|
364 */ |
|
365 |
|
366 int |
|
367 CHOL::init (const Matrix& a) |
|
368 { |
|
369 if (a.nr != a.nc) |
227
|
370 { |
|
371 (*current_liboctave_error_handler) ("CHOL requires square matrix"); |
|
372 return -1; |
|
373 } |
182
|
374 |
|
375 char uplo = 'U'; |
|
376 |
|
377 int n = a.nc; |
|
378 int info; |
|
379 |
|
380 double *h = dup (a.data, a.len); |
|
381 |
|
382 F77_FCN (dpotrf) (&uplo, &n, h, &n, &info, 1L); |
|
383 |
|
384 chol_mat = Matrix (h, n, n); |
|
385 |
|
386 // If someone thinks of a more graceful way of doing this (or faster for |
|
387 // that matter :-)), please let me know! |
|
388 |
|
389 if (n > 1) |
|
390 for (int j = 0; j < a.nc; j++) |
|
391 for (int i = j+1; i < a.nr; i++) |
|
392 chol_mat.elem (i, j) = 0.0; |
|
393 |
|
394 |
|
395 return info; |
|
396 } |
|
397 |
|
398 |
|
399 int |
|
400 ComplexCHOL::init (const ComplexMatrix& a) |
|
401 { |
|
402 if (a.nr != a.nc) |
227
|
403 { |
|
404 (*current_liboctave_error_handler) |
|
405 ("ComplexCHOL requires square matrix"); |
|
406 return -1; |
|
407 } |
182
|
408 |
|
409 char uplo = 'U'; |
|
410 |
|
411 int n = a.nc; |
|
412 int info; |
|
413 |
|
414 Complex *h = dup (a.data, a.len); |
|
415 |
|
416 F77_FCN (zpotrf) (&uplo, &n, h, &n, &info, 1L); |
|
417 |
|
418 chol_mat = ComplexMatrix (h, n, n); |
|
419 |
|
420 // If someone thinks of a more graceful way of doing this (or faster for |
|
421 // that matter :-)), please let me know! |
|
422 |
|
423 if (n > 1) |
|
424 for (int j = 0; j < a.nc; j++) |
|
425 for (int i = j+1; i < a.nr; i++) |
|
426 chol_mat.elem (i, j) = 0.0; |
|
427 |
|
428 return info; |
|
429 } |
|
430 |
|
431 |
|
432 /* |
3
|
433 * HESS stuff |
|
434 */ |
|
435 |
|
436 int |
|
437 HESS::init (const Matrix& a) |
|
438 { |
|
439 if (a.nr != a.nc) |
227
|
440 { |
|
441 (*current_liboctave_error_handler) ("HESS requires square matrix"); |
|
442 return -1; |
|
443 } |
3
|
444 |
22
|
445 char jobbal = 'N'; |
3
|
446 char side = 'R'; |
|
447 |
|
448 int n = a.nc; |
|
449 int lwork = 32 * n; |
|
450 int info; |
|
451 int ilo; |
|
452 int ihi; |
|
453 |
|
454 double *h = dup(a.data, a.len); |
|
455 |
|
456 double *tau = new double [n+1]; |
|
457 double *scale = new double [n]; |
|
458 double *z = new double [n*n]; |
|
459 double *work = new double [lwork]; |
|
460 |
|
461 F77_FCN (dgebal) (&jobbal, &n, h, &n, &ilo, &ihi, scale, &info, |
|
462 1L, 1L); |
|
463 |
|
464 F77_FCN (dgehrd) (&n, &ilo, &ihi, h, &n, tau, work, &lwork, &info, |
|
465 1L, 1L); |
|
466 |
|
467 copy(z,h,n*n); |
|
468 |
|
469 F77_FCN (dorghr) (&n, &ilo, &ihi, z, &n, tau, work, &lwork, &info, |
|
470 1L, 1L); |
|
471 |
|
472 F77_FCN (dgebak) (&jobbal, &side, &n, &ilo, &ihi, scale, &n, z, &n, |
|
473 &info, 1L, 1L); |
|
474 |
|
475 // We need to clear out all of the area below the sub-diagonal which was used |
|
476 // to store the unitary matrix. |
|
477 |
|
478 hess_mat = Matrix(h,n,n); |
|
479 unitary_hess_mat = Matrix(z,n,n); |
|
480 |
|
481 // If someone thinks of a more graceful way of doing this (or faster for |
|
482 // that matter :-)), please let me know! |
|
483 |
|
484 if (n > 2) |
|
485 for (int j = 0; j < a.nc; j++) |
|
486 for (int i = j+2; i < a.nr; i++) |
|
487 hess_mat.elem(i,j) = 0; |
|
488 |
|
489 delete [] tau; |
|
490 delete [] work; |
|
491 delete [] scale; |
|
492 |
|
493 return info; |
|
494 } |
|
495 |
|
496 |
|
497 int |
|
498 ComplexHESS::init (const ComplexMatrix& a) |
|
499 { |
|
500 if (a.nr != a.nc) |
227
|
501 { |
|
502 (*current_liboctave_error_handler) |
|
503 ("ComplexHESS requires square matrix"); |
|
504 return -1; |
|
505 } |
3
|
506 |
22
|
507 char job = 'N'; |
3
|
508 char side = 'R'; |
|
509 |
|
510 int n = a.nc; |
|
511 int lwork = 32 * n; |
|
512 int info; |
|
513 int ilo; |
|
514 int ihi; |
|
515 |
|
516 Complex *h = dup(a.data,a.len); |
|
517 |
|
518 double *scale = new double [n]; |
|
519 Complex *tau = new Complex [n-1]; |
|
520 Complex *work = new Complex [lwork]; |
|
521 Complex *z = new Complex [n*n]; |
|
522 |
|
523 F77_FCN (zgebal) (&job, &n, h, &n, &ilo, &ihi, scale, &info, 1L, 1L); |
|
524 |
|
525 F77_FCN (zgehrd) (&n, &ilo, &ihi, h, &n, tau, work, &lwork, &info, 1L, |
233
|
526 1L); |
3
|
527 |
|
528 copy(z,h,n*n); |
|
529 |
|
530 F77_FCN (zunghr) (&n, &ilo, &ihi, z, &n, tau, work, &lwork, &info, 1L, |
233
|
531 1L); |
3
|
532 |
|
533 F77_FCN (zgebak) (&job, &side, &n, &ilo, &ihi, scale, &n, z, &n, &info, |
233
|
534 1L, 1L); |
3
|
535 |
|
536 hess_mat = ComplexMatrix (h,n,n); |
|
537 unitary_hess_mat = ComplexMatrix (z,n,n); |
|
538 |
|
539 // If someone thinks of a more graceful way of doing this (or faster for |
|
540 // that matter :-)), please let me know! |
|
541 |
|
542 if (n > 2) |
|
543 for (int j = 0; j < a.nc; j++) |
|
544 for (int i = j+2; i < a.nr; i++) |
|
545 hess_mat.elem(i,j) = 0; |
|
546 |
|
547 delete [] work; |
|
548 delete [] tau; |
|
549 delete [] scale; |
|
550 |
|
551 return info; |
|
552 } |
|
553 |
|
554 /* |
|
555 * SCHUR stuff |
|
556 */ |
|
557 |
|
558 static int |
|
559 select_ana (double *a, double *b) |
|
560 { |
|
561 return (*a < 0.0); |
|
562 } |
|
563 |
|
564 static int |
|
565 select_dig (double *a, double *b) |
|
566 { |
|
567 return (hypot (*a, *b) < 1.0); |
|
568 } |
|
569 |
|
570 // GAG. |
|
571 extern "C" { static int (*dummy_select)(); } |
|
572 |
|
573 int |
|
574 SCHUR::init (const Matrix& a, const char *ord) |
|
575 { |
|
576 if (a.nr != a.nc) |
227
|
577 { |
|
578 (*current_liboctave_error_handler) ("SCHUR requires square matrix"); |
|
579 return -1; |
|
580 } |
3
|
581 |
|
582 char jobvs = 'V'; |
|
583 char sort; |
|
584 |
|
585 if (*ord == 'A' || *ord == 'D' || *ord == 'a' || *ord == 'd') |
|
586 sort = 'S'; |
|
587 else |
|
588 sort = 'N'; |
|
589 |
|
590 char sense = 'N'; |
|
591 |
|
592 int n = a.nc; |
|
593 int lwork = 8 * n; |
|
594 int liwork = 1; |
|
595 int info; |
|
596 int sdim; |
|
597 double rconde; |
|
598 double rcondv; |
|
599 |
|
600 double *s = dup(a.data,a.len); |
|
601 |
|
602 double *wr = new double [n]; |
|
603 double *wi = new double [n]; |
|
604 double *q = new double [n*n]; |
|
605 double *work = new double [lwork]; |
|
606 |
|
607 // These are not referenced for the non-ordered Schur routine. |
|
608 |
|
609 int *iwork = (int *) NULL; |
|
610 int *bwork = (int *) NULL; |
|
611 if (*ord == 'A' || *ord == 'D' || *ord == 'a' || *ord == 'd') |
|
612 { |
|
613 iwork = new int [liwork]; |
|
614 bwork = new int [n]; |
|
615 } |
|
616 |
|
617 if (*ord == 'A' || *ord == 'a') |
|
618 { |
|
619 F77_FCN (dgeesx) (&jobvs, &sort, select_ana, &sense, &n, s, &n, |
|
620 &sdim, wr, wi, q, &n, &rconde, &rcondv, work, |
|
621 &lwork, iwork, &liwork, bwork, &info, 1L, 1L); |
|
622 } |
|
623 else if (*ord == 'D' || *ord == 'd') |
|
624 { |
|
625 F77_FCN (dgeesx) (&jobvs, &sort, select_dig, &sense, &n, s, &n, |
|
626 &sdim, wr, wi, q, &n, &rconde, &rcondv, work, |
|
627 &lwork, iwork, &liwork, bwork, &info, 1L, 1L); |
|
628 |
|
629 } |
|
630 else |
|
631 { |
|
632 F77_FCN (dgeesx) (&jobvs, &sort, dummy_select, &sense, &n, s, |
|
633 &n, &sdim, wr, wi, q, &n, &rconde, &rcondv, |
|
634 work, &lwork, iwork, &liwork, bwork, &info, |
|
635 1L, 1L); |
|
636 } |
|
637 |
|
638 |
|
639 schur_mat = Matrix (s, n, n); |
|
640 unitary_mat = Matrix (q, n, n); |
|
641 |
|
642 delete [] wr; |
|
643 delete [] wi; |
|
644 delete [] work; |
|
645 delete [] iwork; |
|
646 delete [] bwork; |
|
647 |
|
648 return info; |
|
649 } |
|
650 |
|
651 static int |
|
652 complex_select_ana (Complex *a) |
|
653 { |
|
654 return (real (*a) < 0.0); |
|
655 } |
|
656 |
|
657 static int |
|
658 complex_select_dig (Complex *a) |
|
659 { |
|
660 return (abs (*a) < 1.0); |
|
661 } |
|
662 |
|
663 int |
|
664 ComplexSCHUR::init (const ComplexMatrix& a, const char *ord) |
|
665 { |
|
666 if (a.nr != a.nc) |
227
|
667 { |
|
668 (*current_liboctave_error_handler) |
|
669 ("ComplexSCHUR requires square matrix"); |
|
670 return -1; |
|
671 } |
3
|
672 |
|
673 char jobvs = 'V'; |
|
674 char sort; |
|
675 if (*ord == 'A' || *ord == 'D' || *ord == 'a' || *ord == 'd') |
|
676 sort = 'S'; |
|
677 else |
|
678 sort = 'N'; |
|
679 |
|
680 char sense = 'N'; |
|
681 |
|
682 int n = a.nc; |
|
683 int lwork = 8 * n; |
|
684 int info; |
|
685 int sdim; |
|
686 double rconde; |
|
687 double rcondv; |
|
688 |
|
689 double *rwork = new double [n]; |
|
690 |
|
691 // bwork is not referenced for non-ordered Schur. |
|
692 |
|
693 int *bwork = (int *) NULL; |
|
694 if (*ord == 'A' || *ord == 'D' || *ord == 'a' || *ord == 'd') |
|
695 bwork = new int [n]; |
|
696 |
|
697 Complex *s = dup(a.data,a.len); |
|
698 |
|
699 Complex *work = new Complex [lwork]; |
|
700 Complex *q = new Complex [n*n]; |
|
701 Complex *w = new Complex [n]; |
|
702 |
|
703 if (*ord == 'A' || *ord == 'a') |
|
704 { |
|
705 F77_FCN (zgeesx) (&jobvs, &sort, complex_select_ana, &sense, |
|
706 &n, s, &n, &sdim, w, q, &n, &rconde, &rcondv, |
|
707 work, &lwork, rwork, bwork, &info, 1L, 1L); |
|
708 } |
|
709 else if (*ord == 'D' || *ord == 'd') |
|
710 { |
|
711 F77_FCN (zgeesx) (&jobvs, &sort, complex_select_dig, &sense, |
|
712 &n, s, &n, &sdim, w, q, &n, &rconde, &rcondv, |
|
713 work, &lwork, rwork, bwork, &info, 1L, 1L); |
|
714 } |
|
715 else |
|
716 { |
|
717 F77_FCN (zgeesx) (&jobvs, &sort, dummy_select, &sense, &n, s, |
|
718 &n, &sdim, w, q, &n, &rconde, &rcondv, work, |
|
719 &lwork, rwork, bwork, &info, 1L, 1L); |
|
720 } |
|
721 |
|
722 schur_mat = ComplexMatrix (s,n,n); |
|
723 unitary_mat = ComplexMatrix (q,n,n); |
|
724 |
|
725 delete [] w; |
|
726 delete [] work; |
|
727 delete [] rwork; |
|
728 delete [] bwork; |
|
729 |
|
730 return info; |
|
731 } |
|
732 |
|
733 ostream& |
|
734 operator << (ostream& os, const SCHUR& a) |
|
735 { |
|
736 os << a.schur_matrix () << "\n"; |
|
737 os << a.unitary_matrix () << "\n"; |
|
738 |
|
739 return os; |
|
740 } |
|
741 |
|
742 /* |
|
743 * SVD stuff |
|
744 */ |
|
745 |
|
746 int |
|
747 SVD::init (const Matrix& a) |
|
748 { |
|
749 int info; |
|
750 |
|
751 int m = a.nr; |
|
752 int n = a.nc; |
|
753 |
|
754 char jobu = 'A'; |
|
755 char jobv = 'A'; |
|
756 |
|
757 double *tmp_data = dup (a.data, a.len); |
|
758 |
|
759 int min_mn = m < n ? m : n; |
|
760 int max_mn = m > n ? m : n; |
|
761 |
|
762 double *u = new double[m*m]; |
|
763 double *s_vec = new double[min_mn]; |
|
764 double *vt = new double[n*n]; |
|
765 |
|
766 int tmp1 = 3*min_mn + max_mn; |
|
767 int tmp2 = 5*min_mn - 4; |
|
768 int lwork = tmp1 > tmp2 ? tmp1 : tmp2; |
|
769 double *work = new double[lwork]; |
|
770 |
|
771 F77_FCN (dgesvd) (&jobu, &jobv, &m, &n, tmp_data, &m, s_vec, u, &m, |
|
772 vt, &n, work, &lwork, &info, 1L, 1L); |
|
773 |
|
774 left_sm = Matrix (u, m, m); |
|
775 sigma = DiagMatrix (s_vec, m, n); |
|
776 Matrix vt_m (vt, n, n); |
|
777 right_sm = Matrix (vt_m.transpose ()); |
|
778 |
|
779 delete [] tmp_data; |
|
780 delete [] work; |
|
781 |
|
782 return info; |
|
783 } |
|
784 |
|
785 ostream& |
|
786 operator << (ostream& os, const SVD& a) |
|
787 { |
|
788 os << a.left_singular_matrix () << "\n"; |
|
789 os << a.singular_values () << "\n"; |
|
790 os << a.right_singular_matrix () << "\n"; |
|
791 |
|
792 return os; |
|
793 } |
|
794 |
|
795 int |
|
796 ComplexSVD::init (const ComplexMatrix& a) |
|
797 { |
|
798 int info; |
|
799 |
|
800 int m = a.nr; |
|
801 int n = a.nc; |
|
802 |
|
803 char jobu = 'A'; |
|
804 char jobv = 'A'; |
|
805 |
|
806 Complex *tmp_data = dup (a.data, a.len); |
|
807 |
|
808 int min_mn = m < n ? m : n; |
|
809 int max_mn = m > n ? m : n; |
|
810 |
|
811 Complex *u = new Complex[m*m]; |
|
812 double *s_vec = new double[min_mn]; |
|
813 Complex *vt = new Complex[n*n]; |
|
814 |
|
815 int lwork = 2*min_mn + max_mn; |
|
816 Complex *work = new Complex[lwork]; |
|
817 |
|
818 int lrwork = 5*max_mn; |
|
819 double *rwork = new double[lrwork]; |
|
820 |
|
821 F77_FCN (zgesvd) (&jobu, &jobv, &m, &n, tmp_data, &m, s_vec, u, &m, |
|
822 vt, &n, work, &lwork, rwork, &info, 1L, 1L); |
|
823 |
|
824 left_sm = ComplexMatrix (u, m, m); |
|
825 sigma = DiagMatrix (s_vec, m, n); |
|
826 ComplexMatrix vt_m (vt, n, n); |
|
827 right_sm = ComplexMatrix (vt_m.hermitian ()); |
|
828 |
|
829 delete [] tmp_data; |
|
830 delete [] work; |
|
831 |
|
832 return info; |
|
833 } |
|
834 |
|
835 /* |
|
836 * EIG stuff. |
|
837 */ |
|
838 |
|
839 int |
|
840 EIG::init (const Matrix& a) |
|
841 { |
|
842 if (a.nr != a.nc) |
227
|
843 { |
|
844 (*current_liboctave_error_handler) ("EIG requires square matrix"); |
|
845 return -1; |
|
846 } |
3
|
847 |
|
848 int n = a.nr; |
|
849 |
|
850 int info; |
|
851 |
|
852 char jobvl = 'N'; |
|
853 char jobvr = 'V'; |
|
854 |
|
855 double *tmp_data = dup (a.data, a.len); |
|
856 double *wr = new double[n]; |
|
857 double *wi = new double[n]; |
|
858 Matrix vr (n, n); |
|
859 double *pvr = vr.fortran_vec (); |
|
860 int lwork = 8*n; |
|
861 double *work = new double[lwork]; |
|
862 |
|
863 double dummy; |
|
864 int idummy = 1; |
|
865 |
|
866 F77_FCN (dgeev) (&jobvl, &jobvr, &n, tmp_data, &n, wr, wi, &dummy, |
|
867 &idummy, pvr, &n, work, &lwork, &info, 1L, 1L); |
|
868 |
|
869 lambda.resize (n); |
|
870 v.resize (n, n); |
|
871 |
|
872 for (int j = 0; j < n; j++) |
|
873 { |
|
874 if (wi[j] == 0.0) |
|
875 { |
|
876 lambda.elem (j) = Complex (wr[j]); |
|
877 for (int i = 0; i < n; i++) |
|
878 v.elem (i, j) = vr.elem (i, j); |
|
879 } |
|
880 else |
|
881 { |
|
882 if (j+1 >= n) |
227
|
883 { |
|
884 (*current_liboctave_error_handler) ("EIG: internal error"); |
|
885 return -1; |
|
886 } |
3
|
887 |
|
888 for (int i = 0; i < n; i++) |
|
889 { |
|
890 lambda.elem (j) = Complex (wr[j], wi[j]); |
|
891 lambda.elem (j+1) = Complex (wr[j+1], wi[j+1]); |
|
892 double real_part = vr.elem (i, j); |
|
893 double imag_part = vr.elem (i, j+1); |
|
894 v.elem (i, j) = Complex (real_part, imag_part); |
|
895 v.elem (i, j+1) = Complex (real_part, -imag_part); |
|
896 } |
|
897 j++; |
|
898 } |
|
899 } |
|
900 |
|
901 delete [] tmp_data; |
|
902 delete [] wr; |
|
903 delete [] wi; |
|
904 delete [] work; |
|
905 |
|
906 return info; |
|
907 } |
|
908 |
|
909 int |
|
910 EIG::init (const ComplexMatrix& a) |
|
911 { |
|
912 |
|
913 if (a.nr != a.nc) |
227
|
914 { |
|
915 (*current_liboctave_error_handler) ("EIG requires square matrix"); |
|
916 return -1; |
|
917 } |
3
|
918 |
|
919 int n = a.nr; |
|
920 |
|
921 int info; |
|
922 |
|
923 char jobvl = 'N'; |
|
924 char jobvr = 'V'; |
|
925 |
|
926 lambda.resize (n); |
|
927 v.resize (n, n); |
|
928 |
|
929 Complex *pw = lambda.fortran_vec (); |
|
930 Complex *pvr = v.fortran_vec (); |
|
931 |
|
932 Complex *tmp_data = dup (a.data, a.len); |
|
933 |
|
934 int lwork = 8*n; |
|
935 Complex *work = new Complex[lwork]; |
|
936 double *rwork = new double[4*n]; |
|
937 |
|
938 Complex dummy; |
|
939 int idummy = 1; |
|
940 |
|
941 F77_FCN (zgeev) (&jobvl, &jobvr, &n, tmp_data, &n, pw, &dummy, |
|
942 &idummy, pvr, &n, work, &lwork, rwork, &info, 1L, |
|
943 1L); |
|
944 |
|
945 delete [] tmp_data; |
|
946 delete [] work; |
|
947 delete [] rwork; |
|
948 |
|
949 return info; |
|
950 } |
|
951 |
|
952 /* |
|
953 * LU stuff. |
|
954 */ |
|
955 |
|
956 LU::LU (const Matrix& a) |
|
957 { |
|
958 if (a.nr == 0 || a.nc == 0 || a.nr != a.nc) |
227
|
959 { |
|
960 (*current_liboctave_error_handler) ("LU requires square matrix"); |
|
961 return; |
|
962 } |
3
|
963 |
|
964 int n = a.nr; |
|
965 |
|
966 int *ipvt = new int [n]; |
|
967 int *pvt = new int [n]; |
|
968 double *tmp_data = dup (a.data, a.len); |
|
969 int info = 0; |
|
970 int zero = 0; |
|
971 double b; |
|
972 |
|
973 F77_FCN (dgesv) (&n, &zero, tmp_data, &n, ipvt, &b, &n, &info); |
|
974 |
|
975 Matrix A_fact (tmp_data, n, n); |
|
976 |
|
977 int i; |
|
978 |
|
979 for (i = 0; i < n; i++) |
|
980 { |
|
981 ipvt[i] -= 1; |
|
982 pvt[i] = i; |
|
983 } |
|
984 |
|
985 for (i = 0; i < n - 1; i++) |
|
986 { |
|
987 int k = ipvt[i]; |
|
988 if (k != i) |
|
989 { |
|
990 int tmp = pvt[k]; |
|
991 pvt[k] = pvt[i]; |
|
992 pvt[i] = tmp; |
|
993 } |
|
994 } |
|
995 |
|
996 l.resize (n, n, 0.0); |
|
997 u.resize (n, n, 0.0); |
|
998 p.resize (n, n, 0.0); |
|
999 |
|
1000 for (i = 0; i < n; i++) |
|
1001 { |
|
1002 p.elem (i, pvt[i]) = 1.0; |
|
1003 |
|
1004 int j; |
|
1005 |
|
1006 l.elem (i, i) = 1.0; |
|
1007 for (j = 0; j < i; j++) |
|
1008 l.elem (i, j) = A_fact.elem (i, j); |
|
1009 |
|
1010 for (j = i; j < n; j++) |
|
1011 u.elem (i, j) = A_fact.elem (i, j); |
|
1012 } |
|
1013 |
|
1014 delete [] ipvt; |
|
1015 delete [] pvt; |
|
1016 } |
|
1017 |
|
1018 ComplexLU::ComplexLU (const ComplexMatrix& a) |
|
1019 { |
|
1020 if (a.nr == 0 || a.nc == 0 || a.nr != a.nc) |
227
|
1021 { |
|
1022 (*current_liboctave_error_handler) ("ComplexLU requires square matrix"); |
|
1023 return; |
|
1024 } |
3
|
1025 |
|
1026 int n = a.nr; |
|
1027 |
|
1028 int *ipvt = new int [n]; |
|
1029 int *pvt = new int [n]; |
|
1030 Complex *tmp_data = dup (a.data, a.len); |
|
1031 int info = 0; |
|
1032 int zero = 0; |
|
1033 Complex b; |
|
1034 |
|
1035 F77_FCN (zgesv) (&n, &zero, tmp_data, &n, ipvt, &b, &n, &info); |
|
1036 |
|
1037 ComplexMatrix A_fact (tmp_data, n, n); |
|
1038 |
|
1039 int i; |
|
1040 |
|
1041 for (i = 0; i < n; i++) |
|
1042 { |
|
1043 ipvt[i] -= 1; |
|
1044 pvt[i] = i; |
|
1045 } |
|
1046 |
|
1047 for (i = 0; i < n - 1; i++) |
|
1048 { |
|
1049 int k = ipvt[i]; |
|
1050 if (k != i) |
|
1051 { |
|
1052 int tmp = pvt[k]; |
|
1053 pvt[k] = pvt[i]; |
|
1054 pvt[i] = tmp; |
|
1055 } |
|
1056 } |
|
1057 |
|
1058 l.resize (n, n, 0.0); |
|
1059 u.resize (n, n, 0.0); |
|
1060 p.resize (n, n, 0.0); |
|
1061 |
|
1062 for (i = 0; i < n; i++) |
|
1063 { |
|
1064 p.elem (i, pvt[i]) = 1.0; |
|
1065 |
|
1066 int j; |
|
1067 |
|
1068 l.elem (i, i) = 1.0; |
|
1069 for (j = 0; j < i; j++) |
|
1070 l.elem (i, j) = A_fact.elem (i, j); |
|
1071 |
|
1072 for (j = i; j < n; j++) |
|
1073 u.elem (i, j) = A_fact.elem (i, j); |
|
1074 } |
|
1075 |
|
1076 delete [] ipvt; |
|
1077 delete [] pvt; |
|
1078 } |
|
1079 |
|
1080 /* |
|
1081 * QR stuff. |
|
1082 */ |
|
1083 |
|
1084 QR::QR (const Matrix& a) |
|
1085 { |
|
1086 int m = a.nr; |
|
1087 int n = a.nc; |
|
1088 |
|
1089 if (m == 0 || n == 0) |
227
|
1090 { |
|
1091 (*current_liboctave_error_handler) ("QR must have non-empty matrix"); |
|
1092 return; |
|
1093 } |
3
|
1094 |
|
1095 double *tmp_data; |
|
1096 int min_mn = m < n ? m : n; |
|
1097 double *tau = new double[min_mn]; |
|
1098 int lwork = 32*n; |
|
1099 double *work = new double[lwork]; |
|
1100 int info = 0; |
|
1101 |
|
1102 if (m > n) |
|
1103 { |
|
1104 tmp_data = new double [m*m]; |
|
1105 copy (tmp_data, a.data, a.len); |
|
1106 } |
|
1107 else |
|
1108 tmp_data = dup (a.data, a.len); |
|
1109 |
|
1110 F77_FCN (dgeqrf) (&m, &n, tmp_data, &m, tau, work, &lwork, &info); |
|
1111 |
|
1112 delete [] work; |
|
1113 |
|
1114 r.resize (m, n, 0.0); |
|
1115 for (int j = 0; j < n; j++) |
|
1116 { |
|
1117 int limit = j < min_mn-1 ? j : min_mn-1; |
|
1118 for (int i = 0; i <= limit; i++) |
|
1119 r.elem (i, j) = tmp_data[m*j+i]; |
|
1120 } |
|
1121 |
|
1122 lwork = 32*m; |
|
1123 work = new double[lwork]; |
|
1124 |
|
1125 F77_FCN (dorgqr) (&m, &m, &min_mn, tmp_data, &m, tau, work, &lwork, &info); |
|
1126 |
|
1127 q = Matrix (tmp_data, m, m); |
|
1128 |
|
1129 delete [] tau; |
|
1130 delete [] work; |
|
1131 } |
|
1132 |
|
1133 ComplexQR::ComplexQR (const ComplexMatrix& a) |
|
1134 { |
|
1135 int m = a.nr; |
|
1136 int n = a.nc; |
|
1137 |
|
1138 if (m == 0 || n == 0) |
227
|
1139 { |
|
1140 (*current_liboctave_error_handler) |
|
1141 ("ComplexQR must have non-empty matrix"); |
|
1142 return; |
|
1143 } |
3
|
1144 |
|
1145 Complex *tmp_data; |
|
1146 int min_mn = m < n ? m : n; |
|
1147 Complex *tau = new Complex[min_mn]; |
|
1148 int lwork = 32*n; |
|
1149 Complex *work = new Complex[lwork]; |
|
1150 int info = 0; |
|
1151 |
|
1152 if (m > n) |
|
1153 { |
|
1154 tmp_data = new Complex [m*m]; |
|
1155 copy (tmp_data, a.data, a.len); |
|
1156 } |
|
1157 else |
|
1158 tmp_data = dup (a.data, a.len); |
|
1159 |
|
1160 F77_FCN (zgeqrf) (&m, &n, tmp_data, &m, tau, work, &lwork, &info); |
|
1161 |
|
1162 delete [] work; |
|
1163 |
|
1164 r.resize (m, n, 0.0); |
|
1165 for (int j = 0; j < n; j++) |
|
1166 { |
|
1167 int limit = j < min_mn-1 ? j : min_mn-1; |
|
1168 for (int i = 0; i <= limit; i++) |
|
1169 r.elem (i, j) = tmp_data[m*j+i]; |
|
1170 } |
|
1171 |
|
1172 lwork = 32*m; |
|
1173 work = new Complex[lwork]; |
|
1174 |
|
1175 F77_FCN (zungqr) (&m, &m, &min_mn, tmp_data, &m, tau, work, &lwork, &info); |
|
1176 |
|
1177 q = ComplexMatrix (tmp_data, m, m); |
|
1178 |
|
1179 delete [] tau; |
|
1180 delete [] work; |
|
1181 } |
|
1182 |
|
1183 /* |
|
1184 ;;; Local Variables: *** |
|
1185 ;;; mode: C++ *** |
|
1186 ;;; page-delimiter: "^/\\*" *** |
|
1187 ;;; End: *** |
|
1188 */ |