Mercurial > octave-nkf
annotate scripts/sparse/ilu.m @ 19088: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 |
---|---|
19087
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 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
5 ## 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
|
6 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
7 ## 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
|
8 ## 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
|
9 ## 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
|
10 ## version. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
11 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
12 ## 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
|
13 ## 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
|
14 ## 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
|
15 ## details. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
16 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
17 ## 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
|
18 ## 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
|
19 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
20 ## -*- texinfo -*- |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
21 ## @deftypefn {Function File} ilu (@var{A}, @var{setup}) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
22 ## @deftypefnx {Function File} {[@var{L}, @var{U}] =} ilu (@var{A}, @var{setup}) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
23 ## @deftypefnx {Function File} {[@var{L}, @var{U}, @var{P}] =} ilu (@var{A}, @var{setup}) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
24 ## ilu produces a unit lower triangular matrix, an upper triangular matrix, and |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
25 ## a permutation matrix. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
26 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
27 ## These incomplete factorizations may be useful as preconditioners for a system |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
28 ## of linear equations being solved by iterative methods such as BICG |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
29 ## (BiConjugate Gradients), GMRES (Generalized Minimum Residual Method). |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
30 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
31 ## @code{ilu (@var{A}, @var{setup})} computes the incomplete LU factorization |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
32 ## of @var{A}. @var{setup} is an input structure with up to five setup options. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
33 ## The fields must be named exactly as shown below. You can include any number |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
34 ## of these fields in the structure and define them in any order. Any |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
35 ## additional fields are ignored. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
36 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
37 ## @table @asis |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
38 ## @item type |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
39 ## Type of factorization. Values for type include: |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
40 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
41 ## @table @asis |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
42 ## @item @samp{nofill} |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
43 ## Performs ILU factorization with 0 level of fill in, known as ILU(0). With |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
44 ## type set to @samp{nofill}, only the milu setup option is used; all other |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
45 ## fields are ignored. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
46 ## @item @samp{crout} |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
47 ## Performs the Crout version of ILU factorization, known as ILUC. With type |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
48 ## set to @samp{crout}, only the droptol and milu setup options are used; all |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
49 ## other fields are ignored. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
50 ## @item @samp{ilutp} |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
51 ## (default) Performs ILU factorization with threshold and pivoting. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
52 ## @end table |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
53 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
54 ## If type is not specified, the ILU factorization with pivoting ILUTP is |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
55 ## performed. Pivoting is never performed with type set to @samp{nofill} or |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
56 ## @samp{crout}. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
57 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
58 ## @item droptol |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
59 ## Drop tolerance of the incomplete LU factorization. droptol is a non-negative |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
60 ## scalar. The default value is 0, which produces the complete LU factorization. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
61 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
62 ## The nonzero entries of U satisfy |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
63 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
64 ## @code{abs (@var{U}(i,j)) >= droptol * norm ((@var{A}:,j))} |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
65 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
66 ## with the exception of the diagonal entries, which are retained regardless of |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
67 ## satisfying the criterion. The entries of @var{L} are tested against the |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
68 ## local drop tolerance before being scaled by the pivot, so for nonzeros in |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
69 ## @var{L} |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
70 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
71 ## @code{abs(@var{L}(i,j)) >= droptol * norm(@var{A}(:,j))/@var{U}(j,j)}. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
72 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
73 ## @item milu |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
74 ## Modified incomplete LU factorization. Values for milu |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
75 ## include: |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
76 ## @table @asis |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
77 ## @item @samp{row} |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
78 ## Produces the row-sum modified incomplete LU factorization. Entries from the |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
79 ## newly-formed column of the factors are subtracted from the diagonal of the |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
80 ## upper triangular factor, @var{U}, preserving column sums. That is, |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
81 ## @code{@var{A} * e = @var{L} * @var{U} * e}, where e is the vector of ones. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
82 ## @item @samp{col} |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
83 ## Produces the column-sum modified incomplete LU factorization. Entries from |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
84 ## the newly-formed column of the factors are subtracted from the diagonal of |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
85 ## the upper triangular factor, @var{U}, preserving column sums. That is, |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
86 ## @code{e'*@var{A} = e'*@var{L}*@var{U}}. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
87 ## @item @samp{off} |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
88 ## (default) No modified incomplete LU factorization is produced. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
89 ## @end table |
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 ## @item udiag |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
92 ## If udiag is 1, any zeros on the diagonal of the upper |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
93 ## triangular factor are replaced by the local drop tolerance. The default is 0. |
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 ## @item thresh |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
96 ## Pivot threshold between 0 (forces diagonal pivoting) and 1, |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
97 ## the default, which always chooses the maximum magnitude entry in the column |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
98 ## to be the pivot. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
99 ## @end table |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
100 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
101 ## @code{ilu (@var{A},@var{setup})} returns |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
102 ## @code{@var{L} + @var{U} - speye (size (@var{A}))}, where @var{L} is a unit |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
103 ## lower triangular matrix and @var{U} is an upper triangular matrix. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
104 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
105 ## @code{[@var{L}, @var{U}] = ilu (@var{A},@var{setup})} returns a unit lower |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
106 ## triangular matrix in @var{L} and an upper triangular matrix in @var{U}. When |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
107 ## SETUP.type = 'ilutp', the role of @var{P} is determined by the value of |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
108 ## SETUP.milu. For SETUP.type == 'ilutp', one of the factors is permuted |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
109 ## based on the value of SETUP.milu. When SETUP.milu == 'row', U is a column |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
110 ## permuted upper triangular factor. Otherwise, L is a row-permuted unit lower |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
111 ## triangular factor. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
112 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
113 ## @code{[@var{L}, @var{U}, @var{P}] = ilu (@var{A},@var{setup})} returns a |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
114 ## unit lower triangular matrix in @var{L}, an upper triangular matrix in |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
115 ## @var{U}, and a permutation matrix in @var{P}. When SETUP.milu ~= 'row', @var{P} |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
116 ## is returned such that @var{L} and @var{U} are incomplete factors of @var{P}*@var{A}. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
117 ## When SETUP.milu == 'row', @var{P} is returned such that and @var{U} are |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
118 ## incomplete factors of A*P. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
119 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
120 ## @strong{NOTE}: ilu works on sparse square matrices only. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
121 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
122 ## EXAMPLES |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
123 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
124 ## @example |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
125 ## A = gallery('neumann', 1600) + speye(1600); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
126 ## setup.type = 'nofill'; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
127 ## nnz(A) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
128 ## ans = 7840 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
129 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
130 ## nnz(lu(A)) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
131 ## ans = 126478 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
132 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
133 ## nnz(ilu(A,setup)) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
134 ## ans = 7840 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
135 ## @end example |
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 ## This shows that @var{A} has 7840 nonzeros, the complete LU factorization has |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
138 ## 126478 nonzeros, and the incomplete LU factorization, with 0 level of |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
139 ## fill-in, has 7840 nonzeros, the same amount as @var{A}. Taken from: |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
140 ## http://www.mathworks.com/help/matlab/ref/ilu.html |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
141 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
142 ## @example |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
143 ## A = gallery ('wathen', 10, 10); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
144 ## b = sum (A,2); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
145 ## tol = 1e-8; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
146 ## maxit = 50; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
147 ## opts.type = 'crout'; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
148 ## opts.droptol = 1e-4; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
149 ## [L, U] = ilu (A, opts); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
150 ## x = bicg (A, b, tol, maxit, L, U); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
151 ## norm(A * x - b, inf) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
152 ## @end example |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
153 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
154 ## This example uses ILU as preconditioner for a random FEM-Matrix, which has a |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
155 ## bad condition. Without @var{L} and @var{U} BICG would not converge. |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
156 ## |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
157 ## @end deftypefn |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
158 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
159 function [L, U, P] = ilu (A, setup) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
160 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
161 if ((nargin > 2) || (nargin < 1) || (nargout > 3)) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
162 print_usage (); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
163 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
164 |
19088
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:
19087
diff
changeset
|
165 |
19087
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
166 % Check input matrix |
19088
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:
19087
diff
changeset
|
167 if (~issparse (A) || ~issquare (A)) |
19087
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
168 error ("ilu: Input A must be a 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
|
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 |
19088
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:
19087
diff
changeset
|
171 % If A is empty and sparse then return empty L, U and P |
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:
19087
diff
changeset
|
172 % 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:
19087
diff
changeset
|
173 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:
19087
diff
changeset
|
174 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:
19087
diff
changeset
|
175 U = 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:
19087
diff
changeset
|
176 P = 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:
19087
diff
changeset
|
177 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:
19087
diff
changeset
|
178 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:
19087
diff
changeset
|
179 |
19087
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
180 % 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
|
181 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
|
182 if (~isstruct (setup)) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
183 error ("ilu: Input 'setup' 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
|
184 endif |
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 setup = struct (); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
187 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
188 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
189 if (~isfield (setup, "type")) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
190 setup.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
|
191 else |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
192 type = tolower (getfield (setup, "type")); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
193 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
|
194 && (strcmp (type, "crout") == 0) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
195 && (strcmp (type, "ilutp") == 0)) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
196 error ("ilu: 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
|
197 else |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
198 setup.type = type; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
199 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
200 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
201 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
202 if (~isfield (setup, "droptol")) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
203 setup.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
|
204 else |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
205 if (~isscalar (setup.droptol) || (setup.droptol < 0)) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
206 error ("ilu: 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
|
207 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
208 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
209 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
210 if (~isfield (setup, "milu")) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
211 setup.milu = "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
|
212 else |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
213 milu = tolower (getfield (setup, "milu")); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
214 if ((strcmp (milu, "off") == 0) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
215 && (strcmp (milu, "col") == 0) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
216 && (strcmp (milu, "row") == 0)) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
217 error ("ilu: Invalid field \"milu\" 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
|
218 else |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
219 setup.milu = milu; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
220 endif |
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 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
223 if (~isfield (setup, "udiag")) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
224 setup.udiag = 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
|
225 else |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
226 if (~isscalar (setup.udiag) || ((setup.udiag ~= 0) && (setup.udiag ~= 1))) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
227 error ("ilu: Invalid field \"udiag\" 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
|
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 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
230 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
231 if (~isfield (setup, "thresh")) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
232 setup.thresh = 1; % set default |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
233 else |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
234 if (~isscalar (setup.thresh) || (setup.thresh < 0) || (setup.thresh > 1)) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
235 error ("ilu: Invalid field \"thresh\" 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
|
236 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
237 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
238 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
239 n = length (A); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
240 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
241 % Delegate to specialized ILU |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
242 switch (setup.type) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
243 case "nofill" |
19088
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:
19087
diff
changeset
|
244 [L, U] = __ilu0__ (A, setup.milu); |
19087
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
245 if (nargout == 3) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
246 P = speye (length (A)); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
247 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
248 case "crout" |
19088
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:
19087
diff
changeset
|
249 [L, U] = __iluc__ (A, setup.droptol, setup.milu); |
19087
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
250 if (nargout == 3) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
251 P = speye (length (A)); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
252 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
253 case "ilutp" |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
254 if (nargout == 2) |
19088
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:
19087
diff
changeset
|
255 [L, U] = __ilutp__ (A, setup.droptol, setup.thresh, setup.milu, setup.udiag); |
19087
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
256 elseif (nargout == 3) |
19088
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:
19087
diff
changeset
|
257 [L, U, P] = __ilutp__ (A, setup.droptol, setup.thresh, setup.milu, setup.udiag); |
19087
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
258 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
259 otherwise |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
260 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
|
261 endswitch |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
262 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
263 if (nargout == 1) |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
264 L = L + U - speye (n); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
265 endif |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
266 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
267 endfunction |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
268 |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
269 %!shared n, dtol, A |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
270 %! n = 1600; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
271 %! dtol = 0.1; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
272 %! A = gallery ('neumann', n) + speye (n); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
273 %!test |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
274 %! setup.type = 'nofill'; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
275 %! assert (nnz (ilu (A, setup)), 7840); |
19088
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:
19087
diff
changeset
|
276 %! # This test is taken from the mathworks and should work for full support. |
19087
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
277 %!test |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
278 %! setup.type = 'crout'; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
279 %! setup.milu = 'row'; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
280 %! setup.droptol = dtol; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
281 %! [L, U] = ilu (A, setup); |
19088
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:
19087
diff
changeset
|
282 %! e = ones (size (A, 2),1); |
19087
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
283 %! assert (norm (A*e - L*U*e), 1e-14, 1e-14); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
284 %!test |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
285 %! setup.type = 'crout'; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
286 %! setup.droptol = dtol; |
19088
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:
19087
diff
changeset
|
287 %! [L, U] = ilu(A, setup); |
19087
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
288 %! assert (norm (A - L * U, 'fro') / norm (A, 'fro'), 0.05, 1e-2); |
19088
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:
19087
diff
changeset
|
289 |
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:
19087
diff
changeset
|
290 %! # Tests for input validation |
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:
19087
diff
changeset
|
291 %!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:
19087
diff
changeset
|
292 %! [L, U] = ilu (sparse ([])); |
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:
19087
diff
changeset
|
293 %! assert (isempty (L), logical (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:
19087
diff
changeset
|
294 %! assert (isempty (U), logical (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:
19087
diff
changeset
|
295 %! setup.type = 'crout'; |
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:
19087
diff
changeset
|
296 %! [L, U] = ilu (sparse ([]), 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:
19087
diff
changeset
|
297 %! assert (isempty (L), logical (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:
19087
diff
changeset
|
298 %! assert (isempty (U), logical (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:
19087
diff
changeset
|
299 %! setup.type = 'ilutp'; |
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:
19087
diff
changeset
|
300 %! [L, U] = ilu (sparse ([]), 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:
19087
diff
changeset
|
301 %! assert (isempty (L), logical (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:
19087
diff
changeset
|
302 %! assert (isempty (U), logical (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:
19087
diff
changeset
|
303 %!error [L, U] = ilu (0); |
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:
19087
diff
changeset
|
304 %!error [L, U] = ilu ([]); |
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:
19087
diff
changeset
|
305 %!error [L, U] = ilu (sparse (0)); |
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:
19087
diff
changeset
|
306 %!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:
19087
diff
changeset
|
307 %! setup.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:
19087
diff
changeset
|
308 %!error [L, U] = ilu (A_tiny, 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:
19087
diff
changeset
|
309 %! setup.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:
19087
diff
changeset
|
310 %!error [L, U] = ilu (A_tiny, 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:
19087
diff
changeset
|
311 %! setup.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:
19087
diff
changeset
|
312 %!error [L, U] = ilu (A_tiny, 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:
19087
diff
changeset
|
313 %!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:
19087
diff
changeset
|
314 %! setup.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:
19087
diff
changeset
|
315 %!error [L, U] = ilu (A_tiny, 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:
19087
diff
changeset
|
316 %! setup.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:
19087
diff
changeset
|
317 %!error [L, U] = ilu (A_tiny, 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:
19087
diff
changeset
|
318 %! setup.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:
19087
diff
changeset
|
319 %!error [L, U] = ilu (A_tiny, 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:
19087
diff
changeset
|
320 %!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:
19087
diff
changeset
|
321 %! setup.thresh= -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:
19087
diff
changeset
|
322 %!error [L, U] = ilu (A_tiny, 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:
19087
diff
changeset
|
323 %! setup.thresh = 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:
19087
diff
changeset
|
324 %!error [L, U] = ilu (A_tiny, 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:
19087
diff
changeset
|
325 %! setup.thresh = []; |
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:
19087
diff
changeset
|
326 %!error [L, U] = ilu (A_tiny, 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:
19087
diff
changeset
|
327 %! setup.thresh = 2; |
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:
19087
diff
changeset
|
328 %!error [L, U] = ilu (A_tiny, 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:
19087
diff
changeset
|
329 %!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:
19087
diff
changeset
|
330 %! setup.diag = 0.5; |
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:
19087
diff
changeset
|
331 %!error [L, U] = ilu (A_tiny, 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:
19087
diff
changeset
|
332 %! setup.diag = []; |
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:
19087
diff
changeset
|
333 %!error [L, U] = ilu (A_tiny, 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:
19087
diff
changeset
|
334 %! setup.diag = -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:
19087
diff
changeset
|
335 %!error [L, U] = ilu (A_tiny, 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:
19087
diff
changeset
|
336 %!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:
19087
diff
changeset
|
337 %! setup.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:
19087
diff
changeset
|
338 %!error [L, U] = ilu (A_tiny, 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:
19087
diff
changeset
|
339 %! setup.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:
19087
diff
changeset
|
340 %!error [L, U] = ilu (A_tiny, 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:
19087
diff
changeset
|
341 %! setup.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:
19087
diff
changeset
|
342 %!error [L, U] = ilu (A_tiny, 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:
19087
diff
changeset
|
343 |
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:
19087
diff
changeset
|
344 %! # Check if the elements in U satisfy the non-dropping condition. |
19087
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
345 %!test |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
346 %! setup.type = 'crout'; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
347 %! setup.droptol = dtol; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
348 %! [L, U] = ilu (A, setup); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
349 %! for j = 1:n |
19088
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:
19087
diff
changeset
|
350 %! cmp_value = dtol * norm (A(:, j)); |
19087
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
351 %! non_zeros = nonzeros (U(:, j)); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
352 %! for i = 1:length (non_zeros); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
353 %! assert (abs (non_zeros (i)) >= cmp_value, logical (1)); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
354 %! endfor |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
355 %! endfor |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
356 %!test |
19088
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:
19087
diff
changeset
|
357 %! setup.type = 'ilutp'; |
19087
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
358 %! setup.droptol = dtol; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
359 %! [L, U] = ilu (A, setup); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
360 %! for j = 1:n |
19088
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:
19087
diff
changeset
|
361 %! cmp_value = dtol * norm (A(:, j)); |
19087
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
362 %! non_zeros = nonzeros (U(:, j)); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
363 %! for i = 1:length (non_zeros); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
364 %! assert (abs (non_zeros (i)) >= cmp_value, logical (1)); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
365 %! endfor |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
366 %! endfor |
19088
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:
19087
diff
changeset
|
367 |
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:
19087
diff
changeset
|
368 %! # Check that the complete LU factorisation with crout and ilutp algorithms |
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:
19087
diff
changeset
|
369 %! # output the same result. |
19087
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
370 %!test |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
371 %! setup.type = 'crout'; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
372 %! setup.droptol = 0; |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
373 %! [L1, U1] = ilu (A, setup); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
374 %! setup.type = 'ilutp'; |
19088
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:
19087
diff
changeset
|
375 %! setup.thresh = 0; |
19087
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
376 %! [L2, U2] = ilu (A, setup); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
377 %! assert (norm (L1 - L2, 'fro') / norm (L1, 'fro'), 0, eps); |
168c0aa9bb05
Added all the files related with ilu.m and ichol.m functions.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff
changeset
|
378 %! assert (norm (U1 - U2, 'fro') / norm (U1, 'fro'), 0, eps); |
19088
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:
19087
diff
changeset
|
379 |
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:
19087
diff
changeset
|
380 %! # Tests for real matrices of different sizes for ilu0, iluc and ilutp. |
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:
19087
diff
changeset
|
381 %! # The difference A - L*U should be not greater than eps because with 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:
19087
diff
changeset
|
382 %! # equaling 0, the LU complete factorization is performed. |
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:
19087
diff
changeset
|
383 %!shared n_tiny, n_small, n_medium, n_large, A_tiny, A_small, A_medium, A_large |
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:
19087
diff
changeset
|
384 %! n_tiny = 5; |
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:
19087
diff
changeset
|
385 %! n_small = 40; |
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:
19087
diff
changeset
|
386 %! n_medium = 600; |
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:
19087
diff
changeset
|
387 %! n_large = 10000; |
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:
19087
diff
changeset
|
388 %! A_tiny = spconvert ([1 4 2 3 3 4 2 5; 1 1 2 3 4 4 5 5; 1 2 3 4 5 6 7 8]'); |
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:
19087
diff
changeset
|
389 %! A_small = sprand (n_small, n_small, 1/n_small) + speye (n_small); |
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:
19087
diff
changeset
|
390 %! A_medium = sprand (n_medium, n_medium, 1/n_medium) + speye (n_medium); |
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:
19087
diff
changeset
|
391 %! A_large = sprand (n_large, n_large, 1/n_large/10) + speye (n_large); |
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:
19087
diff
changeset
|
392 %! |
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:
19087
diff
changeset
|
393 %!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:
19087
diff
changeset
|
394 %! setup.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:
19087
diff
changeset
|
395 %! [L, U] = ilu (A_tiny); |
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:
19087
diff
changeset
|
396 %! assert (norm (A_tiny - L * U, "fro") / norm (A_tiny, "fro"), 0, n_tiny * 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:
19087
diff
changeset
|
397 %!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:
19087
diff
changeset
|
398 %! setup.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:
19087
diff
changeset
|
399 %! [L, U] = ilu (A_small); |
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:
19087
diff
changeset
|
400 %! assert (norm (A_small - L * U, "fro") / norm (A_small, "fro"), 0, 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:
19087
diff
changeset
|
401 %!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:
19087
diff
changeset
|
402 %! setup.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:
19087
diff
changeset
|
403 %! [L, U] = ilu (A_medium); |
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:
19087
diff
changeset
|
404 %! assert (norm (A_medium - L * U, "fro") / norm (A_medium, "fro"), 0, 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:
19087
diff
changeset
|
405 %!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:
19087
diff
changeset
|
406 %! setup.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:
19087
diff
changeset
|
407 %! [L, U] = ilu (A_large); |
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:
19087
diff
changeset
|
408 %! assert (norm (A_large - L * U, "fro") / norm (A_large, "fro"), 0, 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:
19087
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:
19087
diff
changeset
|
410 %!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:
19087
diff
changeset
|
411 %! setup.type = "crout"; |
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:
19087
diff
changeset
|
412 %! [L, U] = ilu (A_tiny, 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:
19087
diff
changeset
|
413 %! assert (norm (A_tiny - L * U, "fro") / norm (A_tiny, "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:
19087
diff
changeset
|
414 %!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:
19087
diff
changeset
|
415 %! setup.type = "crout"; |
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:
19087
diff
changeset
|
416 %! [L, U] = ilu (A_small, 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:
19087
diff
changeset
|
417 %! assert (norm (A_small - L * U, "fro") / norm (A_small, "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:
19087
diff
changeset
|
418 %!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:
19087
diff
changeset
|
419 %! setup.type = "crout"; |
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:
19087
diff
changeset
|
420 %! [L, U] = ilu (A_medium, 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:
19087
diff
changeset
|
421 %! assert (norm (A_medium - L * U, "fro") / norm (A_medium, "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:
19087
diff
changeset
|
422 %!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:
19087
diff
changeset
|
423 %! setup.type = "crout"; |
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:
19087
diff
changeset
|
424 %! [L, U] = ilu (A_large, 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:
19087
diff
changeset
|
425 %! assert (norm (A_large - L * U, "fro") / norm (A_large, "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:
19087
diff
changeset
|
426 %! |
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:
19087
diff
changeset
|
427 %!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:
19087
diff
changeset
|
428 %! setup.type = "ilutp"; |
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:
19087
diff
changeset
|
429 %! setup.droptol = 0; |
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:
19087
diff
changeset
|
430 %! setup.thresh = 0; |
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:
19087
diff
changeset
|
431 %! [L, U] = ilu (A_tiny, 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:
19087
diff
changeset
|
432 %! assert (norm (A_tiny - L * U, "fro") / norm (A_tiny, "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:
19087
diff
changeset
|
433 %!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:
19087
diff
changeset
|
434 %! setup.type = "ilutp"; |
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:
19087
diff
changeset
|
435 %! setup.droptol = 0; |
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:
19087
diff
changeset
|
436 %! setup.thresh = 0; |
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:
19087
diff
changeset
|
437 %! [L, U] = ilu (A_small, 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:
19087
diff
changeset
|
438 %! assert (norm (A_small - L * U, "fro") / norm (A_small, "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:
19087
diff
changeset
|
439 %!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:
19087
diff
changeset
|
440 %! setup.type = "ilutp"; |
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:
19087
diff
changeset
|
441 %! setup.droptol = 0; |
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:
19087
diff
changeset
|
442 %! setup.thresh = 0; |
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:
19087
diff
changeset
|
443 %! [L, U] = ilu (A_medium, 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:
19087
diff
changeset
|
444 %! assert (norm (A_medium - L * U, "fro") / norm (A_medium, "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:
19087
diff
changeset
|
445 %!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:
19087
diff
changeset
|
446 %! setup.type = "ilutp"; |
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:
19087
diff
changeset
|
447 %! setup.droptol = 0; |
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:
19087
diff
changeset
|
448 %! setup.thresh = 0; |
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:
19087
diff
changeset
|
449 %! [L, U] = ilu (A_large, 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:
19087
diff
changeset
|
450 %! assert (norm (A_large - L * U, "fro") / norm (A_large, "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:
19087
diff
changeset
|
451 %! |
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:
19087
diff
changeset
|
452 |
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:
19087
diff
changeset
|
453 %! # Tests for complex matrices of different sizes for ilu0, iluc and ilutp. |
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:
19087
diff
changeset
|
454 %!shared n_tiny, n_small, n_medium, n_large, A_tiny, A_small, A_medium, A_large |
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:
19087
diff
changeset
|
455 %! n_tiny = 5; |
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:
19087
diff
changeset
|
456 %! n_small = 40; |
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:
19087
diff
changeset
|
457 %! n_medium = 600; |
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:
19087
diff
changeset
|
458 %! n_large = 10000; |
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:
19087
diff
changeset
|
459 %! A_tiny = spconvert ([1 4 2 3 3 4 2 5; 1 1 2 3 4 4 5 5; 1 2 3 4 5 6 7 8]'); |
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:
19087
diff
changeset
|
460 %! A_tiny(1,1) += 1i; |
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:
19087
diff
changeset
|
461 %! A_small = sprand(n_small, n_small, 1/n_small) + ... |
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:
19087
diff
changeset
|
462 %! i * sprand(n_small, n_small, 1/n_small) + speye (n_small); |
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:
19087
diff
changeset
|
463 %! A_medium = sprand(n_medium, n_medium, 1/n_medium) + ... |
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:
19087
diff
changeset
|
464 %! i * sprand(n_medium, n_medium, 1/n_medium) + speye (n_medium); |
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:
19087
diff
changeset
|
465 %! A_large = sprand(n_large, n_large, 1/n_large/10) + ... |
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:
19087
diff
changeset
|
466 %! i * sprand(n_large, n_large, 1/n_large/10) + speye (n_large); |
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:
19087
diff
changeset
|
467 %! |
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:
19087
diff
changeset
|
468 %!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:
19087
diff
changeset
|
469 %! setup.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:
19087
diff
changeset
|
470 %! [L, U] = ilu (A_tiny); |
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:
19087
diff
changeset
|
471 %! assert (norm (A_tiny - L * U, "fro") / norm (A_tiny, "fro"), 0, n_tiny * 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:
19087
diff
changeset
|
472 %!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:
19087
diff
changeset
|
473 %! setup.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:
19087
diff
changeset
|
474 %! [L, U] = ilu (A_small); |
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:
19087
diff
changeset
|
475 %! assert (norm (A_small - L * U, "fro") / norm (A_small, "fro"), 0, 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:
19087
diff
changeset
|
476 %!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:
19087
diff
changeset
|
477 %! setup.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:
19087
diff
changeset
|
478 %! [L, U] = ilu (A_medium); |
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:
19087
diff
changeset
|
479 %! assert (norm (A_medium - L * U, "fro") / norm (A_medium, "fro"), 0, 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:
19087
diff
changeset
|
480 %!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:
19087
diff
changeset
|
481 %! setup.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:
19087
diff
changeset
|
482 %! [L, U] = ilu (A_large); |
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:
19087
diff
changeset
|
483 %! assert (norm (A_large - L * U, "fro") / norm (A_large, "fro"), 0, 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:
19087
diff
changeset
|
484 %! |
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:
19087
diff
changeset
|
485 %!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:
19087
diff
changeset
|
486 %! setup.type = "crout"; |
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:
19087
diff
changeset
|
487 %! [L, U] = ilu (A_tiny, 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:
19087
diff
changeset
|
488 %! assert (norm (A_tiny - L * U, "fro") / norm (A_tiny, "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:
19087
diff
changeset
|
489 %!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:
19087
diff
changeset
|
490 %! setup.type = "crout"; |
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:
19087
diff
changeset
|
491 %! [L, U] = ilu (A_small, 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:
19087
diff
changeset
|
492 %! assert (norm (A_small - L * U, "fro") / norm (A_small, "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:
19087
diff
changeset
|
493 %!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:
19087
diff
changeset
|
494 %! setup.type = "crout"; |
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:
19087
diff
changeset
|
495 %! [L, U] = ilu (A_medium, 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:
19087
diff
changeset
|
496 %! assert (norm (A_medium - L * U, "fro") / norm (A_medium, "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:
19087
diff
changeset
|
497 %!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:
19087
diff
changeset
|
498 %! setup.type = "crout"; |
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:
19087
diff
changeset
|
499 %! [L, U] = ilu (A_large, 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:
19087
diff
changeset
|
500 %! assert (norm (A_large - L * U, "fro") / norm (A_large, "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:
19087
diff
changeset
|
501 %! |
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:
19087
diff
changeset
|
502 %!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:
19087
diff
changeset
|
503 %! setup.type = "ilutp"; |
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:
19087
diff
changeset
|
504 %! setup.droptol = 0; |
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:
19087
diff
changeset
|
505 %! setup.thresh = 0; |
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:
19087
diff
changeset
|
506 %! [L, U] = ilu (A_tiny, 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:
19087
diff
changeset
|
507 %! assert (norm (A_tiny - L * U, "fro") / norm (A_tiny, "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:
19087
diff
changeset
|
508 %!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:
19087
diff
changeset
|
509 %! setup.type = "ilutp"; |
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:
19087
diff
changeset
|
510 %! setup.droptol = 0; |
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:
19087
diff
changeset
|
511 %! setup.thresh = 0; |
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:
19087
diff
changeset
|
512 %! [L, U] = ilu (A_small, 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:
19087
diff
changeset
|
513 %! assert (norm (A_small - L * U, "fro") / norm (A_small, "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:
19087
diff
changeset
|
514 %!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:
19087
diff
changeset
|
515 %! setup.type = "ilutp"; |
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:
19087
diff
changeset
|
516 %! setup.droptol = 0; |
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:
19087
diff
changeset
|
517 %! setup.thresh = 0; |
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:
19087
diff
changeset
|
518 %! [L, U] = ilu (A_medium, 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:
19087
diff
changeset
|
519 %! assert (norm (A_medium - L * U, "fro") / norm (A_medium, "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:
19087
diff
changeset
|
520 %!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:
19087
diff
changeset
|
521 %! setup.type = "ilutp"; |
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:
19087
diff
changeset
|
522 %! setup.droptol = 0; |
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:
19087
diff
changeset
|
523 %! setup.thresh = 0; |
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:
19087
diff
changeset
|
524 %! [L, U] = ilu (A_large, 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:
19087
diff
changeset
|
525 %! assert (norm (A_large - L * U, "fro") / norm (A_large, "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:
19087
diff
changeset
|
526 |
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:
19087
diff
changeset
|
527 %! #Specific tests for ilutp |
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:
19087
diff
changeset
|
528 %!shared a1, a2 |
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:
19087
diff
changeset
|
529 %! a1 = sparse ([0 0 4 3 1; 5 1 2.3 2 4.5; 0 0 0 2 1;0 0 8 0 2.2; 0 0 9 9 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:
19087
diff
changeset
|
530 %! a2 = sparse ([3 1 0 0 4; 3 1 0 0 -2;0 0 8 0 0; 0 4 0 4 -4.5; 0 -1 0 0 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:
19087
diff
changeset
|
531 %!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:
19087
diff
changeset
|
532 %! setup.udiag = 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:
19087
diff
changeset
|
533 %! setup.type = "ilutp"; |
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:
19087
diff
changeset
|
534 %! setup.droptol = 0.2; |
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:
19087
diff
changeset
|
535 %! [L, U, P] = ilu (a1, 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:
19087
diff
changeset
|
536 %! assert (norm (U, "fro"), 17.4577, 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:
19087
diff
changeset
|
537 %! assert (norm (L, "fro"), 2.4192, 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:
19087
diff
changeset
|
538 %! setup.udiag = 0; |
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:
19087
diff
changeset
|
539 %!error [L, U, P] = ilu (a1, 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:
19087
diff
changeset
|
540 %! |
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:
19087
diff
changeset
|
541 %!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:
19087
diff
changeset
|
542 %! setup.type = "ilutp"; |
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:
19087
diff
changeset
|
543 %! setup.droptol = 0; |
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:
19087
diff
changeset
|
544 %! setup.thresh = 0; |
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:
19087
diff
changeset
|
545 %! setup.milu = "row"; |
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:
19087
diff
changeset
|
546 %!error [L, U] = ilu (a2, 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:
19087
diff
changeset
|
547 %! |