Mercurial > hg > octave-nkf
annotate scripts/sparse/ichol.m @ 19246:df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
* libinterp/dldfcn/__ichol__.cc: File created now contains __ichol0__ and __icholt__ functions.
* libinterp/dldfcn/__ilu__.cc: File created now contains __ilu0__ __iluc__ and __ilutp__ functions.
* scripts/sparse/ichol.m: Tests have been moved from .cc files to this one.
* changed scripts/sparse/ilu.m: Tests have been moved from .cc files to this one.
author | Eduardo Ramos (edu159) <eduradical951@gmail.com> |
---|---|
date | Mon, 18 Aug 2014 12:32:16 +0100 |
parents | 168c0aa9bb05 |
children | 38937efbee21 |
rev | line source |
---|---|
19245
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
1 ## Copyright (C) 2013 Kai T. Ohlhus <k.ohlhus@gmail.com> |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
2 ## Copyright (C) 2014 Eduardo Ramos Fernández <eduradical951@gmail.com> |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
3 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
4 ## This file is part of Octave. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
5 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
6 ## Octave is free software; you can redistribute it and/or modify it under the |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
7 ## terms of the GNU General Public License as published by the Free Software |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
8 ## Foundation; either version 3 of the License, or (at your option) any later |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
9 ## version. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
10 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
11 ## Octave is distributed in the hope that it will be useful, but WITHOUT ANY |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
12 ## WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
13 ## FOR A PARTICULAR PURPOSE. See the GNU General Public License for more |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
14 ## details. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
15 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
16 ## You should have received a copy of the GNU General Public License along with |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
17 ## Octave; see the file COPYING. If not, see <http://www.gnu.org/licenses/>. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
18 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
19 ## -*- texinfo -*- |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
20 ## @deftypefn {Function File} ichol (@var{A}, @var{opts}) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
21 ## @deftypefnx {Function File} {@var{L} =} ichol (@var{A}, @var{opts}) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
22 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
23 ## @code{@var{L} = ichol (@var{A})} performs the incomplete Cholesky |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
24 ## factorization of A with zero-fill. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
25 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
26 ## @code{@var{L} = ichol (@var{A}, @var{opts})} performs the incomplete Cholesky |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
27 ## factorization of A with options specified by opts. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
28 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
29 ## By default, ichol references the lower triangle of A and produces lower |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
30 ## triangular factors. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
31 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
32 ## The factor given by this routine may be useful as a preconditioner for a |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
33 ## system of linear equations being solved by iterative methods such as |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
34 ## PCG (Preconditioned conjugate gradient). |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
35 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
36 ## ichol works only for sparse square matrices. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
37 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
38 ## The fields of @var{opts} must be named exactly as shown below. You can |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
39 ## include any number of these fields in the structure and define them in any |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
40 ## order. Any additional fields are ignored. Names and specifiers are case |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
41 ## sensitive. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
42 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
43 ## @table @asis |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
44 ## @item type |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
45 ## Type of factorization. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
46 ## String indicating which flavor of incomplete Cholesky to perform. Valid |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
47 ## values of this field are @samp{nofill} and @samp{ict}. The |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
48 ## @samp{nofill} variant performs incomplete Cholesky with zero-fill [IC(0)]. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
49 ## The @samp{ict} variant performs incomplete Cholesky with threshold dropping |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
50 ## [ICT]. The default value is @samp{nofill}. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
51 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
52 ## @item droptol |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
53 ## Drop tolerance when type is @samp{ict}. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
54 ## Nonnegative scalar used as a drop tolerance when performing ICT. Elements |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
55 ## which are smaller in magnitude than a local drop tolerance are dropped from |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
56 ## the resulting factor except for the diagonal element which is never dropped. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
57 ## The local drop tolerance at step j of the factorization is |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
58 ## @code{norm (@var{A}(j:end, j), 1) * droptol}. @samp{droptol} is ignored if |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
59 ## @samp{type} is @samp{nofill}. The default value is 0. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
60 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
61 ## @item michol |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
62 ## Indicates whether to perform modified incomplete Cholesky. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
63 ## Indicates whether or not modified incomplete Cholesky [MIC] is performed. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
64 ## The field may be @samp{on} or @samp{off}. When performing MIC, the diagonal |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
65 ## is compensated for dropped elements to enforce the relationship |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
66 ## @code{@var{A} * @var{e} = @var{L} * @var{L}' * @var{e}} where |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
67 ## @code{@var{e} = ones (size (@var{A}, 2), 1))}. The default value is |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
68 ## @samp{off}. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
69 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
70 ## @item diagcomp |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
71 ## Perform compensated incomplete Cholesky with the specified coefficient. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
72 ## Real nonnegative scalar used as a global diagonal shift @code{@var{alpha}} |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
73 ## in forming the incomplete Cholesky factor. That is, instead of performing |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
74 ## incomplete Cholesky on @code{@var{A}}, the factorization of |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
75 ## @code{@var{A} + @var{alpha} * diag (diag (@var{A}))} is formed. The default |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
76 ## value is 0. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
77 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
78 ## @item shape |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
79 ## Determines which triangle is referenced and returned. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
80 ## Valid values are @samp{upper} and @samp{lower}. If @samp{upper} is specified, |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
81 ## only the upper triangle of @code{@var{A}} is referenced and @code{@var{R}} |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
82 ## is constructed such that @code{@var{A}} is approximated by |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
83 ## @code{@var{R}' * @var{R}}. If @samp{lower} is specified, only the lower |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
84 ## triangle of @code{@var{A}} is referenced and @code{@var{L}} is constructed |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
85 ## such that @code{@var{A}} is approximated by @code{@var{L} * @var{L}'}. The |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
86 ## default value is @samp{lower}. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
87 ## @end table |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
88 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
89 ## EXAMPLES |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
90 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
91 ## The following problem demonstrates how to factorize a sample symmetric |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
92 ## positive definite matrix with the full Cholesky decomposition and with the |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
93 ## incomplete one. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
94 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
95 ## @example |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
96 ## A = [ 0.37, -0.05, -0.05, -0.07; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
97 ## -0.05, 0.116, 0.0, -0.05; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
98 ## -0.05, 0.0, 0.116, -0.05; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
99 ## -0.07, -0.05, -0.05, 0.202]; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
100 ## A = sparse(A); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
101 ## nnz(tril (A)) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
102 ## ans = 9 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
103 ## L = chol(A, "lower"); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
104 ## nnz (L) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
105 ## ans = 10 |
19246
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
106 ## norm (A - L * L', "fro") / norm (A, "fro") |
19245
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
107 ## ans = 1.1993e-16 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
108 ## opts.type = 'nofill'; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
109 ## L = ichol(A,opts); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
110 ## nnz (L) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
111 ## ans = 9 |
19246
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
112 ## norm (A - L * L', "fro") / norm (A, "fro") |
19245
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
113 ## ans = 0.019736 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
114 ## @end example |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
115 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
116 ## Another example for decomposition is finite difference matrix to solve a |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
117 ## boundary value problem on the unit square. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
118 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
119 ## @example |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
120 ## nx = 400; ny = 200; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
121 ## hx = 1 / (nx + 1); hy = 1 / (ny + 1); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
122 ## Dxx = spdiags ([ones(nx, 1), -2 * ones(nx, 1), ones(nx, 1)], [-1 0 1 ], nx, nx) / (hx ^ 2); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
123 ## Dyy = spdiags ([ones(ny, 1), -2 * ones(ny, 1), ones(ny, 1)], [-1 0 1 ], ny, ny) / (hy ^ 2); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
124 ## A = -kron (Dxx, speye (ny)) - kron (speye (nx), Dyy); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
125 ## nnz (tril (A)) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
126 ## ans = 239400 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
127 ## opts.type = 'nofill'; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
128 ## L = ichol (A, opts); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
129 ## nnz (tril (A)) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
130 ## ans = 239400 |
19246
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
131 ## norm (A - L * L', "fro") / norm (A, "fro") |
19245
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
132 ## ans = 0.062327 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
133 ## @end example |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
134 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
135 ## References for the implemented algorithms: |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
136 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
137 ## [1] Saad, Yousef. "Preconditioning Techniques." Iterative Methods for Sparse Linear |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
138 ## Systems. PWS Publishing Company, 1996. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
139 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
140 ## [2] Jones, Mark T. and Plassmann, Paul E.: An Improved Incomplete Cholesky |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
141 ## Factorization (1992). |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
142 ## @end deftypefn |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
143 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
144 function [L] = ichol (A, opts) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
145 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
146 if ((nargin > 2) || (nargin < 1) || (nargout > 1)) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
147 print_usage (); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
148 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
149 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
150 % Check input matrix |
19246
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
151 if (~issparse(A) || ~issquare (A)) |
19245
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
152 error ("ichol: Input A must be a non-empty sparse square matrix"); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
153 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
154 |
19246
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
155 % If A is empty and sparse then return empty L |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
156 % Compatibility with Matlab |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
157 if (isempty (A)) |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
158 L = A; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
159 return; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
160 endif |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
161 |
19245
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
162 % Check input structure, otherwise set default values |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
163 if (nargin == 2) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
164 if (~isstruct (opts)) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
165 error ("ichol: Input \"opts\" must be a valid structure."); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
166 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
167 else |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
168 opts = struct (); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
169 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
170 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
171 if (~isfield (opts, "type")) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
172 opts.type = "nofill"; % set default |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
173 else |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
174 type = tolower (getfield (opts, "type")); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
175 if ((strcmp (type, "nofill") == 0) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
176 && (strcmp (type, "ict") == 0)) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
177 error ("ichol: Invalid field \"type\" in input structure."); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
178 else |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
179 opts.type = type; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
180 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
181 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
182 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
183 if (~isfield (opts, "droptol")) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
184 opts.droptol = 0; % set default |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
185 else |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
186 if (~isscalar (opts.droptol) || (opts.droptol < 0)) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
187 error ("ichol: Invalid field \"droptol\" in input structure."); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
188 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
189 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
190 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
191 michol = ""; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
192 if (~isfield (opts, "michol")) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
193 opts.michol = "off"; % set default |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
194 else |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
195 michol = tolower (getfield (opts, "michol")); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
196 if ((strcmp (michol, "off") == 0) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
197 && (strcmp (michol, "on") == 0)) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
198 error ("ichol: Invalid field \"michol\" in input structure."); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
199 else |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
200 opts.michol = michol; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
201 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
202 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
203 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
204 if (~isfield (opts, "diagcomp")) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
205 opts.diagcomp = 0; % set default |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
206 else |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
207 if (~isscalar (opts.diagcomp) || (opts.diagcomp < 0)) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
208 error ("ichol: Invalid field \"diagcomp\" in input structure."); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
209 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
210 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
211 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
212 if (~isfield (opts, "shape")) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
213 opts.shape = "lower"; % set default |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
214 else |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
215 shape = tolower (getfield (opts, "shape")); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
216 if ((strcmp (shape, "lower") == 0) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
217 && (strcmp (shape, "upper") == 0)) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
218 error ("ichol: Invalid field \"shape\" in input structure."); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
219 else |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
220 opts.shape = shape; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
221 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
222 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
223 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
224 % Prepare input for specialized ICHOL |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
225 A_in = []; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
226 if (opts.diagcomp > 0) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
227 A += opts.diagcomp * diag (diag (A)); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
228 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
229 if (strcmp (opts.shape, "upper") == 1) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
230 A_in = triu (A); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
231 A_in = A_in'; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
232 else |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
233 A_in = tril (A); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
234 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
235 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
236 % Delegate to specialized ICHOL |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
237 switch (opts.type) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
238 case "nofill" |
19246
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
239 L = __ichol0__ (A_in, opts.michol); |
19245
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
240 case "ict" |
19246
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
241 L = __icholt__ (A_in, opts.droptol, opts.michol); |
19245
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
242 otherwise |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
243 printf ("The input structure is invalid.\n"); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
244 endswitch |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
245 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
246 if (strcmp (opts.shape, "upper") == 1) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
247 L = L'; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
248 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
249 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
250 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
251 endfunction |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
252 |
19246
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
253 %!shared A1, A2, A3, A4, A5, A6, A7 |
19245
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
254 %! A1 = [ 0.37, -0.05, -0.05, -0.07; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
255 %! -0.05, 0.116, 0.0, -0.05; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
256 %! -0.05, 0.0, 0.116, -0.05; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
257 %! -0.07, -0.05, -0.05, 0.202]; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
258 %! A1 = sparse(A1); |
19246
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
259 %! A2 = gallery ('poisson', 30); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
260 %! A3 = gallery ('tridiag', 50); |
19245
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
261 %! nx = 400; ny = 200; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
262 %! hx = 1 / (nx + 1); hy = 1 / (ny + 1); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
263 %! Dxx = spdiags ([ones(nx, 1), -2 * ones(nx, 1), ones(nx, 1)], [-1 0 1 ], nx, nx) / (hx ^ 2); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
264 %! Dyy = spdiags ([ones(ny, 1), -2 * ones(ny, 1), ones(ny, 1)], [-1 0 1 ], ny, ny) / (hy ^ 2); |
19246
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
265 %! A4 = -kron (Dxx, speye (ny)) - kron (speye (nx), Dyy); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
266 %! A5 = [ 0.37, -0.05, -0.05, -0.07; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
267 %! -0.05, 0.116, 0.0, -0.05 + 0.05i; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
268 %! -0.05, 0.0, 0.116, -0.05; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
269 %! -0.07, -0.05 - 0.05i, -0.05, 0.202]; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
270 %! A5 = sparse(A5); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
271 %! A6 = [ 0.37, -0.05 - i, -0.05, -0.07; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
272 %! -0.05 + i, 0.116, 0.0, -0.05; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
273 %! -0.05, 0.0, 0.116, -0.05; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
274 %! -0.07, -0.05, -0.05, 0.202]; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
275 %! A6 = sparse(A6); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
276 %! A7 = A5; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
277 %! A7(1) = 2i; |
19245
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
278 %! |
19246
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
279 |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
280 %!# Input validation tests |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
281 |
19245
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
282 %!test |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
283 %!error ichol ([]); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
284 %!error ichol (0); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
285 %!error ichol (-0); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
286 %!error ichol (1); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
287 %!error ichol (-1); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
288 %!error ichol (i); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
289 %!error ichol (-i); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
290 %!error ichol (1 + 1i); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
291 %!error ichol (1 - 1i); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
292 %!error ichol (sparse (0)); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
293 %!error ichol (sparse (-0)); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
294 %!error ichol (sparse (-1)); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
295 %!error ichol (sparse (-1)); |
19246
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
296 %!test |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
297 %! opts.milu = 'foo'; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
298 %!error L = ichol (A1, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
299 %! opts.milu = 1; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
300 %!error L = ichol (A1, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
301 %! opts.milu = []; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
302 %!error L = ichol (A1, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
303 %!test |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
304 %! opts.droptol = -1; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
305 %!error L = ichol (A1, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
306 %! opts.droptol = 0.5i; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
307 %!error L = ichol (A1, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
308 %! opts.droptol = []; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
309 %!error L = ichol (A1, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
310 %!test |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
311 %! opts.type = 'foo'; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
312 %!error L = ichol (A1, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
313 %! opts.type = 1; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
314 %!error L = ichol (A1, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
315 %! opts.type = []; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
316 %!error L = ichol (A1, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
317 %!test |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
318 %! opts.shape = 'foo'; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
319 %!error L = ichol (A1, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
320 %! opts.shape = 1; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
321 %!error L = ichol (A1, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
322 %! opts.shape = []; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
323 %!error L = ichol (A1, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
324 %!test |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
325 %! opts.diagcomp = -1; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
326 %!error L = ichol (A1, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
327 %! opts.diagcomp = 0.5i; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
328 %!error L = ichol (A1, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
329 %! opts.diagcomp = []; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
330 %!error L = ichol (A1, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
331 |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
332 %!# ICHOL0 tests |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
333 |
19245
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
334 %!test |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
335 %! opts.type = "nofill"; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
336 %! opts.michol = "off"; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
337 %! assert (nnz (tril (A1)), nnz (ichol (A1, opts))); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
338 %! assert (nnz (tril (A2)), nnz (ichol (A2, opts))); |
19246
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
339 %! assert (nnz (tril (A3)), nnz (ichol (A3, opts))); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
340 %! assert (nnz (tril (A4)), nnz (ichol (A4, opts))); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
341 %! assert (nnz (tril (A5)), nnz (ichol (A5, opts))); |
19245
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
342 %! |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
343 %!test |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
344 %! opts.type = "nofill"; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
345 %! opts.michol = "off"; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
346 %! L = ichol (A1, opts); |
19246
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
347 %! assert (norm (A1 - L * L', "fro") / norm (A1, "fro"), 0.0197, 1e-4); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
348 %! opts.shape = "upper"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
349 %! U = ichol (A1, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
350 %! assert (norm (A1 - U' * U, "fro") / norm (A1, "fro"), 0.0197, 1e-4); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
351 %! opts.shape = "lower"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
352 %! L = ichol (A1, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
353 %! assert (norm (A1 - L * L', "fro") / norm (A1, "fro"), 0.0197, 1e-4); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
354 %! |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
355 %!test |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
356 %! opts.michol = "on"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
357 %! opts.shape = "lower"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
358 %! opts.type = "nofill"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
359 %! L = ichol (A1, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
360 %! assert (norm (A1 - L * L', "fro") / norm (A1, "fro"), 0.0279, 1e-4); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
361 %! opts.shape = "upper"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
362 %! U = ichol (A1, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
363 %! assert (norm (A1 - U' * U, "fro") / norm (A1, "fro"), 0.0279, 1e-4); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
364 %! opts.shape = "lower"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
365 %! opts.diagcomp = 3e-3; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
366 %! L = ichol (A1, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
367 %! assert (norm (A1 - L * L', "fro") / norm (A1, "fro"), 0.0272, 1e-4); |
19245
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
368 %! |
19246
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
369 %!test |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
370 %! opts.type = "nofill"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
371 %! opts.michol = "off"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
372 %! L = ichol (A2, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
373 %! assert (norm (A2 - L*L', "fro") / norm (A2, "fro"), 0.0893, 1e-4) |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
374 %! opts.michol = "on"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
375 %! L = ichol (A2, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
376 %! assert (norm (A2 - L*L', "fro") / norm (A2, "fro"), 0.2377, 1e-4) |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
377 %! |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
378 %!test |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
379 %! opts.type = "nofill"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
380 %! opts.michol = "off"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
381 %! L = ichol (A3, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
382 %! assert (norm (A3 - L*L', "fro") / norm (A3, "fro"), eps, eps); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
383 %! opts.michol = "on"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
384 %! L = ichol (A3, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
385 %! assert (norm (A3 - L*L', "fro") / norm (A3, "fro"), eps, eps); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
386 %! |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
387 %!test |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
388 %! opts.type = "nofill"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
389 %! opts.michol = "off"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
390 %! L = ichol (A4, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
391 %! assert (norm (A4 - L*L', "fro") / norm (A4, "fro"), 0.0623, 1e-4); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
392 %! opts.michol = "on"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
393 %! L = ichol (A4, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
394 %! assert (norm (A4 - L*L', "fro") / norm (A4, "fro"), 0.1664, 1e-4); |
19245
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
395 %! |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
396 %!test |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
397 %! opts.type = "nofill"; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
398 %! opts.michol = "off"; |
19246
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
399 %! L = ichol (A5, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
400 %! assert (norm (A5 - L*L', "fro") / norm (A5, "fro"), 0.0195, 1e-4); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
401 %! opts.michol = "on"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
402 %! L = ichol (A5, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
403 %! assert (norm (A5 - L*L', "fro") / norm (A5, "fro"), 0.0276, 1e-4); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
404 %!test |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
405 %% Negative pivot |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
406 %!error ichol (A6); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
407 %% Complex entry in the diagonal |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
408 %!error ichol (A7); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
409 |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
410 %%!ICHOLT tests |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
411 |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
412 %%!test |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
413 %! opts.type = "ict"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
414 %! opts.droptol = 1e-1; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
415 %! opts.michol = "off"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
416 %! L = ichol (A1, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
417 %! assert (norm (A1 - L * L', "fro") / norm (A1, "fro"), 0.2065, 1e-4); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
418 %! opts.shape = "upper"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
419 %! U = ichol (A1, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
420 %! assert (norm (A1 - U' * U, "fro") / norm (A1, "fro"), 0.2065, 1e-4); |
19245
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
421 %! opts.shape = "lower"; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
422 %! L = ichol (A1, opts); |
19246
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
423 %! assert (norm (A1 - L * L', "fro") / norm (A1, "fro"), 0.2065, 1e-4); |
19245
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
424 %! |
19246
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
425 %%!test |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
426 %! opts.type = "ict"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
427 %! opts.droptol = 1e-1; |
19245
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
428 %! opts.michol = "on"; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
429 %! L = ichol (A1, opts); |
19246
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
430 %! assert (norm (A1 - L * L', "fro") / norm (A1, "fro"), 0.3266, 1e-4); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
431 %! opts.shape = "upper"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
432 %! U = ichol (A1, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
433 %! assert (norm (A1 - U' * U, "fro") / norm (A1, "fro"), 0.3266, 1e-4); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
434 %! opts.shape = "lower"; |
19245
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
435 %! opts.diagcomp = 3e-3; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
436 %! L = ichol (A1, opts); |
19246
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
437 %! assert (norm (A1 - L * L', "fro") / norm (A1, "fro"), 0.3266, 1e-4); |
19245
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
438 %! |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
439 %!test |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
440 %! opts.type = "ict"; |
19246
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
441 %! opts.droptol = 1e-1; |
19245
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
442 %! opts.michol = "off"; |
19246
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
443 %! L = ichol (A2, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
444 %! assert (norm (A2 - L*L', "fro") / norm (A2, "fro"), 0.0893, 1e-4) |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
445 %! opts.michol = "on"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
446 %! L = ichol (A2, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
447 %! assert (norm (A2 - L*L', "fro") / norm (A2, "fro"), 0.2377, 1e-4) |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
448 %! |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
449 %!test |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
450 %! opts.type = "ict"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
451 %! opts.droptol = 1e-1; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
452 %! opts.michol = "off"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
453 %! L = ichol (A3, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
454 %! assert (norm (A3 - L*L', "fro") / norm (A3, "fro"), eps, eps); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
455 %! opts.michol = "on"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
456 %! L = ichol (A3, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
457 %! assert (norm (A3 - L*L', "fro") / norm (A3, "fro"), eps, eps); |
19245
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
458 %! |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
459 %!test |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
460 %! opts.type = "ict"; |
19246
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
461 %! opts.droptol = 1e-1; |
19245
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
462 %! opts.michol = "off"; |
19246
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
463 %! L = ichol (A4, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
464 %! assert (norm (A4 - L*L', "fro") / norm (A4, "fro"), 0.1224, 1e-4); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
465 %! opts.michol = "on"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
466 %! L = ichol (A4, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
467 %! assert (norm (A4 - L*L', "fro") / norm (A4, "fro"), 0.2118, 1e-4); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
468 %! |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
469 %!test |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
470 %! opts.type = "ict"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
471 %! opts.droptol = 1e-1; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
472 %! opts.michol = "off"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
473 %! L = ichol (A5, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
474 %! assert (norm (A5 - L*L', "fro") / norm (A5, "fro"), 0.2044, 1e-4); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
475 %! opts.michol = "on"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
476 %! L = ichol (A5, opts); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
477 %! assert (norm (A5 - L*L', "fro") / norm (A5, "fro"), 0.3231, 1e-4); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
478 %!test |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
479 %% Negative pivot |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
480 %! opts.type = "ict"; |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
481 %!error ichol (A6, setup); |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
482 %% Complex entry in the diagonal |
df64071e538c
Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
19245
diff
changeset
|
483 %!error ichol (A7, setup); |