annotate libinterp/dldfcn/__ichol__.cc @ 19269:65554f5847ac

don't include oct-locbuf.h in header files unnecessarily * oct-stream.h: Don't include oct-locbuf.h. * Sparse-op-decls.h: New file. * module.mk (OPERATORS_INC): Include it in the list. * Sparse-op-defs.h: Move decls to Sparse-op-decls.h. * CSparse.h, boolSparse.h, dSparse.h: Include Sparse-op-decls.h instead of Sparse-op-defs.h. * help.cc, oct-map.cc, oct-stream.cc, quadcc.cc, strfind.cc, sub2ind.cc, utils.cc, xpow.cc, __delaunayn__.cc, __ichol__.cc, __ilu__.cc, __voronoi__.cc, convhulln.cc, ov-base-mat.cc, op-dm-scm.cc, op-dm-template.cc, op-fcm-fs.cc, op-fcs-fm.cc, op-i16-i16.cc, op-i32-i32.cc, op-i64-i64.cc, op-i8-i8.cc, op-pm-template.cc, op-ui16-ui16.cc, op-ui32-ui32.cc, op-ui64-ui64.cc, op-ui8-ui8.cc, pt-mat.cc, CMatrix.cc, CNDArray.cc, CSparse.cc, MatrixType.cc, boolMatrix.cc, boolNDArray.cc, boolSparse.cc, dMatrix.cc, dSparse.cc, fCMatrix.cc, fCNDArray.cc, fMatrix.cc, eigs-base.cc, oct-norm.cc, sparse-base-lu.cc, * mx-op-defs.h: Update list of include files accordingly. * sparse-mk-ops.awk: Update emitted list of include files accordingly.
author John W. Eaton <jwe@octave.org>
date Sun, 19 Oct 2014 08:42:58 -0400
parents af9c22e20a57
children 9f8ec58b5c74
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
1 /*
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
2
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
3 Copyright (C) 2014 Eduardo Ramos Fernández <eduradical951@gmail.com>
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
4 Copyright (C) 2013 Kai T. Ohlhus <k.ohlhus@gmail.com>
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
5
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
6 This file is part of Octave.
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
7
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
8 Octave is free software; you can redistribute it and/or modify it
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
9 under the terms of the GNU General Public License as published by the
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
10 Free Software Foundation; either version 3 of the License, or (at your
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
11 option) any later version.
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
12
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
13 Octave is distributed in the hope that it will be useful, but WITHOUT
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
14 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
15 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
16 for more details.
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
17
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
18 You should have received a copy of the GNU General Public License
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
19 along with Octave; see the file COPYING. If not, see
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
20 <http://www.gnu.org/licenses/>.
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
21
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
22 */
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
23
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
24 #ifdef HAVE_CONFIG_H
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
25 #include <config.h>
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
26 #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:
diff changeset
27
19269
65554f5847ac don't include oct-locbuf.h in header files unnecessarily
John W. Eaton <jwe@octave.org>
parents: 19080
diff changeset
28 #include "oct-locbuf.h"
65554f5847ac don't include oct-locbuf.h in header files unnecessarily
John W. Eaton <jwe@octave.org>
parents: 19080
diff changeset
29
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
30 #include "defun-dld.h"
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
31 #include "parse.h"
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
32
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
33 // Secondary functions for complex and real case used in ichol algorithms.
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
34 Complex ichol_mult_complex (Complex a, Complex b)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
35 {
19080
af9c22e20a57 Apply std::complex feature tests to ichol helper function
Mike Miller <mtmiller@ieee.org>
parents: 19055
diff changeset
36 #if defined (HAVE_CXX_COMPLEX_SETTERS)
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
37 b.imag (-std::imag (b));
19080
af9c22e20a57 Apply std::complex feature tests to ichol helper function
Mike Miller <mtmiller@ieee.org>
parents: 19055
diff changeset
38 #elif defined (HAVE_CXX_COMPLEX_REFERENCE_ACCESSORS)
af9c22e20a57 Apply std::complex feature tests to ichol helper function
Mike Miller <mtmiller@ieee.org>
parents: 19055
diff changeset
39 b.imag () = -std::imag (b);
af9c22e20a57 Apply std::complex feature tests to ichol helper function
Mike Miller <mtmiller@ieee.org>
parents: 19055
diff changeset
40 #else
af9c22e20a57 Apply std::complex feature tests to ichol helper function
Mike Miller <mtmiller@ieee.org>
parents: 19055
diff changeset
41 b = std::conj (b);
af9c22e20a57 Apply std::complex feature tests to ichol helper function
Mike Miller <mtmiller@ieee.org>
parents: 19055
diff changeset
42 #endif
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
43 return a * b;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
44 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
45
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
46 double ichol_mult_real (double a, double b)
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
47 {
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
48 return a * b;
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
49 }
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
50
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
51 bool ichol_checkpivot_complex (Complex pivot)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
52 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
53 if (pivot.imag () != 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:
diff changeset
54 {
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
55 error ("ichol: non-real pivot encountered. The matrix must be hermitian.");
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
56 return false;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
57 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
58 else if (pivot.real () < 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:
diff changeset
59 {
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
60 error ("ichol: negative pivot encountered");
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
61 return false;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
62 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
63 return true;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
64 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
65
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
66 bool ichol_checkpivot_real (double pivot)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
67 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
68 if (pivot < 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:
diff changeset
69 {
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
70 error ("ichol: negative pivot encountered");
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
71 return false;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
72 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
73 return true;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
74 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
75
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
76 template <typename octave_matrix_t, typename T, T (*ichol_mult) (T, T),
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
77 bool (*ichol_checkpivot) (T)>
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
78 void ichol_0 (octave_matrix_t& sm, const std::string michol = "off")
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
79 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
80
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
81 const octave_idx_type n = sm.cols ();
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
82 octave_idx_type j1, jend, j2, jrow, jjrow, j, jw, i, k, jj, Llist_len, r;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
83 T tl;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
84 char opt;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
85 enum {OFF, ON};
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
86 if (michol == "on")
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
87 opt = ON;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
88 else
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
89 opt = OFF;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
90
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
91 // Input matrix pointers
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
92 octave_idx_type* cidx = sm.cidx ();
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
93 octave_idx_type* ridx = sm.ridx ();
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
94 T* data = sm.data ();
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
95
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
96 // Working arrays
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
97 OCTAVE_LOCAL_BUFFER (octave_idx_type, Lfirst, n);
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
98 OCTAVE_LOCAL_BUFFER (octave_idx_type, Llist, n);
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
99 OCTAVE_LOCAL_BUFFER (octave_idx_type, iw, n);
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
100 OCTAVE_LOCAL_BUFFER (T, dropsums, n);
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
101
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
102 // Initialize working arrays
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
103 for (i = 0; i < n; i++)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
104 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
105 iw[i] = -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:
diff changeset
106 Llist[i] = -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:
diff changeset
107 Lfirst[i] = -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:
diff changeset
108 dropsums[i] = 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:
diff changeset
109 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
110
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
111 // Main loop
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
112 for (k = 0; k < n; k++)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
113 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
114 j1 = cidx[k];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
115 j2 = cidx[k+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:
diff changeset
116 for (j = j1; j < j2; j++)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
117 iw[ridx[j]] = j;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
118
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
119 jrow = Llist [k];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
120 // Iterate over each non-zero element in the actual row.
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
121 while (jrow != -1)
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
122 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
123 jjrow = Lfirst[jrow];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
124 jend = cidx[jrow+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:
diff changeset
125 for (jj = jjrow; jj < jend; jj++)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
126 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
127 r = ridx[jj];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
128 jw = iw[r];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
129 tl = ichol_mult (data[jj], data[jjrow]);
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
130 if (jw != -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:
diff changeset
131 data[jw] -= tl;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
132 else
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
133 // Because of the symmetry of the matrix, we know
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
134 // the drops in the column r are also in the column k.
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
135 if (opt == ON)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
136 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
137 dropsums[r] -= tl;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
138 dropsums[k] -= tl;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
139 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
140 }
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
141 // Update the linked list and the first entry of the actual column.
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
142 if ((jjrow + 1) < jend)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
143 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
144 Lfirst[jrow]++;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
145 j = jrow;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
146 jrow = Llist[jrow];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
147 Llist[j] = Llist[ridx[Lfirst[j]]];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
148 Llist[ridx[Lfirst[j]]] = j;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
149 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
150 else
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
151 jrow = Llist[jrow];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
152 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
153
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
154 if (opt == ON)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
155 data[j1] += dropsums[k];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
156
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
157 if (ridx[j1] != k)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
158 {
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
159 error ("ichol: encountered a pivot equal to 0");
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
160 break;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
161 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
162
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
163 if (! ichol_checkpivot (data[j1]))
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
164 break;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
165
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
166 data[cidx[k]] = std::sqrt (data[j1]);
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
167
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
168 // Update Llist and Lfirst with the k-column information. Also,
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
169 // scale the column elements by the pivot and reset the working array iw.
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
170 if (k < (n - 1))
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
171 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
172 iw[ridx[j1]] = -1;
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
173 for (i = j1 + 1; i < j2; i++)
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
174 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
175 iw[ridx[i]] = -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:
diff changeset
176 data[i] /= data[j1];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
177 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
178 Lfirst[k] = j1;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
179 if ((Lfirst[k] + 1) < j2)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
180 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
181 Lfirst[k]++;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
182 jjrow = ridx[Lfirst[k]];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
183 Llist[k] = Llist[jjrow];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
184 Llist[jjrow] = k;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
185 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
186 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
187 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
188 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
189
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
190 DEFUN_DLD (__ichol0__, args, nargout, "-*- texinfo -*-\n\
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
191 @deftypefn {Loadable Function} {@var{L} =} __ichol0__ (@var{A})\n\
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
192 @deftypefnx {Loadable Function} {@var{L} =} __ichol0__ (@var{A}, @var{michol})\n\
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
193 Undocumented internal function.\n\
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
194 @end deftypefn")
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
195
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
196 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
197 octave_value_list retval;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
198
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
199 int nargin = args.length ();
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
200 std::string michol = "off";
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
201
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
202 if (nargout > 1 || nargin < 1 || nargin > 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:
diff changeset
203 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
204 print_usage ();
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
205 return retval;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
206 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
207
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
208 if (nargin == 2)
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
209 michol = args(1).string_value ();
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
210
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
211 // In ICHOL0 algorithm the zero-pattern of the input matrix is preserved
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
212 // so it's structure does not change during the algorithm. The same input
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
213 // matrix is used to build the output matrix due to that fact.
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
214 octave_value_list param_list;
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
215 if (!args(0).is_complex_type ())
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
216 {
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
217 SparseMatrix sm = args(0).sparse_matrix_value ();
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
218 param_list.append (sm);
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
219 sm = feval ("tril", param_list)(0).sparse_matrix_value ();
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
220 ichol_0 <SparseMatrix, double, ichol_mult_real,
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
221 ichol_checkpivot_real> (sm, michol);
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
222 if (! error_state)
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
223 retval(0) = sm;
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
224 }
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
225 else
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
226 {
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
227 SparseComplexMatrix sm = args(0).sparse_complex_matrix_value ();
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
228 param_list.append (sm);
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
229 sm = feval ("tril", param_list)(0).sparse_complex_matrix_value ();
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
230 ichol_0 <SparseComplexMatrix, Complex, ichol_mult_complex,
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
231 ichol_checkpivot_complex> (sm, michol);
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
232 if (! error_state)
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
233 retval(0) = sm;
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
234 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
235
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
236 return retval;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
237 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
238
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
239 template <typename octave_matrix_t, typename T, T (*ichol_mult) (T, T),
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
240 bool (*ichol_checkpivot) (T)>
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
241 void ichol_t (const octave_matrix_t& sm, octave_matrix_t& L, const T* cols_norm,
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
242 const T droptol, const std::string michol = "off")
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
243
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
244 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
245
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
246 const octave_idx_type n = sm.cols ();
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
247 octave_idx_type j, jrow, jend, jjrow, jw, i, k, jj, Llist_len, total_len,
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
248 w_len, max_len, ind;
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
249 char opt;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
250 enum {OFF, ON};
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
251 if (michol == "on")
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
252 opt = ON;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
253 else
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
254 opt = OFF;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
255
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
256 // Input matrix pointers
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
257 octave_idx_type* cidx = sm.cidx ();
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
258 octave_idx_type* ridx = sm.ridx ();
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
259 T* data = sm.data ();
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
260
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
261 // Output matrix data structures. Because the final zero pattern pattern of
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
262 // the output matrix is not known due to fill-in elements, a heuristic
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
263 // approach has been adopted for memory allocation. The size of ridx_out_l
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
264 // and data_out_l is incremented 10% of their actual size (nnz (A) in the
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
265 // beginning). If that amount is less than n, their size is just incremented
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
266 // in n elements. This way the number of reallocations decreases throughout
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
267 // the process, obtaining a good performance.
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
268 max_len = sm.nnz ();
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
269 max_len += (0.1 * max_len) > n ? 0.1 * max_len : n;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
270 Array <octave_idx_type> cidx_out_l (dim_vector (n + 1, 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:
diff changeset
271 octave_idx_type* cidx_l = cidx_out_l.fortran_vec ();
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
272 Array <octave_idx_type> ridx_out_l (dim_vector (max_len ,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:
diff changeset
273 octave_idx_type* ridx_l = ridx_out_l.fortran_vec ();
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
274 Array <T> data_out_l (dim_vector (max_len, 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:
diff changeset
275 T* data_l = data_out_l.fortran_vec ();
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
276
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
277 // Working arrays
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
278 OCTAVE_LOCAL_BUFFER (T, w_data, n);
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
279 OCTAVE_LOCAL_BUFFER (octave_idx_type, Lfirst, n);
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
280 OCTAVE_LOCAL_BUFFER (octave_idx_type, Llist, n);
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
281 OCTAVE_LOCAL_BUFFER (T, col_drops, n);
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
282 std::vector <octave_idx_type> vec;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
283 vec.resize (n);
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
284
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
285 T zero = T (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:
diff changeset
286 cidx_l[0] = cidx[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:
diff changeset
287 for (i = 0; i < n; i++)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
288 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
289 Llist[i] = -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:
diff changeset
290 Lfirst[i] = -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:
diff changeset
291 w_data[i] = 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:
diff changeset
292 col_drops[i] = zero;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
293 vec[i] = 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:
diff changeset
294 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
295
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
296 total_len = 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:
diff changeset
297 for (k = 0; k < n; k++)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
298 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
299 ind = 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:
diff changeset
300 for (j = cidx[k]; j < cidx[k+1]; j++)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
301 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
302 w_data[ridx[j]] = data[j];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
303 if (ridx[j] != k)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
304 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
305 vec[ind] = ridx[j];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
306 ind++;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
307 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
308 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
309 jrow = Llist[k];
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
310 while (jrow != -1)
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
311 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
312 jjrow = Lfirst[jrow];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
313 jend = cidx_l[jrow+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:
diff changeset
314 for (jj = jjrow; jj < jend; jj++)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
315 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
316 j = ridx_l[jj];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
317 // If the element in the j position of the row is zero,
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
318 // then it will become non-zero, so we add it to the
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
319 // vector that tracks non-zero elements in the working row.
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
320 if (w_data[j] == zero)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
321 {
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
322 vec[ind] = j;
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
323 ind++;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
324 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
325 w_data[j] -= ichol_mult (data_l[jj], data_l[jjrow]);
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
326 }
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
327 // Update the actual column first element and
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
328 // update the linked list of the jrow row.
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
329 if ((jjrow + 1) < jend)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
330 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
331 Lfirst[jrow]++;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
332 j = jrow;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
333 jrow = Llist[jrow];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
334 Llist[j] = Llist[ridx_l[Lfirst[j]]];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
335 Llist[ridx_l[Lfirst[j]]] = j;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
336 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
337 else
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
338 jrow = Llist[jrow];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
339 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
340
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
341 // Resizing output arrays
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
342 if ((max_len - total_len) < n)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
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:
diff changeset
344 max_len += (0.1 * max_len) > n ? 0.1 * max_len : n;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
345 data_out_l.resize (dim_vector (max_len, 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:
diff changeset
346 data_l = data_out_l.fortran_vec ();
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
347 ridx_out_l.resize (dim_vector (max_len, 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:
diff changeset
348 ridx_l = ridx_out_l.fortran_vec ();
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
349 }
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
350
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
351 // The sorting of the non-zero elements of the working column can be
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
352 // handled in a couple of ways. The most efficient two I found, are
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
353 // keeping the elements in an ordered binary search tree dynamically or
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
354 // keep them unsorted in a vector and at the end of the outer iteration
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
355 // order them. The last approach exhibits lower execution times.
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
356 std::sort (vec.begin (), vec.begin () + ind);
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
357
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
358 data_l[total_len] = w_data[k];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
359 ridx_l[total_len] = k;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
360 w_len = 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:
diff changeset
361
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
362 // Extract the non-zero elements of working column and
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
363 // drop the elements that are lower than droptol * cols_norm[k].
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
364 for (i = 0; i < ind ; i++)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
365 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
366 jrow = vec[i];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
367 if (w_data[jrow] != zero)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
368 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
369 if (std::abs (w_data[jrow]) < (droptol * cols_norm[k]))
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
370 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
371 if (opt == ON)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
372 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
373 col_drops[k] += w_data[jrow];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
374 col_drops[jrow] += w_data[jrow];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
375 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
376 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
377 else
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
378 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
379 data_l[total_len + w_len] = w_data[jrow];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
380 ridx_l[total_len + w_len] = jrow;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
381 w_len++;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
382 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
383 vec[i] = 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:
diff changeset
384 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
385 w_data[jrow] = zero;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
386 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
387
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
388 // Compensate column sums --> michol option
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
389 if (opt == ON)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
390 data_l[total_len] += col_drops[k];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
391
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
392 if (data_l[total_len] == zero)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
393 {
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
394 error ("ichol: encountered a pivot equal to 0");
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
395 break;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
396 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
397 else if (! ichol_checkpivot (data_l[total_len]))
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
398 break;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
399
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
400 // Once elements are dropped and compensation of column sums are done,
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
401 // scale the elements by the pivot.
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
402 data_l[total_len] = std::sqrt (data_l[total_len]);
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
403 for (jj = total_len + 1; jj < (total_len + w_len); jj++)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
404 data_l[jj] /= data_l[total_len];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
405 total_len += w_len;
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
406 // Check if there are too many elements to be indexed with
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
407 // octave_idx_type type due to fill-in during the process.
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
408 if (total_len < 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:
diff changeset
409 {
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
410 error ("ichol: integer overflow. Too many fill-in elements in L");
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
411 break;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
412 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
413 cidx_l[k+1] = cidx_l[k] - cidx_l[0] + w_len;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
414
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
415 // Update Llist and Lfirst with the k-column information.
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
416 if (k < (n - 1))
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
417 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
418 Lfirst[k] = cidx_l[k];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
419 if ((Lfirst[k] + 1) < cidx_l[k+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:
diff changeset
420 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
421 Lfirst[k]++;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
422 jjrow = ridx_l[Lfirst[k]];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
423 Llist[k] = Llist[jjrow];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
424 Llist[jjrow] = k;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
425 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
426 }
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
427 }
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
428
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
429 if (! error_state)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
430 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
431 // Build the output matrices
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
432 L = octave_matrix_t (n, n, total_len);
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
433 for (i = 0; i <= n; i++)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
434 L.cidx (i) = cidx_l[i];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
435 for (i = 0; i < total_len; i++)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
436 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
437 L.ridx (i) = ridx_l[i];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
438 L.data (i) = data_l[i];
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
439 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
440 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
441 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
442
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
443 DEFUN_DLD (__icholt__, args, nargout, "-*- texinfo -*-\n\
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
444 @deftypefn {Loadable Function} {@var{L} =} __icholt__ (@var{A})\n\
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
445 @deftypefnx {Loadable Function} {@var{L} =} __icholt__ (@var{A}, @var{droptol})\n\
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
446 @deftypefnx {Loadable Function} {@var{L} =} __icholt__ (@var{A}, @var{droptol}, @var{michol})\n\
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
447 Undocumented internal function.\n\
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
448 @end deftypefn")
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
449 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
450 octave_value_list retval;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
451 int nargin = args.length ();
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
452 // Default values of parameters
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
453 std::string michol = "off";
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
454 double 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:
diff changeset
455
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
456 if (nargout > 1 || nargin < 1 || nargin > 3)
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
457 {
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
458 print_usage ();
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
459 return retval;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
460 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
461
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
462 // Don't repeat input validation of arguments done in ichol.m
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
463
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
464 if (nargin >= 2)
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
465 droptol = args(1).double_value ();
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
466
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
467 if (nargin == 3)
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
468 michol = args(2).string_value ();
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
469
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
470 octave_value_list param_list;
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
471 if (! args(0).is_complex_type ())
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
472 {
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
473 Array <double> cols_norm;
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
474 SparseMatrix L;
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
475 param_list.append (args(0).sparse_matrix_value ());
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
476 SparseMatrix sm_l =
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
477 feval ("tril", param_list)(0).sparse_matrix_value ();
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
478 param_list(0) = sm_l;
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
479 param_list(1) = 1;
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
480 param_list(2) = "cols";
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
481 cols_norm = feval ("norm", param_list)(0).vector_value ();
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
482 param_list.clear ();
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
483 ichol_t <SparseMatrix,
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
484 double, ichol_mult_real, ichol_checkpivot_real>
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
485 (sm_l, L, cols_norm.fortran_vec (), droptol, michol);
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
486 if (! error_state)
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
487 retval(0) = L;
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
488 }
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
489 else
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
490 {
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
491 Array <Complex> cols_norm;
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
492 SparseComplexMatrix L;
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
493 param_list.append (args(0).sparse_complex_matrix_value ());
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
494 SparseComplexMatrix sm_l =
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
495 feval ("tril", param_list)(0).sparse_complex_matrix_value ();
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
496 param_list(0) = sm_l;
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
497 param_list(1) = 1;
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
498 param_list(2) = "cols";
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
499 cols_norm = feval ("norm", param_list)(0).complex_vector_value ();
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
500 param_list.clear ();
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
501 ichol_t <SparseComplexMatrix,
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
502 Complex, ichol_mult_complex, ichol_checkpivot_complex>
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
503 (sm_l, L, cols_norm.fortran_vec (),
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
504 Complex (droptol), michol);
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
505 if (! error_state)
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
506 retval(0) = L;
19054
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
507 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
508
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
509 return retval;
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
510 }
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
511
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
512 /*
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
513 ## No test needed for internal helper function.
df64071e538c Removed ichol0.cc, icholt.cc, ilu0.cc, iluc.cc, ilutp.cc. Created __ichol__.cc and __ilu__.cc. Minor bugs fixed.
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
diff changeset
514 %!assert (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:
diff changeset
515 */
19055
38937efbee21 Incorporate new functions ichol and ilu into Octave.
Rik <rik@octave.org>
parents: 19054
diff changeset
516