Mercurial > hg > octave-jordi
annotate src/DLD-FUNCTIONS/find.cc @ 9019:12ca81f1fa99
compatibility fix for find called for empty arguments
author | John W. Eaton <jwe@octave.org> |
---|---|
date | Wed, 25 Mar 2009 19:02:15 -0400 |
parents | 6d3fcbf89267 |
children | e67dc11ed6e8 97aa01a85ea4 |
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 |
9019
12ca81f1fa99
compatibility fix for find called for empty arguments
John W. Eaton <jwe@octave.org>
parents:
8955
diff
changeset
|
92 bool column_vector_arg = false; |
5130 | 93 bool scalar_arg = false; |
94 | |
9019
12ca81f1fa99
compatibility fix for find called for empty arguments
John W. Eaton <jwe@octave.org>
parents:
8955
diff
changeset
|
95 if (nda.ndims () == 2) |
4826 | 96 { |
9019
12ca81f1fa99
compatibility fix for find called for empty arguments
John W. Eaton <jwe@octave.org>
parents:
8955
diff
changeset
|
97 octave_idx_type nr = nda.rows (); |
12ca81f1fa99
compatibility fix for find called for empty arguments
John W. Eaton <jwe@octave.org>
parents:
8955
diff
changeset
|
98 octave_idx_type nc = nda.columns (); |
5130 | 99 |
9019
12ca81f1fa99
compatibility fix for find called for empty arguments
John W. Eaton <jwe@octave.org>
parents:
8955
diff
changeset
|
100 if (nr == 1) |
12ca81f1fa99
compatibility fix for find called for empty arguments
John W. Eaton <jwe@octave.org>
parents:
8955
diff
changeset
|
101 { |
12ca81f1fa99
compatibility fix for find called for empty arguments
John W. Eaton <jwe@octave.org>
parents:
8955
diff
changeset
|
102 result_nr = 1; |
12ca81f1fa99
compatibility fix for find called for empty arguments
John W. Eaton <jwe@octave.org>
parents:
8955
diff
changeset
|
103 result_nc = count; |
12ca81f1fa99
compatibility fix for find called for empty arguments
John W. Eaton <jwe@octave.org>
parents:
8955
diff
changeset
|
104 |
12ca81f1fa99
compatibility fix for find called for empty arguments
John W. Eaton <jwe@octave.org>
parents:
8955
diff
changeset
|
105 scalar_arg = (nc == 1); |
12ca81f1fa99
compatibility fix for find called for empty arguments
John W. Eaton <jwe@octave.org>
parents:
8955
diff
changeset
|
106 } |
12ca81f1fa99
compatibility fix for find called for empty arguments
John W. Eaton <jwe@octave.org>
parents:
8955
diff
changeset
|
107 else if (nc == 1) |
12ca81f1fa99
compatibility fix for find called for empty arguments
John W. Eaton <jwe@octave.org>
parents:
8955
diff
changeset
|
108 column_vector_arg = true; |
4826 | 109 } |
3418 | 110 |
4826 | 111 Matrix idx (result_nr, result_nc); |
4678 | 112 |
4826 | 113 Matrix i_idx (result_nr, result_nc); |
114 Matrix j_idx (result_nr, result_nc); | |
4678 | 115 |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
116 ArrayN<T> val (dim_vector (result_nr, result_nc)); |
4826 | 117 |
118 if (count > 0) | |
4678 | 119 { |
4826 | 120 count = 0; |
121 | |
5275 | 122 octave_idx_type nr = nda.rows (); |
4678 | 123 |
5275 | 124 octave_idx_type i = 0; |
4826 | 125 |
6002 | 126 // Search for elements to return. Only search the region where |
127 // there are elements to be found using the count that we want | |
128 // to find. | |
6005 | 129 |
130 // For compatibility, all N-d arrays are handled as if they are | |
131 // 2-d, with the number of columns equal to "prod (dims (2:end))". | |
132 | |
6002 | 133 for (k = start_el; k < end_el; k++) |
4678 | 134 { |
4826 | 135 OCTAVE_QUIT; |
4678 | 136 |
8810
c9e1db15035b
eliminate unnecessary casts
John W. Eaton <jwe@octave.org>
parents:
7814
diff
changeset
|
137 if (nda(k) != T ()) |
4826 | 138 { |
139 idx(count) = k + 1; | |
4678 | 140 |
6005 | 141 octave_idx_type xr = k % nr; |
142 i_idx(count) = xr + 1; | |
143 j_idx(count) = (k - xr) / nr + 1; | |
4678 | 144 |
4826 | 145 val(count) = nda(k); |
146 | |
147 count++; | |
148 } | |
4678 | 149 |
4826 | 150 i++; |
4678 | 151 } |
152 } | |
9019
12ca81f1fa99
compatibility fix for find called for empty arguments
John W. Eaton <jwe@octave.org>
parents:
8955
diff
changeset
|
153 else if (scalar_arg || (nda.rows () == 0 && ! column_vector_arg)) |
5131 | 154 { |
155 idx.resize (0, 0); | |
156 | |
157 i_idx.resize (0, 0); | |
158 j_idx.resize (0, 0); | |
159 | |
160 val.resize (dim_vector (0, 0)); | |
161 } | |
2928 | 162 |
163 switch (nargout) | |
164 { | |
6254 | 165 default: |
2928 | 166 case 3: |
167 retval(2) = val; | |
168 // Fall through! | |
169 | |
170 case 2: | |
3418 | 171 retval(1) = j_idx; |
172 retval(0) = i_idx; | |
2928 | 173 break; |
174 | |
6254 | 175 case 1: |
176 case 0: | |
177 retval(0) = idx; | |
2928 | 178 break; |
179 } | |
180 | |
181 return retval; | |
182 } | |
183 | |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
184 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
|
185 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
|
186 |
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 octave_value_list find_nonzero_elem_idx (const Array<Complex>&, int, |
6002 | 188 octave_idx_type, int); |
2928 | 189 |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
190 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
|
191 octave_idx_type, int); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
192 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
193 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
|
194 int, octave_idx_type, int); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
195 |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
196 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
|
197 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
|
198 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
|
199 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
|
200 { |
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_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
|
202 |
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 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
|
205 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
|
206 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
|
207 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
208 // 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
|
209 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
|
210 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
|
211 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
|
212 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
213 // 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
|
214 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
|
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 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
|
217 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
|
218 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
|
219 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
|
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 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
|
222 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
223 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
|
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 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
|
226 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
|
227 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
|
228 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
|
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 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
|
231 break; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
232 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
233 } |
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 else |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
236 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
237 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
|
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 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
|
240 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
|
241 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
|
242 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
|
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 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
|
245 break; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
246 } |
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 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
249 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
250 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
|
251 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
|
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 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
|
254 // 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
|
255 // case... |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
256 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
257 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
|
258 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
|
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 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
|
261 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
262 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
|
263 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
264 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
|
265 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
|
266 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
267 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
|
268 } |
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 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
|
271 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
272 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
|
273 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
|
274 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
275 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
|
276 |
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 (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
|
278 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
279 // 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
|
280 // 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
|
281 // 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
|
282 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
|
283 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
|
284 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
285 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
|
286 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
|
287 continue; |
8810
c9e1db15035b
eliminate unnecessary casts
John W. Eaton <jwe@octave.org>
parents:
7814
diff
changeset
|
288 i_idx(cx) = static_cast<double> (v.ridx(i) + 1); |
c9e1db15035b
eliminate unnecessary casts
John W. Eaton <jwe@octave.org>
parents:
7814
diff
changeset
|
289 j_idx(cx) = static_cast<double> (j + 1); |
c9e1db15035b
eliminate unnecessary casts
John W. Eaton <jwe@octave.org>
parents:
7814
diff
changeset
|
290 idx(cx) = j * nr + v.ridx(i) + 1; |
c9e1db15035b
eliminate unnecessary casts
John W. Eaton <jwe@octave.org>
parents:
7814
diff
changeset
|
291 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
|
292 cx++; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
293 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
|
294 break; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
295 } |
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 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
|
298 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
299 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
|
300 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
301 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
|
302 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
|
303 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
304 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
|
305 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
306 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
307 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
|
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 0: |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
310 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
|
311 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
|
312 break; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
313 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
314 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
|
315 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
|
316 // 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
|
317 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
318 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
|
319 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
|
320 // 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
|
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 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
|
323 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
|
324 // 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
|
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 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
|
327 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
|
328 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
|
329 break; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
330 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
331 default: |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
332 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
|
333 break; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
334 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
335 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
336 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
|
337 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
338 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
339 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
|
340 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
|
341 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
342 template octave_value_list find_nonzero_elem_idx (const Sparse<Complex>&, int, |
6002 | 343 octave_idx_type, int); |
2928 | 344 |
8955
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
345 octave_value_list |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
346 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
|
347 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
|
348 { |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
349 // 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
|
350 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
|
351 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
352 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
|
353 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
|
354 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
355 // 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
|
356 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
|
357 { |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
358 start_nc = 0; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
359 end_nc = nc; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
360 n_to_find = nc; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
361 count = nc; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
362 } |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
363 else if (direction > 0) |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
364 { |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
365 start_nc = 0; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
366 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
|
367 count = n_to_find; |
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 else |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
370 { |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
371 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
|
372 end_nc = nc; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
373 count = n_to_find; |
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 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
376 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
|
377 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
378 Matrix idx (count, 1); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
379 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
|
380 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
|
381 // Every value is 1. |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
382 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
|
383 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
384 if (count > 0) |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
385 { |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
386 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
|
387 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
|
388 { |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
389 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
|
390 { |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
391 OCTAVE_QUIT; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
392 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
|
393 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
|
394 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
|
395 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
|
396 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
|
397 } |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
398 } |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
399 else |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
400 { |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
401 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
|
402 { |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
403 OCTAVE_QUIT; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
404 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
|
405 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
|
406 // 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
|
407 // 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
|
408 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
|
409 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
|
410 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
|
411 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
|
412 } |
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 else if (scalar_arg) |
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 // 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
|
418 idx.resize (0, 0); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
419 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
|
420 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
|
421 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
|
422 } |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
423 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
424 switch (nargout) |
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 0: |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
427 case 1: |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
428 retval(0) = idx; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
429 break; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
430 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
431 case 5: |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
432 retval(4) = nc; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
433 // Fall through |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
434 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
435 case 4: |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
436 retval(3) = nc; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
437 // Fall through |
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 case 3: |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
440 retval(2) = val; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
441 // Fall through! |
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 case 2: |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
444 retval(1) = j_idx; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
445 retval(0) = i_idx; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
446 break; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
447 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
448 default: |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
449 panic_impossible (); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
450 break; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
451 } |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
452 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
453 return retval; |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
454 } |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
455 |
2928 | 456 DEFUN_DLD (find, args, nargout, |
3369 | 457 "-*- texinfo -*-\n\ |
458 @deftypefn {Loadable Function} {} find (@var{x})\n\ | |
6002 | 459 @deftypefnx {Loadable Function} {} find (@var{x}, @var{n})\n\ |
460 @deftypefnx {Loadable Function} {} find (@var{x}, @var{n}, @var{direction})\n\ | |
6524 | 461 Return a vector of indices of nonzero elements of a matrix, as a row if\n\ |
462 @var{x} is a row or as a column otherwise. To obtain a single index for\n\ | |
463 each matrix element, Octave pretends that the columns of a matrix form one\n\ | |
464 long vector (like Fortran arrays are stored). For example,\n\ | |
3369 | 465 \n\ |
466 @example\n\ | |
467 @group\n\ | |
468 find (eye (2))\n\ | |
469 @result{} [ 1; 4 ]\n\ | |
470 @end group\n\ | |
471 @end example\n\ | |
472 \n\ | |
473 If two outputs are requested, @code{find} returns the row and column\n\ | |
474 indices of nonzero elements of a matrix. For example,\n\ | |
475 \n\ | |
476 @example\n\ | |
477 @group\n\ | |
478 [i, j] = find (2 * eye (2))\n\ | |
479 @result{} i = [ 1; 2 ]\n\ | |
480 @result{} j = [ 1; 2 ]\n\ | |
481 @end group\n\ | |
482 @end example\n\ | |
483 \n\ | |
484 If three outputs are requested, @code{find} also returns a vector\n\ | |
485 containing the nonzero values. For example,\n\ | |
486 \n\ | |
487 @example\n\ | |
488 @group\n\ | |
489 [i, j, v] = find (3 * eye (2))\n\ | |
490 @result{} i = [ 1; 2 ]\n\ | |
491 @result{} j = [ 1; 2 ]\n\ | |
492 @result{} v = [ 3; 3 ]\n\ | |
493 @end group\n\ | |
494 @end example\n\ | |
6002 | 495 \n\ |
496 If two inputs are given, @var{n} indicates the number of elements to\n\ | |
497 find from the beginning of the matrix or vector.\n\ | |
498 \n\ | |
499 If three inputs are given, @var{direction} should be one of \"first\" or\n\ | |
500 \"last\" indicating that it should start counting found elements from the\n\ | |
501 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
|
502 \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 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
|
504 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
|
505 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
|
506 \n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
507 @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
|
508 @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
|
509 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
|
510 [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
|
511 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
|
512 @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
|
513 @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
|
514 @seealso{sparse}\n\ |
3369 | 515 @end deftypefn") |
2928 | 516 { |
517 octave_value_list retval; | |
518 | |
519 int nargin = args.length (); | |
520 | |
6002 | 521 if (nargin > 3 || nargin < 1) |
2928 | 522 { |
5823 | 523 print_usage (); |
2928 | 524 return retval; |
525 } | |
526 | |
6002 | 527 // Setup the default options. |
528 octave_idx_type n_to_find = -1; | |
529 if (nargin > 1) | |
530 { | |
531 n_to_find = args(1).int_value (); | |
532 if (error_state) | |
533 { | |
534 error ("find: expecting second argument to be an integer"); | |
535 return retval; | |
536 } | |
537 } | |
538 | |
539 // Direction to do the searching (1 == forward, -1 == reverse). | |
540 int direction = 1; | |
541 if (nargin > 2) | |
542 { | |
543 direction = 0; | |
544 | |
545 std::string s_arg = args(2).string_value (); | |
546 | |
547 if (! error_state) | |
548 { | |
549 if (s_arg == "first") | |
550 direction = 1; | |
551 else if (s_arg == "last") | |
552 direction = -1; | |
553 } | |
554 | |
555 if (direction == 0) | |
556 { | |
557 error ("find: expecting third argument to be \"first\" or \"last\""); | |
558 return retval; | |
559 } | |
560 } | |
561 | |
2928 | 562 octave_value arg = args(0); |
563 | |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
564 if (arg.is_sparse_type ()) |
2928 | 565 { |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
566 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
|
567 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
568 SparseMatrix v = arg.sparse_matrix_value (); |
2928 | 569 |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
570 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
|
571 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
|
572 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
|
573 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
574 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
|
575 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
576 SparseComplexMatrix v = arg.sparse_complex_matrix_value (); |
5107 | 577 |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
578 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
|
579 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
|
580 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
|
581 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
582 else |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
583 gripe_wrong_type_arg ("find", arg); |
5107 | 584 } |
8955
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
585 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
|
586 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
|
587 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
588 if (! error_state) |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
589 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
|
590 } |
2928 | 591 else |
592 { | |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
593 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
|
594 { |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
595 if (arg.is_real_type ()) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
596 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
597 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
|
598 |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
599 if (! error_state) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
600 retval = find_nonzero_elem_idx (nda, nargout, |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
601 n_to_find, direction); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
602 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
603 else if (arg.is_complex_type ()) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
604 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
605 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
|
606 |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
607 if (! error_state) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
608 retval = find_nonzero_elem_idx (cnda, nargout, |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
609 n_to_find, direction); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
610 } |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
611 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
612 else |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
613 { |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
614 if (arg.is_real_type ()) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
615 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
616 NDArray nda = arg.array_value (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
617 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
618 if (! error_state) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
619 retval = find_nonzero_elem_idx (nda, nargout, |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
620 n_to_find, direction); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
621 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
622 else if (arg.is_complex_type ()) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
623 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
624 ComplexNDArray cnda = arg.complex_array_value (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
625 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
626 if (! error_state) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
627 retval = find_nonzero_elem_idx (cnda, nargout, |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
628 n_to_find, direction); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
629 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
630 else if (arg.is_string ()) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
631 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
632 charNDArray cnda = arg.char_array_value (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
633 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
634 if (! error_state) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
635 retval = find_nonzero_elem_idx (cnda, nargout, |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
636 n_to_find, direction); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
637 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
638 else |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
639 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
640 gripe_wrong_type_arg ("find", arg); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7505
diff
changeset
|
641 } |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
642 } |
2928 | 643 } |
644 | |
645 return retval; | |
646 } | |
647 | |
648 /* | |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
649 %!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
|
650 %!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
|
651 %!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
|
652 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
653 %!test |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
654 %! [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
|
655 %! |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
656 %! 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
|
657 %! 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
|
658 %! 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
|
659 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
660 %!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
|
661 %!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
|
662 %!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
|
663 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
664 %!test |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
665 %! [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
|
666 %! |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
667 %! 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
|
668 %! 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
|
669 %! 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
|
670 |
8955
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 %! 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
|
673 %! P = eye (5) (:, pcol); |
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 |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
680 %!test |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
681 %! 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
|
682 %! P = eye (5) (prow, :); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
683 %! [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
|
684 %! [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
|
685 %! assert (i, ifull); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
686 %! assert (j, jfull); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
687 %! assert (all (v == 1)); |
6d3fcbf89267
Add an override to Octave's find() for permutation matrices.
Jason Riedy <jason@acm.org>
parents:
8920
diff
changeset
|
688 |
7814
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
689 %!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
|
690 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
691 */ |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
692 |
87865ed7405f
Second set of single precision test code and fix of resulting bugs
David Bateman <dbateman@free.fr>
parents:
7789
diff
changeset
|
693 /* |
2928 | 694 ;;; Local Variables: *** |
695 ;;; mode: C++ *** | |
696 ;;; End: *** | |
697 */ |