GNU Octave 7.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
dMatrix.cc
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////
2//
3// Copyright (C) 1994-2022 The Octave Project Developers
4//
5// See the file COPYRIGHT.md in the top-level directory of this
6// distribution or <https://octave.org/copyright/>.
7//
8// This file is part of Octave.
9//
10// Octave is free software: you can redistribute it and/or modify it
11// under the terms of the GNU General Public License as published by
12// the Free Software Foundation, either version 3 of the License, or
13// (at your option) any later version.
14//
15// Octave is distributed in the hope that it will be useful, but
16// WITHOUT ANY WARRANTY; without even the implied warranty of
17// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18// GNU General Public License for more details.
19//
20// You should have received a copy of the GNU General Public License
21// along with Octave; see the file COPYING. If not, see
22// <https://www.gnu.org/licenses/>.
23//
24////////////////////////////////////////////////////////////////////////
25
26#if defined (HAVE_CONFIG_H)
27# include "config.h"
28#endif
29
30#include <algorithm>
31#include <istream>
32#include <limits>
33#include <ostream>
34
35#include "Array-util.h"
36#include "CColVector.h"
37#include "CMatrix.h"
38#include "DET.h"
39#include "PermMatrix.h"
40#include "boolMatrix.h"
41#include "byte-swap.h"
42#include "chMatrix.h"
43#include "chol.h"
44#include "dColVector.h"
45#include "dDiagMatrix.h"
46#include "dMatrix.h"
47#include "dRowVector.h"
48#include "lo-blas-proto.h"
49#include "lo-error.h"
50#include "lo-ieee.h"
51#include "lo-lapack-proto.h"
52#include "lo-mappers.h"
53#include "lo-utils.h"
54#include "mx-dm-m.h"
55#include "mx-inlines.cc"
56#include "mx-m-dm.h"
57#include "mx-op-defs.h"
58#include "oct-cmplx.h"
59#include "oct-fftw.h"
60#include "oct-locbuf.h"
61#include "oct-norm.h"
62#include "quit.h"
63#include "schur.h"
64#include "svd.h"
65
66// Matrix class.
67
69 : NDArray (rv)
70{ }
71
73 : NDArray (cv)
74{ }
75
77 : NDArray (a.dims (), 0.0)
78{
79 for (octave_idx_type i = 0; i < a.length (); i++)
80 elem (i, i) = a.elem (i, i);
81}
82
84 : NDArray (a.dims (), 0.0)
85{
86 for (octave_idx_type i = 0; i < a.length (); i++)
87 elem (i, i) = a.elem (i, i);
88}
89
91 : NDArray (a.dims (), 0.0)
92{
93 for (octave_idx_type i = 0; i < a.length (); i++)
94 elem (i, i) = a.elem (i, i);
95}
96
98 : NDArray (a.dims (), 0.0)
99{
100 const Array<octave_idx_type> ia (a.col_perm_vec ());
101 octave_idx_type len = a.rows ();
102 for (octave_idx_type i = 0; i < len; i++)
103 elem (ia(i), i) = 1.0;
104}
105
106// FIXME: could we use a templated mixed-type copy function here?
107
109 : NDArray (a)
110{ }
111
113 : NDArray (a.dims ())
114{
115 for (octave_idx_type i = 0; i < a.rows (); i++)
116 for (octave_idx_type j = 0; j < a.cols (); j++)
117 elem (i, j) = static_cast<unsigned char> (a.elem (i, j));
118}
119
120bool
122{
123 if (rows () != a.rows () || cols () != a.cols ())
124 return false;
125
126 return mx_inline_equal (numel (), data (), a.data ());
127}
128
129bool
131{
132 return !(*this == a);
133}
135bool
137{
138 if (issquare () && rows () > 0)
139 {
140 for (octave_idx_type i = 0; i < rows (); i++)
141 for (octave_idx_type j = i+1; j < cols (); j++)
142 if (elem (i, j) != elem (j, i))
143 return false;
144
145 return true;
146 }
147
148 return false;
149}
150
151Matrix&
153{
154 Array<double>::insert (a, r, c);
155 return *this;
156}
157
158Matrix&
160{
161 octave_idx_type a_len = a.numel ();
162
163 if (r < 0 || r >= rows () || c < 0 || c + a_len > cols ())
164 (*current_liboctave_error_handler) ("range error for insert");
165
166 if (a_len > 0)
167 {
168 make_unique ();
169
170 for (octave_idx_type i = 0; i < a_len; i++)
171 xelem (r, c+i) = a.elem (i);
172 }
173
174 return *this;
175}
176
177Matrix&
179{
180 octave_idx_type a_len = a.numel ();
181
182 if (r < 0 || r + a_len > rows () || c < 0 || c >= cols ())
183 (*current_liboctave_error_handler) ("range error for insert");
184
185 if (a_len > 0)
186 {
187 make_unique ();
188
189 for (octave_idx_type i = 0; i < a_len; i++)
190 xelem (r+i, c) = a.elem (i);
191 }
192
193 return *this;
194}
195
196Matrix&
198{
199 octave_idx_type a_nr = a.rows ();
200 octave_idx_type a_nc = a.cols ();
201
202 if (r < 0 || r + a_nr > rows () || c < 0 || c + a_nc > cols ())
203 (*current_liboctave_error_handler) ("range error for insert");
204
205 fill (0.0, r, c, r + a_nr - 1, c + a_nc - 1);
206
207 octave_idx_type a_len = a.length ();
208
209 if (a_len > 0)
210 {
211 make_unique ();
212
213 for (octave_idx_type i = 0; i < a_len; i++)
214 xelem (r+i, c+i) = a.elem (i, i);
215 }
216
217 return *this;
218}
219
220Matrix&
221Matrix::fill (double val)
222{
223 octave_idx_type nr = rows ();
224 octave_idx_type nc = cols ();
225
226 if (nr > 0 && nc > 0)
227 {
228 make_unique ();
229
230 for (octave_idx_type j = 0; j < nc; j++)
231 for (octave_idx_type i = 0; i < nr; i++)
232 xelem (i, j) = val;
233 }
234
235 return *this;
236}
237
238Matrix&
241{
242 octave_idx_type nr = rows ();
243 octave_idx_type nc = cols ();
244
245 if (r1 < 0 || r2 < 0 || c1 < 0 || c2 < 0
246 || r1 >= nr || r2 >= nr || c1 >= nc || c2 >= nc)
247 (*current_liboctave_error_handler) ("range error for fill");
248
249 if (r1 > r2) { std::swap (r1, r2); }
250 if (c1 > c2) { std::swap (c1, c2); }
251
252 if (r2 >= r1 && c2 >= c1)
253 {
254 make_unique ();
255
256 for (octave_idx_type j = c1; j <= c2; j++)
257 for (octave_idx_type i = r1; i <= r2; i++)
258 xelem (i, j) = val;
259 }
260
261 return *this;
262}
263
264Matrix
265Matrix::append (const Matrix& a) const
266{
267 octave_idx_type nr = rows ();
268 octave_idx_type nc = cols ();
269 if (nr != a.rows ())
270 (*current_liboctave_error_handler) ("row dimension mismatch for append");
271
272 octave_idx_type nc_insert = nc;
273 Matrix retval (nr, nc + a.cols ());
274 retval.insert (*this, 0, 0);
275 retval.insert (a, 0, nc_insert);
276 return retval;
277}
278
279Matrix
280Matrix::append (const RowVector& a) const
281{
282 octave_idx_type nr = rows ();
283 octave_idx_type nc = cols ();
284 if (nr != 1)
285 (*current_liboctave_error_handler) ("row dimension mismatch for append");
286
287 octave_idx_type nc_insert = nc;
288 Matrix retval (nr, nc + a.numel ());
289 retval.insert (*this, 0, 0);
290 retval.insert (a, 0, nc_insert);
291 return retval;
292}
293
294Matrix
296{
297 octave_idx_type nr = rows ();
298 octave_idx_type nc = cols ();
299 if (nr != a.numel ())
300 (*current_liboctave_error_handler) ("row dimension mismatch for append");
301
302 octave_idx_type nc_insert = nc;
303 Matrix retval (nr, nc + 1);
304 retval.insert (*this, 0, 0);
305 retval.insert (a, 0, nc_insert);
306 return retval;
307}
308
309Matrix
311{
312 octave_idx_type nr = rows ();
313 octave_idx_type nc = cols ();
314 if (nr != a.rows ())
315 (*current_liboctave_error_handler) ("row dimension mismatch for append");
316
317 octave_idx_type nc_insert = nc;
318 Matrix retval (nr, nc + a.cols ());
319 retval.insert (*this, 0, 0);
320 retval.insert (a, 0, nc_insert);
321 return retval;
322}
323
324Matrix
325Matrix::stack (const Matrix& a) const
326{
327 octave_idx_type nr = rows ();
328 octave_idx_type nc = cols ();
329 if (nc != a.cols ())
330 (*current_liboctave_error_handler) ("column dimension mismatch for stack");
331
332 octave_idx_type nr_insert = nr;
333 Matrix retval (nr + a.rows (), nc);
334 retval.insert (*this, 0, 0);
335 retval.insert (a, nr_insert, 0);
336 return retval;
337}
338
339Matrix
340Matrix::stack (const RowVector& a) const
341{
342 octave_idx_type nr = rows ();
343 octave_idx_type nc = cols ();
344 if (nc != a.numel ())
345 (*current_liboctave_error_handler) ("column dimension mismatch for stack");
346
347 octave_idx_type nr_insert = nr;
348 Matrix retval (nr + 1, nc);
349 retval.insert (*this, 0, 0);
350 retval.insert (a, nr_insert, 0);
351 return retval;
352}
353
354Matrix
356{
357 octave_idx_type nr = rows ();
358 octave_idx_type nc = cols ();
359 if (nc != 1)
360 (*current_liboctave_error_handler) ("column dimension mismatch for stack");
361
362 octave_idx_type nr_insert = nr;
363 Matrix retval (nr + a.numel (), nc);
364 retval.insert (*this, 0, 0);
365 retval.insert (a, nr_insert, 0);
366 return retval;
367}
368
369Matrix
370Matrix::stack (const DiagMatrix& a) const
371{
372 octave_idx_type nr = rows ();
373 octave_idx_type nc = cols ();
374 if (nc != a.cols ())
375 (*current_liboctave_error_handler) ("column dimension mismatch for stack");
376
377 octave_idx_type nr_insert = nr;
378 Matrix retval (nr + a.rows (), nc);
379 retval.insert (*this, 0, 0);
380 retval.insert (a, nr_insert, 0);
381 return retval;
382}
383
384Matrix
386{
387 return do_mx_unary_op<double, Complex> (a, mx_inline_real);
388}
389
390Matrix
392{
393 return do_mx_unary_op<double, Complex> (a, mx_inline_imag);
394}
395
396Matrix
398 octave_idx_type r2, octave_idx_type c2) const
399{
400 if (r1 > r2) { std::swap (r1, r2); }
401 if (c1 > c2) { std::swap (c1, c2); }
402
403 return index (octave::idx_vector (r1, r2+1), octave::idx_vector (c1, c2+1));
404}
405
406Matrix
408 octave_idx_type nc) const
409{
410 return index (octave::idx_vector (r1, r1 + nr), octave::idx_vector (c1, c1 + nc));
411}
412
413// extract row or column i.
414
417{
419}
420
423{
425}
426
427// Local function to calculate the 1-norm.
428static
429double
430norm1 (const Matrix& a)
431{
432 double anorm = 0.0;
433 RowVector colsum = a.abs ().sum ().row (0);
434
435 for (octave_idx_type i = 0; i < colsum.numel (); i++)
436 {
437 double sum = colsum.xelem (i);
438 if (octave::math::isinf (sum) || octave::math::isnan (sum))
439 {
440 anorm = sum; // Pass Inf or NaN to output
441 break;
442 }
443 else
444 anorm = std::max (anorm, sum);
445 }
446
447 return anorm;
448}
449
450Matrix
451Matrix::inverse (void) const
452{
453 octave_idx_type info;
454 double rcon;
455 MatrixType mattype (*this);
456 return inverse (mattype, info, rcon, 0, 0);
457}
458
459Matrix
461{
462 double rcon;
463 MatrixType mattype (*this);
464 return inverse (mattype, info, rcon, 0, 0);
465}
466
467Matrix
468Matrix::inverse (octave_idx_type& info, double& rcon, bool force,
469 bool calc_cond) const
470{
471 MatrixType mattype (*this);
472 return inverse (mattype, info, rcon, force, calc_cond);
473}
474
475Matrix
477{
478 octave_idx_type info;
479 double rcon;
480 return inverse (mattype, info, rcon, 0, 0);
481}
482
483Matrix
485{
486 double rcon;
487 return inverse (mattype, info, rcon, 0, 0);
488}
489
490Matrix
491Matrix::tinverse (MatrixType& mattype, octave_idx_type& info, double& rcon,
492 bool force, bool calc_cond) const
493{
494 Matrix retval;
495
496 F77_INT nr = octave::to_f77_int (rows ());
497 F77_INT nc = octave::to_f77_int (cols ());
498
499 if (nr != nc || nr == 0 || nc == 0)
500 (*current_liboctave_error_handler) ("inverse requires square matrix");
501
502 int typ = mattype.type ();
503 char uplo = (typ == MatrixType::Lower ? 'L' : 'U');
504 char udiag = 'N';
505 retval = *this;
506 double *tmp_data = retval.fortran_vec ();
507
508 F77_INT tmp_info = 0;
509
510 F77_XFCN (dtrtri, DTRTRI, (F77_CONST_CHAR_ARG2 (&uplo, 1),
511 F77_CONST_CHAR_ARG2 (&udiag, 1),
512 nr, tmp_data, nr, tmp_info
513 F77_CHAR_ARG_LEN (1)
514 F77_CHAR_ARG_LEN (1)));
515
516 info = tmp_info;
517
518 // Throw away extra info LAPACK gives so as to not change output.
519 rcon = 0.0;
520 if (info != 0)
521 info = -1;
522 else if (calc_cond)
523 {
524 F77_INT dtrcon_info = 0;
525 char job = '1';
526
527 OCTAVE_LOCAL_BUFFER (double, work, 3 * nr);
528 OCTAVE_LOCAL_BUFFER (F77_INT, iwork, nr);
529
530 F77_XFCN (dtrcon, DTRCON, (F77_CONST_CHAR_ARG2 (&job, 1),
531 F77_CONST_CHAR_ARG2 (&uplo, 1),
532 F77_CONST_CHAR_ARG2 (&udiag, 1),
533 nr, tmp_data, nr, rcon,
534 work, iwork, dtrcon_info
535 F77_CHAR_ARG_LEN (1)
536 F77_CHAR_ARG_LEN (1)
537 F77_CHAR_ARG_LEN (1)));
538
539 if (dtrcon_info != 0)
540 info = -1;
541 }
542
543 if (info == -1 && ! force)
544 retval = *this; // Restore matrix contents.
545
546 return retval;
547}
548
549Matrix
550Matrix::finverse (MatrixType& mattype, octave_idx_type& info, double& rcon,
551 bool force, bool calc_cond) const
552{
553 Matrix retval;
554
555 F77_INT nr = octave::to_f77_int (rows ());
556 F77_INT nc = octave::to_f77_int (cols ());
557
558 if (nr != nc || nr == 0 || nc == 0)
559 (*current_liboctave_error_handler) ("inverse requires square matrix");
560
561 Array<F77_INT> ipvt (dim_vector (nr, 1));
562 F77_INT *pipvt = ipvt.fortran_vec ();
563
564 retval = *this;
565 double *tmp_data = retval.fortran_vec ();
566
567 Array<double> z (dim_vector (1, 1));
568 F77_INT lwork = -1;
569
570 F77_INT tmp_info = 0;
571
572 // Query the optimum work array size.
573 F77_XFCN (dgetri, DGETRI, (nc, tmp_data, nr, pipvt,
574 z.fortran_vec (), lwork, tmp_info));
575
576 lwork = static_cast<F77_INT> (z(0));
577 lwork = (lwork < 4 * nc ? 4 * nc : lwork);
578 z.resize (dim_vector (lwork, 1));
579 double *pz = z.fortran_vec ();
580
581 info = 0;
582 tmp_info = 0;
583
584 // Calculate the norm of the matrix for later use when determining rcon.
585 double anorm;
586 if (calc_cond)
587 anorm = norm1 (retval);
588
589 F77_XFCN (dgetrf, DGETRF, (nc, nc, tmp_data, nr, pipvt, tmp_info));
590
591 info = tmp_info;
592
593 // Throw away extra info LAPACK gives so as to not change output.
594 rcon = 0.0;
595 if (info != 0)
596 info = -1;
597 else if (calc_cond)
598 {
599 F77_INT dgecon_info = 0;
600
601 // Now calculate the condition number for non-singular matrix.
602 char job = '1';
603 Array<F77_INT> iz (dim_vector (nc, 1));
604 F77_INT *piz = iz.fortran_vec ();
605 F77_XFCN (dgecon, DGECON, (F77_CONST_CHAR_ARG2 (&job, 1),
606 nc, tmp_data, nr, anorm,
607 rcon, pz, piz, dgecon_info
608 F77_CHAR_ARG_LEN (1)));
609
610 if (dgecon_info != 0)
611 info = -1;
612 }
613
614 if (info == -1 && ! force)
615 retval = *this; // Restore matrix contents.
616 else
617 {
618 F77_INT dgetri_info = 0;
619
620 F77_XFCN (dgetri, DGETRI, (nc, tmp_data, nr, pipvt,
621 pz, lwork, dgetri_info));
622
623 if (dgetri_info != 0)
624 info = -1;
625 }
626
627 if (info != 0)
628 mattype.mark_as_rectangular ();
629
630 return retval;
631}
632
633Matrix
634Matrix::inverse (MatrixType& mattype, octave_idx_type& info, double& rcon,
635 bool force, bool calc_cond) const
636{
637 int typ = mattype.type (false);
638 Matrix ret;
639
640 if (typ == MatrixType::Unknown)
641 typ = mattype.type (*this);
642
643 if (typ == MatrixType::Diagonal) // a scalar is classified as Diagonal.
644 {
645 ret = 1 / (*this);
646 if (calc_cond)
647 {
648 double scalar = this->elem (0);
649 if (octave::math::isfinite (scalar) && scalar != 0)
650 rcon = 1.0;
651 else if (octave::math::isinf (scalar) || scalar == 0)
652 rcon = 0.0;
653 else
655 }
656 }
657 else if (typ == MatrixType::Upper || typ == MatrixType::Lower)
658 ret = tinverse (mattype, info, rcon, force, calc_cond);
659 else
660 {
661 if (mattype.ishermitian ())
662 {
663 octave::math::chol<Matrix> chol (*this, info, true, calc_cond);
664 if (info == 0)
665 {
666 if (calc_cond)
667 rcon = chol.rcond ();
668 else
669 rcon = 1.0;
670 ret = chol.inverse ();
671 }
672 else
673 mattype.mark_as_unsymmetric ();
674 }
675
676 if (! mattype.ishermitian ())
677 ret = finverse (mattype, info, rcon, force, calc_cond);
678
679 if ((calc_cond || mattype.ishermitian ()) && rcon == 0.0)
680 ret = Matrix (rows (), columns (),
682 }
683
684 return ret;
685}
686
687Matrix
688Matrix::pseudo_inverse (double tol) const
689{
690 octave::math::svd<Matrix> result (*this,
692
693 DiagMatrix S = result.singular_values ();
694 Matrix U = result.left_singular_matrix ();
695 Matrix V = result.right_singular_matrix ();
696
697 ColumnVector sigma = S.extract_diag ();
698
699 octave_idx_type r = sigma.numel () - 1;
700 octave_idx_type nr = rows ();
701 octave_idx_type nc = cols ();
702
703 if (tol <= 0.0)
704 {
705 tol = std::max (nr, nc) * sigma.elem (0)
706 * std::numeric_limits<double>::epsilon ();
707
708 if (tol == 0)
710 }
711
712 while (r >= 0 && sigma.elem (r) < tol)
713 r--;
714
715 if (r < 0)
716 return Matrix (nc, nr, 0.0);
717 else
718 {
719 Matrix Ur = U.extract (0, 0, nr-1, r);
720 DiagMatrix D = DiagMatrix (sigma.extract (0, r)).inverse ();
721 Matrix Vr = V.extract (0, 0, nc-1, r);
722 return Vr * D * Ur.transpose ();
723 }
724}
725
726#if defined (HAVE_FFTW)
727
729Matrix::fourier (void) const
730{
731 std::size_t nr = rows ();
732 std::size_t nc = cols ();
733
734 ComplexMatrix retval (nr, nc);
735
736 std::size_t npts, nsamples;
737
738 if (nr == 1 || nc == 1)
739 {
740 npts = (nr > nc ? nr : nc);
741 nsamples = 1;
742 }
743 else
744 {
745 npts = nr;
746 nsamples = nc;
747 }
748
749 const double *in (data ());
750 Complex *out (retval.fortran_vec ());
751
752 octave::fftw::fft (in, out, npts, nsamples);
753
754 return retval;
755}
756
759{
760 std::size_t nr = rows ();
761 std::size_t nc = cols ();
762
763 ComplexMatrix retval (nr, nc);
764
765 std::size_t npts, nsamples;
766
767 if (nr == 1 || nc == 1)
768 {
769 npts = (nr > nc ? nr : nc);
770 nsamples = 1;
771 }
772 else
773 {
774 npts = nr;
775 nsamples = nc;
776 }
777
778 ComplexMatrix tmp (*this);
779 const Complex *in (tmp.data ());
780 Complex *out (retval.fortran_vec ());
781
782 octave::fftw::ifft (in, out, npts, nsamples);
783
784 return retval;
785}
786
789{
790 dim_vector dv (rows (), cols ());
791
792 const double *in = data ();
793 ComplexMatrix retval (rows (), cols ());
794 octave::fftw::fftNd (in, retval.fortran_vec (), 2, dv);
795
796 return retval;
797}
798
801{
802 dim_vector dv (rows (), cols ());
803
804 ComplexMatrix retval (*this);
805 Complex *out (retval.fortran_vec ());
806
807 octave::fftw::ifftNd (out, out, 2, dv);
808
809 return retval;
810}
811
812#else
813
815Matrix::fourier (void) const
816{
817 (*current_liboctave_error_handler)
818 ("support for FFTW was unavailable or disabled when liboctave was built");
819
820 return ComplexMatrix ();
821}
822
824Matrix::ifourier (void) const
825{
826 (*current_liboctave_error_handler)
827 ("support for FFTW was unavailable or disabled when liboctave was built");
828
829 return ComplexMatrix ();
830}
831
833Matrix::fourier2d (void) const
834{
835 (*current_liboctave_error_handler)
836 ("support for FFTW was unavailable or disabled when liboctave was built");
837
838 return ComplexMatrix ();
839}
840
842Matrix::ifourier2d (void) const
843{
844 (*current_liboctave_error_handler)
845 ("support for FFTW was unavailable or disabled when liboctave was built");
846
847 return ComplexMatrix ();
848}
849
850#endif
851
852DET
854{
855 octave_idx_type info;
856 double rcon;
857 return determinant (info, rcon, 0);
858}
859
860DET
862{
863 double rcon;
864 return determinant (info, rcon, 0);
865}
866
867DET
868Matrix::determinant (octave_idx_type& info, double& rcon, bool calc_cond) const
869{
870 MatrixType mattype (*this);
871 return determinant (mattype, info, rcon, calc_cond);
872}
873
874DET
876 octave_idx_type& info, double& rcon, bool calc_cond) const
877{
878 DET retval (1.0);
879
880 info = 0;
881 rcon = 0.0;
882
883 F77_INT nr = octave::to_f77_int (rows ());
884 F77_INT nc = octave::to_f77_int (cols ());
885
886 if (nr != nc)
887 (*current_liboctave_error_handler) ("matrix must be square");
888
889 volatile int typ = mattype.type ();
890
891 // Even though the matrix is marked as singular (Rectangular), we may still
892 // get a useful number from the LU factorization, because it always completes.
893
894 if (typ == MatrixType::Unknown)
895 typ = mattype.type (*this);
896 else if (typ == MatrixType::Rectangular)
897 typ = MatrixType::Full;
898
899 if (typ == MatrixType::Lower || typ == MatrixType::Upper)
900 {
901 for (F77_INT i = 0; i < nc; i++)
902 retval *= elem (i, i);
903 }
904 else if (typ == MatrixType::Hermitian)
905 {
906 Matrix atmp = *this;
907 double *tmp_data = atmp.fortran_vec ();
908
909 // Calculate the norm of the matrix for later use when determining rcon.
910 double anorm;
911 if (calc_cond)
912 anorm = norm1 (*this);
913
914 F77_INT tmp_info = 0;
915
916 char job = 'L';
917 F77_XFCN (dpotrf, DPOTRF, (F77_CONST_CHAR_ARG2 (&job, 1), nr,
918 tmp_data, nr, tmp_info
919 F77_CHAR_ARG_LEN (1)));
920
921 info = tmp_info;
922
923 if (info != 0)
924 {
925 rcon = 0.0;
926 mattype.mark_as_unsymmetric ();
927 typ = MatrixType::Full;
928 }
929 else
930 {
931 if (calc_cond)
932 {
933 Array<double> z (dim_vector (3 * nc, 1));
934 double *pz = z.fortran_vec ();
935 Array<F77_INT> iz (dim_vector (nc, 1));
936 F77_INT *piz = iz.fortran_vec ();
937
938 F77_XFCN (dpocon, DPOCON, (F77_CONST_CHAR_ARG2 (&job, 1),
939 nr, tmp_data, nr, anorm,
940 rcon, pz, piz, tmp_info
941 F77_CHAR_ARG_LEN (1)));
942
943 info = tmp_info;
944
945 if (info != 0)
946 rcon = 0.0;
947 }
948
949 for (F77_INT i = 0; i < nc; i++)
950 retval *= atmp(i, i);
951
952 retval = retval.square ();
953 }
954 }
955 else if (typ != MatrixType::Full)
956 (*current_liboctave_error_handler) ("det: invalid dense matrix type");
957
958 if (typ == MatrixType::Full)
959 {
960 Array<F77_INT> ipvt (dim_vector (nr, 1));
961 F77_INT *pipvt = ipvt.fortran_vec ();
962
963 Matrix atmp = *this;
964 double *tmp_data = atmp.fortran_vec ();
965
966 info = 0;
967 F77_INT tmp_info = 0;
968
969 // Calculate the norm of the matrix for later use when determining rcon.
970 double anorm;
971 if (calc_cond)
972 anorm = norm1 (*this);
973
974 F77_XFCN (dgetrf, DGETRF, (nr, nr, tmp_data, nr, pipvt, tmp_info));
975
976 info = tmp_info;
977
978 // Throw away extra info LAPACK gives so as to not change output.
979 rcon = 0.0;
980 if (info != 0)
981 {
982 info = -1;
983 retval = DET ();
984 }
985 else
986 {
987 if (calc_cond)
988 {
989 // Now calc the condition number for non-singular matrix.
990 char job = '1';
991 Array<double> z (dim_vector (4 * nc, 1));
992 double *pz = z.fortran_vec ();
993 Array<F77_INT> iz (dim_vector (nc, 1));
994 F77_INT *piz = iz.fortran_vec ();
995
996 F77_XFCN (dgecon, DGECON, (F77_CONST_CHAR_ARG2 (&job, 1),
997 nc, tmp_data, nr, anorm,
998 rcon, pz, piz, tmp_info
999 F77_CHAR_ARG_LEN (1)));
1000
1001 info = tmp_info;
1002 }
1003
1004 if (info != 0)
1005 {
1006 info = -1;
1007 retval = DET ();
1008 }
1009 else
1010 {
1011 for (F77_INT i = 0; i < nc; i++)
1012 {
1013 double c = atmp(i, i);
1014 retval *= (ipvt(i) != (i+1)) ? -c : c;
1015 }
1016 }
1017 }
1018 }
1019
1020 return retval;
1021}
1022
1023double
1024Matrix::rcond (void) const
1025{
1026 MatrixType mattype (*this);
1027 return rcond (mattype);
1028}
1029
1030double
1032{
1033 double rcon = octave::numeric_limits<double>::NaN ();
1034 F77_INT nr = octave::to_f77_int (rows ());
1035 F77_INT nc = octave::to_f77_int (cols ());
1036
1037 if (nr != nc)
1038 (*current_liboctave_error_handler) ("matrix must be square");
1039
1040 if (nr == 0 || nc == 0)
1042 else
1043 {
1044 volatile int typ = mattype.type ();
1045
1046 if (typ == MatrixType::Unknown)
1047 typ = mattype.type (*this);
1048
1049 // Only calculate the condition number for LU/Cholesky
1050 if (typ == MatrixType::Upper)
1051 {
1052 const double *tmp_data = data ();
1053 F77_INT info = 0;
1054 char norm = '1';
1055 char uplo = 'U';
1056 char dia = 'N';
1057
1058 Array<double> z (dim_vector (3 * nc, 1));
1059 double *pz = z.fortran_vec ();
1060 Array<F77_INT> iz (dim_vector (nc, 1));
1061 F77_INT *piz = iz.fortran_vec ();
1062
1063 F77_XFCN (dtrcon, DTRCON, (F77_CONST_CHAR_ARG2 (&norm, 1),
1064 F77_CONST_CHAR_ARG2 (&uplo, 1),
1065 F77_CONST_CHAR_ARG2 (&dia, 1),
1066 nr, tmp_data, nr, rcon,
1067 pz, piz, info
1068 F77_CHAR_ARG_LEN (1)
1069 F77_CHAR_ARG_LEN (1)
1070 F77_CHAR_ARG_LEN (1)));
1071
1072 if (info != 0)
1073 rcon = 0.0;
1074 }
1075 else if (typ == MatrixType::Permuted_Upper)
1076 (*current_liboctave_error_handler)
1077 ("permuted triangular matrix not implemented");
1078 else if (typ == MatrixType::Lower)
1079 {
1080 const double *tmp_data = data ();
1081 F77_INT info = 0;
1082 char norm = '1';
1083 char uplo = 'L';
1084 char dia = 'N';
1085
1086 Array<double> z (dim_vector (3 * nc, 1));
1087 double *pz = z.fortran_vec ();
1088 Array<F77_INT> iz (dim_vector (nc, 1));
1089 F77_INT *piz = iz.fortran_vec ();
1090
1091 F77_XFCN (dtrcon, DTRCON, (F77_CONST_CHAR_ARG2 (&norm, 1),
1092 F77_CONST_CHAR_ARG2 (&uplo, 1),
1093 F77_CONST_CHAR_ARG2 (&dia, 1),
1094 nr, tmp_data, nr, rcon,
1095 pz, piz, info
1096 F77_CHAR_ARG_LEN (1)
1097 F77_CHAR_ARG_LEN (1)
1098 F77_CHAR_ARG_LEN (1)));
1099
1100 if (info != 0)
1101 rcon = 0.0;
1102 }
1103 else if (typ == MatrixType::Permuted_Lower)
1104 (*current_liboctave_error_handler)
1105 ("permuted triangular matrix not implemented");
1106 else if (typ == MatrixType::Full || typ == MatrixType::Hermitian)
1107 {
1108 double anorm = -1.0;
1109
1110 if (typ == MatrixType::Hermitian)
1111 {
1112 F77_INT info = 0;
1113 char job = 'L';
1114
1115 Matrix atmp = *this;
1116 double *tmp_data = atmp.fortran_vec ();
1117
1118 anorm = norm1 (atmp);
1119
1120 F77_XFCN (dpotrf, DPOTRF, (F77_CONST_CHAR_ARG2 (&job, 1), nr,
1121 tmp_data, nr, info
1122 F77_CHAR_ARG_LEN (1)));
1123
1124 if (info != 0)
1125 {
1126 rcon = 0.0;
1127 mattype.mark_as_unsymmetric ();
1128 typ = MatrixType::Full;
1129 }
1130 else
1131 {
1132 Array<double> z (dim_vector (3 * nc, 1));
1133 double *pz = z.fortran_vec ();
1134 Array<F77_INT> iz (dim_vector (nc, 1));
1135 F77_INT *piz = iz.fortran_vec ();
1136
1137 F77_XFCN (dpocon, DPOCON, (F77_CONST_CHAR_ARG2 (&job, 1),
1138 nr, tmp_data, nr, anorm,
1139 rcon, pz, piz, info
1140 F77_CHAR_ARG_LEN (1)));
1141
1142 if (info != 0)
1143 rcon = 0.0;
1144 }
1145 }
1146
1147 if (typ == MatrixType::Full)
1148 {
1149 F77_INT info = 0;
1150
1151 Matrix atmp = *this;
1152 double *tmp_data = atmp.fortran_vec ();
1153
1154 Array<F77_INT> ipvt (dim_vector (nr, 1));
1155 F77_INT *pipvt = ipvt.fortran_vec ();
1156
1157 if (anorm < 0.0)
1158 anorm = norm1 (atmp);
1159
1160 Array<double> z (dim_vector (4 * nc, 1));
1161 double *pz = z.fortran_vec ();
1162 Array<F77_INT> iz (dim_vector (nc, 1));
1163 F77_INT *piz = iz.fortran_vec ();
1164
1165 F77_XFCN (dgetrf, DGETRF, (nr, nr, tmp_data, nr, pipvt, info));
1166
1167 if (info != 0)
1168 {
1169 rcon = 0.0;
1170 mattype.mark_as_rectangular ();
1171 }
1172 else
1173 {
1174 char job = '1';
1175 F77_XFCN (dgecon, DGECON, (F77_CONST_CHAR_ARG2 (&job, 1),
1176 nc, tmp_data, nr, anorm,
1177 rcon, pz, piz, info
1178 F77_CHAR_ARG_LEN (1)));
1179
1180 if (info != 0)
1181 rcon = 0.0;
1182 }
1183 }
1184 }
1185 else
1186 rcon = 0.0;
1187 }
1188
1189 return rcon;
1190}
1191
1192Matrix
1194 double& rcon, solve_singularity_handler sing_handler,
1195 bool calc_cond, blas_trans_type transt) const
1196{
1197 Matrix retval;
1198
1199 F77_INT nr = octave::to_f77_int (rows ());
1200 F77_INT nc = octave::to_f77_int (cols ());
1201
1202 F77_INT b_nr = octave::to_f77_int (b.rows ());
1203 F77_INT b_nc = octave::to_f77_int (b.cols ());
1204
1205 if (nr != b_nr)
1206 (*current_liboctave_error_handler)
1207 ("matrix dimension mismatch solution of linear equations");
1208
1209 if (nr == 0 || nc == 0 || b_nc == 0)
1210 retval = Matrix (nc, b_nc, 0.0);
1211 else
1212 {
1213 volatile int typ = mattype.type ();
1214
1215 if (typ != MatrixType::Permuted_Upper && typ != MatrixType::Upper)
1216 (*current_liboctave_error_handler) ("incorrect matrix type");
1217
1218 rcon = 1.0;
1219 info = 0;
1220
1221 if (typ == MatrixType::Permuted_Upper)
1222 (*current_liboctave_error_handler)
1223 ("permuted triangular matrix not implemented");
1224
1225 const double *tmp_data = data ();
1226
1227 retval = b;
1228 double *result = retval.fortran_vec ();
1229
1230 char uplo = 'U';
1231 char trans = get_blas_char (transt);
1232 char dia = 'N';
1233
1234 F77_INT tmp_info = 0;
1235
1236 F77_XFCN (dtrtrs, DTRTRS, (F77_CONST_CHAR_ARG2 (&uplo, 1),
1237 F77_CONST_CHAR_ARG2 (&trans, 1),
1238 F77_CONST_CHAR_ARG2 (&dia, 1),
1239 nr, b_nc, tmp_data, nr,
1240 result, nr, tmp_info
1241 F77_CHAR_ARG_LEN (1)
1242 F77_CHAR_ARG_LEN (1)
1243 F77_CHAR_ARG_LEN (1)));
1244
1245 info = tmp_info;
1246
1247 if (calc_cond)
1248 {
1249 char norm = '1';
1250 uplo = 'U';
1251 dia = 'N';
1252
1253 Array<double> z (dim_vector (3 * nc, 1));
1254 double *pz = z.fortran_vec ();
1255 Array<F77_INT> iz (dim_vector (nc, 1));
1256 F77_INT *piz = iz.fortran_vec ();
1257
1258 F77_XFCN (dtrcon, DTRCON, (F77_CONST_CHAR_ARG2 (&norm, 1),
1259 F77_CONST_CHAR_ARG2 (&uplo, 1),
1260 F77_CONST_CHAR_ARG2 (&dia, 1),
1261 nr, tmp_data, nr, rcon,
1262 pz, piz, tmp_info
1263 F77_CHAR_ARG_LEN (1)
1264 F77_CHAR_ARG_LEN (1)
1265 F77_CHAR_ARG_LEN (1)));
1266
1267 info = tmp_info;
1268
1269 if (info != 0)
1270 info = -2;
1271
1272 // FIXME: Why calculate this, rather than just compare to 0.0?
1273 volatile double rcond_plus_one = rcon + 1.0;
1274
1275 if (rcond_plus_one == 1.0 || octave::math::isnan (rcon))
1276 {
1277 info = -2;
1278
1279 if (sing_handler)
1280 sing_handler (rcon);
1281 else
1283 }
1284 }
1285 }
1286
1287 return retval;
1288}
1289
1290Matrix
1292 double& rcon, solve_singularity_handler sing_handler,
1293 bool calc_cond, blas_trans_type transt) const
1294{
1295 Matrix retval;
1296
1297 F77_INT nr = octave::to_f77_int (rows ());
1298 F77_INT nc = octave::to_f77_int (cols ());
1299
1300 F77_INT b_nr = octave::to_f77_int (b.rows ());
1301 F77_INT b_nc = octave::to_f77_int (b.cols ());
1302
1303 if (nr != b_nr)
1304 (*current_liboctave_error_handler)
1305 ("matrix dimension mismatch solution of linear equations");
1306
1307 if (nr == 0 || nc == 0 || b_nc == 0)
1308 retval = Matrix (nc, b_nc, 0.0);
1309 else
1310 {
1311 volatile int typ = mattype.type ();
1312
1313 if (typ != MatrixType::Permuted_Lower && typ != MatrixType::Lower)
1314 (*current_liboctave_error_handler) ("incorrect matrix type");
1315
1316 rcon = 1.0;
1317 info = 0;
1318
1319 if (typ == MatrixType::Permuted_Lower)
1320 (*current_liboctave_error_handler)
1321 ("permuted triangular matrix not implemented");
1322
1323 const double *tmp_data = data ();
1324
1325 retval = b;
1326 double *result = retval.fortran_vec ();
1327
1328 char uplo = 'L';
1329 char trans = get_blas_char (transt);
1330 char dia = 'N';
1331
1332 F77_INT tmp_info = 0;
1333
1334 F77_XFCN (dtrtrs, DTRTRS, (F77_CONST_CHAR_ARG2 (&uplo, 1),
1335 F77_CONST_CHAR_ARG2 (&trans, 1),
1336 F77_CONST_CHAR_ARG2 (&dia, 1),
1337 nr, b_nc, tmp_data, nr,
1338 result, nr, tmp_info
1339 F77_CHAR_ARG_LEN (1)
1340 F77_CHAR_ARG_LEN (1)
1341 F77_CHAR_ARG_LEN (1)));
1342
1343 info = tmp_info;
1344
1345 if (calc_cond)
1346 {
1347 char norm = '1';
1348 uplo = 'L';
1349 dia = 'N';
1350
1351 Array<double> z (dim_vector (3 * nc, 1));
1352 double *pz = z.fortran_vec ();
1353 Array<F77_INT> iz (dim_vector (nc, 1));
1354 F77_INT *piz = iz.fortran_vec ();
1355
1356 F77_XFCN (dtrcon, DTRCON, (F77_CONST_CHAR_ARG2 (&norm, 1),
1357 F77_CONST_CHAR_ARG2 (&uplo, 1),
1358 F77_CONST_CHAR_ARG2 (&dia, 1),
1359 nr, tmp_data, nr, rcon,
1360 pz, piz, tmp_info
1361 F77_CHAR_ARG_LEN (1)
1362 F77_CHAR_ARG_LEN (1)
1363 F77_CHAR_ARG_LEN (1)));
1364
1365 info = tmp_info;
1366
1367 if (info != 0)
1368 info = -2;
1369
1370 volatile double rcond_plus_one = rcon + 1.0;
1371
1372 if (rcond_plus_one == 1.0 || octave::math::isnan (rcon))
1373 {
1374 info = -2;
1375
1376 if (sing_handler)
1377 sing_handler (rcon);
1378 else
1380 }
1381 }
1382 }
1383
1384 return retval;
1385}
1386
1387Matrix
1389 double& rcon, solve_singularity_handler sing_handler,
1390 bool calc_cond) const
1391{
1392 Matrix retval;
1393
1394 F77_INT nr = octave::to_f77_int (rows ());
1395 F77_INT nc = octave::to_f77_int (cols ());
1396
1397 if (nr != nc || nr != b.rows ())
1398 (*current_liboctave_error_handler)
1399 ("matrix dimension mismatch solution of linear equations");
1400
1401 if (nr == 0 || b.cols () == 0)
1402 retval = Matrix (nc, b.cols (), 0.0);
1403 else
1404 {
1405 volatile int typ = mattype.type ();
1406
1407 // Calculate the norm of the matrix for later use when determining rcon.
1408 double anorm = -1.0;
1409
1410 if (typ == MatrixType::Hermitian)
1411 {
1412 info = 0;
1413 char job = 'L';
1414
1415 Matrix atmp = *this;
1416 double *tmp_data = atmp.fortran_vec ();
1417
1418 // The norm of the matrix for later use when determining rcon.
1419 if (calc_cond)
1420 anorm = norm1 (atmp);
1421
1422 F77_INT tmp_info = 0;
1423
1424 F77_XFCN (dpotrf, DPOTRF, (F77_CONST_CHAR_ARG2 (&job, 1), nr,
1425 tmp_data, nr, tmp_info
1426 F77_CHAR_ARG_LEN (1)));
1427
1428 info = tmp_info;
1429
1430 // Throw away extra info LAPACK gives so as to not change output.
1431 rcon = 0.0;
1432 if (info != 0)
1433 {
1434 info = -2;
1435
1436 mattype.mark_as_unsymmetric ();
1437 typ = MatrixType::Full;
1438 }
1439 else
1440 {
1441 if (calc_cond)
1442 {
1443 Array<double> z (dim_vector (3 * nc, 1));
1444 double *pz = z.fortran_vec ();
1445 Array<F77_INT> iz (dim_vector (nc, 1));
1446 F77_INT *piz = iz.fortran_vec ();
1447
1448 F77_XFCN (dpocon, DPOCON, (F77_CONST_CHAR_ARG2 (&job, 1),
1449 nr, tmp_data, nr, anorm,
1450 rcon, pz, piz, tmp_info
1451 F77_CHAR_ARG_LEN (1)));
1452
1453 info = tmp_info;
1454
1455 if (info != 0)
1456 info = -2;
1457
1458 volatile double rcond_plus_one = rcon + 1.0;
1459
1460 if (rcond_plus_one == 1.0 || octave::math::isnan (rcon))
1461 {
1462 info = -2;
1463
1464 if (sing_handler)
1465 sing_handler (rcon);
1466 else
1468 }
1469 }
1470
1471 if (info == 0)
1472 {
1473 retval = b;
1474 double *result = retval.fortran_vec ();
1475
1476 F77_INT b_nr = octave::to_f77_int (b.rows ());
1477 F77_INT b_nc = octave::to_f77_int (b.cols ());
1478
1479 F77_XFCN (dpotrs, DPOTRS, (F77_CONST_CHAR_ARG2 (&job, 1),
1480 nr, b_nc, tmp_data, nr,
1481 result, b_nr, tmp_info
1482 F77_CHAR_ARG_LEN (1)));
1483
1484 info = tmp_info;
1485 }
1486 else
1487 {
1488 mattype.mark_as_unsymmetric ();
1489 typ = MatrixType::Full;
1490 }
1491 }
1492 }
1493
1494 if (typ == MatrixType::Full)
1495 {
1496 info = 0;
1497
1498 Array<F77_INT> ipvt (dim_vector (nr, 1));
1499 F77_INT *pipvt = ipvt.fortran_vec ();
1500
1501 Matrix atmp = *this;
1502 double *tmp_data = atmp.fortran_vec ();
1503
1504 if (calc_cond && anorm < 0.0)
1505 anorm = norm1 (atmp);
1506
1507 Array<double> z (dim_vector (4 * nc, 1));
1508 double *pz = z.fortran_vec ();
1509 Array<F77_INT> iz (dim_vector (nc, 1));
1510 F77_INT *piz = iz.fortran_vec ();
1511
1512 F77_INT tmp_info = 0;
1513
1514 F77_XFCN (dgetrf, DGETRF, (nr, nr, tmp_data, nr, pipvt, tmp_info));
1515
1516 info = tmp_info;
1517
1518 // Throw away extra info LAPACK gives so as to not change output.
1519 rcon = 0.0;
1520 if (info != 0)
1521 {
1522 info = -2;
1523
1524 if (sing_handler)
1525 sing_handler (rcon);
1526 else
1528
1529 mattype.mark_as_rectangular ();
1530 }
1531 else
1532 {
1533 if (calc_cond)
1534 {
1535 // Calculate the condition number for non-singular matrix.
1536 char job = '1';
1537 F77_XFCN (dgecon, DGECON, (F77_CONST_CHAR_ARG2 (&job, 1),
1538 nc, tmp_data, nr, anorm,
1539 rcon, pz, piz, tmp_info
1540 F77_CHAR_ARG_LEN (1)));
1541
1542 info = tmp_info;
1543
1544 if (info != 0)
1545 info = -2;
1546
1547 volatile double rcond_plus_one = rcon + 1.0;
1548
1549 if (rcond_plus_one == 1.0 || octave::math::isnan (rcon))
1550 {
1551 if (sing_handler)
1552 sing_handler (rcon);
1553 else
1555 }
1556 }
1557
1558 if (info == 0)
1559 {
1560 retval = b;
1561 double *result = retval.fortran_vec ();
1562
1563 F77_INT b_nr = octave::to_f77_int (b.rows ());
1564 F77_INT b_nc = octave::to_f77_int (b.cols ());
1565
1566 char job = 'N';
1567 F77_XFCN (dgetrs, DGETRS, (F77_CONST_CHAR_ARG2 (&job, 1),
1568 nr, b_nc, tmp_data, nr,
1569 pipvt, result, b_nr, tmp_info
1570 F77_CHAR_ARG_LEN (1)));
1571
1572 info = tmp_info;
1573 }
1574 else
1575 mattype.mark_as_rectangular ();
1576 }
1577 }
1578 else if (typ != MatrixType::Hermitian)
1579 (*current_liboctave_error_handler) ("incorrect matrix type");
1580 }
1581
1582 return retval;
1583}
1584
1585Matrix
1586Matrix::solve (MatrixType& mattype, const Matrix& b) const
1587{
1588 octave_idx_type info;
1589 double rcon;
1590 return solve (mattype, b, info, rcon, nullptr);
1591}
1592
1593Matrix
1594Matrix::solve (MatrixType& mattype, const Matrix& b,
1595 octave_idx_type& info) const
1596{
1597 double rcon;
1598 return solve (mattype, b, info, rcon, nullptr);
1599}
1600
1601Matrix
1603 double& rcon) const
1604{
1605 return solve (mattype, b, info, rcon, nullptr);
1606}
1607
1608Matrix
1610 double& rcon, solve_singularity_handler sing_handler,
1611 bool singular_fallback, blas_trans_type transt) const
1612{
1613 Matrix retval;
1614 int typ = mattype.type ();
1615
1616 if (typ == MatrixType::Unknown)
1617 typ = mattype.type (*this);
1618
1619 // Only calculate the condition number for LU/Cholesky
1620 if (typ == MatrixType::Upper || typ == MatrixType::Permuted_Upper)
1621 retval = utsolve (mattype, b, info, rcon, sing_handler, true, transt);
1622 else if (typ == MatrixType::Lower || typ == MatrixType::Permuted_Lower)
1623 retval = ltsolve (mattype, b, info, rcon, sing_handler, true, transt);
1624 else if (transt == blas_trans || transt == blas_conj_trans)
1625 return transpose ().solve (mattype, b, info, rcon, sing_handler,
1626 singular_fallback);
1627 else if (typ == MatrixType::Full || typ == MatrixType::Hermitian)
1628 retval = fsolve (mattype, b, info, rcon, sing_handler, true);
1629 else if (typ != MatrixType::Rectangular)
1630 (*current_liboctave_error_handler) ("unknown matrix type");
1631
1632 // Rectangular or one of the above solvers flags a singular matrix
1633 if (singular_fallback && mattype.type () == MatrixType::Rectangular)
1634 {
1635 octave_idx_type rank;
1636 retval = lssolve (b, info, rank, rcon);
1637 }
1638
1639 return retval;
1640}
1641
1643Matrix::solve (MatrixType& mattype, const ComplexMatrix& b) const
1644{
1645 octave_idx_type info;
1646 double rcon;
1647 return solve (mattype, b, info, rcon, nullptr);
1648}
1649
1652 octave_idx_type& info) const
1653{
1654 double rcon;
1655 return solve (mattype, b, info, rcon, nullptr);
1656}
1657
1660 octave_idx_type& info, double& rcon) const
1661{
1662 return solve (mattype, b, info, rcon, nullptr);
1663}
1664
1665static Matrix
1667{
1668 octave_idx_type m = cm.rows ();
1669 octave_idx_type n = cm.cols ();
1670 octave_idx_type nel = m*n;
1671 Matrix retval (m, 2*n);
1672 const Complex *cmd = cm.data ();
1673 double *rd = retval.fortran_vec ();
1674 for (octave_idx_type i = 0; i < nel; i++)
1675 {
1676 rd[i] = std::real (cmd[i]);
1677 rd[nel+i] = std::imag (cmd[i]);
1678 }
1679 return retval;
1680}
1681
1682static ComplexMatrix
1684{
1685 octave_idx_type m = sm.rows ();
1686 octave_idx_type n = sm.cols () / 2;
1687 octave_idx_type nel = m*n;
1688 ComplexMatrix retval (m, n);
1689 const double *smd = sm.data ();
1690 Complex *rd = retval.fortran_vec ();
1691 for (octave_idx_type i = 0; i < nel; i++)
1692 rd[i] = Complex (smd[i], smd[nel+i]);
1693 return retval;
1694}
1695
1698 octave_idx_type& info, double& rcon,
1699 solve_singularity_handler sing_handler,
1700 bool singular_fallback, blas_trans_type transt) const
1701{
1702 Matrix tmp = stack_complex_matrix (b);
1703 tmp = solve (mattype, tmp, info, rcon, sing_handler, singular_fallback,
1704 transt);
1705 return unstack_complex_matrix (tmp);
1706}
1707
1709Matrix::solve (MatrixType& mattype, const ColumnVector& b) const
1710{
1711 octave_idx_type info; double rcon;
1712 return solve (mattype, b, info, rcon);
1713}
1714
1717 octave_idx_type& info) const
1718{
1719 double rcon;
1720 return solve (mattype, b, info, rcon);
1721}
1722
1725 octave_idx_type& info, double& rcon) const
1726{
1727 return solve (mattype, b, info, rcon, nullptr);
1728}
1729
1732 octave_idx_type& info, double& rcon,
1733 solve_singularity_handler sing_handler,
1734 blas_trans_type transt) const
1735{
1736 Matrix tmp (b);
1737 tmp = solve (mattype, tmp, info, rcon, sing_handler, true, transt);
1738 return tmp.column (static_cast<octave_idx_type> (0));
1739}
1740
1743{
1744 ComplexMatrix tmp (*this);
1745 return tmp.solve (mattype, b);
1746}
1747
1750 octave_idx_type& info) const
1751{
1752 ComplexMatrix tmp (*this);
1753 return tmp.solve (mattype, b, info);
1754}
1755
1758 octave_idx_type& info, double& rcon) const
1759{
1760 ComplexMatrix tmp (*this);
1761 return tmp.solve (mattype, b, info, rcon);
1762}
1763
1766 octave_idx_type& info, double& rcon,
1767 solve_singularity_handler sing_handler,
1768 blas_trans_type transt) const
1769{
1770 ComplexMatrix tmp (*this);
1771 return tmp.solve (mattype, b, info, rcon, sing_handler, transt);
1772}
1773
1774Matrix
1775Matrix::solve (const Matrix& b) const
1776{
1777 octave_idx_type info;
1778 double rcon;
1779 return solve (b, info, rcon, nullptr);
1780}
1781
1782Matrix
1783Matrix::solve (const Matrix& b, octave_idx_type& info) const
1784{
1785 double rcon;
1786 return solve (b, info, rcon, nullptr);
1787}
1788
1789Matrix
1790Matrix::solve (const Matrix& b, octave_idx_type& info, double& rcon) const
1791{
1792 return solve (b, info, rcon, nullptr);
1793}
1794
1795Matrix
1797 double& rcon, solve_singularity_handler sing_handler,
1798 blas_trans_type transt) const
1799{
1800 MatrixType mattype (*this);
1801 return solve (mattype, b, info, rcon, sing_handler, true, transt);
1802}
1803
1806{
1807 ComplexMatrix tmp (*this);
1808 return tmp.solve (b);
1809}
1810
1813{
1814 ComplexMatrix tmp (*this);
1815 return tmp.solve (b, info);
1816}
1817
1820 double& rcon) const
1821{
1822 ComplexMatrix tmp (*this);
1823 return tmp.solve (b, info, rcon);
1824}
1825
1827Matrix::solve (const ComplexMatrix& b, octave_idx_type& info, double& rcon,
1828 solve_singularity_handler sing_handler,
1829 blas_trans_type transt) const
1830{
1831 ComplexMatrix tmp (*this);
1832 return tmp.solve (b, info, rcon, sing_handler, transt);
1833}
1834
1837{
1838 octave_idx_type info; double rcon;
1839 return solve (b, info, rcon);
1840}
1841
1844{
1845 double rcon;
1846 return solve (b, info, rcon);
1847}
1848
1850Matrix::solve (const ColumnVector& b, octave_idx_type& info, double& rcon) const
1851{
1852 return solve (b, info, rcon, nullptr);
1853}
1854
1856Matrix::solve (const ColumnVector& b, octave_idx_type& info, double& rcon,
1857 solve_singularity_handler sing_handler,
1858 blas_trans_type transt) const
1859{
1860 MatrixType mattype (*this);
1861 return solve (mattype, b, info, rcon, sing_handler, transt);
1862}
1863
1866{
1867 ComplexMatrix tmp (*this);
1868 return tmp.solve (b);
1869}
1870
1873{
1874 ComplexMatrix tmp (*this);
1875 return tmp.solve (b, info);
1876}
1877
1880 double& rcon) const
1881{
1882 ComplexMatrix tmp (*this);
1883 return tmp.solve (b, info, rcon);
1884}
1885
1888 double& rcon,
1889 solve_singularity_handler sing_handler,
1890 blas_trans_type transt) const
1891{
1892 ComplexMatrix tmp (*this);
1893 return tmp.solve (b, info, rcon, sing_handler, transt);
1894}
1895
1896Matrix
1897Matrix::lssolve (const Matrix& b) const
1898{
1899 octave_idx_type info;
1900 octave_idx_type rank;
1901 double rcon;
1902 return lssolve (b, info, rank, rcon);
1903}
1904
1905Matrix
1907{
1908 octave_idx_type rank;
1909 double rcon;
1910 return lssolve (b, info, rank, rcon);
1911}
1912
1913Matrix
1915 octave_idx_type& rank) const
1916{
1917 double rcon;
1918 return lssolve (b, info, rank, rcon);
1919}
1920
1921Matrix
1923 octave_idx_type& rank, double& rcon) const
1924{
1925 Matrix retval;
1926
1927 F77_INT m = octave::to_f77_int (rows ());
1928 F77_INT n = octave::to_f77_int (cols ());
1929
1930 F77_INT b_nr = octave::to_f77_int (b.rows ());
1931 F77_INT b_nc = octave::to_f77_int (b.cols ());
1932 F77_INT nrhs = b_nc; // alias for code readability
1933
1934 if (m != b_nr)
1935 (*current_liboctave_error_handler)
1936 ("matrix dimension mismatch solution of linear equations");
1937
1938 if (m == 0 || n == 0 || b_nc == 0)
1939 retval = Matrix (n, b_nc, 0.0);
1940 else
1941 {
1942 volatile F77_INT minmn = (m < n ? m : n);
1943 F77_INT maxmn = (m > n ? m : n);
1944 rcon = -1.0;
1945 if (m != n)
1946 {
1947 retval = Matrix (maxmn, nrhs, 0.0);
1948
1949 for (F77_INT j = 0; j < nrhs; j++)
1950 for (F77_INT i = 0; i < m; i++)
1951 retval.elem (i, j) = b.elem (i, j);
1952 }
1953 else
1954 retval = b;
1955
1956 Matrix atmp = *this;
1957 double *tmp_data = atmp.fortran_vec ();
1958
1959 double *pretval = retval.fortran_vec ();
1960 Array<double> s (dim_vector (minmn, 1));
1961 double *ps = s.fortran_vec ();
1962
1963 // Ask DGELSD what the dimension of WORK should be.
1964 F77_INT lwork = -1;
1965
1966 Array<double> work (dim_vector (1, 1));
1967
1968 F77_INT smlsiz;
1969 F77_FUNC (xilaenv, XILAENV) (9, F77_CONST_CHAR_ARG2 ("DGELSD", 6),
1970 F77_CONST_CHAR_ARG2 (" ", 1),
1971 0, 0, 0, 0, smlsiz
1972 F77_CHAR_ARG_LEN (6)
1973 F77_CHAR_ARG_LEN (1));
1974
1975 F77_INT mnthr;
1976 F77_FUNC (xilaenv, XILAENV) (6, F77_CONST_CHAR_ARG2 ("DGELSD", 6),
1977 F77_CONST_CHAR_ARG2 (" ", 1),
1978 m, n, nrhs, -1, mnthr
1979 F77_CHAR_ARG_LEN (6)
1980 F77_CHAR_ARG_LEN (1));
1981
1982 // We compute the size of iwork because DGELSD in older versions
1983 // of LAPACK does not return it on a query call.
1984 double dminmn = static_cast<double> (minmn);
1985 double dsmlsizp1 = static_cast<double> (smlsiz+1);
1986 double tmp = octave::math::log2 (dminmn / dsmlsizp1);
1987
1988 F77_INT nlvl = static_cast<F77_INT> (tmp) + 1;
1989 if (nlvl < 0)
1990 nlvl = 0;
1991
1992 F77_INT liwork = 3 * minmn * nlvl + 11 * minmn;
1993 if (liwork < 1)
1994 liwork = 1;
1995 Array<F77_INT> iwork (dim_vector (liwork, 1));
1996 F77_INT *piwork = iwork.fortran_vec ();
1997
1998 F77_INT tmp_info = 0;
1999 F77_INT tmp_rank = 0;
2000
2001 F77_XFCN (dgelsd, DGELSD, (m, n, nrhs, tmp_data, m, pretval, maxmn,
2002 ps, rcon, tmp_rank, work.fortran_vec (),
2003 lwork, piwork, tmp_info));
2004
2005 info = tmp_info;
2006 rank = tmp_rank;
2007
2008 // The workspace query is broken in at least LAPACK 3.0.0
2009 // through 3.1.1 when n >= mnthr. The obtuse formula below
2010 // should provide sufficient workspace for DGELSD to operate
2011 // efficiently.
2012 if (n > m && n >= mnthr)
2013 {
2014 const F77_INT wlalsd
2015 = 9*m + 2*m*smlsiz + 8*m*nlvl + m*nrhs + (smlsiz+1)*(smlsiz+1);
2016
2017 F77_INT addend = m;
2018
2019 if (2*m-4 > addend)
2020 addend = 2*m-4;
2021
2022 if (nrhs > addend)
2023 addend = nrhs;
2024
2025 if (n-3*m > addend)
2026 addend = n-3*m;
2027
2028 if (wlalsd > addend)
2029 addend = wlalsd;
2030
2031 const F77_INT lworkaround = 4*m + m*m + addend;
2032
2033 if (work(0) < lworkaround)
2034 work(0) = lworkaround;
2035 }
2036 else if (m >= n)
2037 {
2038 F77_INT lworkaround
2039 = 12*n + 2*n*smlsiz + 8*n*nlvl + n*nrhs + (smlsiz+1)*(smlsiz+1);
2040
2041 if (work(0) < lworkaround)
2042 work(0) = lworkaround;
2043 }
2044
2045 lwork = static_cast<F77_INT> (work(0));
2046 work.resize (dim_vector (lwork, 1));
2047
2048 double anorm = norm1 (*this);
2049
2050 if (octave::math::isinf (anorm))
2051 {
2052 rcon = 0.0;
2053 retval = Matrix (n, b_nc, 0.0);
2054 }
2055 else if (octave::math::isnan (anorm))
2056 {
2058 retval = Matrix (n, b_nc, octave::numeric_limits<double>::NaN ());
2059 }
2060 else
2061 {
2062 F77_XFCN (dgelsd, DGELSD, (m, n, nrhs, tmp_data, m, pretval,
2063 maxmn, ps, rcon, tmp_rank,
2064 work.fortran_vec (), lwork,
2065 piwork, tmp_info));
2066
2067 info = tmp_info;
2068 rank = tmp_rank;
2069
2070 if (s.elem (0) == 0.0)
2071 rcon = 0.0;
2072 else
2073 rcon = s.elem (minmn - 1) / s.elem (0);
2074
2075 retval.resize (n, nrhs);
2076 }
2077 }
2078
2079 return retval;
2080}
2081
2084{
2085 ComplexMatrix tmp (*this);
2086 octave_idx_type info;
2087 octave_idx_type rank;
2088 double rcon;
2089 return tmp.lssolve (b, info, rank, rcon);
2090}
2091
2094{
2095 ComplexMatrix tmp (*this);
2096 octave_idx_type rank;
2097 double rcon;
2098 return tmp.lssolve (b, info, rank, rcon);
2099}
2100
2103 octave_idx_type& rank) const
2104{
2105 ComplexMatrix tmp (*this);
2106 double rcon;
2107 return tmp.lssolve (b, info, rank, rcon);
2108}
2109
2112 octave_idx_type& rank, double& rcon) const
2113{
2114 ComplexMatrix tmp (*this);
2115 return tmp.lssolve (b, info, rank, rcon);
2116}
2117
2120{
2121 octave_idx_type info;
2122 octave_idx_type rank;
2123 double rcon;
2124 return lssolve (b, info, rank, rcon);
2125}
2126
2129{
2130 octave_idx_type rank;
2131 double rcon;
2132 return lssolve (b, info, rank, rcon);
2133}
2134
2137 octave_idx_type& rank) const
2138{
2139 double rcon;
2140 return lssolve (b, info, rank, rcon);
2141}
2142
2145 octave_idx_type& rank, double& rcon) const
2146{
2147 ColumnVector retval;
2148
2149 F77_INT nrhs = 1;
2150
2151 F77_INT m = octave::to_f77_int (rows ());
2152 F77_INT n = octave::to_f77_int (cols ());
2153
2154 F77_INT b_nel = octave::to_f77_int (b.numel ());
2155
2156 if (m != b_nel)
2157 (*current_liboctave_error_handler)
2158 ("matrix dimension mismatch solution of linear equations");
2159
2160 if (m == 0 || n == 0)
2161 retval = ColumnVector (n, 0.0);
2162 else
2163 {
2164 volatile F77_INT minmn = (m < n ? m : n);
2165 F77_INT maxmn = (m > n ? m : n);
2166 rcon = -1.0;
2167
2168 if (m != n)
2169 {
2170 retval = ColumnVector (maxmn, 0.0);
2171
2172 for (F77_INT i = 0; i < m; i++)
2173 retval.elem (i) = b.elem (i);
2174 }
2175 else
2176 retval = b;
2177
2178 Matrix atmp = *this;
2179 double *tmp_data = atmp.fortran_vec ();
2180
2181 double *pretval = retval.fortran_vec ();
2182 Array<double> s (dim_vector (minmn, 1));
2183 double *ps = s.fortran_vec ();
2184
2185 // Ask DGELSD what the dimension of WORK should be.
2186 F77_INT lwork = -1;
2187
2188 Array<double> work (dim_vector (1, 1));
2189
2190 F77_INT smlsiz;
2191 F77_FUNC (xilaenv, XILAENV) (9, F77_CONST_CHAR_ARG2 ("DGELSD", 6),
2192 F77_CONST_CHAR_ARG2 (" ", 1),
2193 0, 0, 0, 0, smlsiz
2194 F77_CHAR_ARG_LEN (6)
2195 F77_CHAR_ARG_LEN (1));
2196
2197 // We compute the size of iwork because DGELSD in older versions
2198 // of LAPACK does not return it on a query call.
2199 double dminmn = static_cast<double> (minmn);
2200 double dsmlsizp1 = static_cast<double> (smlsiz+1);
2201 double tmp = octave::math::log2 (dminmn / dsmlsizp1);
2202
2203 F77_INT nlvl = static_cast<F77_INT> (tmp) + 1;
2204 if (nlvl < 0)
2205 nlvl = 0;
2206
2207 F77_INT liwork = 3 * minmn * nlvl + 11 * minmn;
2208 if (liwork < 1)
2209 liwork = 1;
2210 Array<F77_INT> iwork (dim_vector (liwork, 1));
2211 F77_INT *piwork = iwork.fortran_vec ();
2212
2213 F77_INT tmp_info = 0;
2214 F77_INT tmp_rank = 0;
2215
2216 F77_XFCN (dgelsd, DGELSD, (m, n, nrhs, tmp_data, m, pretval, maxmn,
2217 ps, rcon, tmp_rank, work.fortran_vec (),
2218 lwork, piwork, tmp_info));
2219
2220 info = tmp_info;
2221 rank = tmp_rank;
2222
2223 lwork = static_cast<F77_INT> (work(0));
2224 work.resize (dim_vector (lwork, 1));
2225
2226 F77_XFCN (dgelsd, DGELSD, (m, n, nrhs, tmp_data, m, pretval,
2227 maxmn, ps, rcon, tmp_rank,
2228 work.fortran_vec (), lwork,
2229 piwork, tmp_info));
2230
2231 info = tmp_info;
2232 rank = tmp_rank;
2233
2234 if (rank < minmn)
2235 {
2236 if (s.elem (0) == 0.0)
2237 rcon = 0.0;
2238 else
2239 rcon = s.elem (minmn - 1) / s.elem (0);
2240 }
2241
2242 retval.resize (n);
2243 }
2244
2245 return retval;
2246}
2247
2250{
2251 ComplexMatrix tmp (*this);
2252 octave_idx_type info;
2253 octave_idx_type rank;
2254 double rcon;
2255 return tmp.lssolve (b, info, rank, rcon);
2256}
2257
2260{
2261 ComplexMatrix tmp (*this);
2262 octave_idx_type rank;
2263 double rcon;
2264 return tmp.lssolve (b, info, rank, rcon);
2265}
2266
2269 octave_idx_type& rank) const
2270{
2271 ComplexMatrix tmp (*this);
2272 double rcon;
2273 return tmp.lssolve (b, info, rank, rcon);
2274}
2275
2278 octave_idx_type& rank, double& rcon) const
2279{
2280 ComplexMatrix tmp (*this);
2281 return tmp.lssolve (b, info, rank, rcon);
2282}
2283
2284Matrix&
2286{
2287 octave_idx_type nr = rows ();
2288 octave_idx_type nc = cols ();
2289
2290 octave_idx_type a_nr = a.rows ();
2291 octave_idx_type a_nc = a.cols ();
2292
2293 if (nr != a_nr || nc != a_nc)
2294 octave::err_nonconformant ("operator +=", nr, nc, a_nr, a_nc);
2295
2296 for (octave_idx_type i = 0; i < a.length (); i++)
2297 elem (i, i) += a.elem (i, i);
2298
2299 return *this;
2300}
2301
2302Matrix&
2304{
2305 octave_idx_type nr = rows ();
2306 octave_idx_type nc = cols ();
2307
2308 octave_idx_type a_nr = a.rows ();
2309 octave_idx_type a_nc = a.cols ();
2310
2311 if (nr != a_nr || nc != a_nc)
2312 octave::err_nonconformant ("operator -=", nr, nc, a_nr, a_nc);
2313
2314 for (octave_idx_type i = 0; i < a.length (); i++)
2315 elem (i, i) -= a.elem (i, i);
2316
2317 return *this;
2318}
2319
2320// unary operations
2321
2322// column vector by row vector -> matrix operations
2323
2324Matrix
2326{
2327 Matrix retval;
2328
2329 F77_INT len = octave::to_f77_int (v.numel ());
2330
2331 if (len != 0)
2332 {
2333 F77_INT a_len = octave::to_f77_int (a.numel ());
2334
2335 retval = Matrix (len, a_len);
2336 double *c = retval.fortran_vec ();
2337
2338 F77_XFCN (dgemm, DGEMM, (F77_CONST_CHAR_ARG2 ("N", 1),
2339 F77_CONST_CHAR_ARG2 ("N", 1),
2340 len, a_len, 1, 1.0, v.data (), len,
2341 a.data (), 1, 0.0, c, len
2342 F77_CHAR_ARG_LEN (1)
2343 F77_CHAR_ARG_LEN (1)));
2344 }
2345
2346 return retval;
2347}
2348
2349// other operations.
2350
2351// FIXME: Do these really belong here? Maybe they should be in a base class?
2352
2354Matrix::all (int dim) const
2355{
2356 return NDArray::all (dim);
2357}
2358
2360Matrix::any (int dim) const
2361{
2362 return NDArray::any (dim);
2363}
2364
2365Matrix
2366Matrix::cumprod (int dim) const
2367{
2368 return NDArray::cumprod (dim);
2369}
2370
2371Matrix
2372Matrix::cumsum (int dim) const
2373{
2374 return NDArray::cumsum (dim);
2375}
2376
2377Matrix
2378Matrix::prod (int dim) const
2379{
2380 return NDArray::prod (dim);
2381}
2382
2383Matrix
2384Matrix::sum (int dim) const
2385{
2386 return NDArray::sum (dim);
2387}
2388
2389Matrix
2390Matrix::sumsq (int dim) const
2391{
2392 return NDArray::sumsq (dim);
2393}
2394
2395Matrix
2396Matrix::abs (void) const
2397{
2398 return NDArray::abs ();
2399}
2400
2401Matrix
2403{
2404 return NDArray::diag (k);
2405}
2406
2409{
2410 DiagMatrix retval;
2411
2412 octave_idx_type nr = rows ();
2413 octave_idx_type nc = cols ();
2414
2415 if (nr == 1 || nc == 1)
2416 retval = DiagMatrix (*this, m, n);
2417 else
2418 (*current_liboctave_error_handler) ("diag: expecting vector argument");
2419
2420 return retval;
2421}
2422
2425{
2426 Array<octave_idx_type> dummy_idx;
2427 return row_min (dummy_idx);
2428}
2429
2432{
2433 ColumnVector result;
2434
2435 octave_idx_type nr = rows ();
2436 octave_idx_type nc = cols ();
2437
2438 if (nr > 0 && nc > 0)
2439 {
2440 result.resize (nr);
2441 idx_arg.resize (dim_vector (nr, 1));
2442
2443 for (octave_idx_type i = 0; i < nr; i++)
2444 {
2445 octave_idx_type idx_j;
2446
2447 double tmp_min = octave::numeric_limits<double>::NaN ();
2448
2449 for (idx_j = 0; idx_j < nc; idx_j++)
2450 {
2451 tmp_min = elem (i, idx_j);
2452
2453 if (! octave::math::isnan (tmp_min))
2454 break;
2455 }
2456
2457 for (octave_idx_type j = idx_j+1; j < nc; j++)
2458 {
2459 double tmp = elem (i, j);
2460
2461 if (octave::math::isnan (tmp))
2462 continue;
2463 else if (tmp < tmp_min)
2464 {
2465 idx_j = j;
2466 tmp_min = tmp;
2467 }
2468 }
2469
2470 result.elem (i) = tmp_min;
2471 idx_arg.elem (i) = (octave::math::isnan (tmp_min) ? 0 : idx_j);
2472 }
2473 }
2474
2475 return result;
2476}
2477
2480{
2481 Array<octave_idx_type> dummy_idx;
2482 return row_max (dummy_idx);
2483}
2484
2487{
2488 ColumnVector result;
2489
2490 octave_idx_type nr = rows ();
2491 octave_idx_type nc = cols ();
2492
2493 if (nr > 0 && nc > 0)
2494 {
2495 result.resize (nr);
2496 idx_arg.resize (dim_vector (nr, 1));
2497
2498 for (octave_idx_type i = 0; i < nr; i++)
2499 {
2500 octave_idx_type idx_j;
2501
2502 double tmp_max = octave::numeric_limits<double>::NaN ();
2503
2504 for (idx_j = 0; idx_j < nc; idx_j++)
2505 {
2506 tmp_max = elem (i, idx_j);
2507
2508 if (! octave::math::isnan (tmp_max))
2509 break;
2510 }
2511
2512 for (octave_idx_type j = idx_j+1; j < nc; j++)
2513 {
2514 double tmp = elem (i, j);
2515
2516 if (octave::math::isnan (tmp))
2517 continue;
2518 else if (tmp > tmp_max)
2519 {
2520 idx_j = j;
2521 tmp_max = tmp;
2522 }
2523 }
2524
2525 result.elem (i) = tmp_max;
2526 idx_arg.elem (i) = (octave::math::isnan (tmp_max) ? 0 : idx_j);
2527 }
2528 }
2529
2530 return result;
2531}
2532
2535{
2536 Array<octave_idx_type> dummy_idx;
2537 return column_min (dummy_idx);
2538}
2539
2542{
2543 RowVector result;
2544
2545 octave_idx_type nr = rows ();
2546 octave_idx_type nc = cols ();
2547
2548 if (nr > 0 && nc > 0)
2549 {
2550 result.resize (nc);
2551 idx_arg.resize (dim_vector (1, nc));
2552
2553 for (octave_idx_type j = 0; j < nc; j++)
2554 {
2555 octave_idx_type idx_i;
2556
2557 double tmp_min = octave::numeric_limits<double>::NaN ();
2558
2559 for (idx_i = 0; idx_i < nr; idx_i++)
2560 {
2561 tmp_min = elem (idx_i, j);
2562
2563 if (! octave::math::isnan (tmp_min))
2564 break;
2565 }
2566
2567 for (octave_idx_type i = idx_i+1; i < nr; i++)
2568 {
2569 double tmp = elem (i, j);
2570
2571 if (octave::math::isnan (tmp))
2572 continue;
2573 else if (tmp < tmp_min)
2574 {
2575 idx_i = i;
2576 tmp_min = tmp;
2577 }
2578 }
2579
2580 result.elem (j) = tmp_min;
2581 idx_arg.elem (j) = (octave::math::isnan (tmp_min) ? 0 : idx_i);
2582 }
2583 }
2584
2585 return result;
2586}
2587
2590{
2591 Array<octave_idx_type> dummy_idx;
2592 return column_max (dummy_idx);
2593}
2594
2597{
2598 RowVector result;
2599
2600 octave_idx_type nr = rows ();
2601 octave_idx_type nc = cols ();
2602
2603 if (nr > 0 && nc > 0)
2604 {
2605 result.resize (nc);
2606 idx_arg.resize (dim_vector (1, nc));
2607
2608 for (octave_idx_type j = 0; j < nc; j++)
2609 {
2610 octave_idx_type idx_i;
2611
2612 double tmp_max = octave::numeric_limits<double>::NaN ();
2613
2614 for (idx_i = 0; idx_i < nr; idx_i++)
2615 {
2616 tmp_max = elem (idx_i, j);
2617
2618 if (! octave::math::isnan (tmp_max))
2619 break;
2620 }
2621
2622 for (octave_idx_type i = idx_i+1; i < nr; i++)
2623 {
2624 double tmp = elem (i, j);
2625
2626 if (octave::math::isnan (tmp))
2627 continue;
2628 else if (tmp > tmp_max)
2629 {
2630 idx_i = i;
2631 tmp_max = tmp;
2632 }
2633 }
2634
2635 result.elem (j) = tmp_max;
2636 idx_arg.elem (j) = (octave::math::isnan (tmp_max) ? 0 : idx_i);
2637 }
2638 }
2639
2640 return result;
2641}
2642
2643std::ostream&
2644operator << (std::ostream& os, const Matrix& a)
2645{
2646 for (octave_idx_type i = 0; i < a.rows (); i++)
2647 {
2648 for (octave_idx_type j = 0; j < a.cols (); j++)
2649 {
2650 os << ' ';
2651 octave::write_value<double> (os, a.elem (i, j));
2652 }
2653 os << "\n";
2654 }
2655 return os;
2656}
2657
2658std::istream&
2659operator >> (std::istream& is, Matrix& a)
2660{
2661 octave_idx_type nr = a.rows ();
2662 octave_idx_type nc = a.cols ();
2663
2664 if (nr > 0 && nc > 0)
2665 {
2666 double tmp;
2667 for (octave_idx_type i = 0; i < nr; i++)
2668 for (octave_idx_type j = 0; j < nc; j++)
2669 {
2670 tmp = octave::read_value<double> (is);
2671 if (is)
2672 a.elem (i, j) = tmp;
2673 else
2674 return is;
2675 }
2676 }
2677
2678 return is;
2679}
2680
2681Matrix
2682Givens (double x, double y)
2683{
2684 double cc, s, temp_r;
2685
2686 F77_FUNC (dlartg, DLARTG) (x, y, cc, s, temp_r);
2687
2688 Matrix g (2, 2);
2689
2690 g.elem (0, 0) = cc;
2691 g.elem (1, 1) = cc;
2692 g.elem (0, 1) = s;
2693 g.elem (1, 0) = -s;
2694
2695 return g;
2696}
2697
2698Matrix
2699Sylvester (const Matrix& a, const Matrix& b, const Matrix& c)
2700{
2701 Matrix retval;
2702
2703 // FIXME: need to check that a, b, and c are all the same size.
2704
2705 // Compute Schur decompositions.
2706
2707 octave::math::schur<Matrix> as (a, "U");
2708 octave::math::schur<Matrix> bs (b, "U");
2709
2710 // Transform c to new coordinates.
2711
2712 Matrix ua = as.unitary_schur_matrix ();
2713 Matrix sch_a = as.schur_matrix ();
2714
2715 Matrix ub = bs.unitary_schur_matrix ();
2716 Matrix sch_b = bs.schur_matrix ();
2717
2718 Matrix cx = ua.transpose () * c * ub;
2719
2720 // Solve the sylvester equation, back-transform, and return the solution.
2721
2722 F77_INT a_nr = octave::to_f77_int (a.rows ());
2723 F77_INT b_nr = octave::to_f77_int (b.rows ());
2724
2725 double scale;
2726 F77_INT info;
2727
2728 double *pa = sch_a.fortran_vec ();
2729 double *pb = sch_b.fortran_vec ();
2730 double *px = cx.fortran_vec ();
2731
2732 F77_XFCN (dtrsyl, DTRSYL, (F77_CONST_CHAR_ARG2 ("N", 1),
2733 F77_CONST_CHAR_ARG2 ("N", 1),
2734 1, a_nr, b_nr, pa, a_nr, pb,
2735 b_nr, px, a_nr, scale, info
2736 F77_CHAR_ARG_LEN (1)
2737 F77_CHAR_ARG_LEN (1)));
2738
2739 // FIXME: check info?
2740
2741 retval = ua*cx*ub.transpose ();
2742
2743 return retval;
2744}
2745
2746// matrix by matrix -> matrix operations
2747
2748/*
2749
2750## Simple Dot Product, Matrix-Vector and Matrix-Matrix Unit tests
2751%!assert ([1 2 3] * [ 4 ; 5 ; 6], 32, 1e-14)
2752%!assert ([1 2 ; 3 4 ] * [5 ; 6], [17 ; 39 ], 1e-14)
2753%!assert ([1 2 ; 3 4 ] * [5 6 ; 7 8], [19 22; 43 50], 1e-14)
2754
2755## Test some simple identities
2756%!shared M, cv, rv, Mt, rvt
2757%! M = randn (10,10) + 100*eye (10,10);
2758%! Mt = M';
2759%! cv = randn (10,1);
2760%! rv = randn (1,10);
2761%! rvt = rv';
2762%!assert ([M*cv,M*cv], M*[cv,cv], 2e-13)
2763%!assert ([M'*cv,M'*cv], M'*[cv,cv], 2e-13)
2764%!assert ([rv*M;rv*M], [rv;rv]*M, 2e-13)
2765%!assert ([rv*M';rv*M'], [rv;rv]*M', 2e-13)
2766%!assert (2*rv*cv, [rv,rv]*[cv;cv], 2e-13)
2767%!assert (M'\cv, Mt\cv, 1e-14)
2768%!assert (M'\rv', Mt\rvt, 1e-14)
2769
2770*/
2771
2772static inline char
2774{
2775 return trans ? 'T' : 'N';
2776}
2777
2778// the general GEMM operation
2779
2780Matrix
2781xgemm (const Matrix& a, const Matrix& b,
2782 blas_trans_type transa, blas_trans_type transb)
2783{
2784 Matrix retval;
2785
2786 bool tra = transa != blas_no_trans;
2787 bool trb = transb != blas_no_trans;
2788
2789 F77_INT a_nr = octave::to_f77_int (tra ? a.cols () : a.rows ());
2790 F77_INT a_nc = octave::to_f77_int (tra ? a.rows () : a.cols ());
2791
2792 F77_INT b_nr = octave::to_f77_int (trb ? b.cols () : b.rows ());
2793 F77_INT b_nc = octave::to_f77_int (trb ? b.rows () : b.cols ());
2794
2795 if (a_nc != b_nr)
2796 octave::err_nonconformant ("operator *", a_nr, a_nc, b_nr, b_nc);
2797
2798 if (a_nr == 0 || a_nc == 0 || b_nc == 0)
2799 retval = Matrix (a_nr, b_nc, 0.0);
2800 else if (a.data () == b.data () && a_nr == b_nc && tra != trb)
2801 {
2802 F77_INT lda = octave::to_f77_int (a.rows ());
2803
2804 retval = Matrix (a_nr, b_nc);
2805 double *c = retval.fortran_vec ();
2806
2807 const char ctra = get_blas_trans_arg (tra);
2808 F77_XFCN (dsyrk, DSYRK, (F77_CONST_CHAR_ARG2 ("U", 1),
2809 F77_CONST_CHAR_ARG2 (&ctra, 1),
2810 a_nr, a_nc, 1.0,
2811 a.data (), lda, 0.0, c, a_nr
2812 F77_CHAR_ARG_LEN (1)
2813 F77_CHAR_ARG_LEN (1)));
2814 for (int j = 0; j < a_nr; j++)
2815 for (int i = 0; i < j; i++)
2816 retval.xelem (j, i) = retval.xelem (i, j);
2817
2818 }
2819 else
2820 {
2821 F77_INT lda = octave::to_f77_int (a.rows ());
2822 F77_INT tda = octave::to_f77_int (a.cols ());
2823 F77_INT ldb = octave::to_f77_int (b.rows ());
2824 F77_INT tdb = octave::to_f77_int (b.cols ());
2825
2826 retval = Matrix (a_nr, b_nc);
2827 double *c = retval.fortran_vec ();
2828
2829 if (b_nc == 1)
2830 {
2831 if (a_nr == 1)
2832 F77_FUNC (xddot, XDDOT) (a_nc, a.data (), 1, b.data (), 1, *c);
2833 else
2834 {
2835 const char ctra = get_blas_trans_arg (tra);
2836 F77_XFCN (dgemv, DGEMV, (F77_CONST_CHAR_ARG2 (&ctra, 1),
2837 lda, tda, 1.0, a.data (), lda,
2838 b.data (), 1, 0.0, c, 1
2839 F77_CHAR_ARG_LEN (1)));
2840 }
2841 }
2842 else if (a_nr == 1)
2843 {
2844 const char crevtrb = get_blas_trans_arg (! trb);
2845 F77_XFCN (dgemv, DGEMV, (F77_CONST_CHAR_ARG2 (&crevtrb, 1),
2846 ldb, tdb, 1.0, b.data (), ldb,
2847 a.data (), 1, 0.0, c, 1
2848 F77_CHAR_ARG_LEN (1)));
2849 }
2850 else
2851 {
2852 const char ctra = get_blas_trans_arg (tra);
2853 const char ctrb = get_blas_trans_arg (trb);
2854 F77_XFCN (dgemm, DGEMM, (F77_CONST_CHAR_ARG2 (&ctra, 1),
2855 F77_CONST_CHAR_ARG2 (&ctrb, 1),
2856 a_nr, b_nc, a_nc, 1.0, a.data (),
2857 lda, b.data (), ldb, 0.0, c, a_nr
2858 F77_CHAR_ARG_LEN (1)
2859 F77_CHAR_ARG_LEN (1)));
2860 }
2861 }
2862
2863 return retval;
2864}
2865
2866Matrix
2867operator * (const Matrix& a, const Matrix& b)
2868{
2869 return xgemm (a, b);
2870}
2871
2872// FIXME: it would be nice to share code among the min/max functions below.
2873
2874#define EMPTY_RETURN_CHECK(T) \
2875 if (nr == 0 || nc == 0) \
2876 return T (nr, nc);
2877
2878Matrix
2879min (double d, const Matrix& m)
2880{
2881 octave_idx_type nr = m.rows ();
2882 octave_idx_type nc = m.columns ();
2883
2885
2886 Matrix result (nr, nc);
2887
2888 for (octave_idx_type j = 0; j < nc; j++)
2889 for (octave_idx_type i = 0; i < nr; i++)
2890 {
2891 octave_quit ();
2892 result(i, j) = octave::math::min (d, m(i, j));
2893 }
2894
2895 return result;
2896}
2897
2898Matrix
2899min (const Matrix& m, double d)
2900{
2901 octave_idx_type nr = m.rows ();
2902 octave_idx_type nc = m.columns ();
2903
2905
2906 Matrix result (nr, nc);
2907
2908 for (octave_idx_type j = 0; j < nc; j++)
2909 for (octave_idx_type i = 0; i < nr; i++)
2910 {
2911 octave_quit ();
2912 result(i, j) = octave::math::min (m(i, j), d);
2913 }
2914
2915 return result;
2916}
2917
2918Matrix
2919min (const Matrix& a, const Matrix& b)
2920{
2921 octave_idx_type nr = a.rows ();
2922 octave_idx_type nc = a.columns ();
2923
2924 if (nr != b.rows () || nc != b.columns ())
2925 (*current_liboctave_error_handler)
2926 ("two-arg min requires same size arguments");
2927
2929
2930 Matrix result (nr, nc);
2931
2932 for (octave_idx_type j = 0; j < nc; j++)
2933 for (octave_idx_type i = 0; i < nr; i++)
2934 {
2935 octave_quit ();
2936 result(i, j) = octave::math::min (a(i, j), b(i, j));
2937 }
2938
2939 return result;
2940}
2941
2942Matrix
2943max (double d, const Matrix& m)
2944{
2945 octave_idx_type nr = m.rows ();
2946 octave_idx_type nc = m.columns ();
2947
2949
2950 Matrix result (nr, nc);
2951
2952 for (octave_idx_type j = 0; j < nc; j++)
2953 for (octave_idx_type i = 0; i < nr; i++)
2954 {
2955 octave_quit ();
2956 result(i, j) = octave::math::max (d, m(i, j));
2957 }
2958
2959 return result;
2960}
2961
2962Matrix
2963max (const Matrix& m, double d)
2964{
2965 octave_idx_type nr = m.rows ();
2966 octave_idx_type nc = m.columns ();
2967
2969
2970 Matrix result (nr, nc);
2971
2972 for (octave_idx_type j = 0; j < nc; j++)
2973 for (octave_idx_type i = 0; i < nr; i++)
2974 {
2975 octave_quit ();
2976 result(i, j) = octave::math::max (m(i, j), d);
2977 }
2978
2979 return result;
2980}
2981
2982Matrix
2983max (const Matrix& a, const Matrix& b)
2984{
2985 octave_idx_type nr = a.rows ();
2986 octave_idx_type nc = a.columns ();
2987
2988 if (nr != b.rows () || nc != b.columns ())
2989 (*current_liboctave_error_handler)
2990 ("two-arg max requires same size arguments");
2991
2993
2994 Matrix result (nr, nc);
2995
2996 for (octave_idx_type j = 0; j < nc; j++)
2997 for (octave_idx_type i = 0; i < nr; i++)
2998 {
2999 octave_quit ();
3000 result(i, j) = octave::math::max (a(i, j), b(i, j));
3001 }
3002
3003 return result;
3004}
3005
3007 const ColumnVector& x2,
3009
3010{
3011 octave_idx_type m = x1.numel ();
3012
3013 if (x2.numel () != m)
3014 (*current_liboctave_error_handler)
3015 ("linspace: vectors must be of equal length");
3016
3017 Matrix retval;
3018
3019 if (n < 1)
3020 {
3021 retval.clear (m, 0);
3022 return retval;
3023 }
3024
3025 retval.clear (m, n);
3026 for (octave_idx_type i = 0; i < m; i++)
3027 retval.xelem (i, 0) = x1(i);
3028
3029 // The last column is unused so temporarily store delta there
3030 double *delta = &retval.xelem (0, n-1);
3031 for (octave_idx_type i = 0; i < m; i++)
3032 delta[i] = (x1(i) == x2(i)) ? 0 : (x2(i) - x1(i)) / (n - 1);
3033
3034 for (octave_idx_type j = 1; j < n-1; j++)
3035 for (octave_idx_type i = 0; i < m; i++)
3036 retval.xelem (i, j) = x1(i) + j*delta[i];
3037
3038 for (octave_idx_type i = 0; i < m; i++)
3039 retval.xelem (i, n-1) = x2(i);
3040
3041 return retval;
3042}
3043
3046
3049
base_det< double > DET
Definition: DET.h:93
#define Inf
Definition: Faddeeva.cc:260
#define NaN
Definition: Faddeeva.cc:261
bool issquare(void) const
Definition: Array.h:605
double & xelem(octave_idx_type n)
Definition: Array.h:504
void make_unique(void)
Definition: Array.h:215
OCTARRAY_API void clear(void)
Definition: Array.cc:87
OCTARRAY_API Array< double, Alloc > index(const octave::idx_vector &i) const
Indexing without resizing.
Definition: Array.cc:697
octave_idx_type numel(void) const
Number of elements in the array.
Definition: Array.h:411
octave_idx_type cols(void) const
Definition: Array.h:457
double & elem(octave_idx_type n)
Definition: Array.h:534
octave_idx_type rows(void) const
Definition: Array.h:449
OCTARRAY_API void resize(const dim_vector &dv, const T &rfv)
Size of the specified dimension.
Definition: Array.cc:1010
OCTARRAY_API Array< T, Alloc > & insert(const Array< T, Alloc > &a, const Array< octave_idx_type > &idx)
Insert an array into another at a specified position.
Definition: Array.cc:1588
octave_idx_type columns(void) const
Definition: Array.h:458
const double * data(void) const
Definition: Array.h:616
OCTARRAY_API T * fortran_vec(void)
Size of the specified dimension.
Definition: Array.cc:1744
void resize(octave_idx_type n, const double &rfv=0)
Definition: dColVector.h:112
OCTAVE_API ColumnVector extract(octave_idx_type r1, octave_idx_type r2) const
Definition: dColVector.cc:151
OCTAVE_API ComplexMatrix lssolve(const Matrix &b) const
Definition: CMatrix.cc:2215
OCTAVE_API ComplexMatrix solve(MatrixType &mattype, const Matrix &b) const
Definition: CMatrix.cc:1926
T elem(octave_idx_type r, octave_idx_type c) const
Definition: DiagArray2.h:117
octave_idx_type length(void) const
Definition: DiagArray2.h:95
octave_idx_type cols(void) const
Definition: DiagArray2.h:90
octave_idx_type rows(void) const
Definition: DiagArray2.h:89
ColumnVector extract_diag(octave_idx_type k=0) const
Definition: dDiagMatrix.h:107
bool ishermitian(void) const
Definition: MatrixType.h:122
void mark_as_rectangular(void)
Definition: MatrixType.h:155
@ Permuted_Lower
Definition: MatrixType.h:48
@ Permuted_Upper
Definition: MatrixType.h:47
OCTAVE_API void mark_as_unsymmetric(void)
Definition: MatrixType.cc:920
OCTAVE_API int type(bool quiet=true)
Definition: MatrixType.cc:656
Definition: dMatrix.h:42
OCTAVE_API Matrix pseudo_inverse(double tol=0.0) const
Definition: dMatrix.cc:688
OCTAVE_API Matrix & fill(double val)
Definition: dMatrix.cc:221
OCTAVE_API Matrix diag(octave_idx_type k=0) const
Definition: dMatrix.cc:2402
Matrix ltsolve(MatrixType &mattype, const Matrix &b, octave_idx_type &info, double &rcon, solve_singularity_handler sing_handler, bool calc_cond=false, blas_trans_type transt=blas_no_trans) const
Definition: dMatrix.cc:1291
OCTAVE_API ComplexMatrix fourier(void) const
Definition: dMatrix.cc:729
OCTAVE_API ComplexMatrix fourier2d(void) const
Definition: dMatrix.cc:788
OCTAVE_API RowVector column_min(void) const
Definition: dMatrix.cc:2534
OCTAVE_API boolMatrix all(int dim=-1) const
Definition: dMatrix.cc:2354
OCTAVE_API ColumnVector row_max(void) const
Definition: dMatrix.cc:2479
OCTAVE_API Matrix sum(int dim=-1) const
Definition: dMatrix.cc:2384
OCTAVE_API double rcond(void) const
Definition: dMatrix.cc:1024
OCTAVE_API Matrix sumsq(int dim=-1) const
Definition: dMatrix.cc:2390
OCTAVE_API Matrix & insert(const Matrix &a, octave_idx_type r, octave_idx_type c)
Definition: dMatrix.cc:152
OCTAVE_API Matrix prod(int dim=-1) const
Definition: dMatrix.cc:2378
OCTAVE_API RowVector row(octave_idx_type i) const
Definition: dMatrix.cc:416
OCTAVE_API Matrix stack(const Matrix &a) const
Definition: dMatrix.cc:325
OCTAVE_API Matrix append(const Matrix &a) const
Definition: dMatrix.cc:265
Matrix fsolve(MatrixType &mattype, const Matrix &b, octave_idx_type &info, double &rcon, solve_singularity_handler sing_handler, bool calc_cond=false) const
Definition: dMatrix.cc:1388
OCTAVE_API Matrix extract(octave_idx_type r1, octave_idx_type c1, octave_idx_type r2, octave_idx_type c2) const
Definition: dMatrix.cc:397
OCTAVE_API RowVector column_max(void) const
Definition: dMatrix.cc:2589
Matrix utsolve(MatrixType &mattype, const Matrix &b, octave_idx_type &info, double &rcon, solve_singularity_handler sing_handler, bool calc_cond=false, blas_trans_type transt=blas_no_trans) const
Definition: dMatrix.cc:1193
OCTAVE_API Matrix cumprod(int dim=-1) const
Definition: dMatrix.cc:2366
OCTAVE_API Matrix inverse(void) const
Definition: dMatrix.cc:451
Matrix tinverse(MatrixType &mattype, octave_idx_type &info, double &rcon, bool force, bool calc_cond) const
Definition: dMatrix.cc:491
OCTAVE_API Matrix extract_n(octave_idx_type r1, octave_idx_type c1, octave_idx_type nr, octave_idx_type nc) const
Definition: dMatrix.cc:407
OCTAVE_API Matrix abs(void) const
Definition: dMatrix.cc:2396
OCTAVE_API Matrix & operator+=(const DiagMatrix &a)
Definition: dMatrix.cc:2285
OCTAVE_API bool issymmetric(void) const
Definition: dMatrix.cc:136
OCTAVE_API boolMatrix any(int dim=-1) const
Definition: dMatrix.cc:2360
OCTAVE_API ColumnVector row_min(void) const
Definition: dMatrix.cc:2424
OCTAVE_API Matrix & operator-=(const DiagMatrix &a)
Definition: dMatrix.cc:2303
OCTAVE_API ComplexMatrix ifourier2d(void) const
Definition: dMatrix.cc:800
friend class ComplexMatrix
Definition: dMatrix.h:137
Matrix(void)=default
Matrix transpose(void) const
Definition: dMatrix.h:140
OCTAVE_API ComplexMatrix ifourier(void) const
Definition: dMatrix.cc:758
OCTAVE_API DET determinant(void) const
Definition: dMatrix.cc:853
Matrix finverse(MatrixType &mattype, octave_idx_type &info, double &rcon, bool force, bool calc_cond) const
Definition: dMatrix.cc:550
OCTAVE_API bool operator==(const Matrix &a) const
Definition: dMatrix.cc:121
OCTAVE_API Matrix solve(MatrixType &mattype, const Matrix &b) const
Definition: dMatrix.cc:1586
OCTAVE_API Matrix cumsum(int dim=-1) const
Definition: dMatrix.cc:2372
OCTAVE_API bool operator!=(const Matrix &a) const
Definition: dMatrix.cc:130
OCTAVE_API Matrix lssolve(const Matrix &b) const
Definition: dMatrix.cc:1897
void resize(octave_idx_type nr, octave_idx_type nc, double rfv=0)
Definition: dMatrix.h:158
OCTAVE_API ColumnVector column(octave_idx_type i) const
Definition: dMatrix.cc:422
OCTAVE_API NDArray diag(octave_idx_type k=0) const
Definition: dNDArray.cc:609
OCTAVE_API NDArray abs(void) const
Definition: dNDArray.cc:570
OCTAVE_API NDArray cumsum(int dim=-1) const
Definition: dNDArray.cc:413
OCTAVE_API NDArray sumsq(int dim=-1) const
Definition: dNDArray.cc:437
OCTAVE_API NDArray prod(int dim=-1) const
Definition: dNDArray.cc:419
OCTAVE_API NDArray cumprod(int dim=-1) const
Definition: dNDArray.cc:407
OCTAVE_API boolNDArray all(int dim=-1) const
Definition: dNDArray.cc:395
OCTAVE_API boolNDArray any(int dim=-1) const
Definition: dNDArray.cc:401
OCTAVE_API NDArray sum(int dim=-1) const
Definition: dNDArray.cc:425
octave_idx_type rows(void) const
Definition: PermMatrix.h:62
const Array< octave_idx_type > & col_perm_vec(void) const
Definition: PermMatrix.h:83
void resize(octave_idx_type n, const double &rfv=0)
Definition: dRowVector.h:100
Definition: DET.h:39
base_det square() const
Definition: DET.h:75
Definition: mx-defs.h:39
Vector representing the dimensions (size) of an Array.
Definition: dim-vector.h:94
static int ifft(const Complex *in, Complex *out, std::size_t npts, std::size_t nsamples=1, octave_idx_type stride=1, octave_idx_type dist=-1)
Definition: oct-fftw.cc:900
static int fftNd(const double *, Complex *, const int, const dim_vector &)
Definition: oct-fftw.cc:924
static int ifftNd(const Complex *, Complex *, const int, const dim_vector &)
Definition: oct-fftw.cc:970
static int fft(const double *in, Complex *out, std::size_t npts, std::size_t nsamples=1, octave_idx_type stride=1, octave_idx_type dist=-1)
Definition: oct-fftw.cc:858
static const idx_vector colon
Definition: idx-vector.h:483
T unitary_schur_matrix(void) const
Definition: schur.h:91
T schur_matrix(void) const
Definition: schur.h:89
DM_T singular_values(void) const
Definition: svd.h:90
T right_singular_matrix(void) const
Definition: svd.cc:321
T left_singular_matrix(void) const
Definition: svd.cc:310
Matrix imag(const ComplexMatrix &a)
Definition: dMatrix.cc:391
Matrix xgemm(const Matrix &a, const Matrix &b, blas_trans_type transa, blas_trans_type transb)
Definition: dMatrix.cc:2781
static ComplexMatrix unstack_complex_matrix(const Matrix &sm)
Definition: dMatrix.cc:1683
static Matrix stack_complex_matrix(const ComplexMatrix &cm)
Definition: dMatrix.cc:1666
Matrix Givens(double x, double y)
Definition: dMatrix.cc:2682
static double norm1(const Matrix &a)
Definition: dMatrix.cc:430
Matrix max(double d, const Matrix &m)
Definition: dMatrix.cc:2943
std::istream & operator>>(std::istream &is, Matrix &a)
Definition: dMatrix.cc:2659
Matrix Sylvester(const Matrix &a, const Matrix &b, const Matrix &c)
Definition: dMatrix.cc:2699
Matrix operator*(const ColumnVector &v, const RowVector &a)
Definition: dMatrix.cc:2325
#define EMPTY_RETURN_CHECK(T)
Definition: dMatrix.cc:2874
Matrix linspace(const ColumnVector &x1, const ColumnVector &x2, octave_idx_type n)
Definition: dMatrix.cc:3006
Matrix min(double d, const Matrix &m)
Definition: dMatrix.cc:2879
static char get_blas_trans_arg(bool trans)
Definition: dMatrix.cc:2773
std::ostream & operator<<(std::ostream &os, const Matrix &a)
Definition: dMatrix.cc:2644
Matrix real(const ComplexMatrix &a)
Definition: dMatrix.cc:385
#define F77_XFCN(f, F, args)
Definition: f77-fcn.h:45
octave_f77_int_type F77_INT
Definition: f77-fcn.h:306
double norm(const ColumnVector &v)
Definition: graphics.cc:5940
void scale(Matrix &m, double x, double y, double z)
Definition: graphics.cc:5893
OCTAVE_NORETURN liboctave_error_handler current_liboctave_error_handler
Definition: lo-error.c:41
F77_RET_T const F77_INT const F77_INT const F77_INT const F77_DBLE const F77_DBLE F77_INT F77_DBLE * V
F77_RET_T const F77_DBLE const F77_DBLE F77_DBLE * d
F77_RET_T const F77_DBLE * x
blas_trans_type
Definition: mx-defs.h:80
@ blas_no_trans
Definition: mx-defs.h:81
@ blas_conj_trans
Definition: mx-defs.h:83
@ blas_trans
Definition: mx-defs.h:82
char get_blas_char(blas_trans_type transt)
Definition: mx-defs.h:87
class OCTAVE_API Matrix
Definition: mx-fwd.h:31
class OCTAVE_API DiagMatrix
Definition: mx-fwd.h:59
class OCTAVE_API ColumnVector
Definition: mx-fwd.h:45
bool mx_inline_equal(std::size_t n, const T1 *x, const T2 *y)
Definition: mx-inlines.cc:570
void mx_inline_real(std::size_t n, T *r, const std::complex< T > *x)
Definition: mx-inlines.cc:321
void mx_inline_imag(std::size_t n, T *r, const std::complex< T > *x)
Definition: mx-inlines.cc:328
#define MM_BOOL_OPS(M1, M2)
Definition: mx-op-defs.h:213
#define MM_CMP_OPS(M1, M2)
Definition: mx-op-defs.h:196
#define MS_CMP_OPS(M, S)
Definition: mx-op-defs.h:110
#define SM_CMP_OPS(S, M)
Definition: mx-op-defs.h:153
#define SM_BOOL_OPS(S, M)
Definition: mx-op-defs.h:170
#define MS_BOOL_OPS(M, S)
Definition: mx-op-defs.h:127
T max(T x, T y)
Definition: lo-mappers.h:368
bool isfinite(double x)
Definition: lo-mappers.h:192
bool isnan(bool)
Definition: lo-mappers.h:178
bool isinf(double x)
Definition: lo-mappers.h:203
Complex log2(const Complex &x)
Definition: lo-mappers.cc:139
T min(T x, T y)
Definition: lo-mappers.h:361
void warn_singular_matrix(double rcond)
void err_nonconformant(const char *op, octave_idx_type op1_len, octave_idx_type op2_len)
std::complex< double > Complex
Definition: oct-cmplx.h:33
#define OCTAVE_LOCAL_BUFFER(T, buf, size)
Definition: oct-locbuf.h:44
static bool scalar(const dim_vector &dims)
Definition: ov-struct.cc:679
subroutine xddot(n, dx, incx, dy, incy, retval)
Definition: xddot.f:2
F77_RET_T F77_FUNC(xerbla, XERBLA)(F77_CONST_CHAR_ARG_DEF(s_arg
F77_RET_T len
Definition: xerbla.cc:61
subroutine xilaenv(ispec, name, opts, n1, n2, n3, n4, retval)
Definition: xilaenv.f:2