Mercurial > octave-nkf
annotate src/DLD-FUNCTIONS/inv.cc @ 7911:3b46230f7a4d
DLD-FUNCTIONS/inv.cc (Finv): Avoid GCC warning
author | John W. Eaton <jwe@octave.org> |
---|---|
date | Wed, 09 Jul 2008 12:28:26 -0400 |
parents | 87865ed7405f |
children | 8b1a2555c4e2 |
rev | line source |
---|---|
2928 | 1 /* |
2 | |
7017 | 3 Copyright (C) 1996, 1997, 1999, 2000, 2001, 2002, 2004, 2005, 2006, |
4 2007 John W. Eaton | |
2928 | 5 |
6 This file is part of Octave. | |
7 | |
8 Octave is free software; you can redistribute it and/or modify it | |
9 under the terms of the GNU General Public License as published by the | |
7016 | 10 Free Software Foundation; either version 3 of the License, or (at your |
11 option) any later version. | |
2928 | 12 |
13 Octave is distributed in the hope that it will be useful, but WITHOUT | |
14 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or | |
15 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License | |
16 for more details. | |
17 | |
18 You should have received a copy of the GNU General Public License | |
7016 | 19 along with Octave; see the file COPYING. If not, see |
20 <http://www.gnu.org/licenses/>. | |
2928 | 21 |
22 */ | |
23 | |
24 #ifdef HAVE_CONFIG_H | |
25 #include <config.h> | |
26 #endif | |
27 | |
28 #include "defun-dld.h" | |
29 #include "error.h" | |
30 #include "gripes.h" | |
31 #include "oct-obj.h" | |
32 #include "utils.h" | |
33 | |
3808 | 34 DEFUN_DLD (inv, args, nargout, |
3548 | 35 "-*- texinfo -*-\n\ |
7650 | 36 @deftypefn {Loadable Function} {[@var{x}, @var{rcond}] =} inv (@var{a})\n\ |
37 @deftypefnx {Loadable Function} {[@var{x}, @var{rcond}] =} inverse (@var{a})\n\ | |
3808 | 38 Compute the inverse of the square matrix @var{a}. Return an estimate\n\ |
39 of the reciprocal condition number if requested, otherwise warn of an\n\ | |
40 ill-conditioned matrix if the reciprocal condition number is small.\n\ | |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
41 \n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
42 If called with a sparse matrix, then in general @var{x} will be a full\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
43 matrix, and so if possible forming the inverse of a sparse matrix should\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
44 be avoided. It is significantly more accurate and faster to do\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
45 @code{@var{y} = @var{a} \\ @var{b}}, rather than\n\ |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
46 @code{@var{y} = inv (@var{a}) * @var{b}}.\n\ |
3372 | 47 @end deftypefn") |
2928 | 48 { |
49 octave_value_list retval; | |
50 | |
51 int nargin = args.length (); | |
52 | |
53 if (nargin != 1) | |
54 { | |
5823 | 55 print_usage (); |
2928 | 56 return retval; |
57 } | |
58 | |
59 octave_value arg = args(0); | |
60 | |
5275 | 61 octave_idx_type nr = arg.rows (); |
62 octave_idx_type nc = arg.columns (); | |
2928 | 63 |
64 int arg_is_empty = empty_arg ("inverse", nr, nc); | |
65 | |
66 if (arg_is_empty < 0) | |
67 return retval; | |
68 else if (arg_is_empty > 0) | |
4233 | 69 return octave_value (Matrix ()); |
2928 | 70 |
71 if (nr != nc) | |
72 { | |
73 gripe_square_matrix_required ("inverse"); | |
74 return retval; | |
75 } | |
76 | |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
77 octave_value result; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
78 octave_idx_type info; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
79 double rcond = 0.0; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
80 float frcond = 0.0; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
81 bool isfloat = arg.is_single_type (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
82 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
83 if (isfloat) |
2928 | 84 { |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
85 if (arg.is_real_type ()) |
2928 | 86 { |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
87 FloatMatrix m = arg.float_matrix_value (); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
88 if (! error_state) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
89 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
90 MatrixType mattyp = args(0).matrix_type (); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
91 result = m.inverse (mattyp, info, frcond, 1); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
92 args(0).matrix_type (mattyp); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
93 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
94 } |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
95 else if (arg.is_complex_type ()) |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
96 { |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
97 FloatComplexMatrix m = arg.float_complex_matrix_value (); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
98 if (! error_state) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
99 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
100 MatrixType mattyp = args(0).matrix_type (); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
101 result = m.inverse (mattyp, info, frcond, 1); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
102 args(0).matrix_type (mattyp); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
103 } |
2928 | 104 } |
105 } | |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
106 else |
2928 | 107 { |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
108 if (arg.is_real_type ()) |
2928 | 109 { |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
110 if (arg.is_sparse_type ()) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
111 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
112 SparseMatrix m = arg.sparse_matrix_value (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
113 if (! error_state) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
114 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
115 MatrixType mattyp = args(0).matrix_type (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
116 result = m.inverse (mattyp, info, rcond, 1); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
117 args(0).matrix_type (mattyp); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
118 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
119 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
120 else |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
121 { |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
122 Matrix m = arg.matrix_value (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
123 if (! error_state) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
124 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
125 MatrixType mattyp = args(0).matrix_type (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
126 result = m.inverse (mattyp, info, rcond, 1); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
127 args(0).matrix_type (mattyp); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
128 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
129 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
130 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
131 else if (arg.is_complex_type ()) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
132 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
133 if (arg.is_sparse_type ()) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
134 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
135 SparseComplexMatrix m = arg.sparse_complex_matrix_value (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
136 if (! error_state) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
137 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
138 MatrixType mattyp = args(0).matrix_type (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
139 result = m.inverse (mattyp, info, rcond, 1); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
140 args(0).matrix_type (mattyp); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
141 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
142 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
143 else |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
144 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
145 ComplexMatrix m = arg.complex_matrix_value (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
146 if (! error_state) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
147 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
148 MatrixType mattyp = args(0).matrix_type (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
149 result = m.inverse (mattyp, info, rcond, 1); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
150 args(0).matrix_type (mattyp); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
151 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
152 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
153 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
154 else |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7650
diff
changeset
|
155 gripe_wrong_type_arg ("inv", arg); |
2928 | 156 } |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
157 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
158 if (! error_state) |
2928 | 159 { |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
160 if (nargout > 1) |
7911
3b46230f7a4d
DLD-FUNCTIONS/inv.cc (Finv): Avoid GCC warning
John W. Eaton <jwe@octave.org>
parents:
7814
diff
changeset
|
161 retval(1) = isfloat ? octave_value (frcond) : octave_value (rcond); |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
162 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
163 retval(0) = result; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
164 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
165 volatile double xrcond = rcond; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
166 xrcond += 1.0; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
167 if (nargout < 2 && (info == -1 || xrcond == 1.0)) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
168 warning ("inverse: matrix singular to machine precision, rcond = %g", |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
169 rcond); |
2928 | 170 } |
171 | |
172 return retval; | |
173 } | |
174 | |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
175 /* |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
176 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
177 %!assert(inv ([1, 2; 3, 4]), [-2, 1; 1.5, -0.5], sqrt (eps)) |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
178 %!assert(inv (single([1, 2; 3, 4])), single([-2, 1; 1.5, -0.5]), sqrt (eps ('single'))) |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
179 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
180 %!error <Invalid call to inv.*> inv (); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
181 %!error <Invalid call to inv.*> inv ([1, 2; 3, 4], 2); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
182 %!error inv ([1, 2; 3, 4; 5, 6]); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
183 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
184 */ |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
185 |
5775 | 186 // FIXME -- this should really be done with an alias, but |
2928 | 187 // alias_builtin() won't do the right thing if we are actually using |
188 // dynamic linking. | |
189 | |
190 DEFUN_DLD (inverse, args, nargout, | |
3458 | 191 "-*- texinfo -*-\n\ |
192 @deftypefn {Loadable Function} {} inverse (@var{a})\n\ | |
193 See inv.\n\ | |
194 @end deftypefn") | |
2928 | 195 { |
196 return Finv (args, nargout); | |
197 } | |
198 | |
199 /* | |
200 ;;; Local Variables: *** | |
201 ;;; mode: C++ *** | |
202 ;;; End: *** | |
203 */ |