Mercurial > hg > octave-nkf
annotate scripts/linear-algebra/linsolve.m @ 20815:a260a6acb70f
fix test failures introduced by a22d8a2eb0e5
* scripts/ode/private/integrate_adaptive.m: fix stepping backwards, fix
invocation of OutputFcn, fix text of some error messages
* scripts/ode/private/integrate_const.m: remove use of option OutputSave
* scripts/ode/private/integrate_n_steps.m: remove use of option OutputSave
author | Carlo de Falco <carlo.defalco@polimi.it> |
---|---|
date | Sun, 11 Oct 2015 23:09:01 +0200 |
parents | 014e942ac29f |
children |
rev | line source |
---|---|
19898
4197fc428c7d
maint: Update copyright notices for 2015.
John W. Eaton <jwe@octave.org>
parents:
19790
diff
changeset
|
1 ## Copyright (C) 2013-2015 Nir Krakauer |
17499
c3a3532e3d98
linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff
changeset
|
2 ## |
17795
0a8c35ae5ce1
maint: Fix various problems with GPL copyright statements.
Rik <rik@octave.org>
parents:
17500
diff
changeset
|
3 ## This file is part of Octave. |
0a8c35ae5ce1
maint: Fix various problems with GPL copyright statements.
Rik <rik@octave.org>
parents:
17500
diff
changeset
|
4 ## |
0a8c35ae5ce1
maint: Fix various problems with GPL copyright statements.
Rik <rik@octave.org>
parents:
17500
diff
changeset
|
5 ## Octave is free software; you can redistribute it and/or modify |
17499
c3a3532e3d98
linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff
changeset
|
6 ## it under the terms of the GNU General Public License as published by |
c3a3532e3d98
linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff
changeset
|
7 ## the Free Software Foundation; either version 3 of the License, or |
c3a3532e3d98
linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff
changeset
|
8 ## (at your option) any later version. |
c3a3532e3d98
linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff
changeset
|
9 ## |
17795
0a8c35ae5ce1
maint: Fix various problems with GPL copyright statements.
Rik <rik@octave.org>
parents:
17500
diff
changeset
|
10 ## Octave is distributed in the hope that it will be useful, |
17499
c3a3532e3d98
linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff
changeset
|
11 ## but WITHOUT ANY WARRANTY; without even the implied warranty of |
c3a3532e3d98
linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff
changeset
|
12 ## 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
|
13 ## 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
|
14 ## |
c3a3532e3d98
linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff
changeset
|
15 ## You should have received a copy of the GNU General Public License |
17795
0a8c35ae5ce1
maint: Fix various problems with GPL copyright statements.
Rik <rik@octave.org>
parents:
17500
diff
changeset
|
16 ## along with Octave; If not, see <http://www.gnu.org/licenses/>. |
17499
c3a3532e3d98
linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff
changeset
|
17 |
c3a3532e3d98
linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff
changeset
|
18 ## -*- texinfo -*- |
17500
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
19 ## @deftypefn {Function File} {@var{x} =} linsolve (@var{A}, @var{b}) |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
20 ## @deftypefnx {Function File} {@var{x} =} linsolve (@var{A}, @var{b}, @var{opts}) |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
21 ## @deftypefnx {Function File} {[@var{x}, @var{R}] =} linsolve (@dots{}) |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
22 ## 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
|
23 ## |
17500
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
24 ## 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
|
25 ## @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
|
26 ## @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
|
27 ## |
17500
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
28 ## 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
|
29 ## 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
|
30 ## 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
|
31 ## 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
|
32 ## to solving the linear system. |
17499
c3a3532e3d98
linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff
changeset
|
33 ## |
17500
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
34 ## @strong{Warning:} If the matrix @var{A} does not have the properties |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
35 ## listed in the @var{opts} structure then the result will not be accurate |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
36 ## AND no warning will be given. When in doubt, let Octave examine the matrix |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
37 ## and choose the appropriate solver as this step takes little time and the |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
38 ## result is cached so that it is only done once per linear system. |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
39 ## |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
40 ## Possible @var{opts} fields (set value to true/false): |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
41 ## |
17499
c3a3532e3d98
linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff
changeset
|
42 ## @table @asis |
17500
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
43 ## @item LT |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
44 ## @var{A} is lower triangular |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
45 ## |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
46 ## @item UT |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
47 ## @var{A} is upper triangular |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
48 ## |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
49 ## @item UHESS |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
50 ## @var{A} is upper Hessenberg (currently makes no difference) |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
51 ## |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
52 ## @item SYM |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
53 ## @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
|
54 ## |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
55 ## @item POSDEF |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
56 ## @var{A} is positive definite |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
57 ## |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
58 ## @item RECT |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
59 ## @var{A} is general rectangular (currently makes no difference) |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
60 ## |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
61 ## @item TRANSA |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
62 ## Solve @code{A'*x = b} by @code{transpose (A) \ b} |
17499
c3a3532e3d98
linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff
changeset
|
63 ## @end table |
c3a3532e3d98
linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff
changeset
|
64 ## |
17500
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
65 ## 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
|
66 ## @var{A} (zero if matrix is singular). |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
67 ## @seealso{mldivide, matrix_type, rcond} |
17499
c3a3532e3d98
linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff
changeset
|
68 ## @end deftypefn |
c3a3532e3d98
linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff
changeset
|
69 |
c3a3532e3d98
linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff
changeset
|
70 ## Author: Nir Krakauer <nkrakauer@ccny.cuny.edu> |
c3a3532e3d98
linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff
changeset
|
71 |
17500
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
72 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
|
73 |
17500
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
74 if (nargin < 2 || nargin > 3) |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
75 print_usage (); |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
76 endif |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
77 |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
78 if (! (isnumeric (A) && isnumeric (b))) |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
79 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
|
80 endif |
17500
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
81 |
20448
014e942ac29f
linsolve.m: Fix regression when calling linsolve with 2 arguments (bug #45212)
Mike Miller <mtmiller@octave.org>
parents:
20240
diff
changeset
|
82 trans_A = false; |
014e942ac29f
linsolve.m: Fix regression when calling linsolve with 2 arguments (bug #45212)
Mike Miller <mtmiller@octave.org>
parents:
20240
diff
changeset
|
83 |
17500
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
84 ## Process any opts |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
85 if (nargin > 2) |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
86 if (! isstruct (opts)) |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
87 error ("linsolve: OPTS must be a structure"); |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
88 endif |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
89 if (isfield (opts, "TRANSA") && opts.TRANSA) |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
90 trans_A = true; |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
91 endif |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
92 if (isfield (opts, "POSDEF") && opts.POSDEF) |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
93 A = matrix_type (A, "positive definite"); |
19790
446c46af4b42
strip trailing whitespace from most source files
John W. Eaton <jwe@octave.org>
parents:
17795
diff
changeset
|
94 endif |
17500
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
95 if (isfield (opts, "LT") && opts.LT) |
20240
91e1da1d1918
linsolve.m: Calculate solution in a BLAS-optimized way (bug #44722).
Rik <rik@octave.org>
parents:
19898
diff
changeset
|
96 A = matrix_type (A, "lower"); |
91e1da1d1918
linsolve.m: Calculate solution in a BLAS-optimized way (bug #44722).
Rik <rik@octave.org>
parents:
19898
diff
changeset
|
97 elseif (isfield (opts, "UT") && opts.UT) |
91e1da1d1918
linsolve.m: Calculate solution in a BLAS-optimized way (bug #44722).
Rik <rik@octave.org>
parents:
19898
diff
changeset
|
98 A = matrix_type (A, "upper"); |
19790
446c46af4b42
strip trailing whitespace from most source files
John W. Eaton <jwe@octave.org>
parents:
17795
diff
changeset
|
99 endif |
17499
c3a3532e3d98
linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff
changeset
|
100 endif |
17500
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
101 |
20240
91e1da1d1918
linsolve.m: Calculate solution in a BLAS-optimized way (bug #44722).
Rik <rik@octave.org>
parents:
19898
diff
changeset
|
102 ## 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:
19898
diff
changeset
|
103 ## 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:
19898
diff
changeset
|
104 if (trans_A) |
91e1da1d1918
linsolve.m: Calculate solution in a BLAS-optimized way (bug #44722).
Rik <rik@octave.org>
parents:
19898
diff
changeset
|
105 x = A' \ b; |
91e1da1d1918
linsolve.m: Calculate solution in a BLAS-optimized way (bug #44722).
Rik <rik@octave.org>
parents:
19898
diff
changeset
|
106 else |
91e1da1d1918
linsolve.m: Calculate solution in a BLAS-optimized way (bug #44722).
Rik <rik@octave.org>
parents:
19898
diff
changeset
|
107 x = A \ b; |
91e1da1d1918
linsolve.m: Calculate solution in a BLAS-optimized way (bug #44722).
Rik <rik@octave.org>
parents:
19898
diff
changeset
|
108 endif |
17500
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
109 |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
110 if (nargout > 1) |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
111 if (issquare (A)) |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
112 R = rcond (A); |
17499
c3a3532e3d98
linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff
changeset
|
113 else |
17500
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
114 R = 0; |
17499
c3a3532e3d98
linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff
changeset
|
115 endif |
c3a3532e3d98
linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff
changeset
|
116 endif |
20240
91e1da1d1918
linsolve.m: Calculate solution in a BLAS-optimized way (bug #44722).
Rik <rik@octave.org>
parents:
19898
diff
changeset
|
117 |
17500
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
118 endfunction |
17499
c3a3532e3d98
linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff
changeset
|
119 |
c3a3532e3d98
linsolve.m: Add new function for Matlab compatibility.
Nir Krakauer < nkrakauer@ccny.cuny.edu>
parents:
diff
changeset
|
120 |
17500
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
121 %!test |
20240
91e1da1d1918
linsolve.m: Calculate solution in a BLAS-optimized way (bug #44722).
Rik <rik@octave.org>
parents:
19898
diff
changeset
|
122 %! n = 10; |
20448
014e942ac29f
linsolve.m: Fix regression when calling linsolve with 2 arguments (bug #45212)
Mike Miller <mtmiller@octave.org>
parents:
20240
diff
changeset
|
123 %! A = rand (n); |
014e942ac29f
linsolve.m: Fix regression when calling linsolve with 2 arguments (bug #45212)
Mike Miller <mtmiller@octave.org>
parents:
20240
diff
changeset
|
124 %! x = rand (n, 1); |
014e942ac29f
linsolve.m: Fix regression when calling linsolve with 2 arguments (bug #45212)
Mike Miller <mtmiller@octave.org>
parents:
20240
diff
changeset
|
125 %! b = A * x; |
014e942ac29f
linsolve.m: Fix regression when calling linsolve with 2 arguments (bug #45212)
Mike Miller <mtmiller@octave.org>
parents:
20240
diff
changeset
|
126 %! 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:
20240
diff
changeset
|
127 %! 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:
20240
diff
changeset
|
128 |
014e942ac29f
linsolve.m: Fix regression when calling linsolve with 2 arguments (bug #45212)
Mike Miller <mtmiller@octave.org>
parents:
20240
diff
changeset
|
129 %!test |
014e942ac29f
linsolve.m: Fix regression when calling linsolve with 2 arguments (bug #45212)
Mike Miller <mtmiller@octave.org>
parents:
20240
diff
changeset
|
130 %! n = 10; |
20240
91e1da1d1918
linsolve.m: Calculate solution in a BLAS-optimized way (bug #44722).
Rik <rik@octave.org>
parents:
19898
diff
changeset
|
131 %! A = triu (gallery ("condex", n)); |
17500
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
132 %! x = rand (n, 1); |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
133 %! b = A' * x; |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
134 %! opts.UT = true; |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
135 %! opts.TRANSA = true; |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
136 %! 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
|
137 |
17500
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
138 %!error linsolve () |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
139 %!error linsolve (1) |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
140 %!error linsolve (1,2,3) |
b66f068e4468
linsolve.m: Use Octave coding conventions.
Rik <rik@octave.org>
parents:
17499
diff
changeset
|
141 %!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
|
142 %!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
|
143 %!error <OPTS must be a structure> linsolve (1,2,3) |