458
|
1 // Matrix manipulations. -*- C++ -*- |
|
2 /* |
|
3 |
1882
|
4 Copyright (C) 1996 John W. Eaton |
458
|
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 |
1315
|
20 Software Foundation, 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. |
458
|
21 |
|
22 */ |
|
23 |
1296
|
24 #if defined (__GNUG__) |
|
25 #pragma implementation |
|
26 #endif |
|
27 |
458
|
28 #ifdef HAVE_CONFIG_H |
1192
|
29 #include <config.h> |
458
|
30 #endif |
|
31 |
1367
|
32 #include <cfloat> |
|
33 #include <cstdio> |
1472
|
34 #include <cstring> |
1367
|
35 |
458
|
36 #include <iostream.h> |
1367
|
37 |
1472
|
38 #include <sys/types.h> // XXX FIXME XXX |
458
|
39 |
1819
|
40 #include "dbleAEPBAL.h" |
458
|
41 #include "dbleDET.h" |
1819
|
42 #include "dbleSCHUR.h" |
740
|
43 #include "dbleSVD.h" |
1847
|
44 #include "f77-fcn.h" |
458
|
45 #include "lo-error.h" |
1968
|
46 #include "lo-utils.h" |
1367
|
47 #include "mx-base.h" |
|
48 #include "mx-inlines.cc" |
1650
|
49 #include "oct-cmplx.h" |
458
|
50 |
|
51 // Fortran functions we call. |
|
52 |
|
53 extern "C" |
|
54 { |
1253
|
55 int F77_FCN (dgemm, DGEMM) (const char*, const char*, const int&, |
|
56 const int&, const int&, const double&, |
|
57 const double*, const int&, |
|
58 const double*, const int&, |
|
59 const double&, double*, const int&, |
|
60 long, long); |
458
|
61 |
1253
|
62 int F77_FCN (dgeco, DGECO) (double*, const int&, const int&, int*, |
|
63 double&, double*); |
458
|
64 |
1253
|
65 int F77_FCN (dgesl, DGESL) (const double*, const int&, const int&, |
|
66 const int*, double*, const int&); |
458
|
67 |
1253
|
68 int F77_FCN (dgedi, DGEDI) (double*, const int&, const int&, |
|
69 const int*, double*, double*, |
|
70 const int&); |
458
|
71 |
1253
|
72 int F77_FCN (dgelss, DGELSS) (const int&, const int&, const int&, |
|
73 double*, const int&, double*, |
|
74 const int&, double*, double&, int&, |
|
75 double*, const int&, int&); |
458
|
76 |
1360
|
77 // Note that the original complex fft routines were not written for |
|
78 // double complex arguments. They have been modified by adding an |
|
79 // implicit double precision (a-h,o-z) statement at the beginning of |
|
80 // each subroutine. |
458
|
81 |
1253
|
82 int F77_FCN (cffti, CFFTI) (const int&, Complex*); |
458
|
83 |
1253
|
84 int F77_FCN (cfftf, CFFTF) (const int&, Complex*, Complex*); |
458
|
85 |
1253
|
86 int F77_FCN (cfftb, CFFTB) (const int&, Complex*, Complex*); |
1819
|
87 |
|
88 int F77_FCN (dlartg, DLARTG) (const double&, const double&, double&, |
|
89 double&, double&); |
|
90 |
|
91 int F77_FCN (dtrsyl, DTRSYL) (const char*, const char*, const int&, |
|
92 const int&, const int&, const double*, |
|
93 const int&, const double*, const int&, |
|
94 const double*, const int&, double&, |
|
95 int&, long, long); |
|
96 |
|
97 double F77_FCN (dlange, DLANGE) (const char*, const int&, |
|
98 const int&, const double*, |
|
99 const int&, double*); |
1959
|
100 |
|
101 int F77_FCN (qzhes, QZHES) (const int&, const int&, double*, |
|
102 double*, const long&, double*); |
|
103 |
|
104 int F77_FCN (qzit, QZIT) (const int&, const int&, double*, double*, |
|
105 const double&, const long&, double*, |
|
106 int&); |
|
107 |
|
108 int F77_FCN (qzval, QZVAL) (const int&, const int&, double*, |
|
109 double*, double*, double*, double*, |
|
110 const long&, double*); |
458
|
111 } |
|
112 |
1360
|
113 // Matrix class. |
458
|
114 |
|
115 Matrix::Matrix (const DiagMatrix& a) |
1214
|
116 : MArray2<double> (a.rows (), a.cols (), 0.0) |
458
|
117 { |
|
118 for (int i = 0; i < a.length (); i++) |
|
119 elem (i, i) = a.elem (i, i); |
|
120 } |
|
121 |
1574
|
122 // XXX FIXME XXX -- could we use a templated mixed-type copy function |
|
123 // here? |
|
124 |
|
125 Matrix::Matrix (const charMatrix& a) |
|
126 : MArray2<double> (a.rows (), a.cols ()) |
|
127 { |
|
128 for (int i = 0; i < a.rows (); i++) |
|
129 for (int j = 0; j < a.cols (); j++) |
|
130 elem (i, j) = a.elem (i, j); |
|
131 } |
|
132 |
458
|
133 int |
|
134 Matrix::operator == (const Matrix& a) const |
|
135 { |
|
136 if (rows () != a.rows () || cols () != a.cols ()) |
|
137 return 0; |
|
138 |
|
139 return equal (data (), a.data (), length ()); |
|
140 } |
|
141 |
|
142 int |
|
143 Matrix::operator != (const Matrix& a) const |
|
144 { |
|
145 return !(*this == a); |
|
146 } |
|
147 |
|
148 Matrix& |
|
149 Matrix::insert (const Matrix& a, int r, int c) |
|
150 { |
1561
|
151 Array2<double>::insert (a, r, c); |
458
|
152 return *this; |
|
153 } |
|
154 |
|
155 Matrix& |
|
156 Matrix::insert (const RowVector& a, int r, int c) |
|
157 { |
|
158 int a_len = a.length (); |
1698
|
159 if (r < 0 || r >= rows () || c < 0 || c + a_len > cols ()) |
458
|
160 { |
|
161 (*current_liboctave_error_handler) ("range error for insert"); |
|
162 return *this; |
|
163 } |
|
164 |
|
165 for (int i = 0; i < a_len; i++) |
|
166 elem (r, c+i) = a.elem (i); |
|
167 |
|
168 return *this; |
|
169 } |
|
170 |
|
171 Matrix& |
|
172 Matrix::insert (const ColumnVector& a, int r, int c) |
|
173 { |
|
174 int a_len = a.length (); |
1698
|
175 if (r < 0 || r + a_len > rows () || c < 0 || c >= cols ()) |
458
|
176 { |
|
177 (*current_liboctave_error_handler) ("range error for insert"); |
|
178 return *this; |
|
179 } |
|
180 |
|
181 for (int i = 0; i < a_len; i++) |
|
182 elem (r+i, c) = a.elem (i); |
|
183 |
|
184 return *this; |
|
185 } |
|
186 |
|
187 Matrix& |
|
188 Matrix::insert (const DiagMatrix& a, int r, int c) |
|
189 { |
1697
|
190 int a_nr = a.rows (); |
|
191 int a_nc = a.cols (); |
|
192 |
1698
|
193 if (r < 0 || r + a_nr > rows () || c < 0 || c + a_nc > cols ()) |
458
|
194 { |
|
195 (*current_liboctave_error_handler) ("range error for insert"); |
|
196 return *this; |
|
197 } |
|
198 |
1697
|
199 fill (0.0, r, c, r + a_nr - 1, c + a_nc - 1); |
|
200 |
458
|
201 for (int i = 0; i < a.length (); i++) |
|
202 elem (r+i, c+i) = a.elem (i, i); |
|
203 |
|
204 return *this; |
|
205 } |
|
206 |
|
207 Matrix& |
|
208 Matrix::fill (double val) |
|
209 { |
|
210 int nr = rows (); |
|
211 int nc = cols (); |
|
212 if (nr > 0 && nc > 0) |
|
213 for (int j = 0; j < nc; j++) |
|
214 for (int i = 0; i < nr; i++) |
|
215 elem (i, j) = val; |
|
216 |
|
217 return *this; |
|
218 } |
|
219 |
|
220 Matrix& |
|
221 Matrix::fill (double val, int r1, int c1, int r2, int c2) |
|
222 { |
|
223 int nr = rows (); |
|
224 int nc = cols (); |
|
225 if (r1 < 0 || r2 < 0 || c1 < 0 || c2 < 0 |
|
226 || r1 >= nr || r2 >= nr || c1 >= nc || c2 >= nc) |
|
227 { |
|
228 (*current_liboctave_error_handler) ("range error for fill"); |
|
229 return *this; |
|
230 } |
|
231 |
|
232 if (r1 > r2) { int tmp = r1; r1 = r2; r2 = tmp; } |
|
233 if (c1 > c2) { int tmp = c1; c1 = c2; c2 = tmp; } |
|
234 |
|
235 for (int j = c1; j <= c2; j++) |
|
236 for (int i = r1; i <= r2; i++) |
|
237 elem (i, j) = val; |
|
238 |
|
239 return *this; |
|
240 } |
|
241 |
|
242 Matrix |
|
243 Matrix::append (const Matrix& a) const |
|
244 { |
|
245 int nr = rows (); |
|
246 int nc = cols (); |
|
247 if (nr != a.rows ()) |
|
248 { |
|
249 (*current_liboctave_error_handler) ("row dimension mismatch for append"); |
|
250 return Matrix (); |
|
251 } |
|
252 |
|
253 int nc_insert = nc; |
|
254 Matrix retval (nr, nc + a.cols ()); |
|
255 retval.insert (*this, 0, 0); |
|
256 retval.insert (a, 0, nc_insert); |
|
257 return retval; |
|
258 } |
|
259 |
|
260 Matrix |
|
261 Matrix::append (const RowVector& a) const |
|
262 { |
|
263 int nr = rows (); |
|
264 int nc = cols (); |
|
265 if (nr != 1) |
|
266 { |
|
267 (*current_liboctave_error_handler) ("row dimension mismatch for append"); |
|
268 return Matrix (); |
|
269 } |
|
270 |
|
271 int nc_insert = nc; |
|
272 Matrix retval (nr, nc + a.length ()); |
|
273 retval.insert (*this, 0, 0); |
|
274 retval.insert (a, 0, nc_insert); |
|
275 return retval; |
|
276 } |
|
277 |
|
278 Matrix |
|
279 Matrix::append (const ColumnVector& a) const |
|
280 { |
|
281 int nr = rows (); |
|
282 int nc = cols (); |
|
283 if (nr != a.length ()) |
|
284 { |
|
285 (*current_liboctave_error_handler) ("row dimension mismatch for append"); |
|
286 return Matrix (); |
|
287 } |
|
288 |
|
289 int nc_insert = nc; |
|
290 Matrix retval (nr, nc + 1); |
|
291 retval.insert (*this, 0, 0); |
|
292 retval.insert (a, 0, nc_insert); |
|
293 return retval; |
|
294 } |
|
295 |
|
296 Matrix |
|
297 Matrix::append (const DiagMatrix& a) const |
|
298 { |
|
299 int nr = rows (); |
|
300 int nc = cols (); |
|
301 if (nr != a.rows ()) |
|
302 { |
|
303 (*current_liboctave_error_handler) ("row dimension mismatch for append"); |
|
304 return *this; |
|
305 } |
|
306 |
|
307 int nc_insert = nc; |
|
308 Matrix retval (nr, nc + a.cols ()); |
|
309 retval.insert (*this, 0, 0); |
|
310 retval.insert (a, 0, nc_insert); |
|
311 return retval; |
|
312 } |
|
313 |
|
314 Matrix |
|
315 Matrix::stack (const Matrix& a) const |
|
316 { |
|
317 int nr = rows (); |
|
318 int nc = cols (); |
|
319 if (nc != a.cols ()) |
|
320 { |
|
321 (*current_liboctave_error_handler) |
|
322 ("column dimension mismatch for stack"); |
|
323 return Matrix (); |
|
324 } |
|
325 |
|
326 int nr_insert = nr; |
|
327 Matrix retval (nr + a.rows (), nc); |
|
328 retval.insert (*this, 0, 0); |
|
329 retval.insert (a, nr_insert, 0); |
|
330 return retval; |
|
331 } |
|
332 |
|
333 Matrix |
|
334 Matrix::stack (const RowVector& a) const |
|
335 { |
|
336 int nr = rows (); |
|
337 int nc = cols (); |
|
338 if (nc != a.length ()) |
|
339 { |
|
340 (*current_liboctave_error_handler) |
|
341 ("column dimension mismatch for stack"); |
|
342 return Matrix (); |
|
343 } |
|
344 |
|
345 int nr_insert = nr; |
|
346 Matrix retval (nr + 1, nc); |
|
347 retval.insert (*this, 0, 0); |
|
348 retval.insert (a, nr_insert, 0); |
|
349 return retval; |
|
350 } |
|
351 |
|
352 Matrix |
|
353 Matrix::stack (const ColumnVector& a) const |
|
354 { |
|
355 int nr = rows (); |
|
356 int nc = cols (); |
|
357 if (nc != 1) |
|
358 { |
|
359 (*current_liboctave_error_handler) |
|
360 ("column dimension mismatch for stack"); |
|
361 return Matrix (); |
|
362 } |
|
363 |
|
364 int nr_insert = nr; |
|
365 Matrix retval (nr + a.length (), nc); |
|
366 retval.insert (*this, 0, 0); |
|
367 retval.insert (a, nr_insert, 0); |
|
368 return retval; |
|
369 } |
|
370 |
|
371 Matrix |
|
372 Matrix::stack (const DiagMatrix& a) const |
|
373 { |
|
374 int nr = rows (); |
|
375 int nc = cols (); |
|
376 if (nc != a.cols ()) |
|
377 { |
|
378 (*current_liboctave_error_handler) |
|
379 ("column dimension mismatch for stack"); |
|
380 return Matrix (); |
|
381 } |
|
382 |
|
383 int nr_insert = nr; |
|
384 Matrix retval (nr + a.rows (), nc); |
|
385 retval.insert (*this, 0, 0); |
|
386 retval.insert (a, nr_insert, 0); |
|
387 return retval; |
|
388 } |
|
389 |
|
390 Matrix |
|
391 Matrix::transpose (void) const |
|
392 { |
|
393 int nr = rows (); |
|
394 int nc = cols (); |
|
395 Matrix result (nc, nr); |
|
396 if (length () > 0) |
|
397 { |
|
398 for (int j = 0; j < nc; j++) |
|
399 for (int i = 0; i < nr; i++) |
|
400 result.elem (j, i) = elem (i, j); |
|
401 } |
|
402 return result; |
|
403 } |
|
404 |
|
405 Matrix |
1205
|
406 real (const ComplexMatrix& a) |
|
407 { |
|
408 int a_len = a.length (); |
|
409 Matrix retval; |
|
410 if (a_len > 0) |
|
411 retval = Matrix (real_dup (a.data (), a_len), a.rows (), a.cols ()); |
|
412 return retval; |
|
413 } |
|
414 |
|
415 Matrix |
|
416 imag (const ComplexMatrix& a) |
|
417 { |
|
418 int a_len = a.length (); |
|
419 Matrix retval; |
|
420 if (a_len > 0) |
|
421 retval = Matrix (imag_dup (a.data (), a_len), a.rows (), a.cols ()); |
|
422 return retval; |
|
423 } |
|
424 |
|
425 Matrix |
458
|
426 Matrix::extract (int r1, int c1, int r2, int c2) const |
|
427 { |
|
428 if (r1 > r2) { int tmp = r1; r1 = r2; r2 = tmp; } |
|
429 if (c1 > c2) { int tmp = c1; c1 = c2; c2 = tmp; } |
|
430 |
|
431 int new_r = r2 - r1 + 1; |
|
432 int new_c = c2 - c1 + 1; |
|
433 |
|
434 Matrix result (new_r, new_c); |
|
435 |
|
436 for (int j = 0; j < new_c; j++) |
|
437 for (int i = 0; i < new_r; i++) |
|
438 result.elem (i, j) = elem (r1+i, c1+j); |
|
439 |
|
440 return result; |
|
441 } |
|
442 |
|
443 // extract row or column i. |
|
444 |
|
445 RowVector |
|
446 Matrix::row (int i) const |
|
447 { |
|
448 int nc = cols (); |
|
449 if (i < 0 || i >= rows ()) |
|
450 { |
|
451 (*current_liboctave_error_handler) ("invalid row selection"); |
|
452 return RowVector (); |
|
453 } |
|
454 |
|
455 RowVector retval (nc); |
|
456 for (int j = 0; j < nc; j++) |
|
457 retval.elem (j) = elem (i, j); |
|
458 |
|
459 return retval; |
|
460 } |
|
461 |
|
462 RowVector |
|
463 Matrix::row (char *s) const |
|
464 { |
533
|
465 if (! s) |
458
|
466 { |
|
467 (*current_liboctave_error_handler) ("invalid row selection"); |
|
468 return RowVector (); |
|
469 } |
|
470 |
|
471 char c = *s; |
|
472 if (c == 'f' || c == 'F') |
|
473 return row (0); |
|
474 else if (c == 'l' || c == 'L') |
|
475 return row (rows () - 1); |
|
476 else |
|
477 { |
|
478 (*current_liboctave_error_handler) ("invalid row selection"); |
|
479 return RowVector (); |
|
480 } |
|
481 } |
|
482 |
|
483 ColumnVector |
|
484 Matrix::column (int i) const |
|
485 { |
|
486 int nr = rows (); |
|
487 if (i < 0 || i >= cols ()) |
|
488 { |
|
489 (*current_liboctave_error_handler) ("invalid column selection"); |
|
490 return ColumnVector (); |
|
491 } |
|
492 |
|
493 ColumnVector retval (nr); |
|
494 for (int j = 0; j < nr; j++) |
|
495 retval.elem (j) = elem (j, i); |
|
496 |
|
497 return retval; |
|
498 } |
|
499 |
|
500 ColumnVector |
|
501 Matrix::column (char *s) const |
|
502 { |
533
|
503 if (! s) |
458
|
504 { |
|
505 (*current_liboctave_error_handler) ("invalid column selection"); |
|
506 return ColumnVector (); |
|
507 } |
|
508 |
|
509 char c = *s; |
|
510 if (c == 'f' || c == 'F') |
|
511 return column (0); |
|
512 else if (c == 'l' || c == 'L') |
|
513 return column (cols () - 1); |
|
514 else |
|
515 { |
|
516 (*current_liboctave_error_handler) ("invalid column selection"); |
|
517 return ColumnVector (); |
|
518 } |
|
519 } |
|
520 |
|
521 Matrix |
|
522 Matrix::inverse (void) const |
|
523 { |
|
524 int info; |
|
525 double rcond; |
|
526 return inverse (info, rcond); |
|
527 } |
|
528 |
|
529 Matrix |
|
530 Matrix::inverse (int& info) const |
|
531 { |
|
532 double rcond; |
|
533 return inverse (info, rcond); |
|
534 } |
|
535 |
|
536 Matrix |
1656
|
537 Matrix::inverse (int& info, double& rcond, int force) const |
458
|
538 { |
1948
|
539 Matrix retval; |
|
540 |
458
|
541 int nr = rows (); |
|
542 int nc = cols (); |
1948
|
543 |
458
|
544 if (nr != nc || nr == 0 || nc == 0) |
1948
|
545 (*current_liboctave_error_handler) ("inverse requires square matrix"); |
458
|
546 else |
|
547 { |
1948
|
548 info = 0; |
|
549 |
|
550 Array<int> ipvt (nr); |
|
551 int *pipvt = ipvt.fortran_vec (); |
|
552 |
|
553 Array<double> z (nr); |
|
554 double *pz = z.fortran_vec (); |
|
555 |
|
556 retval = *this; |
|
557 double *tmp_data = retval.fortran_vec (); |
|
558 |
|
559 F77_XFCN (dgeco, DGECO, (tmp_data, nr, nc, pipvt, rcond, pz)); |
|
560 |
|
561 if (f77_exception_encountered) |
|
562 (*current_liboctave_error_handler) ("unrecoverable error in dgeco"); |
|
563 else |
|
564 { |
|
565 volatile double rcond_plus_one = rcond + 1.0; |
|
566 |
|
567 if (rcond_plus_one == 1.0) |
|
568 info = -1; |
|
569 |
|
570 if (info == -1 && ! force) |
|
571 retval = *this; // Restore matrix contents. |
|
572 else |
|
573 { |
|
574 double *dummy = 0; |
|
575 |
|
576 F77_XFCN (dgedi, DGEDI, (tmp_data, nr, nc, pipvt, dummy, |
|
577 pz, 1)); |
|
578 |
|
579 if (f77_exception_encountered) |
|
580 (*current_liboctave_error_handler) |
|
581 ("unrecoverable error in dgedi"); |
|
582 } |
|
583 } |
458
|
584 } |
|
585 |
1948
|
586 return retval; |
458
|
587 } |
|
588 |
740
|
589 Matrix |
|
590 Matrix::pseudo_inverse (double tol) |
|
591 { |
|
592 SVD result (*this); |
|
593 |
|
594 DiagMatrix S = result.singular_values (); |
|
595 Matrix U = result.left_singular_matrix (); |
|
596 Matrix V = result.right_singular_matrix (); |
|
597 |
|
598 ColumnVector sigma = S.diag (); |
|
599 |
|
600 int r = sigma.length () - 1; |
|
601 int nr = rows (); |
|
602 int nc = cols (); |
|
603 |
|
604 if (tol <= 0.0) |
|
605 { |
|
606 if (nr > nc) |
|
607 tol = nr * sigma.elem (0) * DBL_EPSILON; |
|
608 else |
|
609 tol = nc * sigma.elem (0) * DBL_EPSILON; |
|
610 } |
|
611 |
|
612 while (r >= 0 && sigma.elem (r) < tol) |
|
613 r--; |
|
614 |
|
615 if (r < 0) |
|
616 return Matrix (nc, nr, 0.0); |
|
617 else |
|
618 { |
|
619 Matrix Ur = U.extract (0, 0, nr-1, r); |
|
620 DiagMatrix D = DiagMatrix (sigma.extract (0, r)) . inverse (); |
|
621 Matrix Vr = V.extract (0, 0, nc-1, r); |
|
622 return Vr * D * Ur.transpose (); |
|
623 } |
|
624 } |
|
625 |
458
|
626 ComplexMatrix |
|
627 Matrix::fourier (void) const |
|
628 { |
1948
|
629 ComplexMatrix retval; |
|
630 |
458
|
631 int nr = rows (); |
|
632 int nc = cols (); |
1948
|
633 |
458
|
634 int npts, nsamples; |
1948
|
635 |
458
|
636 if (nr == 1 || nc == 1) |
|
637 { |
|
638 npts = nr > nc ? nr : nc; |
|
639 nsamples = 1; |
|
640 } |
|
641 else |
|
642 { |
|
643 npts = nr; |
|
644 nsamples = nc; |
|
645 } |
|
646 |
|
647 int nn = 4*npts+15; |
1948
|
648 |
|
649 Array<Complex> wsave (nn); |
|
650 Complex *pwsave = wsave.fortran_vec (); |
|
651 |
|
652 retval = *this; |
|
653 Complex *tmp_data = retval.fortran_vec (); |
|
654 |
|
655 F77_FCN (cffti, CFFTI) (npts, pwsave); |
458
|
656 |
|
657 for (int j = 0; j < nsamples; j++) |
1948
|
658 F77_FCN (cfftf, CFFTF) (npts, &tmp_data[npts*j], pwsave); |
|
659 |
|
660 return retval; |
458
|
661 } |
|
662 |
|
663 ComplexMatrix |
|
664 Matrix::ifourier (void) const |
|
665 { |
1948
|
666 ComplexMatrix retval; |
|
667 |
458
|
668 int nr = rows (); |
|
669 int nc = cols (); |
1948
|
670 |
458
|
671 int npts, nsamples; |
1948
|
672 |
458
|
673 if (nr == 1 || nc == 1) |
|
674 { |
|
675 npts = nr > nc ? nr : nc; |
|
676 nsamples = 1; |
|
677 } |
|
678 else |
|
679 { |
|
680 npts = nr; |
|
681 nsamples = nc; |
|
682 } |
|
683 |
|
684 int nn = 4*npts+15; |
1948
|
685 |
|
686 Array<Complex> wsave (nn); |
|
687 Complex *pwsave = wsave.fortran_vec (); |
|
688 |
|
689 retval = *this; |
|
690 Complex *tmp_data = retval.fortran_vec (); |
|
691 |
|
692 F77_FCN (cffti, CFFTI) (npts, pwsave); |
458
|
693 |
|
694 for (int j = 0; j < nsamples; j++) |
1948
|
695 F77_FCN (cfftb, CFFTB) (npts, &tmp_data[npts*j], pwsave); |
458
|
696 |
1321
|
697 for (int j = 0; j < npts*nsamples; j++) |
458
|
698 tmp_data[j] = tmp_data[j] / (double) npts; |
|
699 |
1948
|
700 return retval; |
458
|
701 } |
|
702 |
677
|
703 ComplexMatrix |
|
704 Matrix::fourier2d (void) const |
|
705 { |
1948
|
706 ComplexMatrix retval; |
|
707 |
677
|
708 int nr = rows (); |
|
709 int nc = cols (); |
1948
|
710 |
677
|
711 int npts, nsamples; |
1948
|
712 |
677
|
713 if (nr == 1 || nc == 1) |
|
714 { |
|
715 npts = nr > nc ? nr : nc; |
|
716 nsamples = 1; |
|
717 } |
|
718 else |
|
719 { |
|
720 npts = nr; |
|
721 nsamples = nc; |
|
722 } |
|
723 |
|
724 int nn = 4*npts+15; |
1948
|
725 |
|
726 Array<Complex> wsave (nn); |
|
727 Complex *pwsave = wsave.fortran_vec (); |
|
728 |
|
729 retval = *this; |
|
730 Complex *tmp_data = retval.fortran_vec (); |
|
731 |
|
732 F77_FCN (cffti, CFFTI) (npts, pwsave); |
677
|
733 |
|
734 for (int j = 0; j < nsamples; j++) |
1948
|
735 F77_FCN (cfftf, CFFTF) (npts, &tmp_data[npts*j], pwsave); |
677
|
736 |
|
737 npts = nc; |
|
738 nsamples = nr; |
|
739 nn = 4*npts+15; |
1948
|
740 |
|
741 wsave.resize (nn); |
|
742 pwsave = wsave.fortran_vec (); |
|
743 |
|
744 Array<Complex> row (npts); |
|
745 Complex *prow = row.fortran_vec (); |
|
746 |
|
747 F77_FCN (cffti, CFFTI) (npts, pwsave); |
677
|
748 |
1321
|
749 for (int j = 0; j < nsamples; j++) |
677
|
750 { |
|
751 for (int i = 0; i < npts; i++) |
1948
|
752 prow[i] = tmp_data[i*nr + j]; |
|
753 |
|
754 F77_FCN (cfftf, CFFTF) (npts, prow, pwsave); |
677
|
755 |
1321
|
756 for (int i = 0; i < npts; i++) |
1948
|
757 tmp_data[i*nr + j] = prow[i]; |
677
|
758 } |
|
759 |
1948
|
760 return retval; |
677
|
761 } |
|
762 |
|
763 ComplexMatrix |
|
764 Matrix::ifourier2d (void) const |
|
765 { |
1948
|
766 ComplexMatrix retval; |
|
767 |
677
|
768 int nr = rows (); |
|
769 int nc = cols (); |
1948
|
770 |
677
|
771 int npts, nsamples; |
1948
|
772 |
677
|
773 if (nr == 1 || nc == 1) |
|
774 { |
|
775 npts = nr > nc ? nr : nc; |
|
776 nsamples = 1; |
|
777 } |
|
778 else |
|
779 { |
|
780 npts = nr; |
|
781 nsamples = nc; |
|
782 } |
|
783 |
|
784 int nn = 4*npts+15; |
1948
|
785 |
|
786 Array<Complex> wsave (nn); |
|
787 Complex *pwsave = wsave.fortran_vec (); |
|
788 |
|
789 retval = *this; |
|
790 Complex *tmp_data = retval.fortran_vec (); |
|
791 |
|
792 F77_FCN (cffti, CFFTI) (npts, pwsave); |
677
|
793 |
|
794 for (int j = 0; j < nsamples; j++) |
1948
|
795 F77_FCN (cfftb, CFFTB) (npts, &tmp_data[npts*j], pwsave); |
677
|
796 |
1321
|
797 for (int j = 0; j < npts*nsamples; j++) |
677
|
798 tmp_data[j] = tmp_data[j] / (double) npts; |
|
799 |
|
800 npts = nc; |
|
801 nsamples = nr; |
|
802 nn = 4*npts+15; |
1948
|
803 |
|
804 wsave.resize (nn); |
|
805 pwsave = wsave.fortran_vec (); |
|
806 |
|
807 Array<Complex> row (npts); |
|
808 Complex *prow = row.fortran_vec (); |
|
809 |
|
810 F77_FCN (cffti, CFFTI) (npts, pwsave); |
677
|
811 |
1321
|
812 for (int j = 0; j < nsamples; j++) |
677
|
813 { |
|
814 for (int i = 0; i < npts; i++) |
1948
|
815 prow[i] = tmp_data[i*nr + j]; |
|
816 |
|
817 F77_FCN (cfftb, CFFTB) (npts, prow, pwsave); |
677
|
818 |
1321
|
819 for (int i = 0; i < npts; i++) |
1948
|
820 tmp_data[i*nr + j] = prow[i] / (double) npts; |
677
|
821 } |
|
822 |
1948
|
823 return retval; |
677
|
824 } |
|
825 |
458
|
826 DET |
|
827 Matrix::determinant (void) const |
|
828 { |
|
829 int info; |
|
830 double rcond; |
|
831 return determinant (info, rcond); |
|
832 } |
|
833 |
|
834 DET |
|
835 Matrix::determinant (int& info) const |
|
836 { |
|
837 double rcond; |
|
838 return determinant (info, rcond); |
|
839 } |
|
840 |
|
841 DET |
532
|
842 Matrix::determinant (int& info, double& rcond) const |
458
|
843 { |
|
844 DET retval; |
|
845 |
|
846 int nr = rows (); |
|
847 int nc = cols (); |
|
848 |
|
849 if (nr == 0 || nc == 0) |
|
850 { |
|
851 double d[2]; |
|
852 d[0] = 1.0; |
|
853 d[1] = 0.0; |
|
854 retval = DET (d); |
|
855 } |
|
856 else |
|
857 { |
|
858 info = 0; |
1948
|
859 |
|
860 Array<int> ipvt (nr); |
|
861 int *pipvt = ipvt.fortran_vec (); |
|
862 |
|
863 Array<double> z (nr); |
|
864 double *pz = z.fortran_vec (); |
|
865 |
|
866 Matrix atmp = *this; |
|
867 double *tmp_data = atmp.fortran_vec (); |
|
868 |
|
869 F77_XFCN (dgeco, DGECO, (tmp_data, nr, nr, pipvt, rcond, pz)); |
|
870 |
|
871 if (f77_exception_encountered) |
|
872 (*current_liboctave_error_handler) ("unrecoverable error in dgeco"); |
458
|
873 else |
|
874 { |
1948
|
875 volatile double rcond_plus_one = rcond + 1.0; |
|
876 |
|
877 if (rcond_plus_one == 1.0) |
|
878 { |
|
879 info = -1; |
|
880 retval = DET (); |
|
881 } |
|
882 else |
|
883 { |
|
884 double d[2]; |
|
885 |
|
886 F77_XFCN (dgedi, DGEDI, (tmp_data, nr, nr, pipvt, d, pz, 10)); |
|
887 |
|
888 if (f77_exception_encountered) |
|
889 (*current_liboctave_error_handler) |
|
890 ("unrecoverable error in dgedi"); |
|
891 else |
|
892 retval = DET (d); |
|
893 } |
458
|
894 } |
|
895 } |
|
896 |
|
897 return retval; |
|
898 } |
|
899 |
|
900 Matrix |
|
901 Matrix::solve (const Matrix& b) const |
|
902 { |
|
903 int info; |
|
904 double rcond; |
|
905 return solve (b, info, rcond); |
|
906 } |
|
907 |
|
908 Matrix |
|
909 Matrix::solve (const Matrix& b, int& info) const |
|
910 { |
|
911 double rcond; |
|
912 return solve (b, info, rcond); |
|
913 } |
|
914 |
|
915 Matrix |
532
|
916 Matrix::solve (const Matrix& b, int& info, double& rcond) const |
458
|
917 { |
|
918 Matrix retval; |
|
919 |
|
920 int nr = rows (); |
|
921 int nc = cols (); |
1948
|
922 |
458
|
923 if (nr == 0 || nc == 0 || nr != nc || nr != b.rows ()) |
1948
|
924 (*current_liboctave_error_handler) |
|
925 ("matrix dimension mismatch solution of linear equations"); |
458
|
926 else |
|
927 { |
1948
|
928 info = 0; |
|
929 |
|
930 Array<int> ipvt (nr); |
|
931 int *pipvt = ipvt.fortran_vec (); |
|
932 |
|
933 Array<double> z (nr); |
|
934 double *pz = z.fortran_vec (); |
|
935 |
|
936 Matrix atmp = *this; |
|
937 double *tmp_data = atmp.fortran_vec (); |
|
938 |
|
939 F77_XFCN (dgeco, DGECO, (tmp_data, nr, nr, pipvt, rcond, pz)); |
|
940 |
|
941 if (f77_exception_encountered) |
|
942 (*current_liboctave_error_handler) ("unrecoverable error in dgeco"); |
|
943 else |
|
944 { |
|
945 volatile double rcond_plus_one = rcond + 1.0; |
|
946 |
|
947 if (rcond_plus_one == 1.0) |
|
948 { |
|
949 info = -2; |
|
950 } |
|
951 else |
|
952 { |
|
953 retval = b; |
|
954 double *result = retval.fortran_vec (); |
|
955 |
|
956 int b_nc = b.cols (); |
|
957 |
|
958 for (volatile int j = 0; j < b_nc; j++) |
|
959 { |
|
960 F77_XFCN (dgesl, DGESL, (tmp_data, nr, nr, pipvt, |
|
961 &result[nr*j], 0)); |
|
962 |
|
963 if (f77_exception_encountered) |
|
964 { |
|
965 (*current_liboctave_error_handler) |
|
966 ("unrecoverable error in dgesl"); |
|
967 |
|
968 break; |
|
969 } |
|
970 } |
|
971 } |
|
972 } |
458
|
973 } |
|
974 |
|
975 return retval; |
|
976 } |
|
977 |
|
978 ComplexMatrix |
|
979 Matrix::solve (const ComplexMatrix& b) const |
|
980 { |
|
981 ComplexMatrix tmp (*this); |
|
982 return tmp.solve (b); |
|
983 } |
|
984 |
|
985 ComplexMatrix |
|
986 Matrix::solve (const ComplexMatrix& b, int& info) const |
|
987 { |
|
988 ComplexMatrix tmp (*this); |
|
989 return tmp.solve (b, info); |
|
990 } |
|
991 |
|
992 ComplexMatrix |
|
993 Matrix::solve (const ComplexMatrix& b, int& info, double& rcond) const |
|
994 { |
|
995 ComplexMatrix tmp (*this); |
|
996 return tmp.solve (b, info, rcond); |
|
997 } |
|
998 |
|
999 ColumnVector |
|
1000 Matrix::solve (const ColumnVector& b) const |
|
1001 { |
|
1002 int info; double rcond; |
|
1003 return solve (b, info, rcond); |
|
1004 } |
|
1005 |
|
1006 ColumnVector |
|
1007 Matrix::solve (const ColumnVector& b, int& info) const |
|
1008 { |
|
1009 double rcond; |
|
1010 return solve (b, info, rcond); |
|
1011 } |
|
1012 |
|
1013 ColumnVector |
532
|
1014 Matrix::solve (const ColumnVector& b, int& info, double& rcond) const |
458
|
1015 { |
|
1016 ColumnVector retval; |
|
1017 |
|
1018 int nr = rows (); |
|
1019 int nc = cols (); |
1948
|
1020 |
458
|
1021 if (nr == 0 || nc == 0 || nr != nc || nr != b.length ()) |
1948
|
1022 (*current_liboctave_error_handler) |
|
1023 ("matrix dimension mismatch solution of linear equations"); |
458
|
1024 else |
|
1025 { |
1948
|
1026 info = 0; |
|
1027 |
|
1028 Array<int> ipvt (nr); |
|
1029 int *pipvt = ipvt.fortran_vec (); |
|
1030 |
|
1031 Array<double> z (nr); |
|
1032 double *pz = z.fortran_vec (); |
|
1033 |
|
1034 Matrix atmp = *this; |
|
1035 double *tmp_data = atmp.fortran_vec (); |
|
1036 |
|
1037 F77_XFCN (dgeco, DGECO, (tmp_data, nr, nr, pipvt, rcond, pz)); |
|
1038 |
|
1039 if (f77_exception_encountered) |
|
1040 (*current_liboctave_error_handler) |
|
1041 ("unrecoverable error in dgeco"); |
|
1042 else |
|
1043 { |
|
1044 volatile double rcond_plus_one = rcond + 1.0; |
|
1045 |
|
1046 if (rcond_plus_one == 1.0) |
|
1047 { |
|
1048 info = -2; |
|
1049 } |
|
1050 else |
|
1051 { |
|
1052 retval = b; |
|
1053 double *result = retval.fortran_vec (); |
|
1054 |
|
1055 F77_XFCN (dgesl, DGESL, (tmp_data, nr, nr, pipvt, result, 0)); |
|
1056 |
|
1057 if (f77_exception_encountered) |
|
1058 (*current_liboctave_error_handler) |
|
1059 ("unrecoverable error in dgesl"); |
|
1060 } |
|
1061 } |
458
|
1062 } |
|
1063 |
|
1064 return retval; |
|
1065 } |
|
1066 |
|
1067 ComplexColumnVector |
|
1068 Matrix::solve (const ComplexColumnVector& b) const |
|
1069 { |
|
1070 ComplexMatrix tmp (*this); |
|
1071 return tmp.solve (b); |
|
1072 } |
|
1073 |
|
1074 ComplexColumnVector |
|
1075 Matrix::solve (const ComplexColumnVector& b, int& info) const |
|
1076 { |
|
1077 ComplexMatrix tmp (*this); |
|
1078 return tmp.solve (b, info); |
|
1079 } |
|
1080 |
|
1081 ComplexColumnVector |
|
1082 Matrix::solve (const ComplexColumnVector& b, int& info, double& rcond) const |
|
1083 { |
|
1084 ComplexMatrix tmp (*this); |
|
1085 return tmp.solve (b, info, rcond); |
|
1086 } |
|
1087 |
|
1088 Matrix |
|
1089 Matrix::lssolve (const Matrix& b) const |
|
1090 { |
|
1091 int info; |
|
1092 int rank; |
|
1093 return lssolve (b, info, rank); |
|
1094 } |
|
1095 |
|
1096 Matrix |
|
1097 Matrix::lssolve (const Matrix& b, int& info) const |
|
1098 { |
|
1099 int rank; |
|
1100 return lssolve (b, info, rank); |
|
1101 } |
|
1102 |
|
1103 Matrix |
|
1104 Matrix::lssolve (const Matrix& b, int& info, int& rank) const |
|
1105 { |
1948
|
1106 Matrix retval; |
|
1107 |
458
|
1108 int nrhs = b.cols (); |
|
1109 |
|
1110 int m = rows (); |
|
1111 int n = cols (); |
|
1112 |
|
1113 if (m == 0 || n == 0 || m != b.rows ()) |
1948
|
1114 (*current_liboctave_error_handler) |
|
1115 ("matrix dimension mismatch in solution of least squares problem"); |
|
1116 else |
458
|
1117 { |
1948
|
1118 Matrix atmp = *this; |
|
1119 double *tmp_data = atmp.fortran_vec (); |
|
1120 |
|
1121 int nrr = m > n ? m : n; |
|
1122 Matrix result (nrr, nrhs); |
|
1123 |
|
1124 for (int j = 0; j < nrhs; j++) |
|
1125 for (int i = 0; i < m; i++) |
|
1126 result.elem (i, j) = b.elem (i, j); |
|
1127 |
|
1128 double *presult = result.fortran_vec (); |
|
1129 |
|
1130 int len_s = m < n ? m : n; |
|
1131 Array<double> s (len_s); |
|
1132 double *ps = s.fortran_vec (); |
|
1133 |
|
1134 double rcond = -1.0; |
|
1135 |
|
1136 int lwork; |
|
1137 if (m < n) |
|
1138 lwork = 3*m + (2*m > nrhs |
|
1139 ? (2*m > n ? 2*m : n) |
|
1140 : (nrhs > n ? nrhs : n)); |
|
1141 else |
|
1142 lwork = 3*n + (2*n > nrhs |
|
1143 ? (2*n > m ? 2*n : m) |
|
1144 : (nrhs > m ? nrhs : m)); |
|
1145 |
|
1146 Array<double> work (lwork); |
|
1147 double *pwork = work.fortran_vec (); |
|
1148 |
|
1149 F77_XFCN (dgelss, DGELSS, (m, n, nrhs, tmp_data, m, presult, nrr, ps, |
|
1150 rcond, rank, pwork, lwork, info)); |
|
1151 |
|
1152 if (f77_exception_encountered) |
|
1153 (*current_liboctave_error_handler) ("unrecoverable error in dgelss"); |
|
1154 else |
|
1155 { |
|
1156 retval.resize (n, nrhs); |
|
1157 for (int j = 0; j < nrhs; j++) |
|
1158 for (int i = 0; i < n; i++) |
|
1159 retval.elem (i, j) = result.elem (i, j); |
|
1160 } |
458
|
1161 } |
|
1162 |
|
1163 return retval; |
|
1164 } |
|
1165 |
|
1166 ComplexMatrix |
|
1167 Matrix::lssolve (const ComplexMatrix& b) const |
|
1168 { |
|
1169 ComplexMatrix tmp (*this); |
1484
|
1170 int info; |
|
1171 int rank; |
|
1172 return tmp.lssolve (b, info, rank); |
458
|
1173 } |
|
1174 |
|
1175 ComplexMatrix |
|
1176 Matrix::lssolve (const ComplexMatrix& b, int& info) const |
|
1177 { |
|
1178 ComplexMatrix tmp (*this); |
1484
|
1179 int rank; |
|
1180 return tmp.lssolve (b, info, rank); |
458
|
1181 } |
|
1182 |
|
1183 ComplexMatrix |
|
1184 Matrix::lssolve (const ComplexMatrix& b, int& info, int& rank) const |
|
1185 { |
|
1186 ComplexMatrix tmp (*this); |
1484
|
1187 return tmp.lssolve (b, info, rank); |
458
|
1188 } |
|
1189 |
|
1190 ColumnVector |
|
1191 Matrix::lssolve (const ColumnVector& b) const |
|
1192 { |
|
1193 int info; |
1484
|
1194 int rank; |
|
1195 return lssolve (b, info, rank); |
458
|
1196 } |
|
1197 |
|
1198 ColumnVector |
|
1199 Matrix::lssolve (const ColumnVector& b, int& info) const |
|
1200 { |
|
1201 int rank; |
|
1202 return lssolve (b, info, rank); |
|
1203 } |
|
1204 |
|
1205 ColumnVector |
|
1206 Matrix::lssolve (const ColumnVector& b, int& info, int& rank) const |
|
1207 { |
1948
|
1208 ColumnVector retval; |
|
1209 |
458
|
1210 int nrhs = 1; |
|
1211 |
|
1212 int m = rows (); |
|
1213 int n = cols (); |
|
1214 |
|
1215 if (m == 0 || n == 0 || m != b.length ()) |
1948
|
1216 (*current_liboctave_error_handler) |
|
1217 ("matrix dimension mismatch in solution of least squares problem"); |
|
1218 else |
458
|
1219 { |
1948
|
1220 Matrix atmp = *this; |
|
1221 double *tmp_data = atmp.fortran_vec (); |
|
1222 |
|
1223 int nrr = m > n ? m : n; |
|
1224 ColumnVector result (nrr); |
|
1225 |
|
1226 for (int i = 0; i < m; i++) |
|
1227 result.elem (i) = b.elem (i); |
|
1228 |
|
1229 double *presult = result.fortran_vec (); |
|
1230 |
|
1231 int len_s = m < n ? m : n; |
|
1232 Array<double> s (len_s); |
|
1233 double *ps = s.fortran_vec (); |
|
1234 |
|
1235 double rcond = -1.0; |
|
1236 |
|
1237 int lwork; |
|
1238 if (m < n) |
|
1239 lwork = 3*m + (2*m > nrhs |
|
1240 ? (2*m > n ? 2*m : n) |
|
1241 : (nrhs > n ? nrhs : n)); |
|
1242 else |
|
1243 lwork = 3*n + (2*n > nrhs |
|
1244 ? (2*n > m ? 2*n : m) |
|
1245 : (nrhs > m ? nrhs : m)); |
|
1246 |
|
1247 Array<double> work (lwork); |
|
1248 double *pwork = work.fortran_vec (); |
|
1249 |
|
1250 F77_XFCN (dgelss, DGELSS, (m, n, nrhs, tmp_data, m, presult, nrr, |
|
1251 ps, rcond, rank, pwork, lwork, info)); |
|
1252 |
|
1253 if (f77_exception_encountered) |
|
1254 (*current_liboctave_error_handler) ("unrecoverable error in dgelss"); |
|
1255 else |
|
1256 { |
|
1257 retval.resize (n); |
|
1258 for (int i = 0; i < n; i++) |
|
1259 retval.elem (i) = result.elem (i); |
|
1260 } |
458
|
1261 } |
|
1262 |
|
1263 return retval; |
|
1264 } |
|
1265 |
|
1266 ComplexColumnVector |
|
1267 Matrix::lssolve (const ComplexColumnVector& b) const |
|
1268 { |
|
1269 ComplexMatrix tmp (*this); |
|
1270 return tmp.lssolve (b); |
|
1271 } |
|
1272 |
|
1273 ComplexColumnVector |
|
1274 Matrix::lssolve (const ComplexColumnVector& b, int& info) const |
|
1275 { |
|
1276 ComplexMatrix tmp (*this); |
|
1277 return tmp.lssolve (b, info); |
|
1278 } |
|
1279 |
|
1280 ComplexColumnVector |
|
1281 Matrix::lssolve (const ComplexColumnVector& b, int& info, int& rank) const |
|
1282 { |
|
1283 ComplexMatrix tmp (*this); |
|
1284 return tmp.lssolve (b, info, rank); |
|
1285 } |
|
1286 |
1819
|
1287 // Constants for matrix exponential calculation. |
|
1288 |
|
1289 static double padec [] = |
|
1290 { |
|
1291 5.0000000000000000e-1, |
|
1292 1.1666666666666667e-1, |
|
1293 1.6666666666666667e-2, |
|
1294 1.6025641025641026e-3, |
|
1295 1.0683760683760684e-4, |
|
1296 4.8562548562548563e-6, |
|
1297 1.3875013875013875e-7, |
|
1298 1.9270852604185938e-9, |
|
1299 }; |
|
1300 |
|
1301 Matrix |
|
1302 Matrix::expm (void) const |
|
1303 { |
|
1304 Matrix retval; |
|
1305 |
|
1306 Matrix m = *this; |
|
1307 |
|
1308 int nc = columns (); |
|
1309 |
|
1310 // trace shift value |
|
1311 double trshift = 0; |
|
1312 |
|
1313 // Preconditioning step 1: trace normalization. |
|
1314 |
|
1315 for (int i = 0; i < nc; i++) |
|
1316 trshift += m.elem (i, i); |
|
1317 |
|
1318 trshift /= nc; |
|
1319 |
|
1320 for (int i = 0; i < nc; i++) |
|
1321 m.elem (i, i) -= trshift; |
|
1322 |
|
1323 // Preconditioning step 2: balancing. |
|
1324 |
|
1325 AEPBALANCE mbal (m, "B"); |
|
1326 m = mbal.balanced_matrix (); |
|
1327 Matrix d = mbal.balancing_matrix (); |
|
1328 |
|
1329 // Preconditioning step 3: scaling. |
|
1330 |
|
1331 ColumnVector work(nc); |
|
1332 double inf_norm |
|
1333 = F77_FCN (dlange, DLANGE) ("I", nc, nc, m.fortran_vec (),nc, |
|
1334 work.fortran_vec ()); |
|
1335 |
|
1336 int sqpow = (int) (inf_norm > 0.0 |
|
1337 ? (1.0 + log (inf_norm) / log (2.0)) |
|
1338 : 0.0); |
|
1339 |
|
1340 // Check whether we need to square at all. |
|
1341 |
|
1342 if (sqpow < 0) |
|
1343 sqpow = 0; |
|
1344 |
|
1345 if (sqpow > 0) |
|
1346 { |
|
1347 double scale_factor = 1.0; |
|
1348 for (int i = 0; i < sqpow; i++) |
|
1349 scale_factor *= 2.0; |
|
1350 |
|
1351 m = m / scale_factor; |
|
1352 } |
|
1353 |
|
1354 // npp, dpp: pade' approx polynomial matrices. |
|
1355 |
|
1356 Matrix npp (nc, nc, 0.0); |
|
1357 Matrix dpp = npp; |
|
1358 |
|
1359 // Now powers a^8 ... a^1. |
|
1360 |
|
1361 int minus_one_j = -1; |
|
1362 for (int j = 7; j >= 0; j--) |
|
1363 { |
|
1364 npp = m * npp + m * padec[j]; |
|
1365 dpp = m * dpp + m * (minus_one_j * padec[j]); |
|
1366 minus_one_j *= -1; |
|
1367 } |
|
1368 |
|
1369 // Zero power. |
|
1370 |
|
1371 dpp = -dpp; |
|
1372 for(int j = 0; j < nc; j++) |
|
1373 { |
|
1374 npp.elem (j, j) += 1.0; |
|
1375 dpp.elem (j, j) += 1.0; |
|
1376 } |
|
1377 |
|
1378 // Compute pade approximation = inverse (dpp) * npp. |
|
1379 |
|
1380 retval = dpp.solve (npp); |
|
1381 |
|
1382 // Reverse preconditioning step 3: repeated squaring. |
|
1383 |
|
1384 while (sqpow) |
|
1385 { |
|
1386 retval = retval * retval; |
|
1387 sqpow--; |
|
1388 } |
|
1389 |
|
1390 // Reverse preconditioning step 2: inverse balancing. |
|
1391 |
|
1392 retval = retval.transpose(); |
|
1393 d = d.transpose (); |
|
1394 retval = retval * d; |
|
1395 retval = d.solve (retval); |
|
1396 retval = retval.transpose (); |
|
1397 |
|
1398 // Reverse preconditioning step 1: fix trace normalization. |
|
1399 |
|
1400 return retval * exp (trshift); |
|
1401 } |
|
1402 |
458
|
1403 Matrix& |
|
1404 Matrix::operator += (const Matrix& a) |
|
1405 { |
|
1406 int nr = rows (); |
|
1407 int nc = cols (); |
|
1408 if (nr != a.rows () || nc != a.cols ()) |
|
1409 { |
|
1410 (*current_liboctave_error_handler) |
|
1411 ("nonconformant matrix += operation attempted"); |
|
1412 return *this; |
|
1413 } |
|
1414 |
|
1415 if (nr == 0 || nc == 0) |
|
1416 return *this; |
|
1417 |
|
1418 double *d = fortran_vec (); // Ensures only one reference to my privates! |
|
1419 |
|
1420 add2 (d, a.data (), length ()); |
|
1421 |
|
1422 return *this; |
|
1423 } |
|
1424 |
|
1425 Matrix& |
|
1426 Matrix::operator -= (const Matrix& a) |
|
1427 { |
|
1428 int nr = rows (); |
|
1429 int nc = cols (); |
|
1430 if (nr != a.rows () || nc != a.cols ()) |
|
1431 { |
|
1432 (*current_liboctave_error_handler) |
|
1433 ("nonconformant matrix -= operation attempted"); |
|
1434 return *this; |
|
1435 } |
|
1436 |
|
1437 if (nr == 0 || nc == 0) |
|
1438 return *this; |
|
1439 |
|
1440 double *d = fortran_vec (); // Ensures only one reference to my privates! |
|
1441 |
|
1442 subtract2 (d, a.data (), length ()); |
|
1443 |
|
1444 return *this; |
|
1445 } |
|
1446 |
|
1447 Matrix& |
|
1448 Matrix::operator += (const DiagMatrix& a) |
|
1449 { |
|
1450 if (rows () != a.rows () || cols () != a.cols ()) |
|
1451 { |
|
1452 (*current_liboctave_error_handler) |
|
1453 ("nonconformant matrix += operation attempted"); |
|
1454 return *this; |
|
1455 } |
|
1456 |
|
1457 for (int i = 0; i < a.length (); i++) |
|
1458 elem (i, i) += a.elem (i, i); |
|
1459 |
|
1460 return *this; |
|
1461 } |
|
1462 |
|
1463 Matrix& |
|
1464 Matrix::operator -= (const DiagMatrix& a) |
|
1465 { |
|
1466 if (rows () != a.rows () || cols () != a.cols ()) |
|
1467 { |
|
1468 (*current_liboctave_error_handler) |
|
1469 ("nonconformant matrix += operation attempted"); |
|
1470 return *this; |
|
1471 } |
|
1472 |
|
1473 for (int i = 0; i < a.length (); i++) |
|
1474 elem (i, i) -= a.elem (i, i); |
|
1475 |
|
1476 return *this; |
|
1477 } |
|
1478 |
|
1479 // unary operations |
|
1480 |
|
1481 Matrix |
|
1482 Matrix::operator ! (void) const |
|
1483 { |
|
1484 int nr = rows (); |
|
1485 int nc = cols (); |
|
1486 |
|
1487 Matrix b (nr, nc); |
|
1488 |
|
1489 for (int j = 0; j < nc; j++) |
|
1490 for (int i = 0; i < nr; i++) |
|
1491 b.elem (i, j) = ! elem (i, j); |
|
1492 |
|
1493 return b; |
|
1494 } |
|
1495 |
1205
|
1496 // column vector by row vector -> matrix operations |
458
|
1497 |
1205
|
1498 Matrix |
|
1499 operator * (const ColumnVector& v, const RowVector& a) |
458
|
1500 { |
1948
|
1501 Matrix retval; |
|
1502 |
1205
|
1503 int len = v.length (); |
|
1504 int a_len = a.length (); |
1948
|
1505 |
1205
|
1506 if (len != a_len) |
1948
|
1507 (*current_liboctave_error_handler) |
|
1508 ("nonconformant vector multiplication attempted"); |
|
1509 else |
1205
|
1510 { |
1948
|
1511 if (len != 0) |
|
1512 { |
|
1513 retval.resize (len, a_len); |
|
1514 double *c = retval.fortran_vec (); |
|
1515 |
|
1516 F77_XFCN (dgemm, DGEMM, ("N", "N", len, a_len, 1, 1.0, |
|
1517 v.data (), len, a.data (), 1, 0.0, |
|
1518 c, len, 1L, 1L)); |
|
1519 |
|
1520 if (f77_exception_encountered) |
|
1521 (*current_liboctave_error_handler) |
|
1522 ("unrecoverable error in dgemm"); |
|
1523 } |
1205
|
1524 } |
458
|
1525 |
1948
|
1526 return retval; |
458
|
1527 } |
|
1528 |
1205
|
1529 // diagonal matrix by scalar -> matrix operations |
|
1530 |
|
1531 Matrix |
|
1532 operator + (const DiagMatrix& a, double s) |
458
|
1533 { |
1205
|
1534 Matrix tmp (a.rows (), a.cols (), s); |
|
1535 return a + tmp; |
458
|
1536 } |
|
1537 |
1205
|
1538 Matrix |
|
1539 operator - (const DiagMatrix& a, double s) |
458
|
1540 { |
1205
|
1541 Matrix tmp (a.rows (), a.cols (), -s); |
|
1542 return a + tmp; |
458
|
1543 } |
|
1544 |
1205
|
1545 // scalar by diagonal matrix -> matrix operations |
|
1546 |
|
1547 Matrix |
|
1548 operator + (double s, const DiagMatrix& a) |
458
|
1549 { |
1205
|
1550 Matrix tmp (a.rows (), a.cols (), s); |
|
1551 return tmp + a; |
|
1552 } |
|
1553 |
|
1554 Matrix |
|
1555 operator - (double s, const DiagMatrix& a) |
|
1556 { |
|
1557 Matrix tmp (a.rows (), a.cols (), s); |
|
1558 return tmp - a; |
458
|
1559 } |
|
1560 |
|
1561 // matrix by diagonal matrix -> matrix operations |
|
1562 |
|
1563 Matrix |
|
1564 operator + (const Matrix& m, const DiagMatrix& a) |
|
1565 { |
|
1566 int nr = m.rows (); |
|
1567 int nc = m.cols (); |
|
1568 if (nr != a.rows () || nc != a.cols ()) |
|
1569 { |
|
1570 (*current_liboctave_error_handler) |
|
1571 ("nonconformant matrix addition attempted"); |
|
1572 return Matrix (); |
|
1573 } |
|
1574 |
|
1575 if (nr == 0 || nc == 0) |
|
1576 return Matrix (nr, nc); |
|
1577 |
|
1578 Matrix result (m); |
|
1579 int a_len = a.length (); |
|
1580 for (int i = 0; i < a_len; i++) |
|
1581 result.elem (i, i) += a.elem (i, i); |
|
1582 |
|
1583 return result; |
|
1584 } |
|
1585 |
|
1586 Matrix |
|
1587 operator - (const Matrix& m, const DiagMatrix& a) |
|
1588 { |
|
1589 int nr = m.rows (); |
|
1590 int nc = m.cols (); |
|
1591 if (nr != a.rows () || nc != a.cols ()) |
|
1592 { |
|
1593 (*current_liboctave_error_handler) |
|
1594 ("nonconformant matrix subtraction attempted"); |
|
1595 return Matrix (); |
|
1596 } |
|
1597 |
|
1598 if (nr == 0 || nc == 0) |
|
1599 return Matrix (nr, nc); |
|
1600 |
|
1601 Matrix result (m); |
|
1602 int a_len = a.length (); |
|
1603 for (int i = 0; i < a_len; i++) |
|
1604 result.elem (i, i) -= a.elem (i, i); |
|
1605 |
|
1606 return result; |
|
1607 } |
|
1608 |
|
1609 Matrix |
|
1610 operator * (const Matrix& m, const DiagMatrix& a) |
|
1611 { |
1948
|
1612 Matrix retval; |
|
1613 |
458
|
1614 int nr = m.rows (); |
|
1615 int nc = m.cols (); |
1948
|
1616 |
458
|
1617 int a_nr = a.rows (); |
|
1618 int a_nc = a.cols (); |
1948
|
1619 |
458
|
1620 if (nc != a_nr) |
1948
|
1621 (*current_liboctave_error_handler) |
|
1622 ("nonconformant matrix multiplication attempted"); |
|
1623 else |
458
|
1624 { |
1948
|
1625 if (nr == 0 || nc == 0 || a_nc == 0) |
|
1626 retval.resize (nr, a_nc, 0.0); |
458
|
1627 else |
|
1628 { |
1948
|
1629 retval.resize (nr, a_nc); |
|
1630 double *c = retval.fortran_vec (); |
|
1631 |
|
1632 double *ctmp = 0; |
|
1633 |
|
1634 int a_len = a.length (); |
|
1635 |
|
1636 for (int j = 0; j < a_len; j++) |
|
1637 { |
|
1638 int idx = j * nr; |
|
1639 ctmp = c + idx; |
|
1640 |
|
1641 if (a.elem (j, j) == 1.0) |
|
1642 { |
|
1643 for (int i = 0; i < nr; i++) |
|
1644 ctmp[i] = m.elem (i, j); |
|
1645 } |
|
1646 else if (a.elem (j, j) == 0.0) |
|
1647 { |
|
1648 for (int i = 0; i < nr; i++) |
|
1649 ctmp[i] = 0.0; |
|
1650 } |
|
1651 else |
|
1652 { |
|
1653 for (int i = 0; i < nr; i++) |
|
1654 ctmp[i] = a.elem (j, j) * m.elem (i, j); |
|
1655 } |
|
1656 } |
|
1657 |
|
1658 if (a_nr < a_nc) |
|
1659 { |
|
1660 for (int i = nr * nc; i < nr * a_nc; i++) |
|
1661 ctmp[i] = 0.0; |
|
1662 } |
458
|
1663 } |
|
1664 } |
|
1665 |
1948
|
1666 return retval; |
458
|
1667 } |
|
1668 |
1205
|
1669 // diagonal matrix by matrix -> matrix operations |
|
1670 |
|
1671 Matrix |
|
1672 operator + (const DiagMatrix& m, const Matrix& a) |
458
|
1673 { |
|
1674 int nr = m.rows (); |
|
1675 int nc = m.cols (); |
|
1676 if (nr != a.rows () || nc != a.cols ()) |
|
1677 { |
|
1678 (*current_liboctave_error_handler) |
|
1679 ("nonconformant matrix addition attempted"); |
1205
|
1680 return Matrix (); |
458
|
1681 } |
|
1682 |
|
1683 if (nr == 0 || nc == 0) |
1205
|
1684 return Matrix (nr, nc); |
458
|
1685 |
1205
|
1686 Matrix result (a); |
|
1687 for (int i = 0; i < m.length (); i++) |
|
1688 result.elem (i, i) += m.elem (i, i); |
458
|
1689 |
|
1690 return result; |
|
1691 } |
|
1692 |
1205
|
1693 Matrix |
|
1694 operator - (const DiagMatrix& m, const Matrix& a) |
458
|
1695 { |
|
1696 int nr = m.rows (); |
|
1697 int nc = m.cols (); |
|
1698 if (nr != a.rows () || nc != a.cols ()) |
|
1699 { |
|
1700 (*current_liboctave_error_handler) |
|
1701 ("nonconformant matrix subtraction attempted"); |
1205
|
1702 return Matrix (); |
458
|
1703 } |
|
1704 |
|
1705 if (nr == 0 || nc == 0) |
1205
|
1706 return Matrix (nr, nc); |
458
|
1707 |
1205
|
1708 Matrix result (-a); |
|
1709 for (int i = 0; i < m.length (); i++) |
|
1710 result.elem (i, i) += m.elem (i, i); |
458
|
1711 |
|
1712 return result; |
|
1713 } |
|
1714 |
1205
|
1715 Matrix |
|
1716 operator * (const DiagMatrix& m, const Matrix& a) |
458
|
1717 { |
|
1718 int nr = m.rows (); |
|
1719 int nc = m.cols (); |
|
1720 int a_nr = a.rows (); |
|
1721 int a_nc = a.cols (); |
|
1722 if (nc != a_nr) |
|
1723 { |
|
1724 (*current_liboctave_error_handler) |
|
1725 ("nonconformant matrix multiplication attempted"); |
1205
|
1726 return Matrix (); |
458
|
1727 } |
|
1728 |
|
1729 if (nr == 0 || nc == 0 || a_nc == 0) |
1205
|
1730 return Matrix (nr, a_nc, 0.0); |
458
|
1731 |
1205
|
1732 Matrix c (nr, a_nc); |
458
|
1733 |
1205
|
1734 for (int i = 0; i < m.length (); i++) |
458
|
1735 { |
1205
|
1736 if (m.elem (i, i) == 1.0) |
458
|
1737 { |
1205
|
1738 for (int j = 0; j < a_nc; j++) |
|
1739 c.elem (i, j) = a.elem (i, j); |
458
|
1740 } |
1205
|
1741 else if (m.elem (i, i) == 0.0) |
458
|
1742 { |
1205
|
1743 for (int j = 0; j < a_nc; j++) |
|
1744 c.elem (i, j) = 0.0; |
458
|
1745 } |
|
1746 else |
|
1747 { |
1205
|
1748 for (int j = 0; j < a_nc; j++) |
|
1749 c.elem (i, j) = m.elem (i, i) * a.elem (i, j); |
458
|
1750 } |
|
1751 } |
|
1752 |
1205
|
1753 if (nr > nc) |
458
|
1754 { |
1205
|
1755 for (int j = 0; j < a_nc; j++) |
|
1756 for (int i = a_nr; i < nr; i++) |
|
1757 c.elem (i, j) = 0.0; |
458
|
1758 } |
|
1759 |
1205
|
1760 return c; |
458
|
1761 } |
|
1762 |
|
1763 // matrix by matrix -> matrix operations |
|
1764 |
|
1765 Matrix |
|
1766 operator * (const Matrix& m, const Matrix& a) |
|
1767 { |
1948
|
1768 Matrix retval; |
|
1769 |
458
|
1770 int nr = m.rows (); |
|
1771 int nc = m.cols (); |
1948
|
1772 |
458
|
1773 int a_nr = a.rows (); |
|
1774 int a_nc = a.cols (); |
1948
|
1775 |
458
|
1776 if (nc != a_nr) |
1948
|
1777 (*current_liboctave_error_handler) |
|
1778 ("nonconformant matrix multiplication attempted"); |
|
1779 else |
458
|
1780 { |
1948
|
1781 if (nr == 0 || nc == 0 || a_nc == 0) |
|
1782 retval.resize (nr, a_nc, 0.0); |
|
1783 else |
|
1784 { |
|
1785 int ld = nr; |
|
1786 int lda = a_nr; |
|
1787 |
|
1788 retval.resize (nr, a_nc); |
|
1789 double *c = retval.fortran_vec (); |
|
1790 |
|
1791 F77_XFCN (dgemm, DGEMM, ("N", "N", nr, a_nc, nc, 1.0, |
|
1792 m.data (), ld, a.data (), lda, 0.0, |
|
1793 c, nr, 1L, 1L)); |
|
1794 |
|
1795 if (f77_exception_encountered) |
|
1796 (*current_liboctave_error_handler) |
|
1797 ("unrecoverable error in dgemm"); |
|
1798 } |
458
|
1799 } |
|
1800 |
1948
|
1801 return retval; |
458
|
1802 } |
|
1803 |
|
1804 // other operations. |
|
1805 |
|
1806 Matrix |
|
1807 map (d_d_Mapper f, const Matrix& a) |
|
1808 { |
|
1809 Matrix b (a); |
|
1810 b.map (f); |
|
1811 return b; |
|
1812 } |
|
1813 |
1205
|
1814 Matrix |
|
1815 map (d_c_Mapper f, const ComplexMatrix& a) |
|
1816 { |
|
1817 int a_nc = a.cols (); |
|
1818 int a_nr = a.rows (); |
|
1819 Matrix b (a_nr, a_nc); |
|
1820 for (int j = 0; j < a_nc; j++) |
|
1821 for (int i = 0; i < a_nr; i++) |
|
1822 b.elem (i, j) = f (a.elem (i, j)); |
|
1823 return b; |
|
1824 } |
|
1825 |
458
|
1826 void |
|
1827 Matrix::map (d_d_Mapper f) |
|
1828 { |
|
1829 double *d = fortran_vec (); // Ensures only one reference to my privates! |
|
1830 |
|
1831 for (int i = 0; i < length (); i++) |
|
1832 d[i] = f (d[i]); |
|
1833 } |
|
1834 |
1968
|
1835 // Return nonzero if any element of M is not an integer. Also extract |
|
1836 // the largest and smallest values and return them in MAX_VAL and MIN_VAL. |
|
1837 |
|
1838 int |
|
1839 Matrix::all_integers (double& max_val, double& min_val) const |
|
1840 { |
|
1841 int nr = rows (); |
|
1842 int nc = cols (); |
|
1843 |
|
1844 if (nr > 0 && nc > 0) |
|
1845 { |
|
1846 max_val = elem (0, 0); |
|
1847 min_val = elem (0, 0); |
|
1848 } |
|
1849 else |
|
1850 return 0; |
|
1851 |
|
1852 for (int j = 0; j < nc; j++) |
|
1853 for (int i = 0; i < nr; i++) |
|
1854 { |
|
1855 double val = elem (i, j); |
|
1856 |
|
1857 if (val > max_val) |
|
1858 max_val = val; |
|
1859 |
|
1860 if (val < min_val) |
|
1861 min_val = val; |
|
1862 |
|
1863 if (D_NINT (val) != val) |
|
1864 return 0; |
|
1865 } |
|
1866 return 1; |
|
1867 } |
|
1868 |
|
1869 int |
|
1870 Matrix::too_large_for_float (void) const |
|
1871 { |
|
1872 int nr = rows (); |
|
1873 int nc = columns (); |
|
1874 |
|
1875 for (int j = 0; j < nc; j++) |
|
1876 for (int i = 0; i < nr; i++) |
|
1877 { |
|
1878 double val = elem (i, j); |
|
1879 |
|
1880 if (val > FLT_MAX || val < FLT_MIN) |
|
1881 return 1; |
|
1882 } |
|
1883 |
|
1884 return 0; |
|
1885 } |
|
1886 |
458
|
1887 // XXX FIXME XXX Do these really belong here? They should maybe be |
|
1888 // cleaned up a bit, no? What about corresponding functions for the |
|
1889 // Vectors? |
|
1890 |
|
1891 Matrix |
|
1892 Matrix::all (void) const |
|
1893 { |
|
1894 int nr = rows (); |
|
1895 int nc = cols (); |
|
1896 Matrix retval; |
|
1897 if (nr > 0 && nc > 0) |
|
1898 { |
|
1899 if (nr == 1) |
|
1900 { |
|
1901 retval.resize (1, 1); |
|
1902 retval.elem (0, 0) = 1.0; |
|
1903 for (int j = 0; j < nc; j++) |
|
1904 { |
|
1905 if (elem (0, j) == 0.0) |
|
1906 { |
|
1907 retval.elem (0, 0) = 0.0; |
|
1908 break; |
|
1909 } |
|
1910 } |
|
1911 } |
|
1912 else if (nc == 1) |
|
1913 { |
|
1914 retval.resize (1, 1); |
|
1915 retval.elem (0, 0) = 1.0; |
|
1916 for (int i = 0; i < nr; i++) |
|
1917 { |
|
1918 if (elem (i, 0) == 0.0) |
|
1919 { |
|
1920 retval.elem (0, 0) = 0.0; |
|
1921 break; |
|
1922 } |
|
1923 } |
|
1924 } |
|
1925 else |
|
1926 { |
|
1927 retval.resize (1, nc); |
|
1928 for (int j = 0; j < nc; j++) |
|
1929 { |
|
1930 retval.elem (0, j) = 1.0; |
|
1931 for (int i = 0; i < nr; i++) |
|
1932 { |
|
1933 if (elem (i, j) == 0.0) |
|
1934 { |
|
1935 retval.elem (0, j) = 0.0; |
|
1936 break; |
|
1937 } |
|
1938 } |
|
1939 } |
|
1940 } |
|
1941 } |
|
1942 return retval; |
|
1943 } |
|
1944 |
|
1945 Matrix |
|
1946 Matrix::any (void) const |
|
1947 { |
|
1948 int nr = rows (); |
|
1949 int nc = cols (); |
|
1950 Matrix retval; |
|
1951 if (nr > 0 && nc > 0) |
|
1952 { |
|
1953 if (nr == 1) |
|
1954 { |
|
1955 retval.resize (1, 1); |
|
1956 retval.elem (0, 0) = 0.0; |
|
1957 for (int j = 0; j < nc; j++) |
|
1958 { |
|
1959 if (elem (0, j) != 0.0) |
|
1960 { |
|
1961 retval.elem (0, 0) = 1.0; |
|
1962 break; |
|
1963 } |
|
1964 } |
|
1965 } |
|
1966 else if (nc == 1) |
|
1967 { |
|
1968 retval.resize (1, 1); |
|
1969 retval.elem (0, 0) = 0.0; |
|
1970 for (int i = 0; i < nr; i++) |
|
1971 { |
|
1972 if (elem (i, 0) != 0.0) |
|
1973 { |
|
1974 retval.elem (0, 0) = 1.0; |
|
1975 break; |
|
1976 } |
|
1977 } |
|
1978 } |
|
1979 else |
|
1980 { |
|
1981 retval.resize (1, nc); |
|
1982 for (int j = 0; j < nc; j++) |
|
1983 { |
|
1984 retval.elem (0, j) = 0.0; |
|
1985 for (int i = 0; i < nr; i++) |
|
1986 { |
|
1987 if (elem (i, j) != 0.0) |
|
1988 { |
|
1989 retval.elem (0, j) = 1.0; |
|
1990 break; |
|
1991 } |
|
1992 } |
|
1993 } |
|
1994 } |
|
1995 } |
|
1996 return retval; |
|
1997 } |
|
1998 |
|
1999 Matrix |
|
2000 Matrix::cumprod (void) const |
|
2001 { |
|
2002 Matrix retval; |
|
2003 |
|
2004 int nr = rows (); |
|
2005 int nc = cols (); |
|
2006 |
|
2007 if (nr == 1) |
|
2008 { |
|
2009 retval.resize (1, nc); |
|
2010 if (nc > 0) |
|
2011 { |
|
2012 double prod = elem (0, 0); |
|
2013 for (int j = 0; j < nc; j++) |
|
2014 { |
|
2015 retval.elem (0, j) = prod; |
|
2016 if (j < nc - 1) |
|
2017 prod *= elem (0, j+1); |
|
2018 } |
|
2019 } |
|
2020 } |
|
2021 else if (nc == 1) |
|
2022 { |
|
2023 retval.resize (nr, 1); |
|
2024 if (nr > 0) |
|
2025 { |
|
2026 double prod = elem (0, 0); |
|
2027 for (int i = 0; i < nr; i++) |
|
2028 { |
|
2029 retval.elem (i, 0) = prod; |
|
2030 if (i < nr - 1) |
|
2031 prod *= elem (i+1, 0); |
|
2032 } |
|
2033 } |
|
2034 } |
|
2035 else |
|
2036 { |
|
2037 retval.resize (nr, nc); |
|
2038 if (nr > 0 && nc > 0) |
|
2039 { |
|
2040 for (int j = 0; j < nc; j++) |
|
2041 { |
|
2042 double prod = elem (0, j); |
|
2043 for (int i = 0; i < nr; i++) |
|
2044 { |
|
2045 retval.elem (i, j) = prod; |
|
2046 if (i < nr - 1) |
|
2047 prod *= elem (i+1, j); |
|
2048 } |
|
2049 } |
|
2050 } |
|
2051 } |
|
2052 return retval; |
|
2053 } |
|
2054 |
|
2055 Matrix |
|
2056 Matrix::cumsum (void) const |
|
2057 { |
|
2058 Matrix retval; |
|
2059 |
|
2060 int nr = rows (); |
|
2061 int nc = cols (); |
|
2062 |
|
2063 if (nr == 1) |
|
2064 { |
|
2065 retval.resize (1, nc); |
|
2066 if (nc > 0) |
|
2067 { |
|
2068 double sum = elem (0, 0); |
|
2069 for (int j = 0; j < nc; j++) |
|
2070 { |
|
2071 retval.elem (0, j) = sum; |
|
2072 if (j < nc - 1) |
|
2073 sum += elem (0, j+1); |
|
2074 } |
|
2075 } |
|
2076 } |
|
2077 else if (nc == 1) |
|
2078 { |
|
2079 retval.resize (nr, 1); |
|
2080 if (nr > 0) |
|
2081 { |
|
2082 double sum = elem (0, 0); |
|
2083 for (int i = 0; i < nr; i++) |
|
2084 { |
|
2085 retval.elem (i, 0) = sum; |
|
2086 if (i < nr - 1) |
|
2087 sum += elem (i+1, 0); |
|
2088 } |
|
2089 } |
|
2090 } |
|
2091 else |
|
2092 { |
|
2093 retval.resize (nr, nc); |
|
2094 if (nr > 0 && nc > 0) |
|
2095 { |
|
2096 for (int j = 0; j < nc; j++) |
|
2097 { |
|
2098 double sum = elem (0, j); |
|
2099 for (int i = 0; i < nr; i++) |
|
2100 { |
|
2101 retval.elem (i, j) = sum; |
|
2102 if (i < nr - 1) |
|
2103 sum += elem (i+1, j); |
|
2104 } |
|
2105 } |
|
2106 } |
|
2107 } |
|
2108 return retval; |
|
2109 } |
|
2110 |
|
2111 Matrix |
|
2112 Matrix::prod (void) const |
|
2113 { |
|
2114 Matrix retval; |
|
2115 |
|
2116 int nr = rows (); |
|
2117 int nc = cols (); |
|
2118 |
|
2119 if (nr == 1) |
|
2120 { |
|
2121 retval.resize (1, 1); |
|
2122 retval.elem (0, 0) = 1.0; |
|
2123 for (int j = 0; j < nc; j++) |
|
2124 retval.elem (0, 0) *= elem (0, j); |
|
2125 } |
|
2126 else if (nc == 1) |
|
2127 { |
|
2128 retval.resize (1, 1); |
|
2129 retval.elem (0, 0) = 1.0; |
|
2130 for (int i = 0; i < nr; i++) |
|
2131 retval.elem (0, 0) *= elem (i, 0); |
|
2132 } |
|
2133 else |
|
2134 { |
|
2135 if (nc == 0) |
|
2136 { |
|
2137 retval.resize (1, 1); |
|
2138 retval.elem (0, 0) = 1.0; |
|
2139 } |
|
2140 else |
|
2141 retval.resize (1, nc); |
|
2142 |
|
2143 for (int j = 0; j < nc; j++) |
|
2144 { |
|
2145 retval.elem (0, j) = 1.0; |
|
2146 for (int i = 0; i < nr; i++) |
|
2147 retval.elem (0, j) *= elem (i, j); |
|
2148 } |
|
2149 } |
|
2150 return retval; |
|
2151 } |
|
2152 |
|
2153 Matrix |
|
2154 Matrix::sum (void) const |
|
2155 { |
|
2156 Matrix retval; |
|
2157 |
|
2158 int nr = rows (); |
|
2159 int nc = cols (); |
|
2160 |
|
2161 if (nr == 1) |
|
2162 { |
|
2163 retval.resize (1, 1); |
|
2164 retval.elem (0, 0) = 0.0; |
|
2165 for (int j = 0; j < nc; j++) |
|
2166 retval.elem (0, 0) += elem (0, j); |
|
2167 } |
|
2168 else if (nc == 1) |
|
2169 { |
|
2170 retval.resize (1, 1); |
|
2171 retval.elem (0, 0) = 0.0; |
|
2172 for (int i = 0; i < nr; i++) |
|
2173 retval.elem (0, 0) += elem (i, 0); |
|
2174 } |
|
2175 else |
|
2176 { |
|
2177 if (nc == 0) |
|
2178 { |
|
2179 retval.resize (1, 1); |
|
2180 retval.elem (0, 0) = 0.0; |
|
2181 } |
|
2182 else |
|
2183 retval.resize (1, nc); |
|
2184 |
|
2185 for (int j = 0; j < nc; j++) |
|
2186 { |
|
2187 retval.elem (0, j) = 0.0; |
|
2188 for (int i = 0; i < nr; i++) |
|
2189 retval.elem (0, j) += elem (i, j); |
|
2190 } |
|
2191 } |
|
2192 return retval; |
|
2193 } |
|
2194 |
|
2195 Matrix |
|
2196 Matrix::sumsq (void) const |
|
2197 { |
|
2198 Matrix retval; |
|
2199 |
|
2200 int nr = rows (); |
|
2201 int nc = cols (); |
|
2202 |
|
2203 if (nr == 1) |
|
2204 { |
|
2205 retval.resize (1, 1); |
|
2206 retval.elem (0, 0) = 0.0; |
|
2207 for (int j = 0; j < nc; j++) |
|
2208 { |
|
2209 double d = elem (0, j); |
|
2210 retval.elem (0, 0) += d * d; |
|
2211 } |
|
2212 } |
|
2213 else if (nc == 1) |
|
2214 { |
|
2215 retval.resize (1, 1); |
|
2216 retval.elem (0, 0) = 0.0; |
|
2217 for (int i = 0; i < nr; i++) |
|
2218 { |
|
2219 double d = elem (i, 0); |
|
2220 retval.elem (0, 0) += d * d; |
|
2221 } |
|
2222 } |
|
2223 else |
|
2224 { |
|
2225 retval.resize (1, nc); |
|
2226 for (int j = 0; j < nc; j++) |
|
2227 { |
|
2228 retval.elem (0, j) = 0.0; |
|
2229 for (int i = 0; i < nr; i++) |
|
2230 { |
|
2231 double d = elem (i, j); |
|
2232 retval.elem (0, j) += d * d; |
|
2233 } |
|
2234 } |
|
2235 } |
|
2236 return retval; |
|
2237 } |
|
2238 |
|
2239 ColumnVector |
|
2240 Matrix::diag (void) const |
|
2241 { |
|
2242 return diag (0); |
|
2243 } |
|
2244 |
|
2245 ColumnVector |
|
2246 Matrix::diag (int k) const |
|
2247 { |
|
2248 int nnr = rows (); |
|
2249 int nnc = cols (); |
|
2250 if (k > 0) |
|
2251 nnc -= k; |
|
2252 else if (k < 0) |
|
2253 nnr += k; |
|
2254 |
|
2255 ColumnVector d; |
|
2256 |
|
2257 if (nnr > 0 && nnc > 0) |
|
2258 { |
|
2259 int ndiag = (nnr < nnc) ? nnr : nnc; |
|
2260 |
|
2261 d.resize (ndiag); |
|
2262 |
|
2263 if (k > 0) |
|
2264 { |
|
2265 for (int i = 0; i < ndiag; i++) |
|
2266 d.elem (i) = elem (i, i+k); |
|
2267 } |
|
2268 else if ( k < 0) |
|
2269 { |
|
2270 for (int i = 0; i < ndiag; i++) |
|
2271 d.elem (i) = elem (i-k, i); |
|
2272 } |
|
2273 else |
|
2274 { |
|
2275 for (int i = 0; i < ndiag; i++) |
|
2276 d.elem (i) = elem (i, i); |
|
2277 } |
|
2278 } |
|
2279 else |
|
2280 cerr << "diag: requested diagonal out of range\n"; |
|
2281 |
|
2282 return d; |
|
2283 } |
|
2284 |
|
2285 ColumnVector |
|
2286 Matrix::row_min (void) const |
|
2287 { |
|
2288 ColumnVector result; |
|
2289 |
|
2290 int nr = rows (); |
|
2291 int nc = cols (); |
|
2292 |
|
2293 if (nr > 0 && nc > 0) |
|
2294 { |
|
2295 result.resize (nr); |
|
2296 |
|
2297 for (int i = 0; i < nr; i++) |
|
2298 { |
|
2299 double res = elem (i, 0); |
|
2300 for (int j = 1; j < nc; j++) |
|
2301 if (elem (i, j) < res) |
|
2302 res = elem (i, j); |
|
2303 result.elem (i) = res; |
|
2304 } |
|
2305 } |
|
2306 |
|
2307 return result; |
|
2308 } |
|
2309 |
|
2310 ColumnVector |
|
2311 Matrix::row_min_loc (void) const |
|
2312 { |
|
2313 ColumnVector result; |
|
2314 |
|
2315 int nr = rows (); |
|
2316 int nc = cols (); |
|
2317 |
|
2318 if (nr > 0 && nc > 0) |
|
2319 { |
|
2320 result.resize (nr); |
|
2321 |
|
2322 for (int i = 0; i < nr; i++) |
|
2323 { |
|
2324 int res = 0; |
|
2325 for (int j = 0; j < nc; j++) |
|
2326 if (elem (i, j) < elem (i, res)) |
|
2327 res = j; |
|
2328 result.elem (i) = (double) (res + 1); |
|
2329 } |
|
2330 } |
|
2331 |
|
2332 return result; |
|
2333 } |
|
2334 |
|
2335 ColumnVector |
|
2336 Matrix::row_max (void) const |
|
2337 { |
|
2338 ColumnVector result; |
|
2339 |
|
2340 int nr = rows (); |
|
2341 int nc = cols (); |
|
2342 |
|
2343 if (nr > 0 && nc > 0) |
|
2344 { |
|
2345 result.resize (nr); |
|
2346 |
|
2347 for (int i = 0; i < nr; i++) |
|
2348 { |
|
2349 double res = elem (i, 0); |
|
2350 for (int j = 1; j < nc; j++) |
|
2351 if (elem (i, j) > res) |
|
2352 res = elem (i, j); |
|
2353 result.elem (i) = res; |
|
2354 } |
|
2355 } |
|
2356 |
|
2357 return result; |
|
2358 } |
|
2359 |
|
2360 ColumnVector |
|
2361 Matrix::row_max_loc (void) const |
|
2362 { |
|
2363 ColumnVector result; |
|
2364 |
|
2365 int nr = rows (); |
|
2366 int nc = cols (); |
|
2367 |
|
2368 if (nr > 0 && nc > 0) |
|
2369 { |
|
2370 result.resize (nr); |
|
2371 |
|
2372 for (int i = 0; i < nr; i++) |
|
2373 { |
|
2374 int res = 0; |
|
2375 for (int j = 0; j < nc; j++) |
|
2376 if (elem (i, j) > elem (i, res)) |
|
2377 res = j; |
|
2378 result.elem (i) = (double) (res + 1); |
|
2379 } |
|
2380 } |
|
2381 |
|
2382 return result; |
|
2383 } |
|
2384 |
|
2385 RowVector |
|
2386 Matrix::column_min (void) const |
|
2387 { |
|
2388 RowVector result; |
|
2389 |
|
2390 int nr = rows (); |
|
2391 int nc = cols (); |
|
2392 |
|
2393 if (nr > 0 && nc > 0) |
|
2394 { |
|
2395 result.resize (nc); |
|
2396 |
|
2397 for (int j = 0; j < nc; j++) |
|
2398 { |
|
2399 double res = elem (0, j); |
|
2400 for (int i = 1; i < nr; i++) |
|
2401 if (elem (i, j) < res) |
|
2402 res = elem (i, j); |
|
2403 result.elem (j) = res; |
|
2404 } |
|
2405 } |
|
2406 |
|
2407 return result; |
|
2408 } |
|
2409 RowVector |
|
2410 Matrix::column_min_loc (void) const |
|
2411 { |
|
2412 RowVector result; |
|
2413 |
|
2414 int nr = rows (); |
|
2415 int nc = cols (); |
|
2416 |
|
2417 if (nr > 0 && nc > 0) |
|
2418 { |
|
2419 result.resize (nc); |
|
2420 |
|
2421 for (int j = 0; j < nc; j++) |
|
2422 { |
|
2423 int res = 0; |
|
2424 for (int i = 0; i < nr; i++) |
|
2425 if (elem (i, j) < elem (res, j)) |
|
2426 res = i; |
|
2427 result.elem (j) = (double) (res + 1); |
|
2428 } |
|
2429 } |
|
2430 |
|
2431 return result; |
|
2432 } |
|
2433 |
|
2434 |
|
2435 RowVector |
|
2436 Matrix::column_max (void) const |
|
2437 { |
|
2438 RowVector result; |
|
2439 |
|
2440 int nr = rows (); |
|
2441 int nc = cols (); |
|
2442 |
|
2443 if (nr > 0 && nc > 0) |
|
2444 { |
|
2445 result.resize (nc); |
|
2446 |
|
2447 for (int j = 0; j < nc; j++) |
|
2448 { |
|
2449 double res = elem (0, j); |
|
2450 for (int i = 1; i < nr; i++) |
|
2451 if (elem (i, j) > res) |
|
2452 res = elem (i, j); |
|
2453 result.elem (j) = res; |
|
2454 } |
|
2455 } |
|
2456 |
|
2457 return result; |
|
2458 } |
|
2459 |
|
2460 RowVector |
|
2461 Matrix::column_max_loc (void) const |
|
2462 { |
|
2463 RowVector result; |
|
2464 |
|
2465 int nr = rows (); |
|
2466 int nc = cols (); |
|
2467 |
|
2468 if (nr > 0 && nc > 0) |
|
2469 { |
|
2470 result.resize (nc); |
|
2471 |
|
2472 for (int j = 0; j < nc; j++) |
|
2473 { |
|
2474 int res = 0; |
|
2475 for (int i = 0; i < nr; i++) |
|
2476 if (elem (i, j) > elem (res, j)) |
|
2477 res = i; |
|
2478 result.elem (j) = (double) (res + 1); |
|
2479 } |
|
2480 } |
|
2481 |
|
2482 return result; |
|
2483 } |
|
2484 |
|
2485 ostream& |
|
2486 operator << (ostream& os, const Matrix& a) |
|
2487 { |
|
2488 // int field_width = os.precision () + 7; |
1360
|
2489 |
458
|
2490 for (int i = 0; i < a.rows (); i++) |
|
2491 { |
|
2492 for (int j = 0; j < a.cols (); j++) |
|
2493 os << " " /* setw (field_width) */ << a.elem (i, j); |
|
2494 os << "\n"; |
|
2495 } |
|
2496 return os; |
|
2497 } |
|
2498 |
|
2499 istream& |
|
2500 operator >> (istream& is, Matrix& a) |
|
2501 { |
|
2502 int nr = a.rows (); |
|
2503 int nc = a.cols (); |
|
2504 |
|
2505 if (nr < 1 || nc < 1) |
|
2506 is.clear (ios::badbit); |
|
2507 else |
|
2508 { |
|
2509 double tmp; |
|
2510 for (int i = 0; i < nr; i++) |
|
2511 for (int j = 0; j < nc; j++) |
|
2512 { |
|
2513 is >> tmp; |
|
2514 if (is) |
|
2515 a.elem (i, j) = tmp; |
|
2516 else |
|
2517 break; |
|
2518 } |
|
2519 } |
|
2520 |
|
2521 return is; |
|
2522 } |
|
2523 |
1365
|
2524 // Read an array of data from a file in binary format. |
1360
|
2525 |
458
|
2526 int |
1365
|
2527 Matrix::read (FILE *fptr, const char *type) |
458
|
2528 { |
1360
|
2529 // Allocate buffer pointers. |
458
|
2530 |
|
2531 union |
|
2532 { |
|
2533 void *vd; |
|
2534 char *ch; |
|
2535 u_char *uc; |
|
2536 short *sh; |
|
2537 u_short *us; |
|
2538 int *in; |
|
2539 u_int *ui; |
|
2540 long *ln; |
|
2541 u_long *ul; |
|
2542 float *fl; |
|
2543 double *db; |
|
2544 } |
|
2545 buf; |
|
2546 |
1360
|
2547 // Convert data to double. |
458
|
2548 |
471
|
2549 if (! type) |
458
|
2550 { |
471
|
2551 (*current_liboctave_error_handler) |
|
2552 ("fread: invalid NULL type parameter"); |
|
2553 return 0; |
|
2554 } |
458
|
2555 |
471
|
2556 int count; |
|
2557 int nitems = length (); |
458
|
2558 |
471
|
2559 double *d = fortran_vec (); // Ensures only one reference to my privates! |
458
|
2560 |
471
|
2561 #define DO_FREAD(TYPE,ELEM) \ |
|
2562 do \ |
|
2563 { \ |
|
2564 size_t size = sizeof (TYPE); \ |
|
2565 buf.ch = new char [size * nitems]; \ |
|
2566 count = fread (buf.ch, size, nitems, fptr); \ |
|
2567 for (int k = 0; k < count; k++) \ |
|
2568 d[k] = buf.ELEM[k]; \ |
|
2569 delete [] buf.ch; \ |
|
2570 } \ |
|
2571 while (0) |
458
|
2572 |
471
|
2573 if (strcasecmp (type, "double") == 0) |
|
2574 DO_FREAD (double, db); |
|
2575 else if (strcasecmp (type, "char") == 0) |
|
2576 DO_FREAD (char, ch); |
|
2577 else if (strcasecmp (type, "uchar") == 0) |
|
2578 DO_FREAD (u_char, uc); |
|
2579 else if (strcasecmp (type, "short") == 0) |
|
2580 DO_FREAD (short, sh); |
|
2581 else if (strcasecmp (type, "ushort") == 0) |
|
2582 DO_FREAD (u_short, us); |
|
2583 else if (strcasecmp (type, "int") == 0) |
|
2584 DO_FREAD (int, in); |
|
2585 else if (strcasecmp (type, "uint") == 0) |
|
2586 DO_FREAD (u_int, ui); |
|
2587 else if (strcasecmp (type, "long") == 0) |
|
2588 DO_FREAD (long, ul); |
|
2589 else if (strcasecmp (type, "float") == 0) |
|
2590 DO_FREAD (float, fl); |
|
2591 else |
|
2592 { |
|
2593 (*current_liboctave_error_handler) |
|
2594 ("fread: invalid NULL type parameter"); |
458
|
2595 return 0; |
|
2596 } |
|
2597 |
|
2598 return count; |
|
2599 } |
|
2600 |
1360
|
2601 // Write the data array to a file in binary format. |
|
2602 |
458
|
2603 int |
1365
|
2604 Matrix::write (FILE *fptr, const char *type) |
458
|
2605 { |
1360
|
2606 // Allocate buffer pointers. |
458
|
2607 |
|
2608 union |
|
2609 { |
|
2610 void *vd; |
|
2611 char *ch; |
|
2612 u_char *uc; |
|
2613 short *sh; |
|
2614 u_short *us; |
|
2615 int *in; |
|
2616 u_int *ui; |
|
2617 long *ln; |
|
2618 u_long *ul; |
|
2619 float *fl; |
|
2620 double *db; |
|
2621 } |
|
2622 buf; |
|
2623 |
471
|
2624 int nitems = length (); |
458
|
2625 |
471
|
2626 double *d = fortran_vec (); |
458
|
2627 |
1360
|
2628 // Convert from double to correct size. |
458
|
2629 |
471
|
2630 if (! type) |
458
|
2631 { |
471
|
2632 (*current_liboctave_error_handler) |
|
2633 ("fwrite: invalid NULL type parameter"); |
|
2634 return 0; |
|
2635 } |
458
|
2636 |
471
|
2637 size_t size; |
|
2638 int count; |
458
|
2639 |
471
|
2640 #define DO_FWRITE(TYPE,ELEM) \ |
|
2641 do \ |
|
2642 { \ |
|
2643 size = sizeof (TYPE); \ |
|
2644 buf.ELEM = new TYPE [nitems]; \ |
|
2645 for (int k = 0; k < nitems; k++) \ |
|
2646 buf.ELEM[k] = (TYPE) d[k]; \ |
|
2647 count = fwrite (buf.ELEM, size, nitems, fptr); \ |
|
2648 delete [] buf.ELEM; \ |
|
2649 } \ |
|
2650 while (0) |
458
|
2651 |
471
|
2652 if (strcasecmp (type, "double") == 0) |
|
2653 DO_FWRITE (double, db); |
|
2654 else if (strcasecmp (type, "char") == 0) |
|
2655 DO_FWRITE (char, ch); |
|
2656 else if (strcasecmp (type, "uchar") == 0) |
|
2657 DO_FWRITE (u_char, uc); |
|
2658 else if (strcasecmp (type, "short") == 0) |
|
2659 DO_FWRITE (short, sh); |
|
2660 else if (strcasecmp (type, "ushort") == 0) |
|
2661 DO_FWRITE (u_short, us); |
|
2662 else if (strcasecmp (type, "int") == 0) |
|
2663 DO_FWRITE (int, in); |
|
2664 else if (strcasecmp (type, "uint") == 0) |
|
2665 DO_FWRITE (u_int, ui); |
|
2666 else if (strcasecmp (type, "long") == 0) |
|
2667 DO_FWRITE (long, ln); |
|
2668 else if (strcasecmp (type, "ulong") == 0) |
|
2669 DO_FWRITE (u_long, ul); |
|
2670 else if (strcasecmp (type, "float") == 0) |
|
2671 DO_FWRITE (float, fl); |
|
2672 else |
|
2673 { |
|
2674 (*current_liboctave_error_handler) |
|
2675 ("fwrite: unrecognized type parameter %s", type); |
458
|
2676 return 0; |
471
|
2677 } |
458
|
2678 |
|
2679 return count; |
|
2680 } |
|
2681 |
1819
|
2682 Matrix |
|
2683 Givens (double x, double y) |
|
2684 { |
|
2685 double cc, s, temp_r; |
|
2686 |
|
2687 F77_FCN (dlartg, DLARTG) (x, y, cc, s, temp_r); |
|
2688 |
|
2689 Matrix g (2, 2); |
|
2690 |
|
2691 g.elem (0, 0) = cc; |
|
2692 g.elem (1, 1) = cc; |
|
2693 g.elem (0, 1) = s; |
|
2694 g.elem (1, 0) = -s; |
|
2695 |
|
2696 return g; |
|
2697 } |
|
2698 |
|
2699 Matrix |
|
2700 Sylvester (const Matrix& a, const Matrix& b, const Matrix& c) |
|
2701 { |
|
2702 Matrix retval; |
|
2703 |
|
2704 // XXX FIXME XXX -- need to check that a, b, and c are all the same |
|
2705 // size. |
|
2706 |
|
2707 // Compute Schur decompositions. |
|
2708 |
|
2709 SCHUR as (a, "U"); |
|
2710 SCHUR bs (b, "U"); |
|
2711 |
|
2712 // Transform c to new coordinates. |
|
2713 |
|
2714 Matrix ua = as.unitary_matrix (); |
|
2715 Matrix sch_a = as.schur_matrix (); |
|
2716 |
|
2717 Matrix ub = bs.unitary_matrix (); |
|
2718 Matrix sch_b = bs.schur_matrix (); |
|
2719 |
|
2720 Matrix cx = ua.transpose () * c * ub; |
|
2721 |
|
2722 // Solve the sylvester equation, back-transform, and return the |
|
2723 // solution. |
|
2724 |
|
2725 int a_nr = a.rows (); |
|
2726 int b_nr = b.rows (); |
|
2727 |
|
2728 double scale; |
|
2729 int info; |
|
2730 |
1950
|
2731 double *pa = sch_a.fortran_vec (); |
|
2732 double *pb = sch_b.fortran_vec (); |
|
2733 double *px = cx.fortran_vec (); |
|
2734 |
|
2735 F77_XFCN (dtrsyl, DTRSYL, ("N", "N", 1, a_nr, b_nr, pa, a_nr, pb, |
|
2736 b_nr, px, a_nr, scale, info, 1L, 1L)); |
|
2737 |
|
2738 |
|
2739 if (f77_exception_encountered) |
|
2740 (*current_liboctave_error_handler) ("unrecoverable error in dtrsyl"); |
|
2741 else |
|
2742 { |
|
2743 // XXX FIXME XXX -- check info? |
1819
|
2744 |
1950
|
2745 retval = -ua*cx*ub.transpose (); |
|
2746 } |
1819
|
2747 |
|
2748 return retval; |
|
2749 } |
|
2750 |
1959
|
2751 ComplexColumnVector |
|
2752 Qzval (const Matrix& a, const Matrix& b) |
|
2753 { |
|
2754 ComplexColumnVector retval; |
|
2755 |
|
2756 int a_nr = a.rows(); |
|
2757 int a_nc = a.cols(); |
|
2758 |
|
2759 int b_nr = b.rows(); |
|
2760 int b_nc = b.cols(); |
|
2761 |
|
2762 if (a_nr == a_nc) |
|
2763 { |
|
2764 if (a_nr == b_nr && a_nc == b_nc) |
|
2765 { |
|
2766 if (a_nr != 0) |
|
2767 { |
|
2768 Matrix jnk (a_nr, a_nr, 0.0); |
|
2769 double *pjnk = jnk.fortran_vec (); |
|
2770 |
|
2771 ColumnVector alfr (a_nr); |
|
2772 double *palfr = alfr.fortran_vec (); |
|
2773 |
|
2774 ColumnVector alfi (a_nr); |
|
2775 double *palfi = alfr.fortran_vec (); |
|
2776 |
|
2777 ColumnVector beta (a_nr); |
|
2778 double *pbeta = alfr.fortran_vec (); |
|
2779 |
|
2780 Matrix atmp = a; |
|
2781 double *pa = atmp.fortran_vec (); |
|
2782 |
|
2783 Matrix btmp = b; |
|
2784 double *pb = btmp.fortran_vec (); |
|
2785 |
|
2786 long matz = 0; |
|
2787 int info; |
|
2788 |
|
2789 // XXX FIXME ??? XXX |
|
2790 double eps = DBL_EPSILON; |
|
2791 |
|
2792 F77_FCN (qzhes, QZHES) (a_nr, a_nr, pa, pb, matz, pjnk); |
|
2793 |
|
2794 F77_FCN (qzit, QZIT) (a_nr, a_nr, pa, pb, eps, matz, pjnk, info); |
|
2795 |
|
2796 if (! info) |
|
2797 { |
|
2798 F77_FCN (qzval, QZVAL) (a_nr, a_nr, pa, pb, palfr, |
|
2799 palfi, pbeta, matz, pjnk); |
|
2800 |
|
2801 // Count and extract finite generalized eigenvalues. |
|
2802 |
|
2803 int cnt = 0; |
|
2804 |
|
2805 for (int i = 0; i < a_nr; i++) |
|
2806 if (beta.elem (i) != 0) |
|
2807 cnt++; |
|
2808 |
|
2809 ComplexColumnVector cx (cnt, 0.0); |
|
2810 |
|
2811 Complex Im (0, 1); |
|
2812 |
|
2813 for (int i = 0; i < a_nr; i++) |
|
2814 { |
|
2815 if (beta.elem (i) != 0) |
|
2816 { |
|
2817 // Finite generalized eigenvalue. |
|
2818 |
|
2819 cnt--; |
|
2820 cx.elem (cnt) = (alfr.elem (i) + Im * alfi.elem (i)) |
|
2821 / beta.elem (i); |
|
2822 } |
|
2823 } |
|
2824 |
|
2825 retval = cx; |
|
2826 } |
|
2827 else |
|
2828 (*current_liboctave_error_handler) |
|
2829 ("qzval: trouble in qzit, info = %d", info); |
|
2830 } |
|
2831 } |
|
2832 else |
|
2833 (*current_liboctave_error_handler) ("nonconformant matrices"); |
|
2834 } |
|
2835 else |
|
2836 (*current_liboctave_error_handler) ("qzval: square matrices required"); |
|
2837 |
|
2838 return retval; |
|
2839 } |
|
2840 |
458
|
2841 /* |
|
2842 ;;; Local Variables: *** |
|
2843 ;;; mode: C++ *** |
|
2844 ;;; page-delimiter: "^/\\*" *** |
|
2845 ;;; End: *** |
|
2846 */ |