Mercurial > hg > octave-nkf
annotate scripts/sparse/cgs.m @ 8596:8833c0b18eb2
enable default settings queries in optim funcs
author | Jaroslav Hajek <highegg@gmail.com> |
---|---|
date | Tue, 27 Jan 2009 08:15:08 +0100 |
parents | cadc73247d65 |
children | eb63fbe60fab |
rev | line source |
---|---|
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
1 ## Copyright (C) 2008 Radek Salac |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
2 ## |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
3 ## This file is part of Octave. |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
4 ## |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
5 ## Octave is free software; you can redistribute it and/or modify it |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
6 ## under the terms of the GNU General Public License as published by |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
7 ## the Free Software Foundation; either version 3 of the License, or (at |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
8 ## your option) any later version. |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
9 ## |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
10 ## Octave is distributed in the hope that it will be useful, but |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
11 ## WITHOUT ANY WARRANTY; without even the implied warranty of |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
12 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
13 ## General Public License for more details. |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
14 ## |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
15 ## You should have received a copy of the GNU General Public License |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
16 ## along with Octave; see the file COPYING. If not, see |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
17 ## <http://www.gnu.org/licenses/>. |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
18 |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
19 ## -*- texinfo -*- |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
20 ## @deftypefn {Function File} {} cgs (@var{A}, @var{b}) |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
21 ## @deftypefnx {Function File} {} cgs (@var{A}, @var{b}, @var{tol}, @var{maxit}, @var{M1}, @var{M2}, @var{x0}) |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
22 ## This procedure attempts to solve a system of linear equations A*x = b for x. |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
23 ## The @var{A} must be square, symmetric and positive definite real matrix N*N. |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
24 ## The @var{b} must be a one column vector with a length of N. |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
25 ## The @var{tol} specifies the tolerance of the method, default value is 1e-6. |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
26 ## The @var{maxit} specifies the maximum number of iteration, default value is MIN(20,N). |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
27 ## The @var{M1} specifies a preconditioner, can also be a function handler which returns M\X. |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
28 ## The @var{M2} combined with @var{M1} defines preconditioner as preconditioner=M1*M2. |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
29 ## The @var{x0} is initial guess, default value is zeros(N,1). |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
30 ## |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
31 ## @end deftypefn |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
32 |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
33 function [x, flag, relres, iter, resvec] = cgs (A, b, tol, maxit, M1, M2, x0) |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
34 |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
35 if (nargin < 2 || nargin > 7 || nargout > 5) |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
36 print_usage (); |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
37 elseif (!isnumeric (A) || rows (A) != columns (A)) |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
38 error ("cgs: first argument must be a n-by-n matrix"); |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
39 elseif (!isvector (b)) |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
40 error ("cgs: b must be a vector"); |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
41 elseif (rows (A) != rows (b)) |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
42 error ("cgs: first and second argument must have the same number of rows"); |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
43 elseif (nargin > 2 && !isscalar (tol)) |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
44 error ("cgs: tol must be a scalar"); |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
45 elseif (nargin > 3 && !isscalar (maxit)) |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
46 error ("cgs: maxit must be a scalar"); |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
47 elseif (nargin > 4 && ismatrix (M1) && (rows (M1) != rows (A) || columns (M1) != columns (A))) |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
48 error ("cgs: M1 must have the same number of rows and columns as A"); |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
49 elseif (nargin > 5 && (!ismatrix (M2) || rows (M2) != rows (A) || columns (M2) != columns (A))) |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
50 error ("cgs: M2 must have the same number of rows and columns as A"); |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
51 elseif (nargin > 6 && !isvector (x0)) |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
52 error ("cgs: x0 must be a vector"); |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
53 elseif (nargin > 6 && rows (x0) != rows (b)) |
8507 | 54 error ("cgs: x0 must have the same number of rows as b"); |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
55 endif |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
56 |
8507 | 57 ## Default tolerance. |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
58 if (nargin < 3) |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
59 tol = 1e-6; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
60 endif |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
61 |
8507 | 62 ## Default maximum number of iteration. |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
63 if (nargin < 4) |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
64 maxit = min (rows (b),20); |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
65 endif |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
66 |
8507 | 67 ## Left preconditioner. |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
68 precon = []; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
69 if (nargin == 5) |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
70 precon = M1; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
71 elseif (nargin > 5) |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
72 if (isparse(M1) && issparse(M2)) |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
73 precon = @(x) M1 * (M2 * x); |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
74 else |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
75 precon = M1 * M2; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
76 endif |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
77 endif |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
78 |
8507 | 79 ## Precon can also be a function. |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
80 if (nargin > 4 && isnumeric (precon)) |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
81 |
8507 | 82 ## We can compute inverse preconditioner and use quicker algorithm. |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
83 if (det (precon) != 0) |
8507 | 84 precon = inv (precon); |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
85 else |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
86 error ("cgs: preconditioner is ill conditioned"); |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
87 endif |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
88 |
8507 | 89 ## We must make test if preconditioner isn't ill conditioned. |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
90 if (isinf (cond (precon))); |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
91 error ("cgs: preconditioner is ill conditioned"); |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
92 endif |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
93 endif |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
94 |
8507 | 95 ## Specifies initial estimate x0. |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
96 if (nargin < 7) |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
97 x = zeros (rows (b), 1); |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
98 else |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
99 x = x0; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
100 endif |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
101 |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
102 relres = b - A * x; |
8507 | 103 ## Vector of the residual norms for each iteration. |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
104 resvec = [norm(relres)]; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
105 ro = 0; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
106 norm_b = norm (b); |
8507 | 107 ## Default behavior we don't reach tolerance tol within maxit iterations. |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
108 flag = 1; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
109 for iter = 1 : maxit |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
110 |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
111 if (nargin > 4 && isnumeric (precon)) |
8507 | 112 ## We have computed inverse matrix so we can use quick algorithm. |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
113 z = precon * relres; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
114 elseif (nargin > 4) |
8507 | 115 ## Our preconditioner is a function. |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
116 z = feval (precon, relres); |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
117 else |
8507 | 118 ## We don't use preconditioning. |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
119 z = relres; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
120 endif |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
121 |
8507 | 122 ## Cache. |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
123 ro_old = ro; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
124 ro = relres' * z; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
125 if (iter == 1) |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
126 p = z; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
127 else |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
128 beta = ro / ro_old; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
129 p = z + beta * p; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
130 endif |
8507 | 131 ## Cache. |
132 q = A * p; | |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
133 alpha = ro / (p' * q); |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
134 x = x + alpha * p; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
135 |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
136 relres = relres - alpha * q; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
137 resvec = [resvec; norm(relres)]; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
138 |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
139 relres_distance = resvec (end) / norm_b; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
140 if (relres_distance <= tol) |
8507 | 141 ## We reach tolerance tol within maxit iterations. |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
142 flag = 0; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
143 break; |
8507 | 144 elseif (resvec (end) == resvec (end - 1)) |
145 ## The method stagnates. | |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
146 flag = 3; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
147 break; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
148 endif |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
149 endfor; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
150 |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
151 relres = relres_distance; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
152 endfunction |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
153 |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
154 |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
155 |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
156 %!demo |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
157 %! % Solve system of A*x=b |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
158 %! A=[5 -1 3;-1 2 -2;3 -2 3] |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
159 %! b=[7;-1;4] |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
160 %! [a,b,c,d,e]=cgs(A,b) |