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