annotate scripts/linear-algebra/linsolve.m @ 33632:fed0dc6fd44c default tip @

remove unused variable from libgui/module.mk * libgui/module.mk: remove empty variable OCTAVE_GUI_EDITOR_MOC
author Torsten Lilge <ttl-octave@mailbox.org>
date Mon, 27 May 2024 19:42:05 +0200
parents 2e484f9f1f18
children
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
27923
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
1 ########################################################################
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
2 ##
32632
2e484f9f1f18 maint: update Octave Project Developers copyright for the new year
John W. Eaton <jwe@octave.org>
parents: 31706
diff changeset
3 ## Copyright (C) 2013-2024 The Octave Project Developers
27918
b442ec6dda5c use centralized file for copyright info for individual contributors
John W. Eaton <jwe@octave.org>
parents: 26376
diff changeset
4 ##
27923
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
5 ## See the file COPYRIGHT.md in the top-level directory of this
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
6 ## distribution or <https://octave.org/copyright/>.
17499
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
7 ##
17795
0a8c35ae5ce1 maint: Fix various problems with GPL copyright statements.
Rik <rik@octave.org>
parents: 17500
diff changeset
8 ## This file is part of Octave.
0a8c35ae5ce1 maint: Fix various problems with GPL copyright statements.
Rik <rik@octave.org>
parents: 17500
diff changeset
9 ##
24534
194eb4bd202b maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents: 24489
diff changeset
10 ## Octave is free software: you can redistribute it and/or modify it
22755
3a2b891d0b33 maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents: 22323
diff changeset
11 ## under the terms of the GNU General Public License as published by
24534
194eb4bd202b maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents: 24489
diff changeset
12 ## the Free Software Foundation, either version 3 of the License, or
17499
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
13 ## (at your option) any later version.
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
14 ##
22755
3a2b891d0b33 maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents: 22323
diff changeset
15 ## Octave is distributed in the hope that it will be useful, but
3a2b891d0b33 maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents: 22323
diff changeset
16 ## WITHOUT ANY WARRANTY; without even the implied warranty of
17499
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
17 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
18 ## GNU General Public License for more details.
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
19 ##
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
20 ## You should have received a copy of the GNU General Public License
22755
3a2b891d0b33 maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents: 22323
diff changeset
21 ## along with Octave; see the file COPYING. If not, see
24534
194eb4bd202b maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents: 24489
diff changeset
22 ## <https://www.gnu.org/licenses/>.
27923
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
23 ##
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
24 ########################################################################
17499
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
25
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
26 ## -*- texinfo -*-
20852
516bb87ea72e 2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents: 20238
diff changeset
27 ## @deftypefn {} {@var{x} =} linsolve (@var{A}, @var{b})
516bb87ea72e 2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents: 20238
diff changeset
28 ## @deftypefnx {} {@var{x} =} linsolve (@var{A}, @var{b}, @var{opts})
516bb87ea72e 2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents: 20238
diff changeset
29 ## @deftypefnx {} {[@var{x}, @var{R}] =} linsolve (@dots{})
17500
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
30 ## Solve the linear system @code{A*x = b}.
17499
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
31 ##
17500
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
32 ## With no options, this function is equivalent to the left division operator
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
33 ## @w{(@code{x = A \ b})} or the matrix-left-divide function
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
34 ## @w{(@code{x = mldivide (A, b)})}.
17499
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
35 ##
17500
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
36 ## Octave ordinarily examines the properties of the matrix @var{A} and chooses
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
37 ## a solver that best matches the matrix. By passing a structure @var{opts}
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
38 ## to @code{linsolve} you can inform Octave directly about the matrix @var{A}.
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
39 ## In this case Octave will skip the matrix examination and proceed directly
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
40 ## to solving the linear system.
17499
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
41 ##
24489
6ea279bbe94a linsolve.m: Clarify documentation about TRANSA input.
Rik <rik@octave.org>
parents: 23220
diff changeset
42 ## @strong{Warning:} If the matrix @var{A} does not have the properties listed
6ea279bbe94a linsolve.m: Clarify documentation about TRANSA input.
Rik <rik@octave.org>
parents: 23220
diff changeset
43 ## in the @var{opts} structure then the result will not be accurate AND no
6ea279bbe94a linsolve.m: Clarify documentation about TRANSA input.
Rik <rik@octave.org>
parents: 23220
diff changeset
44 ## warning will be given. When in doubt, let Octave examine the matrix and
6ea279bbe94a linsolve.m: Clarify documentation about TRANSA input.
Rik <rik@octave.org>
parents: 23220
diff changeset
45 ## choose the appropriate solver as this step takes little time and the result
6ea279bbe94a linsolve.m: Clarify documentation about TRANSA input.
Rik <rik@octave.org>
parents: 23220
diff changeset
46 ## is cached so that it is only done once per linear system.
17500
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
47 ##
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
48 ## Possible @var{opts} fields (set value to true/false):
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
49 ##
17499
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
50 ## @table @asis
17500
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
51 ## @item LT
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
52 ## @var{A} is lower triangular
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
53 ##
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
54 ## @item UT
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
55 ## @var{A} is upper triangular
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
56 ##
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
57 ## @item UHESS
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
58 ## @var{A} is upper Hessenberg (currently makes no difference)
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
59 ##
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
60 ## @item SYM
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
61 ## @var{A} is symmetric or complex Hermitian (currently makes no difference)
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
62 ##
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
63 ## @item POSDEF
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
64 ## @var{A} is positive definite
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
65 ##
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
66 ## @item RECT
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
67 ## @var{A} is general rectangular (currently makes no difference)
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
68 ##
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
69 ## @item TRANSA
24489
6ea279bbe94a linsolve.m: Clarify documentation about TRANSA input.
Rik <rik@octave.org>
parents: 23220
diff changeset
70 ## Solve @code{A'*x = b} if true rather than @code{A*x = b}
17499
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
71 ## @end table
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
72 ##
17500
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
73 ## The optional second output @var{R} is the inverse condition number of
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
74 ## @var{A} (zero if matrix is singular).
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
75 ## @seealso{mldivide, matrix_type, rcond}
17499
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
76 ## @end deftypefn
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
77
17500
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
78 function [x, R] = linsolve (A, b, opts)
17499
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
79
28789
28de41192f3c Eliminate unneeded verification of nargin, nargout in m-files.
Rik <rik@octave.org>
parents: 27978
diff changeset
80 if (nargin < 2)
17500
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
81 print_usage ();
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
82 endif
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
83
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
84 if (! (isnumeric (A) && isnumeric (b)))
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
85 error ("linsolve: A and B must be numeric");
17499
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
86 endif
17500
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
87
20238
014e942ac29f linsolve.m: Fix regression when calling linsolve with 2 arguments (bug #45212)
Mike Miller <mtmiller@octave.org>
parents: 20030
diff changeset
88 trans_A = false;
014e942ac29f linsolve.m: Fix regression when calling linsolve with 2 arguments (bug #45212)
Mike Miller <mtmiller@octave.org>
parents: 20030
diff changeset
89
17500
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
90 ## Process any opts
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
91 if (nargin > 2)
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
92 if (! isstruct (opts))
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
93 error ("linsolve: OPTS must be a structure");
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
94 endif
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
95 if (isfield (opts, "TRANSA") && opts.TRANSA)
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
96 trans_A = true;
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
97 endif
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
98 if (isfield (opts, "POSDEF") && opts.POSDEF)
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
99 A = matrix_type (A, "positive definite");
19593
446c46af4b42 strip trailing whitespace from most source files
John W. Eaton <jwe@octave.org>
parents: 17795
diff changeset
100 endif
17500
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
101 if (isfield (opts, "LT") && opts.LT)
20030
91e1da1d1918 linsolve.m: Calculate solution in a BLAS-optimized way (bug #44722).
Rik <rik@octave.org>
parents: 19697
diff changeset
102 A = matrix_type (A, "lower");
91e1da1d1918 linsolve.m: Calculate solution in a BLAS-optimized way (bug #44722).
Rik <rik@octave.org>
parents: 19697
diff changeset
103 elseif (isfield (opts, "UT") && opts.UT)
91e1da1d1918 linsolve.m: Calculate solution in a BLAS-optimized way (bug #44722).
Rik <rik@octave.org>
parents: 19697
diff changeset
104 A = matrix_type (A, "upper");
19593
446c46af4b42 strip trailing whitespace from most source files
John W. Eaton <jwe@octave.org>
parents: 17795
diff changeset
105 endif
17499
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
106 endif
17500
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
107
20030
91e1da1d1918 linsolve.m: Calculate solution in a BLAS-optimized way (bug #44722).
Rik <rik@octave.org>
parents: 19697
diff changeset
108 ## This way is faster as the transpose is not calculated in Octave,
91e1da1d1918 linsolve.m: Calculate solution in a BLAS-optimized way (bug #44722).
Rik <rik@octave.org>
parents: 19697
diff changeset
109 ## but forwarded as a flag option to BLAS.
91e1da1d1918 linsolve.m: Calculate solution in a BLAS-optimized way (bug #44722).
Rik <rik@octave.org>
parents: 19697
diff changeset
110 if (trans_A)
91e1da1d1918 linsolve.m: Calculate solution in a BLAS-optimized way (bug #44722).
Rik <rik@octave.org>
parents: 19697
diff changeset
111 x = A' \ b;
91e1da1d1918 linsolve.m: Calculate solution in a BLAS-optimized way (bug #44722).
Rik <rik@octave.org>
parents: 19697
diff changeset
112 else
91e1da1d1918 linsolve.m: Calculate solution in a BLAS-optimized way (bug #44722).
Rik <rik@octave.org>
parents: 19697
diff changeset
113 x = A \ b;
91e1da1d1918 linsolve.m: Calculate solution in a BLAS-optimized way (bug #44722).
Rik <rik@octave.org>
parents: 19697
diff changeset
114 endif
17500
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
115
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
116 if (nargout > 1)
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
117 if (issquare (A))
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
118 R = rcond (A);
17499
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
119 else
17500
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
120 R = 0;
17499
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
121 endif
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
122 endif
20030
91e1da1d1918 linsolve.m: Calculate solution in a BLAS-optimized way (bug #44722).
Rik <rik@octave.org>
parents: 19697
diff changeset
123
17500
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
124 endfunction
17499
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
125
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
126
17500
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
127 %!test
20030
91e1da1d1918 linsolve.m: Calculate solution in a BLAS-optimized way (bug #44722).
Rik <rik@octave.org>
parents: 19697
diff changeset
128 %! n = 10;
20238
014e942ac29f linsolve.m: Fix regression when calling linsolve with 2 arguments (bug #45212)
Mike Miller <mtmiller@octave.org>
parents: 20030
diff changeset
129 %! A = rand (n);
014e942ac29f linsolve.m: Fix regression when calling linsolve with 2 arguments (bug #45212)
Mike Miller <mtmiller@octave.org>
parents: 20030
diff changeset
130 %! x = rand (n, 1);
014e942ac29f linsolve.m: Fix regression when calling linsolve with 2 arguments (bug #45212)
Mike Miller <mtmiller@octave.org>
parents: 20030
diff changeset
131 %! b = A * x;
014e942ac29f linsolve.m: Fix regression when calling linsolve with 2 arguments (bug #45212)
Mike Miller <mtmiller@octave.org>
parents: 20030
diff changeset
132 %! assert (linsolve (A, b), A \ b);
014e942ac29f linsolve.m: Fix regression when calling linsolve with 2 arguments (bug #45212)
Mike Miller <mtmiller@octave.org>
parents: 20030
diff changeset
133 %! assert (linsolve (A, b, struct ()), A \ b);
014e942ac29f linsolve.m: Fix regression when calling linsolve with 2 arguments (bug #45212)
Mike Miller <mtmiller@octave.org>
parents: 20030
diff changeset
134
014e942ac29f linsolve.m: Fix regression when calling linsolve with 2 arguments (bug #45212)
Mike Miller <mtmiller@octave.org>
parents: 20030
diff changeset
135 %!test
014e942ac29f linsolve.m: Fix regression when calling linsolve with 2 arguments (bug #45212)
Mike Miller <mtmiller@octave.org>
parents: 20030
diff changeset
136 %! n = 10;
20030
91e1da1d1918 linsolve.m: Calculate solution in a BLAS-optimized way (bug #44722).
Rik <rik@octave.org>
parents: 19697
diff changeset
137 %! A = triu (gallery ("condex", n));
17500
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
138 %! x = rand (n, 1);
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
139 %! b = A' * x;
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
140 %! opts.UT = true;
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
141 %! opts.TRANSA = true;
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
142 %! assert (linsolve (A, b, opts), A' \ b);
17499
c3a3532e3d98 linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff changeset
143
28896
90fea9cc9caa test: Add expected error message <Invalid call> to BIST tests for nargin.
Rik <rik@octave.org>
parents: 28789
diff changeset
144 %!error <Invalid call> linsolve ()
90fea9cc9caa test: Add expected error message <Invalid call> to BIST tests for nargin.
Rik <rik@octave.org>
parents: 28789
diff changeset
145 %!error <Invalid call> linsolve (1)
17500
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
146 %!error linsolve (1,2,3)
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
147 %!error <A and B must be numeric> linsolve ({1},2)
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
148 %!error <A and B must be numeric> linsolve (1,{2})
b66f068e4468 linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents: 17499
diff changeset
149 %!error <OPTS must be a structure> linsolve (1,2,3)