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