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 <cstdlib> |
|
29 #include <string> |
|
30 |
|
31 #include "variables.h" |
|
32 #include "utils.h" |
|
33 #include "pager.h" |
|
34 #include "defun-dld.h" |
|
35 #include "gripes.h" |
|
36 #include "quit.h" |
|
37 |
|
38 #include "ov-re-sparse.h" |
|
39 #include "ov-cx-sparse.h" |
|
40 #include "ov-bool-sparse.h" |
|
41 |
|
42 static bool |
|
43 is_sparse (const octave_value& arg) |
|
44 { |
5631
|
45 return (arg.is_sparse_type ()); |
5164
|
46 } |
|
47 |
|
48 DEFUN_DLD (issparse, args, , |
|
49 "-*- texinfo -*-\n\ |
|
50 @deftypefn {Loadable Function} {} issparse (@var{expr})\n\ |
|
51 Return 1 if the value of the expression @var{expr} is a sparse matrix.\n\ |
|
52 @end deftypefn") |
|
53 { |
|
54 if (args.length() != 1) |
|
55 { |
5823
|
56 print_usage (); |
5164
|
57 return octave_value (); |
|
58 } |
|
59 else |
|
60 return octave_value (is_sparse (args(0))); |
|
61 } |
|
62 |
|
63 DEFUN_DLD (sparse, args, , |
|
64 "-*- texinfo -*-\n\ |
6556
|
65 @deftypefn {Loadable Function} {@var{s} =} sparse (@var{a})\n\ |
|
66 Create a sparse matrix from the full matrix @var{a}.\n\ |
|
67 is forced back to a full matrix is resulting matrix is sparse\n\ |
5164
|
68 \n\ |
6556
|
69 @deftypefnx {Loadable Function} {@var{s} =} sparse (@var{i}, @var{j}, @var{sv}, @var{m}, @var{n}, @var{nzmax})\n\ |
|
70 Create a sparse matrix given integer index vectors @var{i} and @var{j},\n\ |
|
71 a 1-by-@code{nnz} vector of real of complex values @var{sv}, overall\n\ |
|
72 dimensions @var{m} and @var{n} of the sparse matrix. The argument\n\ |
7001
|
73 @code{nzmax} is ignored but accepted for compatibility with @sc{Matlab}.\n\ |
5164
|
74 \n\ |
6556
|
75 @strong{Note}: if multiple values are specified with the same\n\ |
|
76 @var{i}, @var{j} indices, the corresponding values in @var{s} will\n\ |
|
77 be added.\n\ |
|
78 \n\ |
6557
|
79 The following are all equivalent:\n\ |
5164
|
80 \n\ |
6556
|
81 @example\n\ |
|
82 @group\n\ |
|
83 s = sparse (i, j, s, m, n)\n\ |
|
84 s = sparse (i, j, s, m, n, \"summation\")\n\ |
|
85 s = sparse (i, j, s, m, n, \"sum\")\n\ |
|
86 @end group\n\ |
|
87 @end example\n\ |
5164
|
88 \n\ |
6556
|
89 @deftypefnx {Loadable Function} {@var{s} =} sparse (@var{i}, @var{j}, @var{s}, @var{m}, @var{n}, \"unique\")\n\ |
|
90 Same as above, except that if more than two values are specified for the\n\ |
|
91 same @var{i}, @var{j} indices, the last specified value will be used.\n\ |
5164
|
92 \n\ |
6556
|
93 @deftypefnx {Loadable Function} {@var{s} =} sparse (@var{i}, @var{j}, @var{sv})\n\ |
|
94 Uses @code{@var{m} = max (@var{i})}, @code{@var{n} = max (@var{j})}\n\ |
5164
|
95 \n\ |
6556
|
96 @deftypefnx {Loadable Function} {@var{s} =} sparse (@var{m}, @var{n})\n\ |
|
97 Equivalent to @code{sparse ([], [], [], @var{m}, @var{n}, 0)}\n\ |
|
98 \n\ |
|
99 If any of @var{sv}, @var{i} or @var{j} are scalars, they are expanded\n\ |
|
100 to have a common size.\n\ |
5164
|
101 @seealso{full}\n\ |
|
102 @end deftypefn") |
|
103 { |
|
104 octave_value retval; |
|
105 |
|
106 // WARNING: This function should always use constructions like |
|
107 // retval = new octave_sparse_matrix (sm); |
|
108 // To avoid calling the maybe_mutate function. This is the only |
7287
|
109 // function that should not call maybe_mutate |
5164
|
110 |
|
111 int nargin= args.length(); |
|
112 if (nargin < 1 || (nargin == 4 && !args(3).is_string ()) || nargin > 6) |
|
113 { |
5823
|
114 print_usage (); |
5164
|
115 return retval; |
|
116 } |
|
117 |
|
118 bool use_complex = false; |
|
119 bool use_bool = false; |
|
120 if (nargin > 2) |
|
121 { |
|
122 use_complex= args(2).is_complex_type(); |
|
123 use_bool = args(2).is_bool_type (); |
|
124 } |
|
125 else |
|
126 { |
|
127 use_complex= args(0).is_complex_type(); |
|
128 use_bool = args(0).is_bool_type (); |
|
129 } |
|
130 |
7287
|
131 if (nargin == 1) |
5164
|
132 { |
|
133 octave_value arg = args (0); |
|
134 |
|
135 if (is_sparse (arg)) |
|
136 { |
|
137 if (use_complex) |
|
138 { |
5760
|
139 SparseComplexMatrix sm = arg.sparse_complex_matrix_value (); |
5164
|
140 retval = new octave_sparse_complex_matrix (sm); |
|
141 } |
|
142 else if (use_bool) |
|
143 { |
5760
|
144 SparseBoolMatrix sm = arg.sparse_bool_matrix_value (); |
5164
|
145 retval = new octave_sparse_bool_matrix (sm); |
|
146 } |
|
147 else |
|
148 { |
5760
|
149 SparseMatrix sm = arg.sparse_matrix_value (); |
5164
|
150 retval = new octave_sparse_matrix (sm); |
|
151 } |
|
152 } |
|
153 else |
|
154 { |
|
155 if (use_complex) |
|
156 { |
|
157 SparseComplexMatrix sm (args (0).complex_matrix_value ()); |
|
158 if (error_state) |
|
159 return retval; |
|
160 retval = new octave_sparse_complex_matrix (sm); |
|
161 } |
|
162 else if (use_bool) |
|
163 { |
|
164 SparseBoolMatrix sm (args (0).bool_matrix_value ()); |
|
165 if (error_state) |
|
166 return retval; |
|
167 retval = new octave_sparse_bool_matrix (sm); |
|
168 } |
|
169 else |
|
170 { |
|
171 SparseMatrix sm (args (0).matrix_value ()); |
|
172 if (error_state) |
|
173 return retval; |
|
174 retval = new octave_sparse_matrix (sm); |
|
175 } |
|
176 } |
|
177 } |
|
178 else |
|
179 { |
5275
|
180 octave_idx_type m = 1, n = 1; |
5164
|
181 if (nargin == 2) |
|
182 { |
7308
|
183 if (args(0).numel () == 1 && args(1).numel () == 1) |
|
184 { |
|
185 m = args(0).int_value(); |
|
186 n = args(1).int_value(); |
|
187 if (error_state) return retval; |
5164
|
188 |
7308
|
189 if (use_complex) |
|
190 retval = new octave_sparse_complex_matrix |
|
191 (SparseComplexMatrix (m, n)); |
|
192 else if (use_bool) |
|
193 retval = new octave_sparse_bool_matrix |
|
194 (SparseBoolMatrix (m, n)); |
|
195 else |
|
196 retval = new octave_sparse_matrix |
|
197 (SparseMatrix (m, n)); |
|
198 } |
5164
|
199 else |
7308
|
200 error ("sparse: expecting scalar values"); |
5164
|
201 } |
|
202 else |
|
203 { |
|
204 if (args(0).is_empty () || args (1).is_empty () |
|
205 || args(2).is_empty ()) |
|
206 { |
|
207 if (nargin > 4) |
|
208 { |
|
209 m = args(3).int_value(); |
|
210 n = args(4).int_value(); |
|
211 } |
|
212 |
|
213 if (use_bool) |
|
214 retval = new octave_sparse_bool_matrix |
|
215 (SparseBoolMatrix (m, n)); |
|
216 else |
|
217 retval = new octave_sparse_matrix (SparseMatrix (m, n)); |
|
218 } |
|
219 else |
|
220 { |
|
221 // |
|
222 // I use this clumsy construction so that we can use |
|
223 // any orientation of args |
|
224 ColumnVector ridxA = ColumnVector (args(0).vector_value |
|
225 (false, true)); |
|
226 ColumnVector cidxA = ColumnVector (args(1).vector_value |
|
227 (false, true)); |
|
228 ColumnVector coefA; |
|
229 boolNDArray coefAB; |
|
230 ComplexColumnVector coefAC; |
|
231 bool assemble_do_sum = true; // this is the default in matlab6 |
|
232 |
|
233 if (use_complex) |
|
234 { |
|
235 if (args(2).is_empty ()) |
|
236 coefAC = ComplexColumnVector (0); |
|
237 else |
|
238 coefAC = ComplexColumnVector |
|
239 (args(2).complex_vector_value (false, true)); |
|
240 } |
|
241 else if (use_bool) |
|
242 { |
|
243 if (args(2).is_empty ()) |
|
244 coefAB = boolNDArray (dim_vector (1, 0)); |
|
245 else |
|
246 coefAB = args(2).bool_array_value (); |
|
247 dim_vector AB_dims = coefAB.dims (); |
|
248 if (AB_dims.length() > 2 || (AB_dims(0) != 1 && |
|
249 AB_dims(1) != 1)) |
|
250 error ("sparse: vector arguments required"); |
|
251 } |
|
252 else |
|
253 if (args(2).is_empty ()) |
|
254 coefA = ColumnVector (0); |
|
255 else |
|
256 coefA = ColumnVector (args(2).vector_value (false, true)); |
|
257 |
|
258 if (error_state) |
|
259 return retval; |
|
260 |
|
261 // Confirm that i,j,s all have the same number of elements |
5275
|
262 octave_idx_type ns; |
5164
|
263 if (use_complex) |
|
264 ns = coefAC.length(); |
|
265 else if (use_bool) |
|
266 ns = coefAB.length(); |
|
267 else |
|
268 ns = coefA.length(); |
|
269 |
5275
|
270 octave_idx_type ni = ridxA.length(); |
|
271 octave_idx_type nj = cidxA.length(); |
|
272 octave_idx_type nnz = (ni > nj ? ni : nj); |
5164
|
273 if ((ns != 1 && ns != nnz) || |
|
274 (ni != 1 && ni != nnz) || |
|
275 (nj != 1 && nj != nnz)) |
|
276 { |
|
277 error ("sparse i, j and s must have the same length"); |
|
278 return retval; |
|
279 } |
|
280 |
|
281 if (nargin == 3 || nargin == 4) |
|
282 { |
5275
|
283 m = static_cast<octave_idx_type> (ridxA.max()); |
|
284 n = static_cast<octave_idx_type> (cidxA.max()); |
5164
|
285 |
|
286 // if args(3) is not string, then ignore the value |
|
287 // otherwise check for summation or unique |
|
288 if (nargin == 4 && args(3).is_string()) |
|
289 { |
|
290 std::string vv= args(3).string_value(); |
|
291 if (error_state) return retval; |
|
292 |
|
293 if ( vv == "summation" || |
|
294 vv == "sum" ) |
|
295 assemble_do_sum = true; |
|
296 else |
|
297 if ( vv == "unique" ) |
|
298 assemble_do_sum = false; |
|
299 else { |
|
300 error("sparse repeat flag must be 'sum' or 'unique'"); |
|
301 return retval; |
|
302 } |
|
303 } |
|
304 } |
|
305 else |
|
306 { |
|
307 m = args(3).int_value(); |
|
308 n = args(4).int_value(); |
|
309 if (error_state) |
|
310 return retval; |
|
311 |
|
312 // if args(5) is not string, then ignore the value |
|
313 // otherwise check for summation or unique |
|
314 if (nargin >= 6 && args(5).is_string()) |
|
315 { |
|
316 std::string vv= args(5).string_value(); |
|
317 if (error_state) return retval; |
|
318 |
|
319 if ( vv == "summation" || |
|
320 vv == "sum" ) |
|
321 assemble_do_sum = true; |
|
322 else |
|
323 if ( vv == "unique" ) |
|
324 assemble_do_sum = false; |
|
325 else { |
|
326 error("sparse repeat flag must be 'sum' or 'unique'"); |
|
327 return retval; |
|
328 } |
|
329 } |
|
330 |
|
331 } |
|
332 |
|
333 // Convert indexing to zero-indexing used internally |
|
334 ridxA -= 1.; |
|
335 cidxA -= 1.; |
|
336 |
|
337 if (use_complex) |
|
338 retval = new octave_sparse_complex_matrix |
|
339 (SparseComplexMatrix (coefAC, ridxA, cidxA, m, n, |
|
340 assemble_do_sum)); |
|
341 else if (use_bool) |
|
342 retval = new octave_sparse_bool_matrix |
|
343 (SparseBoolMatrix (coefAB, ridxA, cidxA, m, n, |
|
344 assemble_do_sum)); |
|
345 else |
|
346 retval = new octave_sparse_matrix |
|
347 (SparseMatrix (coefA, ridxA, cidxA, m, n, |
|
348 assemble_do_sum)); |
|
349 } |
|
350 } |
|
351 } |
|
352 |
|
353 return retval; |
|
354 } |
|
355 |
|
356 DEFUN_DLD (full, args, , |
|
357 "-*- texinfo -*-\n\ |
|
358 @deftypefn {Loadable Function} {@var{FM} =} full (@var{SM})\n\ |
|
359 returns a full storage matrix from a sparse one\n\ |
|
360 @seealso{sparse}\n\ |
|
361 @end deftypefn") |
|
362 { |
|
363 octave_value retval; |
|
364 |
5760
|
365 if (args.length() < 1) |
|
366 { |
5823
|
367 print_usage (); |
5760
|
368 return retval; |
|
369 } |
5164
|
370 |
5631
|
371 if (args(0).is_sparse_type ()) |
5164
|
372 { |
|
373 if (args(0).type_name () == "sparse matrix") |
|
374 retval = args(0).matrix_value (); |
|
375 else if (args(0).type_name () == "sparse complex matrix") |
|
376 retval = args(0).complex_matrix_value (); |
|
377 else if (args(0).type_name () == "sparse bool matrix") |
|
378 retval = args(0).bool_matrix_value (); |
|
379 } |
|
380 else if (args(0).is_real_type()) |
|
381 retval = args(0).matrix_value(); |
|
382 else if (args(0).is_complex_type()) |
|
383 retval = args(0).complex_matrix_value(); |
|
384 else |
|
385 gripe_wrong_type_arg ("full", args(0)); |
|
386 |
|
387 return retval; |
|
388 } |
|
389 |
|
390 #define SPARSE_DIM_ARG_BODY(NAME, FUNC) \ |
|
391 int nargin = args.length(); \ |
|
392 octave_value retval; \ |
|
393 if ((nargin != 1 ) && (nargin != 2)) \ |
5823
|
394 print_usage (); \ |
5164
|
395 else { \ |
|
396 int dim = (nargin == 1 ? -1 : args(1).int_value(true) - 1); \ |
|
397 if (error_state) return retval; \ |
|
398 if (dim < -1 || dim > 1) { \ |
|
399 error (#NAME ": invalid dimension argument = %d", dim + 1); \ |
|
400 return retval; \ |
|
401 } \ |
|
402 if (args(0).type_id () == \ |
|
403 octave_sparse_matrix::static_type_id () || args(0).type_id () == \ |
|
404 octave_sparse_bool_matrix::static_type_id ()) { \ |
|
405 retval = args(0).sparse_matrix_value () .FUNC (dim); \ |
|
406 } else if (args(0).type_id () == \ |
|
407 octave_sparse_complex_matrix::static_type_id ()) { \ |
|
408 retval = args(0).sparse_complex_matrix_value () .FUNC (dim); \ |
|
409 } else \ |
5823
|
410 print_usage (); \ |
5164
|
411 } \ |
|
412 return retval |
|
413 |
|
414 // PKG_ADD: dispatch ("prod", "spprod", "sparse matrix"); |
|
415 // PKG_ADD: dispatch ("prod", "spprod", "sparse complex matrix"); |
|
416 // PKG_ADD: dispatch ("prod", "spprod", "sparse bool matrix"); |
|
417 DEFUN_DLD (spprod, args, , |
|
418 "-*- texinfo -*-\n\ |
|
419 @deftypefn {Loadable Function} {@var{y} =} spprod (@var{x},@var{dim})\n\ |
|
420 Product of elements along dimension @var{dim}. If @var{dim} is omitted,\n\ |
|
421 it defaults to 1 (column-wise products).\n\ |
5642
|
422 @seealso{spsum, spsumsq}\n\ |
|
423 @end deftypefn") |
5164
|
424 { |
|
425 SPARSE_DIM_ARG_BODY (spprod, prod); |
|
426 } |
|
427 |
|
428 // PKG_ADD: dispatch ("cumprod", "spcumprod", "sparse matrix"); |
|
429 // PKG_ADD: dispatch ("cumprod", "spcumprod", "sparse complex matrix"); |
|
430 // PKG_ADD: dispatch ("cumprod", "spcumprod", "sparse bool matrix"); |
|
431 DEFUN_DLD (spcumprod, args, , |
|
432 "-*- texinfo -*-\n\ |
|
433 @deftypefn {Loadable Function} {@var{y} =} spcumprod (@var{x},@var{dim})\n\ |
|
434 Cumulative product of elements along dimension @var{dim}. If @var{dim}\n\ |
|
435 is omitted, it defaults to 1 (column-wise cumulative products).\n\ |
5642
|
436 @seealso{spcumsum}\n\ |
|
437 @end deftypefn") |
5164
|
438 { |
|
439 SPARSE_DIM_ARG_BODY (spcumprod, cumprod); |
|
440 } |
|
441 |
|
442 // PKG_ADD: dispatch ("sum", "spsum", "sparse matrix"); |
|
443 // PKG_ADD: dispatch ("sum", "spsum", "sparse complex matrix"); |
|
444 // PKG_ADD: dispatch ("sum", "spsum", "sparse bool matrix"); |
|
445 DEFUN_DLD (spsum, args, , |
|
446 "-*- texinfo -*-\n\ |
|
447 @deftypefn {Loadable Function} {@var{y} =} spsum (@var{x},@var{dim})\n\ |
|
448 Sum of elements along dimension @var{dim}. If @var{dim} is omitted, it\n\ |
|
449 defaults to 1 (column-wise sum).\n\ |
5642
|
450 @seealso{spprod, spsumsq}\n\ |
|
451 @end deftypefn") |
5164
|
452 { |
|
453 SPARSE_DIM_ARG_BODY (spsum, sum); |
|
454 } |
|
455 |
|
456 // PKG_ADD: dispatch ("cumsum", "spcumsum", "sparse matrix"); |
|
457 // PKG_ADD: dispatch ("cumsum", "spcumsum", "sparse complex matrix"); |
|
458 // PKG_ADD: dispatch ("cumsum", "spcumsum", "sparse bool matrix"); |
|
459 DEFUN_DLD (spcumsum, args, , |
|
460 "-*- texinfo -*-\n\ |
|
461 @deftypefn {Loadable Function} {@var{y} =} spcumsum (@var{x},@var{dim})\n\ |
|
462 Cumulative sum of elements along dimension @var{dim}. If @var{dim}\n\ |
|
463 is omitted, it defaults to 1 (column-wise cumulative sums).\n\ |
5642
|
464 @seealso{spcumprod}\n\ |
|
465 @end deftypefn") |
5164
|
466 { |
|
467 SPARSE_DIM_ARG_BODY (spcumsum, cumsum); |
|
468 } |
|
469 |
|
470 // PKG_ADD: dispatch ("sumsq", "spsumsq", "sparse matrix"); |
|
471 // PKG_ADD: dispatch ("sumsq", "spsumsq", "sparse complex matrix"); |
|
472 // PKG_ADD: dispatch ("sumsq", "spsumsq", "sparse bool matrix"); |
|
473 DEFUN_DLD (spsumsq, args, , |
|
474 "-*- texinfo -*-\n\ |
|
475 @deftypefn {Loadable Function} {@var{y} =} spsumsq (@var{x},@var{dim})\n\ |
|
476 Sum of squares of elements along dimension @var{dim}. If @var{dim}\n\ |
|
477 is omitted, it defaults to 1 (column-wise sum of squares).\n\ |
|
478 This function is equivalent to computing\n\ |
|
479 @example\n\ |
|
480 spsum (x .* spconj (x), dim)\n\ |
|
481 @end example\n\ |
|
482 but it uses less memory and avoids calling @code{spconj} if @var{x} is\n\ |
|
483 real.\n\ |
5642
|
484 @seealso{spprod, spsum}\n\ |
|
485 @end deftypefn") |
5164
|
486 { |
|
487 SPARSE_DIM_ARG_BODY (spsumsq, sumsq); |
|
488 } |
|
489 |
|
490 #define MINMAX_BODY(FCN) \ |
|
491 \ |
|
492 octave_value_list retval; \ |
|
493 \ |
|
494 int nargin = args.length (); \ |
|
495 \ |
|
496 if (nargin < 1 || nargin > 3 || nargout > 2) \ |
|
497 { \ |
5823
|
498 print_usage (); \ |
5164
|
499 return retval; \ |
|
500 } \ |
|
501 \ |
|
502 octave_value arg1; \ |
|
503 octave_value arg2; \ |
|
504 octave_value arg3; \ |
|
505 \ |
|
506 switch (nargin) \ |
|
507 { \ |
|
508 case 3: \ |
|
509 arg3 = args(2); \ |
|
510 \ |
|
511 case 2: \ |
|
512 arg2 = args(1); \ |
|
513 \ |
|
514 case 1: \ |
|
515 arg1 = args(0); \ |
|
516 break; \ |
|
517 \ |
|
518 default: \ |
|
519 panic_impossible (); \ |
|
520 break; \ |
|
521 } \ |
|
522 \ |
|
523 int dim; \ |
5759
|
524 dim_vector dv = arg1.dims (); \ |
5164
|
525 if (error_state) \ |
|
526 { \ |
|
527 gripe_wrong_type_arg (#FCN, arg1); \ |
|
528 return retval; \ |
|
529 } \ |
|
530 \ |
|
531 if (nargin == 3) \ |
|
532 { \ |
|
533 dim = arg3.nint_value () - 1; \ |
|
534 if (dim < 0 || dim >= dv.length ()) \ |
|
535 { \ |
|
536 error ("%s: invalid dimension", #FCN); \ |
|
537 return retval; \ |
|
538 } \ |
|
539 } \ |
|
540 else \ |
|
541 { \ |
|
542 dim = 0; \ |
|
543 while ((dim < dv.length ()) && (dv (dim) <= 1)) \ |
|
544 dim++; \ |
|
545 if (dim == dv.length ()) \ |
|
546 dim = 0; \ |
|
547 } \ |
|
548 \ |
|
549 bool single_arg = (nargin == 1) || arg2.is_empty(); \ |
|
550 \ |
|
551 if (single_arg && (nargout == 1 || nargout == 0)) \ |
|
552 { \ |
|
553 if (arg1.type_id () == octave_sparse_matrix::static_type_id ()) \ |
|
554 retval(0) = arg1.sparse_matrix_value () .FCN (dim); \ |
|
555 else if (arg1.type_id () == \ |
|
556 octave_sparse_complex_matrix::static_type_id ()) \ |
|
557 retval(0) = arg1.sparse_complex_matrix_value () .FCN (dim); \ |
|
558 else \ |
|
559 gripe_wrong_type_arg (#FCN, arg1); \ |
|
560 } \ |
|
561 else if (single_arg && nargout == 2) \ |
|
562 { \ |
5275
|
563 Array2<octave_idx_type> index; \ |
5164
|
564 \ |
|
565 if (arg1.type_id () == octave_sparse_matrix::static_type_id ()) \ |
|
566 retval(0) = arg1.sparse_matrix_value () .FCN (index, dim); \ |
|
567 else if (arg1.type_id () == \ |
|
568 octave_sparse_complex_matrix::static_type_id ()) \ |
|
569 retval(0) = arg1.sparse_complex_matrix_value () .FCN (index, dim); \ |
|
570 else \ |
|
571 gripe_wrong_type_arg (#FCN, arg1); \ |
|
572 \ |
5275
|
573 octave_idx_type len = index.numel (); \ |
5164
|
574 \ |
|
575 if (len > 0) \ |
|
576 { \ |
|
577 double nan_val = lo_ieee_nan_value (); \ |
|
578 \ |
|
579 NDArray idx (index.dims ()); \ |
|
580 \ |
5275
|
581 for (octave_idx_type i = 0; i < len; i++) \ |
5164
|
582 { \ |
|
583 OCTAVE_QUIT; \ |
5275
|
584 octave_idx_type tmp = index.elem (i) + 1; \ |
5164
|
585 idx.elem (i) = (tmp <= 0) \ |
|
586 ? nan_val : static_cast<double> (tmp); \ |
|
587 } \ |
|
588 \ |
|
589 retval(1) = idx; \ |
|
590 } \ |
|
591 else \ |
|
592 retval(1) = NDArray (); \ |
|
593 } \ |
|
594 else \ |
|
595 { \ |
|
596 int arg1_is_scalar = arg1.is_scalar_type (); \ |
|
597 int arg2_is_scalar = arg2.is_scalar_type (); \ |
|
598 \ |
|
599 int arg1_is_complex = arg1.is_complex_type (); \ |
|
600 int arg2_is_complex = arg2.is_complex_type (); \ |
|
601 \ |
|
602 if (arg1_is_scalar) \ |
|
603 { \ |
|
604 if (arg1_is_complex || arg2_is_complex) \ |
|
605 { \ |
|
606 Complex c1 = arg1.complex_value (); \ |
|
607 \ |
|
608 SparseComplexMatrix m2 = arg2.sparse_complex_matrix_value (); \ |
|
609 \ |
|
610 if (! error_state) \ |
|
611 { \ |
|
612 SparseComplexMatrix result = FCN (c1, m2); \ |
|
613 if (! error_state) \ |
|
614 retval(0) = result; \ |
|
615 } \ |
|
616 } \ |
|
617 else \ |
|
618 { \ |
|
619 double d1 = arg1.double_value (); \ |
|
620 SparseMatrix m2 = arg2.sparse_matrix_value (); \ |
|
621 \ |
|
622 if (! error_state) \ |
|
623 { \ |
|
624 SparseMatrix result = FCN (d1, m2); \ |
|
625 if (! error_state) \ |
|
626 retval(0) = result; \ |
|
627 } \ |
|
628 } \ |
|
629 } \ |
|
630 else if (arg2_is_scalar) \ |
|
631 { \ |
|
632 if (arg1_is_complex || arg2_is_complex) \ |
|
633 { \ |
|
634 SparseComplexMatrix m1 = arg1.sparse_complex_matrix_value (); \ |
|
635 \ |
|
636 if (! error_state) \ |
|
637 { \ |
|
638 Complex c2 = arg2.complex_value (); \ |
|
639 SparseComplexMatrix result = FCN (m1, c2); \ |
|
640 if (! error_state) \ |
|
641 retval(0) = result; \ |
|
642 } \ |
|
643 } \ |
|
644 else \ |
|
645 { \ |
|
646 SparseMatrix m1 = arg1.sparse_matrix_value (); \ |
|
647 \ |
|
648 if (! error_state) \ |
|
649 { \ |
|
650 double d2 = arg2.double_value (); \ |
|
651 SparseMatrix result = FCN (m1, d2); \ |
|
652 if (! error_state) \ |
|
653 retval(0) = result; \ |
|
654 } \ |
|
655 } \ |
|
656 } \ |
|
657 else \ |
|
658 { \ |
|
659 if (arg1_is_complex || arg2_is_complex) \ |
|
660 { \ |
|
661 SparseComplexMatrix m1 = arg1.sparse_complex_matrix_value (); \ |
|
662 \ |
|
663 if (! error_state) \ |
|
664 { \ |
|
665 SparseComplexMatrix m2 = arg2.sparse_complex_matrix_value (); \ |
|
666 \ |
|
667 if (! error_state) \ |
|
668 { \ |
|
669 SparseComplexMatrix result = FCN (m1, m2); \ |
|
670 if (! error_state) \ |
|
671 retval(0) = result; \ |
|
672 } \ |
|
673 } \ |
|
674 } \ |
|
675 else \ |
|
676 { \ |
|
677 SparseMatrix m1 = arg1.sparse_matrix_value (); \ |
|
678 \ |
|
679 if (! error_state) \ |
|
680 { \ |
|
681 SparseMatrix m2 = arg2.sparse_matrix_value (); \ |
|
682 \ |
|
683 if (! error_state) \ |
|
684 { \ |
|
685 SparseMatrix result = FCN (m1, m2); \ |
|
686 if (! error_state) \ |
|
687 retval(0) = result; \ |
|
688 } \ |
|
689 } \ |
|
690 } \ |
|
691 } \ |
|
692 } \ |
|
693 \ |
|
694 return retval |
|
695 |
|
696 // PKG_ADD: dispatch ("min", "spmin", "sparse matrix"); |
|
697 // PKG_ADD: dispatch ("min", "spmin", "sparse complex matrix"); |
|
698 // PKG_ADD: dispatch ("min", "spmin", "sparse bool matrix"); |
|
699 DEFUN_DLD (spmin, args, nargout, |
|
700 "-*- texinfo -*-\n\ |
|
701 @deftypefn {Mapping Function} {} spmin (@var{x}, @var{y}, @var{dim})\n\ |
|
702 @deftypefnx {Mapping Function} {[@var{w}, @var{iw}] =} spmin (@var{x})\n\ |
|
703 @cindex Utility Functions\n\ |
|
704 For a vector argument, return the minimum value. For a matrix\n\ |
|
705 argument, return the minimum value from each column, as a row\n\ |
|
706 vector, or over the dimension @var{dim} if defined. For two matrices\n\ |
|
707 (or a matrix and scalar), return the pair-wise minimum.\n\ |
|
708 Thus,\n\ |
|
709 \n\ |
|
710 @example\n\ |
|
711 min (min (@var{x}))\n\ |
|
712 @end example\n\ |
|
713 \n\ |
|
714 @noindent\n\ |
|
715 returns the smallest element of @var{x}, and\n\ |
|
716 \n\ |
|
717 @example\n\ |
|
718 @group\n\ |
|
719 min (2:5, pi)\n\ |
|
720 @result{} 2.0000 3.0000 3.1416 3.1416\n\ |
|
721 @end group\n\ |
|
722 @end example\n\ |
|
723 @noindent\n\ |
|
724 compares each element of the range @code{2:5} with @code{pi}, and\n\ |
|
725 returns a row vector of the minimum values.\n\ |
|
726 \n\ |
|
727 For complex arguments, the magnitude of the elements are used for\n\ |
|
728 comparison.\n\ |
|
729 \n\ |
|
730 If called with one input and two output arguments,\n\ |
|
731 @code{min} also returns the first index of the\n\ |
|
732 minimum value(s). Thus,\n\ |
|
733 \n\ |
|
734 @example\n\ |
|
735 @group\n\ |
|
736 [x, ix] = min ([1, 3, 0, 2, 5])\n\ |
|
737 @result{} x = 0\n\ |
|
738 ix = 3\n\ |
|
739 @end group\n\ |
|
740 @end example\n\ |
|
741 @end deftypefn") |
|
742 { |
|
743 MINMAX_BODY (min); |
|
744 } |
|
745 |
|
746 // PKG_ADD: dispatch ("max", "spmax", "sparse matrix"); |
|
747 // PKG_ADD: dispatch ("max", "spmax", "sparse complex matrix"); |
|
748 // PKG_ADD: dispatch ("max", "spmax", "sparse bool matrix"); |
|
749 DEFUN_DLD (spmax, args, nargout, |
|
750 "-*- texinfo -*-\n\ |
|
751 @deftypefn {Mapping Function} {} spmax (@var{x}, @var{y}, @var{dim})\n\ |
|
752 @deftypefnx {Mapping Function} {[@var{w}, @var{iw}] =} spmax (@var{x})\n\ |
|
753 @cindex Utility Functions\n\ |
|
754 For a vector argument, return the maximum value. For a matrix\n\ |
|
755 argument, return the maximum value from each column, as a row\n\ |
|
756 vector, or over the dimension @var{dim} if defined. For two matrices\n\ |
|
757 (or a matrix and scalar), return the pair-wise maximum.\n\ |
|
758 Thus,\n\ |
|
759 \n\ |
|
760 @example\n\ |
|
761 max (max (@var{x}))\n\ |
|
762 @end example\n\ |
|
763 \n\ |
|
764 @noindent\n\ |
|
765 returns the largest element of @var{x}, and\n\ |
|
766 \n\ |
|
767 @example\n\ |
|
768 @group\n\ |
|
769 max (2:5, pi)\n\ |
|
770 @result{} 3.1416 3.1416 4.0000 5.0000\n\ |
|
771 @end group\n\ |
|
772 @end example\n\ |
|
773 @noindent\n\ |
|
774 compares each element of the range @code{2:5} with @code{pi}, and\n\ |
|
775 returns a row vector of the maximum values.\n\ |
|
776 \n\ |
|
777 For complex arguments, the magnitude of the elements are used for\n\ |
|
778 comparison.\n\ |
|
779 \n\ |
|
780 If called with one input and two output arguments,\n\ |
|
781 @code{max} also returns the first index of the\n\ |
|
782 maximum value(s). Thus,\n\ |
|
783 \n\ |
|
784 @example\n\ |
|
785 @group\n\ |
|
786 [x, ix] = max ([1, 3, 5, 2, 5])\n\ |
|
787 @result{} x = 5\n\ |
|
788 ix = 3\n\ |
|
789 @end group\n\ |
|
790 @end example\n\ |
|
791 @end deftypefn") |
|
792 { |
|
793 MINMAX_BODY (max); |
|
794 } |
|
795 |
|
796 // PKG_ADD: dispatch ("atan2", "spatan2", "sparse matrix"); |
|
797 // PKG_ADD: dispatch ("atan2", "spatan2", "sparse complex matrix"); |
|
798 // PKG_ADD: dispatch ("atan2", "spatan2", "sparse bool matrix"); |
|
799 DEFUN_DLD (spatan2, args, , |
|
800 "-*- texinfo -*-\n\ |
|
801 @deftypefn {Loadable Function} {} spatan2 (@var{y}, @var{x})\n\ |
|
802 Compute atan (Y / X) for corresponding sparse matrix elements of Y and X.\n\ |
|
803 The result is in range -pi to pi.\n\ |
5646
|
804 @end deftypefn") |
5164
|
805 { |
|
806 octave_value retval; |
|
807 int nargin = args.length (); |
|
808 if (nargin == 2) { |
|
809 SparseMatrix a, b; |
|
810 double da, db; |
|
811 bool is_double_a = false; |
|
812 bool is_double_b = false; |
|
813 |
|
814 if (args(0).is_scalar_type ()) |
|
815 { |
|
816 is_double_a = true; |
|
817 da = args(0).double_value(); |
|
818 } |
|
819 else |
|
820 a = args(0).sparse_matrix_value (); |
|
821 |
|
822 if (args(1).is_scalar_type ()) |
|
823 { |
|
824 is_double_b = true; |
|
825 db = args(1).double_value(); |
|
826 } |
|
827 else |
|
828 b = args(1).sparse_matrix_value (); |
|
829 |
|
830 if (is_double_a && is_double_b) |
|
831 retval = Matrix (1, 1, atan2(da, db)); |
|
832 else if (is_double_a) |
|
833 retval = atan2 (da, b); |
|
834 else if (is_double_b) |
|
835 retval = atan2 (a, db); |
|
836 else |
|
837 retval = atan2 (a, b); |
|
838 |
|
839 } else |
5823
|
840 print_usage (); |
5164
|
841 |
|
842 return retval; |
|
843 } |
|
844 |
|
845 static octave_value |
|
846 make_spdiag (const octave_value& a, const octave_value& b) |
|
847 { |
|
848 octave_value retval; |
|
849 |
|
850 if (a.is_complex_type ()) |
|
851 { |
|
852 SparseComplexMatrix m = a.sparse_complex_matrix_value (); |
5275
|
853 octave_idx_type k = b.nint_value(true); |
5164
|
854 |
|
855 if (error_state) |
|
856 return retval; |
|
857 |
5275
|
858 octave_idx_type nr = m.rows (); |
|
859 octave_idx_type nc = m.columns (); |
5164
|
860 |
|
861 if (nr == 0 || nc == 0) |
|
862 retval = m; |
|
863 else if (nr == 1 || nc == 1) |
|
864 { |
5275
|
865 octave_idx_type roff = 0; |
|
866 octave_idx_type coff = 0; |
5164
|
867 if (k > 0) |
|
868 { |
|
869 roff = 0; |
|
870 coff = k; |
|
871 } |
|
872 else if (k < 0) |
|
873 { |
|
874 k = -k; |
|
875 roff = k; |
|
876 coff = 0; |
|
877 } |
|
878 |
|
879 if (nr == 1) |
|
880 { |
5275
|
881 octave_idx_type n = nc + k; |
5604
|
882 octave_idx_type nz = m.nzmax (); |
5164
|
883 SparseComplexMatrix r (n, n, nz); |
5275
|
884 for (octave_idx_type i = 0; i < coff+1; i++) |
5164
|
885 r.xcidx (i) = 0; |
5275
|
886 for (octave_idx_type j = 0; j < nc; j++) |
5164
|
887 { |
5275
|
888 for (octave_idx_type i = m.cidx(j); i < m.cidx(j+1); i++) |
5164
|
889 { |
|
890 r.xdata (i) = m.data (i); |
|
891 r.xridx (i) = j + roff; |
|
892 } |
|
893 r.xcidx (j+coff+1) = m.cidx(j+1); |
|
894 } |
5275
|
895 for (octave_idx_type i = nc+coff+1; i < n+1; i++) |
5164
|
896 r.xcidx (i) = nz; |
|
897 retval = r; |
|
898 } |
|
899 else |
|
900 { |
5275
|
901 octave_idx_type n = nr + k; |
5604
|
902 octave_idx_type nz = m.nzmax (); |
5275
|
903 octave_idx_type ii = 0; |
|
904 octave_idx_type ir = m.ridx(0); |
5164
|
905 SparseComplexMatrix r (n, n, nz); |
5275
|
906 for (octave_idx_type i = 0; i < coff+1; i++) |
5164
|
907 r.xcidx (i) = 0; |
5275
|
908 for (octave_idx_type i = 0; i < nr; i++) |
5164
|
909 { |
|
910 if (ir == i) |
|
911 { |
|
912 r.xdata (ii) = m.data (ii); |
|
913 r.xridx (ii++) = ir + roff; |
|
914 if (ii != nz) |
|
915 ir = m.ridx (ii); |
|
916 } |
|
917 r.xcidx (i+coff+1) = ii; |
|
918 } |
5275
|
919 for (octave_idx_type i = nr+coff+1; i < n+1; i++) |
5164
|
920 r.xcidx (i) = nz; |
|
921 retval = r; |
|
922 } |
|
923 } |
|
924 else |
|
925 { |
|
926 SparseComplexMatrix r = m.diag (k); |
|
927 // Don't use numel, since it can overflow for very large matrices |
|
928 if (r.rows () > 0 && r.cols () > 0) |
|
929 retval = r; |
|
930 } |
|
931 } |
|
932 else if (a.is_real_type ()) |
|
933 { |
|
934 SparseMatrix m = a.sparse_matrix_value (); |
|
935 |
5275
|
936 octave_idx_type k = b.nint_value(true); |
5164
|
937 |
|
938 if (error_state) |
|
939 return retval; |
|
940 |
5275
|
941 octave_idx_type nr = m.rows (); |
|
942 octave_idx_type nc = m.columns (); |
5164
|
943 |
|
944 if (nr == 0 || nc == 0) |
|
945 retval = m; |
|
946 else if (nr == 1 || nc == 1) |
|
947 { |
5275
|
948 octave_idx_type roff = 0; |
|
949 octave_idx_type coff = 0; |
5164
|
950 if (k > 0) |
|
951 { |
|
952 roff = 0; |
|
953 coff = k; |
|
954 } |
|
955 else if (k < 0) |
|
956 { |
|
957 k = -k; |
|
958 roff = k; |
|
959 coff = 0; |
|
960 } |
|
961 |
|
962 if (nr == 1) |
|
963 { |
5275
|
964 octave_idx_type n = nc + k; |
5604
|
965 octave_idx_type nz = m.nzmax (); |
5164
|
966 SparseMatrix r (n, n, nz); |
|
967 |
5275
|
968 for (octave_idx_type i = 0; i < coff+1; i++) |
5164
|
969 r.xcidx (i) = 0; |
5275
|
970 for (octave_idx_type j = 0; j < nc; j++) |
5164
|
971 { |
5275
|
972 for (octave_idx_type i = m.cidx(j); i < m.cidx(j+1); i++) |
5164
|
973 { |
|
974 r.xdata (i) = m.data (i); |
|
975 r.xridx (i) = j + roff; |
|
976 } |
|
977 r.xcidx (j+coff+1) = m.cidx(j+1); |
|
978 } |
5275
|
979 for (octave_idx_type i = nc+coff+1; i < n+1; i++) |
5164
|
980 r.xcidx (i) = nz; |
|
981 retval = r; |
|
982 } |
|
983 else |
|
984 { |
5275
|
985 octave_idx_type n = nr + k; |
5604
|
986 octave_idx_type nz = m.nzmax (); |
5275
|
987 octave_idx_type ii = 0; |
|
988 octave_idx_type ir = m.ridx(0); |
5164
|
989 SparseMatrix r (n, n, nz); |
5275
|
990 for (octave_idx_type i = 0; i < coff+1; i++) |
5164
|
991 r.xcidx (i) = 0; |
5275
|
992 for (octave_idx_type i = 0; i < nr; i++) |
5164
|
993 { |
|
994 if (ir == i) |
|
995 { |
|
996 r.xdata (ii) = m.data (ii); |
|
997 r.xridx (ii++) = ir + roff; |
|
998 if (ii != nz) |
|
999 ir = m.ridx (ii); |
|
1000 } |
|
1001 r.xcidx (i+coff+1) = ii; |
|
1002 } |
5275
|
1003 for (octave_idx_type i = nr+coff+1; i < n+1; i++) |
5164
|
1004 r.xcidx (i) = nz; |
|
1005 retval = r; |
|
1006 } |
|
1007 } |
|
1008 else |
|
1009 { |
|
1010 SparseMatrix r = m.diag (k); |
|
1011 if (r.rows () > 0 && r.cols () > 0) |
|
1012 retval = r; |
|
1013 } |
|
1014 } |
|
1015 else |
|
1016 gripe_wrong_type_arg ("spdiag", a); |
|
1017 |
|
1018 return retval; |
|
1019 } |
|
1020 |
6771
|
1021 static octave_value |
|
1022 make_spdiag (const octave_value& a) |
|
1023 { |
|
1024 octave_value retval; |
|
1025 octave_idx_type nr = a.rows (); |
|
1026 octave_idx_type nc = a.columns (); |
|
1027 |
|
1028 if (nr == 0 || nc == 0) |
|
1029 retval = SparseMatrix (); |
|
1030 else |
|
1031 retval = make_spdiag (a, octave_value (0.)); |
|
1032 |
|
1033 return retval; |
|
1034 } |
|
1035 |
5164
|
1036 // PKG_ADD: dispatch ("diag", "spdiag", "sparse matrix"); |
|
1037 // PKG_ADD: dispatch ("diag", "spdiag", "sparse complex matrix"); |
|
1038 // PKG_ADD: dispatch ("diag", "spdiag", "sparse bool matrix"); |
|
1039 DEFUN_DLD (spdiag, args, , |
|
1040 "-*- texinfo -*-\n\ |
|
1041 @deftypefn {Loadable Function} {} spdiag (@var{v}, @var{k})\n\ |
|
1042 Return a diagonal matrix with the sparse vector @var{v} on diagonal\n\ |
|
1043 @var{k}. The second argument is optional. If it is positive, the vector is\n\ |
|
1044 placed on the @var{k}-th super-diagonal. If it is negative, it is placed\n\ |
|
1045 on the @var{-k}-th sub-diagonal. The default value of @var{k} is 0, and\n\ |
|
1046 the vector is placed on the main diagonal. For example,\n\ |
|
1047 \n\ |
|
1048 @example\n\ |
6772
|
1049 @group\n\ |
5164
|
1050 spdiag ([1, 2, 3], 1)\n\ |
|
1051 ans =\n\ |
|
1052 \n\ |
|
1053 Compressed Column Sparse (rows=4, cols=4, nnz=3)\n\ |
|
1054 (1 , 2) -> 1\n\ |
|
1055 (2 , 3) -> 2\n\ |
|
1056 (3 , 4) -> 3\n\ |
6772
|
1057 @end group\n\ |
5164
|
1058 @end example\n\ |
6772
|
1059 \n\ |
|
1060 @noindent\n\ |
|
1061 Given a matrix argument, instead of a vector, @code{spdiag} extracts the\n\ |
6774
|
1062 @var{k}-th diagonal of the sparse matrix.\n\ |
5642
|
1063 @seealso{diag}\n\ |
|
1064 @end deftypefn") |
5164
|
1065 { |
|
1066 octave_value retval; |
|
1067 |
|
1068 int nargin = args.length (); |
|
1069 |
|
1070 if (nargin == 1 && args(0).is_defined ()) |
6771
|
1071 retval = make_spdiag (args(0)); |
5164
|
1072 else if (nargin == 2 && args(0).is_defined () && args(1).is_defined ()) |
|
1073 retval = make_spdiag (args(0), args(1)); |
|
1074 else |
5823
|
1075 print_usage (); |
5164
|
1076 |
|
1077 return retval; |
|
1078 } |
|
1079 |
|
1080 /* |
|
1081 ;;; Local Variables: *** |
|
1082 ;;; mode: C++ *** |
|
1083 ;;; End: *** |
|
1084 */ |