Mercurial > hg > octave-nkf
annotate liboctave/boolSparse.cc @ 7515:f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
author | David Bateman <dbateman@free.fr> |
---|---|
date | Fri, 22 Feb 2008 15:50:51 +0100 |
parents | a1dbe9d80eee |
children | 36594d5bbe13 |
rev | line source |
---|---|
5164 | 1 /* |
2 | |
7017 | 3 Copyright (C) 2004, 2005, 2006, 2007 David Bateman |
7016 | 4 Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004 Andy Adler |
5 | |
6 This file is part of Octave. | |
5164 | 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. | |
5164 | 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/>. | |
5164 | 21 |
22 */ | |
23 | |
24 #ifdef HAVE_CONFIG_H | |
25 #include <config.h> | |
26 #endif | |
27 | |
28 #include <iostream> | |
29 #include <vector> | |
30 | |
31 #include "config.h" | |
32 #include "quit.h" | |
33 #include "lo-ieee.h" | |
34 #include "lo-mappers.h" | |
35 | |
36 #include "boolSparse.h" | |
37 | |
38 // SparseBoolMatrix class. | |
39 | |
40 bool | |
41 SparseBoolMatrix::operator == (const SparseBoolMatrix& a) const | |
42 { | |
5275 | 43 octave_idx_type nr = rows (); |
44 octave_idx_type nc = cols (); | |
5604 | 45 octave_idx_type nz = nzmax (); |
5275 | 46 octave_idx_type nr_a = a.rows (); |
47 octave_idx_type nc_a = a.cols (); | |
5604 | 48 octave_idx_type nz_a = a.nzmax (); |
5164 | 49 |
50 if (nr != nr_a || nc != nc_a || nz != nz_a) | |
51 return false; | |
52 | |
5275 | 53 for (octave_idx_type i = 0; i < nc + 1; i++) |
5164 | 54 if (cidx(i) != a.cidx(i)) |
55 return false; | |
56 | |
5275 | 57 for (octave_idx_type i = 0; i < nz; i++) |
5164 | 58 if (data(i) != a.data(i) || ridx(i) != a.ridx(i)) |
59 return false; | |
60 | |
61 return true; | |
62 } | |
63 | |
64 bool | |
65 SparseBoolMatrix::operator != (const SparseBoolMatrix& a) const | |
66 { | |
67 return !(*this == a); | |
68 } | |
69 | |
70 SparseBoolMatrix& | |
5275 | 71 SparseBoolMatrix::insert (const SparseBoolMatrix& a, octave_idx_type r, octave_idx_type c) |
5164 | 72 { |
73 Sparse<bool>::insert (a, r, c); | |
74 return *this; | |
75 } | |
76 | |
6823 | 77 SparseBoolMatrix& |
78 SparseBoolMatrix::insert (const SparseBoolMatrix& a, const Array<octave_idx_type>& indx) | |
79 { | |
80 Sparse<bool>::insert (a, indx); | |
81 return *this; | |
82 } | |
83 | |
5164 | 84 SparseBoolMatrix |
5275 | 85 SparseBoolMatrix::concat (const SparseBoolMatrix& rb, const Array<octave_idx_type>& ra_idx) |
5164 | 86 { |
87 // Don't use numel to avoid all possiblity of an overflow | |
88 if (rb.rows () > 0 && rb.cols () > 0) | |
89 insert (rb, ra_idx(0), ra_idx(1)); | |
90 return *this; | |
91 } | |
92 | |
93 // unary operations | |
94 | |
95 SparseBoolMatrix | |
96 SparseBoolMatrix::operator ! (void) const | |
97 { | |
5275 | 98 octave_idx_type nr = rows (); |
99 octave_idx_type nc = cols (); | |
5604 | 100 octave_idx_type nz1 = nzmax (); |
5275 | 101 octave_idx_type nz2 = nr*nc - nz1; |
5164 | 102 |
103 SparseBoolMatrix r (nr, nc, nz2); | |
104 | |
5275 | 105 octave_idx_type ii = 0; |
106 octave_idx_type jj = 0; | |
6217 | 107 r.cidx (0) = 0; |
5275 | 108 for (octave_idx_type i = 0; i < nc; i++) |
5164 | 109 { |
5275 | 110 for (octave_idx_type j = 0; j < nr; j++) |
5164 | 111 { |
112 if (jj < cidx(i+1) && ridx(jj) == j) | |
113 jj++; | |
114 else | |
115 { | |
116 r.data(ii) = true; | |
117 r.ridx(ii++) = j; | |
118 } | |
119 } | |
6217 | 120 r.cidx (i+1) = ii; |
5164 | 121 } |
122 | |
123 return r; | |
124 } | |
125 | |
126 // other operations | |
127 | |
5775 | 128 // FIXME Do these really belong here? Maybe they should be |
5164 | 129 // in a base class? |
130 | |
131 SparseBoolMatrix | |
132 SparseBoolMatrix::all (int dim) const | |
133 { | |
134 SPARSE_ALL_OP (dim); | |
135 } | |
136 | |
137 SparseBoolMatrix | |
138 SparseBoolMatrix::any (int dim) const | |
139 { | |
140 SPARSE_ANY_OP (dim); | |
141 } | |
142 | |
7515
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
143 SparseBoolMatrix |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
144 SparseBoolMatrix::diag (octave_idx_type k) const |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
145 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
146 octave_idx_type nnr = rows (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
147 octave_idx_type nnc = cols (); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
148 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
149 if (k > 0) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
150 nnc -= k; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
151 else if (k < 0) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
152 nnr += k; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
153 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
154 SparseBoolMatrix d; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
155 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
156 if (nnr > 0 && nnc > 0) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
157 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
158 octave_idx_type ndiag = (nnr < nnc) ? nnr : nnc; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
159 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
160 // Count the number of non-zero elements |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
161 octave_idx_type nel = 0; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
162 if (k > 0) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
163 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
164 for (octave_idx_type i = 0; i < ndiag; i++) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
165 if (elem (i, i+k) != 0.) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
166 nel++; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
167 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
168 else if ( k < 0) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
169 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
170 for (octave_idx_type i = 0; i < ndiag; i++) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
171 if (elem (i-k, i) != 0.) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
172 nel++; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
173 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
174 else |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
175 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
176 for (octave_idx_type i = 0; i < ndiag; i++) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
177 if (elem (i, i) != 0.) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
178 nel++; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
179 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
180 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
181 d = SparseBoolMatrix (ndiag, 1, nel); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
182 d.xcidx (0) = 0; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
183 d.xcidx (1) = nel; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
184 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
185 octave_idx_type ii = 0; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
186 if (k > 0) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
187 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
188 for (octave_idx_type i = 0; i < ndiag; i++) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
189 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
190 bool tmp = elem (i, i+k); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
191 if (tmp != 0.) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
192 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
193 d.xdata (ii) = tmp; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
194 d.xridx (ii++) = i; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
195 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
196 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
197 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
198 else if ( k < 0) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
199 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
200 for (octave_idx_type i = 0; i < ndiag; i++) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
201 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
202 bool tmp = elem (i-k, i); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
203 if (tmp != 0.) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
204 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
205 d.xdata (ii) = tmp; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
206 d.xridx (ii++) = i; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
207 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
208 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
209 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
210 else |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
211 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
212 for (octave_idx_type i = 0; i < ndiag; i++) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
213 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
214 bool tmp = elem (i, i); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
215 if (tmp != 0.) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
216 { |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
217 d.xdata (ii) = tmp; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
218 d.xridx (ii++) = i; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
219 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
220 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
221 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
222 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
223 else |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
224 (*current_liboctave_error_handler) |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
225 ("diag: requested diagonal out of range"); |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
226 |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
227 return d; |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
228 } |
f3c00dc0912b
Eliminate the rest of the dispatched sparse functions
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
229 |
5164 | 230 boolMatrix |
231 SparseBoolMatrix::matrix_value (void) const | |
232 { | |
5275 | 233 octave_idx_type nr = rows (); |
234 octave_idx_type nc = cols (); | |
5164 | 235 |
236 boolMatrix retval (nr, nc, false); | |
5275 | 237 for (octave_idx_type j = 0; j < nc; j++) |
238 for (octave_idx_type i = cidx(j); i < cidx(j+1); i++) | |
5164 | 239 retval.elem (ridx(i), j) = data (i); |
240 | |
241 return retval; | |
242 } | |
243 | |
244 std::ostream& | |
245 operator << (std::ostream& os, const SparseBoolMatrix& a) | |
246 { | |
5275 | 247 octave_idx_type nc = a.cols (); |
5164 | 248 |
249 // add one to the printed indices to go from | |
250 // zero-based to one-based arrays | |
5275 | 251 for (octave_idx_type j = 0; j < nc; j++) |
5164 | 252 { |
253 OCTAVE_QUIT; | |
5275 | 254 for (octave_idx_type i = a.cidx(j); i < a.cidx(j+1); i++) |
5164 | 255 os << a.ridx(i) + 1 << " " << j + 1 << " " << a.data(i) << "\n"; |
256 } | |
257 | |
258 return os; | |
259 } | |
260 | |
261 std::istream& | |
262 operator >> (std::istream& is, SparseBoolMatrix& a) | |
263 { | |
5275 | 264 octave_idx_type nr = a.rows (); |
265 octave_idx_type nc = a.cols (); | |
5604 | 266 octave_idx_type nz = a.nzmax (); |
5164 | 267 |
268 if (nr < 1 || nc < 1) | |
269 is.clear (std::ios::badbit); | |
270 else | |
271 { | |
5275 | 272 octave_idx_type itmp, jtmp, jold = 0; |
5164 | 273 bool tmp; |
5275 | 274 octave_idx_type ii = 0; |
5164 | 275 |
276 a.cidx (0) = 0; | |
5275 | 277 for (octave_idx_type i = 0; i < nz; i++) |
5164 | 278 { |
279 is >> itmp; | |
280 itmp--; | |
281 is >> jtmp; | |
282 jtmp--; | |
283 is >> tmp; | |
284 if (is) | |
285 { | |
286 if (jold != jtmp) | |
287 { | |
5275 | 288 for (octave_idx_type j = jold; j < jtmp; j++) |
5164 | 289 a.cidx(j+1) = ii; |
290 | |
291 jold = jtmp; | |
292 } | |
293 a.data (ii) = tmp; | |
294 a.ridx (ii++) = itmp; | |
295 } | |
296 else | |
297 goto done; | |
298 } | |
299 | |
5275 | 300 for (octave_idx_type j = jold; j < nc; j++) |
5164 | 301 a.cidx(j+1) = ii; |
302 } | |
303 | |
304 done: | |
305 | |
306 return is; | |
307 } | |
308 | |
309 SparseBoolMatrix | |
310 SparseBoolMatrix::squeeze (void) const | |
311 { | |
312 return Sparse<bool>::squeeze (); | |
313 } | |
314 | |
315 SparseBoolMatrix | |
316 SparseBoolMatrix::index (idx_vector& i, int resize_ok) const | |
317 { | |
318 return Sparse<bool>::index (i, resize_ok); | |
319 } | |
320 | |
321 SparseBoolMatrix | |
322 SparseBoolMatrix::index (idx_vector& i, idx_vector& j, int resize_ok) const | |
323 { | |
324 return Sparse<bool>::index (i, j, resize_ok); | |
325 } | |
326 | |
327 SparseBoolMatrix | |
328 SparseBoolMatrix::index (Array<idx_vector>& ra_idx, int resize_ok) const | |
329 { | |
330 return Sparse<bool>::index (ra_idx, resize_ok); | |
331 } | |
332 | |
333 SparseBoolMatrix | |
334 SparseBoolMatrix::reshape (const dim_vector& new_dims) const | |
335 { | |
336 return Sparse<bool>::reshape (new_dims); | |
337 } | |
338 | |
339 SparseBoolMatrix | |
5275 | 340 SparseBoolMatrix::permute (const Array<octave_idx_type>& vec, bool inv) const |
5164 | 341 { |
342 return Sparse<bool>::permute (vec, inv); | |
343 } | |
344 | |
345 SparseBoolMatrix | |
5275 | 346 SparseBoolMatrix::ipermute (const Array<octave_idx_type>& vec) const |
5164 | 347 { |
348 return Sparse<bool>::ipermute (vec); | |
349 } | |
350 | |
351 SPARSE_SMS_EQNE_OPS (SparseBoolMatrix, false, , bool, false, ) | |
352 SPARSE_SMS_BOOL_OPS (SparseBoolMatrix, bool, false) | |
353 | |
354 SPARSE_SSM_EQNE_OPS (bool, false, , SparseBoolMatrix, false, ) | |
355 SPARSE_SSM_BOOL_OPS (bool, SparseBoolMatrix, false) | |
356 | |
357 SPARSE_SMSM_EQNE_OPS (SparseBoolMatrix, false, , SparseBoolMatrix, false, ) | |
358 SPARSE_SMSM_BOOL_OPS (SparseBoolMatrix, SparseBoolMatrix, false) | |
359 | |
360 | |
361 /* | |
362 ;;; Local Variables: *** | |
363 ;;; mode: C++ *** | |
364 ;;; End: *** | |
365 */ |