Mercurial > octave
annotate src/DLD-FUNCTIONS/find.cc @ 8955:6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Because of find()'s count-limiting and direction arguments, this is
slightly more complicated than just copying the permutation vector.
I suspect this is a common operation for people who don't know about
the 'vector' option to lu().
author | Jason Riedy <jason@acm.org> |
---|---|
date | Tue, 10 Mar 2009 21:54:39 -0400 |
parents | eb63fbe60fab |
children | 12ca81f1fa99 853f96e8008f |
rev | line source |
---|---|
2928 | 1 /* |
2 | |
7017 | 3 Copyright (C) 1996, 1997, 1999, 2000, 2002, 2003, 2004, 2005, 2006, |
8920 | 4 2007, 2008, 2009 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 | |
4153 | 28 #include "quit.h" |
29 | |
2928 | 30 #include "defun-dld.h" |
31 #include "error.h" | |
32 #include "gripes.h" | |
33 #include "oct-obj.h" | |
34 | |
6002 | 35 // Find at most N_TO_FIND nonzero elements in NDA. Search forward if |
36 // DIRECTION is 1, backward if it is -1. NARGOUT is the number of | |
37 // output arguments. If N_TO_FIND is -1, find all nonzero elements. | |
4678 | 38 |
39 template <typename T> | |
40 octave_value_list | |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
41 find_nonzero_elem_idx (const Array<T>& nda, int nargout, |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
42 octave_idx_type n_to_find, int direction) |
4678 | 43 { |
6002 | 44 octave_value_list retval ((nargout == 0 ? 1 : nargout), Matrix ()); |
4678 | 45 |
5275 | 46 octave_idx_type count = 0; |
4678 | 47 |
5275 | 48 octave_idx_type nel = nda.nelem (); |
4678 | 49 |
6002 | 50 // Set the starting element to the correct value based on the |
51 // direction to search. | |
52 octave_idx_type k = 0; | |
53 if (direction == -1) | |
54 k = nel - 1; | |
55 | |
56 // Search in the default range. | |
57 octave_idx_type start_el = -1; | |
58 octave_idx_type end_el = -1; | |
59 | |
60 // Search for the number of elements to return. | |
61 while (k < nel && k > -1 && n_to_find != count) | |
4678 | 62 { |
63 OCTAVE_QUIT; | |
64 | |
8810
c9e1db15035b
eliminate unnecessary casts
John W. Eaton <jwe@octave.org>
parents:
7814
diff
changeset
|
65 if (nda(k) != T ()) |
6002 | 66 { |
67 end_el = k; | |
68 if (start_el == -1) | |
69 start_el = k; | |
70 count++; | |
71 } | |
72 k = k + direction; | |
4678 | 73 } |
74 | |
6002 | 75 // Reverse the range if we're looking backward. |
76 if (direction == -1) | |
77 { | |
78 octave_idx_type tmp_el = start_el; | |
79 start_el = end_el; | |
80 end_el = tmp_el; | |
81 } | |
82 // Fix an off by one error. | |
83 end_el++; | |
84 | |
4826 | 85 // If the original argument was a row vector, force a row vector of |
5130 | 86 // the overall indices to be returned. But see below for scalar |
87 // case... | |
4678 | 88 |
5275 | 89 octave_idx_type result_nr = count; |
90 octave_idx_type result_nc = 1; | |
4678 | 91 |
5130 | 92 bool scalar_arg = false; |
93 | |
4826 | 94 if (nda.ndims () == 2 && nda.rows () == 1) |
95 { | |
96 result_nr = 1; | |
97 result_nc = count; | |
5130 | 98 |
99 scalar_arg = (nda.columns () == 1); | |
4826 | 100 } |
3418 | 101 |
4826 | 102 Matrix idx (result_nr, result_nc); |
4678 | 103 |
4826 | 104 Matrix i_idx (result_nr, result_nc); |
105 Matrix j_idx (result_nr, result_nc); | |
4678 | 106 |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
107 ArrayN<T> val (dim_vector (result_nr, result_nc)); |
4826 | 108 |
109 if (count > 0) | |
4678 | 110 { |
4826 | 111 count = 0; |
112 | |
5275 | 113 octave_idx_type nr = nda.rows (); |
4678 | 114 |
5275 | 115 octave_idx_type i = 0; |
4826 | 116 |
6002 | 117 // Search for elements to return. Only search the region where |
118 // there are elements to be found using the count that we want | |
119 // to find. | |
6005 | 120 |
121 // For compatibility, all N-d arrays are handled as if they are | |
122 // 2-d, with the number of columns equal to "prod (dims (2:end))". | |
123 | |
6002 | 124 for (k = start_el; k < end_el; k++) |
4678 | 125 { |
4826 | 126 OCTAVE_QUIT; |
4678 | 127 |
8810
c9e1db15035b
eliminate unnecessary casts
John W. Eaton <jwe@octave.org>
parents:
7814
diff
changeset
|
128 if (nda(k) != T ()) |
4826 | 129 { |
130 idx(count) = k + 1; | |
4678 | 131 |
6005 | 132 octave_idx_type xr = k % nr; |
133 i_idx(count) = xr + 1; | |
134 j_idx(count) = (k - xr) / nr + 1; | |
4678 | 135 |
4826 | 136 val(count) = nda(k); |
137 | |
138 count++; | |
139 } | |
4678 | 140 |
4826 | 141 i++; |
4678 | 142 } |
143 } | |
5130 | 144 else if (scalar_arg) |
5131 | 145 { |
146 idx.resize (0, 0); | |
147 | |
148 i_idx.resize (0, 0); | |
149 j_idx.resize (0, 0); | |
150 | |
151 val.resize (dim_vector (0, 0)); | |
152 } | |
2928 | 153 |
154 switch (nargout) | |
155 { | |
6254 | 156 default: |
2928 | 157 case 3: |
158 retval(2) = val; | |
159 // Fall through! | |
160 | |
161 case 2: | |
3418 | 162 retval(1) = j_idx; |
163 retval(0) = i_idx; | |
2928 | 164 break; |
165 | |
6254 | 166 case 1: |
167 case 0: | |
168 retval(0) = idx; | |
2928 | 169 break; |
170 } | |
171 | |
172 return retval; | |
173 } | |
174 | |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
175 template octave_value_list find_nonzero_elem_idx (const Array<double>&, int, |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
176 octave_idx_type, int); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
177 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
178 template octave_value_list find_nonzero_elem_idx (const Array<Complex>&, int, |
6002 | 179 octave_idx_type, int); |
2928 | 180 |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
181 template octave_value_list find_nonzero_elem_idx (const Array<float>&, int, |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
182 octave_idx_type, int); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
183 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
184 template octave_value_list find_nonzero_elem_idx (const Array<FloatComplex>&, |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
185 int, octave_idx_type, int); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
186 |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
187 template <typename T> |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
188 octave_value_list |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
189 find_nonzero_elem_idx (const Sparse<T>& v, int nargout, |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
190 octave_idx_type n_to_find, int direction) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
191 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
192 octave_value_list retval ((nargout == 0 ? 1 : nargout), Matrix ()); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
193 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
194 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
195 octave_idx_type nc = v.cols(); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
196 octave_idx_type nr = v.rows(); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
197 octave_idx_type nz = v.nnz(); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
198 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
199 // Search in the default range. |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
200 octave_idx_type start_nc = -1; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
201 octave_idx_type end_nc = -1; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
202 octave_idx_type count; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
203 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
204 // Search for the range to search |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
205 if (n_to_find < 0) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
206 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
207 start_nc = 0; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
208 end_nc = nc; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
209 n_to_find = nz; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
210 count = nz; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
211 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
212 else if (direction > 0) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
213 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
214 for (octave_idx_type j = 0; j < nc; j++) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
215 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
216 OCTAVE_QUIT; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
217 if (v.cidx(j) == 0 && v.cidx(j+1) != 0) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
218 start_nc = j; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
219 if (v.cidx(j+1) >= n_to_find) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
220 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
221 end_nc = j + 1; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
222 break; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
223 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
224 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
225 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
226 else |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
227 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
228 for (octave_idx_type j = nc; j > 0; j--) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
229 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
230 OCTAVE_QUIT; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
231 if (v.cidx(j) == nz && v.cidx(j-1) != nz) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
232 end_nc = j; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
233 if (nz - v.cidx(j-1) >= n_to_find) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
234 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
235 start_nc = j - 1; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
236 break; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
237 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
238 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
239 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
240 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
241 count = (n_to_find > v.cidx(end_nc) - v.cidx(start_nc) ? |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
242 v.cidx(end_nc) - v.cidx(start_nc) : n_to_find); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
243 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
244 // If the original argument was a row vector, force a row vector of |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
245 // the overall indices to be returned. But see below for scalar |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
246 // case... |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
247 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
248 octave_idx_type result_nr = count; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
249 octave_idx_type result_nc = 1; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
250 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
251 bool scalar_arg = false; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
252 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
253 if (v.rows () == 1) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
254 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
255 result_nr = 1; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
256 result_nc = count; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
257 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
258 scalar_arg = (v.columns () == 1); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
259 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
260 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
261 Matrix idx (result_nr, result_nc); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
262 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
263 Matrix i_idx (result_nr, result_nc); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
264 Matrix j_idx (result_nr, result_nc); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
265 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
266 ArrayN<T> val (dim_vector (result_nr, result_nc)); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
267 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
268 if (count > 0) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
269 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
270 // Search for elements to return. Only search the region where |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
271 // there are elements to be found using the count that we want |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
272 // to find. |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
273 for (octave_idx_type j = start_nc, cx = 0; j < end_nc; j++) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
274 for (octave_idx_type i = v.cidx(j); i < v.cidx(j+1); i++ ) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
275 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
276 OCTAVE_QUIT; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
277 if (direction < 0 && i < nz - count) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
278 continue; |
8810
c9e1db15035b
eliminate unnecessary casts
John W. Eaton <jwe@octave.org>
parents:
7814
diff
changeset
|
279 i_idx(cx) = static_cast<double> (v.ridx(i) + 1); |
c9e1db15035b
eliminate unnecessary casts
John W. Eaton <jwe@octave.org>
parents:
7814
diff
changeset
|
280 j_idx(cx) = static_cast<double> (j + 1); |
c9e1db15035b
eliminate unnecessary casts
John W. Eaton <jwe@octave.org>
parents:
7814
diff
changeset
|
281 idx(cx) = j * nr + v.ridx(i) + 1; |
c9e1db15035b
eliminate unnecessary casts
John W. Eaton <jwe@octave.org>
parents:
7814
diff
changeset
|
282 val(cx) = v.data(i); |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
283 cx++; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
284 if (cx == count) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
285 break; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
286 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
287 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
288 else if (scalar_arg) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
289 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
290 idx.resize (0, 0); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
291 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
292 i_idx.resize (0, 0); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
293 j_idx.resize (0, 0); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
294 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
295 val.resize (dim_vector (0, 0)); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
296 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
297 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
298 switch (nargout) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
299 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
300 case 0: |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
301 case 1: |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
302 retval(0) = idx; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
303 break; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
304 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
305 case 5: |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
306 retval(4) = nc; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
307 // Fall through |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
308 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
309 case 4: |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
310 retval(3) = nr; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
311 // Fall through |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
312 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
313 case 3: |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
314 retval(2) = val; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
315 // Fall through! |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
316 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
317 case 2: |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
318 retval(1) = j_idx; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
319 retval(0) = i_idx; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
320 break; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
321 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
322 default: |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
323 panic_impossible (); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
324 break; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
325 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
326 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
327 return retval; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
328 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
329 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
330 template octave_value_list find_nonzero_elem_idx (const Sparse<double>&, int, |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
331 octave_idx_type, int); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
332 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
333 template octave_value_list find_nonzero_elem_idx (const Sparse<Complex>&, int, |
6002 | 334 octave_idx_type, int); |
2928 | 335 |
8955
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
336 octave_value_list |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
337 find_nonzero_elem_idx (const PermMatrix& v, int nargout, |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
338 octave_idx_type n_to_find, int direction) |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
339 { |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
340 // There are far fewer special cases to handle for a PermMatrix. |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
341 octave_value_list retval ((nargout == 0 ? 1 : nargout), Matrix ()); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
342 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
343 octave_idx_type nc = v.cols(); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
344 octave_idx_type start_nc, end_nc, count; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
345 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
346 // Determine the range to search. |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
347 if (n_to_find < 0 || n_to_find >= nc) |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
348 { |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
349 start_nc = 0; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
350 end_nc = nc; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
351 n_to_find = nc; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
352 count = nc; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
353 } |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
354 else if (direction > 0) |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
355 { |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
356 start_nc = 0; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
357 end_nc = n_to_find; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
358 count = n_to_find; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
359 } |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
360 else |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
361 { |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
362 start_nc = nc - n_to_find; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
363 end_nc = nc; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
364 count = n_to_find; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
365 } |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
366 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
367 bool scalar_arg = (v.rows () == 1 && v.cols () == 1); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
368 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
369 Matrix idx (count, 1); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
370 Matrix i_idx (count, 1); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
371 Matrix j_idx (count, 1); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
372 // Every value is 1. |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
373 ArrayN<double> val (dim_vector (count, 1), 1.0); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
374 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
375 if (count > 0) |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
376 { |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
377 const octave_idx_type* p = v.data (); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
378 if (v.is_col_perm ()) |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
379 { |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
380 for (octave_idx_type k = 0; k < count; k++) |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
381 { |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
382 OCTAVE_QUIT; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
383 const octave_idx_type j = start_nc + k; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
384 const octave_idx_type i = p[j]; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
385 i_idx(k) = static_cast<double> (1+i); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
386 j_idx(k) = static_cast<double> (1+j); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
387 idx(k) = j * nc + i + 1; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
388 } |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
389 } |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
390 else |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
391 { |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
392 for (octave_idx_type k = 0; k < count; k++) |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
393 { |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
394 OCTAVE_QUIT; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
395 const octave_idx_type i = start_nc + k; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
396 const octave_idx_type j = p[i]; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
397 // Scatter into the index arrays according to |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
398 // j adjusted by the start point. |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
399 const octave_idx_type koff = j - start_nc; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
400 i_idx(koff) = static_cast<double> (1+i); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
401 j_idx(koff) = static_cast<double> (1+j); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
402 idx(koff) = j * nc + i + 1; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
403 } |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
404 } |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
405 } |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
406 else if (scalar_arg) |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
407 { |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
408 // Same odd compatibility case as the other overrides. |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
409 idx.resize (0, 0); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
410 i_idx.resize (0, 0); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
411 j_idx.resize (0, 0); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
412 val.resize (dim_vector (0, 0)); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
413 } |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
414 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
415 switch (nargout) |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
416 { |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
417 case 0: |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
418 case 1: |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
419 retval(0) = idx; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
420 break; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
421 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
422 case 5: |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
423 retval(4) = nc; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
424 // Fall through |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
425 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
426 case 4: |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
427 retval(3) = nc; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
428 // Fall through |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
429 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
430 case 3: |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
431 retval(2) = val; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
432 // Fall through! |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
433 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
434 case 2: |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
435 retval(1) = j_idx; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
436 retval(0) = i_idx; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
437 break; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
438 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
439 default: |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
440 panic_impossible (); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
441 break; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
442 } |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
443 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
444 return retval; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
445 } |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
446 |
2928 | 447 DEFUN_DLD (find, args, nargout, |
3369 | 448 "-*- texinfo -*-\n\ |
449 @deftypefn {Loadable Function} {} find (@var{x})\n\ | |
6002 | 450 @deftypefnx {Loadable Function} {} find (@var{x}, @var{n})\n\ |
451 @deftypefnx {Loadable Function} {} find (@var{x}, @var{n}, @var{direction})\n\ | |
6524 | 452 Return a vector of indices of nonzero elements of a matrix, as a row if\n\ |
453 @var{x} is a row or as a column otherwise. To obtain a single index for\n\ | |
454 each matrix element, Octave pretends that the columns of a matrix form one\n\ | |
455 long vector (like Fortran arrays are stored). For example,\n\ | |
3369 | 456 \n\ |
457 @example\n\ | |
458 @group\n\ | |
459 find (eye (2))\n\ | |
460 @result{} [ 1; 4 ]\n\ | |
461 @end group\n\ | |
462 @end example\n\ | |
463 \n\ | |
464 If two outputs are requested, @code{find} returns the row and column\n\ | |
465 indices of nonzero elements of a matrix. For example,\n\ | |
466 \n\ | |
467 @example\n\ | |
468 @group\n\ | |
469 [i, j] = find (2 * eye (2))\n\ | |
470 @result{} i = [ 1; 2 ]\n\ | |
471 @result{} j = [ 1; 2 ]\n\ | |
472 @end group\n\ | |
473 @end example\n\ | |
474 \n\ | |
475 If three outputs are requested, @code{find} also returns a vector\n\ | |
476 containing the nonzero values. For example,\n\ | |
477 \n\ | |
478 @example\n\ | |
479 @group\n\ | |
480 [i, j, v] = find (3 * eye (2))\n\ | |
481 @result{} i = [ 1; 2 ]\n\ | |
482 @result{} j = [ 1; 2 ]\n\ | |
483 @result{} v = [ 3; 3 ]\n\ | |
484 @end group\n\ | |
485 @end example\n\ | |
6002 | 486 \n\ |
487 If two inputs are given, @var{n} indicates the number of elements to\n\ | |
488 find from the beginning of the matrix or vector.\n\ | |
489 \n\ | |
490 If three inputs are given, @var{direction} should be one of \"first\" or\n\ | |
491 \"last\" indicating that it should start counting found elements from the\n\ | |
492 first or last element.\n\ | |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
493 \n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
494 Note that this function is particularly useful for sparse matrices, as\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
495 it extracts the non-zero elements as vectors, which can then be used to\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
496 create the original matrix. For example,\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
497 \n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
498 @example\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
499 @group\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
500 sz = size(a);\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
501 [i, j, v] = find (a);\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
502 b = sparse(i, j, v, sz(1), sz(2));\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
503 @end group\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
504 @end example\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
505 @seealso{sparse}\n\ |
3369 | 506 @end deftypefn") |
2928 | 507 { |
508 octave_value_list retval; | |
509 | |
510 int nargin = args.length (); | |
511 | |
6002 | 512 if (nargin > 3 || nargin < 1) |
2928 | 513 { |
5823 | 514 print_usage (); |
2928 | 515 return retval; |
516 } | |
517 | |
6002 | 518 // Setup the default options. |
519 octave_idx_type n_to_find = -1; | |
520 if (nargin > 1) | |
521 { | |
522 n_to_find = args(1).int_value (); | |
523 if (error_state) | |
524 { | |
525 error ("find: expecting second argument to be an integer"); | |
526 return retval; | |
527 } | |
528 } | |
529 | |
530 // Direction to do the searching (1 == forward, -1 == reverse). | |
531 int direction = 1; | |
532 if (nargin > 2) | |
533 { | |
534 direction = 0; | |
535 | |
536 std::string s_arg = args(2).string_value (); | |
537 | |
538 if (! error_state) | |
539 { | |
540 if (s_arg == "first") | |
541 direction = 1; | |
542 else if (s_arg == "last") | |
543 direction = -1; | |
544 } | |
545 | |
546 if (direction == 0) | |
547 { | |
548 error ("find: expecting third argument to be \"first\" or \"last\""); | |
549 return retval; | |
550 } | |
551 } | |
552 | |
2928 | 553 octave_value arg = args(0); |
554 | |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
555 if (arg.is_sparse_type ()) |
2928 | 556 { |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
557 if (arg.is_real_type ()) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
558 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
559 SparseMatrix v = arg.sparse_matrix_value (); |
2928 | 560 |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
561 if (! error_state) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
562 retval = find_nonzero_elem_idx (v, nargout, |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
563 n_to_find, direction); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
564 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
565 else if (arg.is_complex_type ()) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
566 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
567 SparseComplexMatrix v = arg.sparse_complex_matrix_value (); |
5107 | 568 |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
569 if (! error_state) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
570 retval = find_nonzero_elem_idx (v, nargout, |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
571 n_to_find, direction); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
572 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
573 else |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
574 gripe_wrong_type_arg ("find", arg); |
5107 | 575 } |
8955
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
576 else if (arg.is_perm_matrix ()) { |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
577 PermMatrix P = arg.perm_matrix_value (); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
578 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
579 if (! error_state) |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
580 retval = find_nonzero_elem_idx (P, nargout, n_to_find, direction); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
581 } |
2928 | 582 else |
583 { | |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
584 if (arg.is_single_type ()) |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
585 { |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
586 if (arg.is_real_type ()) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
587 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
588 FloatNDArray nda = arg.float_array_value (); |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
589 |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
590 if (! error_state) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
591 retval = find_nonzero_elem_idx (nda, nargout, |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
592 n_to_find, direction); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
593 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
594 else if (arg.is_complex_type ()) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
595 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
596 FloatComplexNDArray cnda = arg.float_complex_array_value (); |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
597 |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
598 if (! error_state) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
599 retval = find_nonzero_elem_idx (cnda, nargout, |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
600 n_to_find, direction); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
601 } |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
602 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
603 else |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
604 { |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
605 if (arg.is_real_type ()) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
606 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
607 NDArray nda = arg.array_value (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
608 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
609 if (! error_state) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
610 retval = find_nonzero_elem_idx (nda, nargout, |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
611 n_to_find, direction); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
612 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
613 else if (arg.is_complex_type ()) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
614 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
615 ComplexNDArray cnda = arg.complex_array_value (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
616 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
617 if (! error_state) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
618 retval = find_nonzero_elem_idx (cnda, nargout, |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
619 n_to_find, direction); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
620 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
621 else if (arg.is_string ()) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
622 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
623 charNDArray cnda = arg.char_array_value (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
624 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
625 if (! error_state) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
626 retval = find_nonzero_elem_idx (cnda, nargout, |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
627 n_to_find, direction); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
628 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
629 else |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
630 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
631 gripe_wrong_type_arg ("find", arg); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
632 } |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
633 } |
2928 | 634 } |
635 | |
636 return retval; | |
637 } | |
638 | |
639 /* | |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
640 %!assert(find ([1, 0, 1, 0, 1]), [1, 3, 5]); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
641 %!assert(find ([1; 0; 3; 0; 1]), [1; 3; 5]); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
642 %!assert(find ([0, 0, 2; 0, 3, 0; -1, 0, 0]), [3; 5; 7]); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
643 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
644 %!test |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
645 %! [i, j, v] = find ([0, 0, 2; 0, 3, 0; -1, 0, 0]); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
646 %! |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
647 %! assert(i, [3; 2; 1]); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
648 %! assert(j, [1; 2; 3]); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
649 %! assert(v, [-1; 3; 2]); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
650 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
651 %!assert(find (single([1, 0, 1, 0, 1])), [1, 3, 5]); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
652 %!assert(find (single([1; 0; 3; 0; 1])), [1; 3; 5]); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
653 %!assert(find (single([0, 0, 2; 0, 3, 0; -1, 0, 0])), [3; 5; 7]); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
654 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
655 %!test |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
656 %! [i, j, v] = find (single([0, 0, 2; 0, 3, 0; -1, 0, 0])); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
657 %! |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
658 %! assert(i, [3; 2; 1]); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
659 %! assert(j, [1; 2; 3]); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
660 %! assert(v, single([-1; 3; 2])); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
661 |
8955
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
662 %!test |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
663 %! pcol = [5 1 4 3 2]; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
664 %! P = eye (5) (:, pcol); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
665 %! [i, j, v] = find (P); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
666 %! [ifull, jfull, vfull] = find (full (P)); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
667 %! assert (i, ifull); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
668 %! assert (j, jfull); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
669 %! assert (all (v == 1)); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
670 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
671 %!test |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
672 %! prow = [5 1 4 3 2]; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
673 %! P = eye (5) (prow, :); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
674 %! [i, j, v] = find (P); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
675 %! [ifull, jfull, vfull] = find (full (P)); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
676 %! assert (i, ifull); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
677 %! assert (j, jfull); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
678 %! assert (all (v == 1)); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
679 |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
680 %!error <Invalid call to find.*> find (); |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
681 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
682 */ |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
683 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
684 /* |
2928 | 685 ;;; Local Variables: *** |
686 ;;; mode: C++ *** | |
687 ;;; End: *** | |
688 */ |