Mercurial > hg > octave-lyh
comparison src/data.cc @ 10435:6a271334750c
implement general binary mapping facility
author | Jaroslav Hajek <highegg@gmail.com> |
---|---|
date | Tue, 23 Mar 2010 09:30:46 +0100 |
parents | cc69a17ec801 |
children | 00219bdd2d17 |
comparison
equal
deleted
inserted
replaced
10434:f387c5b3a369 | 10435:6a271334750c |
---|---|
43 #include "lo-math.h" | 43 #include "lo-math.h" |
44 #include "oct-time.h" | 44 #include "oct-time.h" |
45 #include "str-vec.h" | 45 #include "str-vec.h" |
46 #include "quit.h" | 46 #include "quit.h" |
47 #include "mx-base.h" | 47 #include "mx-base.h" |
48 #include "oct-binmap.h" | |
48 | 49 |
49 #include "Cell.h" | 50 #include "Cell.h" |
50 #include "defun.h" | 51 #include "defun.h" |
51 #include "error.h" | 52 #include "error.h" |
52 #include "gripes.h" | 53 #include "gripes.h" |
192 | 193 |
193 */ | 194 */ |
194 | 195 |
195 // These mapping functions may also be useful in other places, eh? | 196 // These mapping functions may also be useful in other places, eh? |
196 | 197 |
197 typedef double (*d_dd_fcn) (double, double); | |
198 typedef float (*f_ff_fcn) (float, float); | |
199 | |
200 static NDArray | |
201 map_d_m (d_dd_fcn f, double x, const NDArray& y) | |
202 { | |
203 NDArray retval (y.dims ()); | |
204 double *r_data = retval.fortran_vec (); | |
205 | |
206 const double *y_data = y.data (); | |
207 | |
208 octave_idx_type nel = y.numel (); | |
209 | |
210 for (octave_idx_type i = 0; i < nel; i++) | |
211 { | |
212 octave_quit (); | |
213 r_data[i] = f (x, y_data[i]); | |
214 } | |
215 | |
216 return retval; | |
217 } | |
218 | |
219 static FloatNDArray | |
220 map_f_fm (f_ff_fcn f, float x, const FloatNDArray& y) | |
221 { | |
222 FloatNDArray retval (y.dims ()); | |
223 float *r_data = retval.fortran_vec (); | |
224 | |
225 const float *y_data = y.data (); | |
226 | |
227 octave_idx_type nel = y.numel (); | |
228 | |
229 for (octave_idx_type i = 0; i < nel; i++) | |
230 { | |
231 octave_quit (); | |
232 r_data[i] = f (x, y_data[i]); | |
233 } | |
234 | |
235 return retval; | |
236 } | |
237 | |
238 static NDArray | |
239 map_m_d (d_dd_fcn f, const NDArray& x, double y) | |
240 { | |
241 NDArray retval (x.dims ()); | |
242 double *r_data = retval.fortran_vec (); | |
243 | |
244 const double *x_data = x.data (); | |
245 | |
246 octave_idx_type nel = x.numel (); | |
247 | |
248 for (octave_idx_type i = 0; i < nel; i++) | |
249 { | |
250 octave_quit (); | |
251 r_data[i] = f (x_data[i], y); | |
252 } | |
253 | |
254 return retval; | |
255 } | |
256 | |
257 static FloatNDArray | |
258 map_fm_f (f_ff_fcn f, const FloatNDArray& x, float y) | |
259 { | |
260 FloatNDArray retval (x.dims ()); | |
261 float *r_data = retval.fortran_vec (); | |
262 | |
263 const float *x_data = x.data (); | |
264 | |
265 octave_idx_type nel = x.numel (); | |
266 | |
267 for (octave_idx_type i = 0; i < nel; i++) | |
268 { | |
269 octave_quit (); | |
270 r_data[i] = f (x_data[i], y); | |
271 } | |
272 | |
273 return retval; | |
274 } | |
275 | |
276 static NDArray | |
277 map_m_m (d_dd_fcn f, const NDArray& x, const NDArray& y) | |
278 { | |
279 assert (x.dims () == y.dims ()); | |
280 | |
281 NDArray retval (x.dims ()); | |
282 double *r_data = retval.fortran_vec (); | |
283 | |
284 const double *x_data = x.data (); | |
285 const double *y_data = y.data (); | |
286 | |
287 octave_idx_type nel = x.numel (); | |
288 | |
289 for (octave_idx_type i = 0; i < nel; i++) | |
290 { | |
291 octave_quit (); | |
292 r_data[i] = f (x_data[i], y_data[i]); | |
293 } | |
294 | |
295 return retval; | |
296 } | |
297 | |
298 static FloatNDArray | |
299 map_fm_fm (f_ff_fcn f, const FloatNDArray& x, const FloatNDArray& y) | |
300 { | |
301 assert (x.dims () == y.dims ()); | |
302 | |
303 FloatNDArray retval (x.dims ()); | |
304 float *r_data = retval.fortran_vec (); | |
305 | |
306 const float *x_data = x.data (); | |
307 const float *y_data = y.data (); | |
308 | |
309 octave_idx_type nel = x.numel (); | |
310 | |
311 for (octave_idx_type i = 0; i < nel; i++) | |
312 { | |
313 octave_quit (); | |
314 r_data[i] = f (x_data[i], y_data[i]); | |
315 } | |
316 | |
317 return retval; | |
318 } | |
319 | |
320 static SparseMatrix | |
321 map_d_s (d_dd_fcn f, double x, const SparseMatrix& y) | |
322 { | |
323 octave_idx_type nr = y.rows (); | |
324 octave_idx_type nc = y.columns (); | |
325 SparseMatrix retval; | |
326 double f_zero = f (x, 0.); | |
327 | |
328 if (f_zero != 0) | |
329 { | |
330 retval = SparseMatrix (nr, nc, f_zero); | |
331 | |
332 for (octave_idx_type j = 0; j < nc; j++) | |
333 for (octave_idx_type i = y.cidx (j); i < y.cidx (j+1); i++) | |
334 { | |
335 octave_quit (); | |
336 retval.data (y.ridx(i) + j * nr) = f (x, y.data (i)); | |
337 } | |
338 | |
339 retval.maybe_compress (true); | |
340 } | |
341 else | |
342 { | |
343 octave_idx_type nz = y.nnz (); | |
344 retval = SparseMatrix (nr, nc, nz); | |
345 octave_idx_type ii = 0; | |
346 retval.cidx (ii) = 0; | |
347 | |
348 for (octave_idx_type j = 0; j < nc; j++) | |
349 { | |
350 for (octave_idx_type i = y.cidx (j); i < y.cidx (j+1); i++) | |
351 { | |
352 octave_quit (); | |
353 double val = f (x, y.data (i)); | |
354 | |
355 if (val != 0.0) | |
356 { | |
357 retval.data (ii) = val; | |
358 retval.ridx (ii++) = y.ridx (i); | |
359 } | |
360 } | |
361 retval.cidx (j + 1) = ii; | |
362 } | |
363 | |
364 retval.maybe_compress (false); | |
365 } | |
366 | |
367 return retval; | |
368 } | |
369 | |
370 static SparseMatrix | |
371 map_s_d (d_dd_fcn f, const SparseMatrix& x, double y) | |
372 { | |
373 octave_idx_type nr = x.rows (); | |
374 octave_idx_type nc = x.columns (); | |
375 SparseMatrix retval; | |
376 double f_zero = f (0., y); | |
377 | |
378 if (f_zero != 0) | |
379 { | |
380 retval = SparseMatrix (nr, nc, f_zero); | |
381 | |
382 for (octave_idx_type j = 0; j < nc; j++) | |
383 for (octave_idx_type i = x.cidx (j); i < x.cidx (j+1); i++) | |
384 { | |
385 octave_quit (); | |
386 retval.data (x.ridx(i) + j * nr) = f (x.data (i), y); | |
387 } | |
388 | |
389 retval.maybe_compress (true); | |
390 } | |
391 else | |
392 { | |
393 octave_idx_type nz = x.nnz (); | |
394 retval = SparseMatrix (nr, nc, nz); | |
395 octave_idx_type ii = 0; | |
396 retval.cidx (ii) = 0; | |
397 | |
398 for (octave_idx_type j = 0; j < nc; j++) | |
399 { | |
400 for (octave_idx_type i = x.cidx (j); i < x.cidx (j+1); i++) | |
401 { | |
402 octave_quit (); | |
403 double val = f (x.data (i), y); | |
404 | |
405 if (val != 0.0) | |
406 { | |
407 retval.data (ii) = val; | |
408 retval.ridx (ii++) = x.ridx (i); | |
409 } | |
410 } | |
411 retval.cidx (j + 1) = ii; | |
412 } | |
413 | |
414 retval.maybe_compress (false); | |
415 } | |
416 | |
417 return retval; | |
418 } | |
419 | |
420 static SparseMatrix | |
421 map_s_s (d_dd_fcn f, const SparseMatrix& x, const SparseMatrix& y) | |
422 { | |
423 octave_idx_type nr = x.rows (); | |
424 octave_idx_type nc = x.columns (); | |
425 | |
426 octave_idx_type y_nr = y.rows (); | |
427 octave_idx_type y_nc = y.columns (); | |
428 | |
429 assert (nr == y_nr && nc == y_nc); | |
430 | |
431 SparseMatrix retval; | |
432 double f_zero = f (0., 0.); | |
433 | |
434 if (f_zero != 0) | |
435 { | |
436 retval = SparseMatrix (nr, nc, f_zero); | |
437 octave_idx_type k1 = 0, k2 = 0; | |
438 | |
439 for (octave_idx_type j = 0; j < nc; j++) | |
440 { | |
441 while (k1 < x.cidx(j+1) && k2 < y.cidx(j+1)) | |
442 { | |
443 octave_quit (); | |
444 if (k1 >= x.cidx(j+1)) | |
445 { | |
446 retval.data (y.ridx(k2) + j * nr) = f (0.0, y.data (k2)); | |
447 k2++; | |
448 } | |
449 else if (k2 >= y.cidx(j+1)) | |
450 { | |
451 retval.data (x.ridx(k1) + j * nr) = f (x.data (k1), 0.0); | |
452 k1++; | |
453 } | |
454 else | |
455 { | |
456 octave_idx_type rx = x.ridx(k1); | |
457 octave_idx_type ry = y.ridx(k2); | |
458 | |
459 if (rx < ry) | |
460 { | |
461 retval.data (rx + j * nr) = f (x.data (k1), 0.0); | |
462 k1++; | |
463 } | |
464 else if (rx > ry) | |
465 { | |
466 retval.data (ry + j * nr) = f (0.0, y.data (k2)); | |
467 k2++; | |
468 } | |
469 else | |
470 { | |
471 retval.data (ry + j * nr) = f (x.data (k1), y.data (k2)); | |
472 k1++; | |
473 k2++; | |
474 } | |
475 } | |
476 } | |
477 } | |
478 | |
479 retval.maybe_compress (true); | |
480 } | |
481 else | |
482 { | |
483 octave_idx_type nz = x.nnz () + y.nnz (); | |
484 retval = SparseMatrix (nr, nc, nz); | |
485 octave_idx_type ii = 0; | |
486 retval.cidx (ii) = 0; | |
487 octave_idx_type k1 = 0, k2 = 0; | |
488 | |
489 for (octave_idx_type j = 0; j < nc; j++) | |
490 { | |
491 while (k1 < x.cidx(j+1) && k2 < y.cidx(j+1)) | |
492 { | |
493 octave_quit (); | |
494 double val; | |
495 octave_idx_type r; | |
496 if (k1 >= x.cidx(j+1)) | |
497 { | |
498 r = y.ridx (k2); | |
499 val = f (0.0, y.data (k2++)); | |
500 } | |
501 else if (k2 >= y.cidx(j+1)) | |
502 { | |
503 r = x.ridx (k1); | |
504 val = f (x.data (k1++), 0.0); | |
505 } | |
506 else | |
507 { | |
508 octave_idx_type rx = x.ridx(k1); | |
509 octave_idx_type ry = y.ridx(k2); | |
510 | |
511 if (rx < ry) | |
512 { | |
513 r = x.ridx (k1); | |
514 val = f (x.data (k1++), 0.0); | |
515 } | |
516 else if (rx > ry) | |
517 { | |
518 r = y.ridx (k2); | |
519 val = f (0.0, y.data (k2++)); | |
520 } | |
521 else | |
522 { | |
523 r = y.ridx (k2); | |
524 val = f (x.data (k1++), y.data (k2++)); | |
525 } | |
526 } | |
527 if (val != 0.0) | |
528 { | |
529 retval.data (ii) = val; | |
530 retval.ridx (ii++) = r; | |
531 } | |
532 } | |
533 retval.cidx (j + 1) = ii; | |
534 } | |
535 | |
536 retval.maybe_compress (false); | |
537 } | |
538 | |
539 return retval; | |
540 } | |
541 | |
542 DEFUN (atan2, args, , | 198 DEFUN (atan2, args, , |
543 "-*- texinfo -*-\n\ | 199 "-*- texinfo -*-\n\ |
544 @deftypefn {Mapping Function} {} atan2 (@var{y}, @var{x})\n\ | 200 @deftypefn {Mapping Function} {} atan2 (@var{y}, @var{x})\n\ |
545 Compute atan (@var{y} / @var{x}) for corresponding elements of @var{y}\n\ | 201 Compute atan (@var{y} / @var{x}) for corresponding elements of @var{y}\n\ |
546 and @var{x}. Signal an error if @var{y} and @var{x} do not match in size\n\ | 202 and @var{x}. Signal an error if @var{y} and @var{x} do not match in size\n\ |
549 { | 205 { |
550 octave_value retval; | 206 octave_value retval; |
551 | 207 |
552 int nargin = args.length (); | 208 int nargin = args.length (); |
553 | 209 |
554 if (nargin == 2 && args(0).is_defined () && args(1).is_defined ()) | 210 if (nargin == 2) |
555 { | 211 { |
556 if (args(0).is_integer_type () || args(1).is_integer_type ()) | 212 if (! args(0).is_numeric_type ()) |
557 error ("atan2: not defined for integer types"); | 213 gripe_wrong_type_arg ("atan2", args(0)); |
214 else if (! args(1).is_numeric_type ()) | |
215 gripe_wrong_type_arg ("atan2", args(1)); | |
216 else if (args(0).is_complex_type () || args(1).is_complex_type ()) | |
217 error ("atan2: not defined for complex numbers"); | |
218 else if (args(0).is_single_type () || args(1).is_single_type ()) | |
219 { | |
220 if (args(0).is_scalar_type () && args(1).is_scalar_type ()) | |
221 retval = atan2f (args(0).float_value (), args(1).float_value ()); | |
222 else | |
223 { | |
224 FloatNDArray a0 = args(0).float_array_value (); | |
225 FloatNDArray a1 = args(1).float_array_value (); | |
226 retval = binmap<float> (a0, a1, ::atan2f, "atan2"); | |
227 } | |
228 } | |
558 else | 229 else |
559 { | 230 { |
560 octave_value arg_y = args(0); | 231 bool a0_scalar = args(0).is_scalar_type (); |
561 octave_value arg_x = args(1); | 232 bool a1_scalar = args(1).is_scalar_type (); |
562 | 233 if (a0_scalar && a1_scalar) |
563 dim_vector y_dims = arg_y.dims (); | 234 retval = atan2 (args(0).scalar_value (), args(1).scalar_value ()); |
564 dim_vector x_dims = arg_x.dims (); | 235 else if ((a0_scalar || args(0).is_sparse_type ()) |
565 | 236 && (a1_scalar || args(1).is_sparse_type ())) |
566 bool y_is_scalar = y_dims.all_ones (); | |
567 bool x_is_scalar = x_dims.all_ones (); | |
568 | |
569 bool is_float = arg_y.is_single_type () || arg_x.is_single_type (); | |
570 | |
571 if (y_is_scalar && x_is_scalar) | |
572 { | 237 { |
573 if (is_float) | 238 SparseMatrix m0 = args(0).sparse_matrix_value (); |
574 { | 239 SparseMatrix m1 = args(1).sparse_matrix_value (); |
575 float y = arg_y.float_value (); | 240 retval = binmap<double> (m0, m1, ::atan2, "atan2"); |
576 | |
577 if (! error_state) | |
578 { | |
579 float x = arg_x.float_value (); | |
580 | |
581 if (! error_state) | |
582 retval = atan2f (y, x); | |
583 } | |
584 } | |
585 else | |
586 { | |
587 double y = arg_y.double_value (); | |
588 | |
589 if (! error_state) | |
590 { | |
591 double x = arg_x.double_value (); | |
592 | |
593 if (! error_state) | |
594 retval = atan2 (y, x); | |
595 } | |
596 } | |
597 } | |
598 else if (y_is_scalar) | |
599 { | |
600 if (is_float) | |
601 { | |
602 float y = arg_y.float_value (); | |
603 | |
604 if (! error_state) | |
605 { | |
606 // Even if x is sparse return a full matrix here | |
607 FloatNDArray x = arg_x.float_array_value (); | |
608 | |
609 if (! error_state) | |
610 retval = map_f_fm (atan2f, y, x); | |
611 } | |
612 } | |
613 else | |
614 { | |
615 double y = arg_y.double_value (); | |
616 | |
617 if (! error_state) | |
618 { | |
619 // Even if x is sparse return a full matrix here | |
620 NDArray x = arg_x.array_value (); | |
621 | |
622 if (! error_state) | |
623 retval = map_d_m (atan2, y, x); | |
624 } | |
625 } | |
626 } | |
627 else if (x_is_scalar) | |
628 { | |
629 if (arg_y.is_sparse_type ()) | |
630 { | |
631 SparseMatrix y = arg_y.sparse_matrix_value (); | |
632 | |
633 if (! error_state) | |
634 { | |
635 double x = arg_x.double_value (); | |
636 | |
637 if (! error_state) | |
638 retval = map_s_d (atan2, y, x); | |
639 } | |
640 } | |
641 else if (is_float) | |
642 { | |
643 FloatNDArray y = arg_y.float_array_value (); | |
644 | |
645 if (! error_state) | |
646 { | |
647 float x = arg_x.float_value (); | |
648 | |
649 if (! error_state) | |
650 retval = map_fm_f (atan2f, y, x); | |
651 } | |
652 } | |
653 else | |
654 { | |
655 NDArray y = arg_y.array_value (); | |
656 | |
657 if (! error_state) | |
658 { | |
659 double x = arg_x.double_value (); | |
660 | |
661 if (! error_state) | |
662 retval = map_m_d (atan2, y, x); | |
663 } | |
664 } | |
665 } | |
666 else if (y_dims == x_dims) | |
667 { | |
668 // Even if y is sparse return a full matrix here | |
669 if (arg_x.is_sparse_type ()) | |
670 { | |
671 SparseMatrix y = arg_y.sparse_matrix_value (); | |
672 | |
673 if (! error_state) | |
674 { | |
675 SparseMatrix x = arg_x.sparse_matrix_value (); | |
676 | |
677 if (! error_state) | |
678 retval = map_s_s (atan2, y, x); | |
679 } | |
680 } | |
681 else if (is_float) | |
682 { | |
683 FloatNDArray y = arg_y.array_value (); | |
684 | |
685 if (! error_state) | |
686 { | |
687 FloatNDArray x = arg_x.array_value (); | |
688 | |
689 if (! error_state) | |
690 retval = map_fm_fm (atan2f, y, x); | |
691 } | |
692 } | |
693 else | |
694 { | |
695 NDArray y = arg_y.array_value (); | |
696 | |
697 if (! error_state) | |
698 { | |
699 NDArray x = arg_x.array_value (); | |
700 | |
701 if (! error_state) | |
702 retval = map_m_m (atan2, y, x); | |
703 } | |
704 } | |
705 } | 241 } |
706 else | 242 else |
707 error ("atan2: nonconformant matrices"); | 243 { |
244 NDArray a0 = args(0).array_value (); | |
245 NDArray a1 = args(1).array_value (); | |
246 retval = binmap<double> (a0, a1, ::atan2, "atan2"); | |
247 } | |
708 } | 248 } |
709 } | 249 } |
710 else | 250 else |
711 print_usage (); | 251 print_usage (); |
712 | 252 |
754 { | 294 { |
755 octave_value retval; | 295 octave_value retval; |
756 | 296 |
757 int nargin = args.length (); | 297 int nargin = args.length (); |
758 | 298 |
759 if (nargin == 2 && args(0).is_defined () && args(1).is_defined ()) | 299 if (nargin == 2) |
760 { | 300 { |
761 if (args(0).is_integer_type () || args(1).is_integer_type ()) | 301 octave_value arg0 = args(0), arg1 = args(1); |
762 error ("hypot: not defined for integer types"); | 302 if (! arg0.is_numeric_type ()) |
303 gripe_wrong_type_arg ("hypot", arg0); | |
304 else if (! arg1.is_numeric_type ()) | |
305 gripe_wrong_type_arg ("hypot", arg1); | |
763 else | 306 else |
764 { | 307 { |
765 octave_value arg_x = args(0); | 308 if (arg0.is_complex_type ()) |
766 octave_value arg_y = args(1); | 309 arg0 = arg0.abs (); |
767 | 310 if (arg1.is_complex_type ()) |
768 dim_vector x_dims = arg_x.dims (); | 311 arg1 = arg1.abs (); |
769 dim_vector y_dims = arg_y.dims (); | 312 |
770 | 313 if (arg0.is_single_type () || arg1.is_single_type ()) |
771 bool x_is_scalar = x_dims.all_ones (); | |
772 bool y_is_scalar = y_dims.all_ones (); | |
773 | |
774 bool is_float = arg_y.is_single_type () || arg_x.is_single_type (); | |
775 | |
776 if (y_is_scalar && x_is_scalar) | |
777 { | 314 { |
778 if (is_float) | 315 if (arg0.is_scalar_type () && arg1.is_scalar_type ()) |
316 retval = hypotf (arg0.float_value (), arg1.float_value ()); | |
317 else | |
779 { | 318 { |
780 float x; | 319 FloatNDArray a0 = arg0.float_array_value (); |
781 if (arg_x.is_complex_type ()) | 320 FloatNDArray a1 = arg1.float_array_value (); |
782 x = abs (arg_x.float_complex_value ()); | 321 retval = binmap<float> (a0, a1, ::hypotf, "hypot"); |
783 else | 322 } |
784 x = arg_x.float_value (); | 323 } |
785 | 324 else |
786 if (! error_state) | 325 { |
787 { | 326 bool a0_scalar = arg0.is_scalar_type (); |
788 float y; | 327 bool a1_scalar = arg1.is_scalar_type (); |
789 if (arg_y.is_complex_type ()) | 328 if (a0_scalar && a1_scalar) |
790 y = abs (arg_y.float_complex_value ()); | 329 retval = hypot (arg0.scalar_value (), arg1.scalar_value ()); |
791 else | 330 else if ((a0_scalar || arg0.is_sparse_type ()) |
792 y = arg_y.float_value (); | 331 && (a1_scalar || arg1.is_sparse_type ())) |
793 | 332 { |
794 if (! error_state) | 333 SparseMatrix m0 = arg0.sparse_matrix_value (); |
795 retval = hypotf (x, y); | 334 SparseMatrix m1 = arg1.sparse_matrix_value (); |
796 } | 335 retval = binmap<double> (m0, m1, ::hypot, "hypot"); |
797 } | 336 } |
798 else | 337 else |
799 { | 338 { |
800 double x; | 339 NDArray a0 = arg0.array_value (); |
801 if (arg_x.is_complex_type ()) | 340 NDArray a1 = arg1.array_value (); |
802 x = abs (arg_x.complex_value ()); | 341 retval = binmap<double> (a0, a1, ::hypot, "hypot"); |
803 else | |
804 x = arg_x.double_value (); | |
805 | |
806 if (! error_state) | |
807 { | |
808 double y; | |
809 if (arg_y.is_complex_type ()) | |
810 y = abs (arg_y.complex_value ()); | |
811 else | |
812 y = arg_y.double_value (); | |
813 | |
814 if (! error_state) | |
815 retval = hypot (x, y); | |
816 } | |
817 } | 342 } |
818 } | 343 } |
819 else if (y_is_scalar) | |
820 { | |
821 if (is_float) | |
822 { | |
823 FloatNDArray x; | |
824 if (arg_x.is_complex_type ()) | |
825 x = arg_x.float_complex_array_value ().abs (); | |
826 else | |
827 x = arg_x.float_array_value (); | |
828 | |
829 if (! error_state) | |
830 { | |
831 float y; | |
832 if (arg_y.is_complex_type ()) | |
833 y = abs (arg_y.float_complex_value ()); | |
834 else | |
835 y = arg_y.float_value (); | |
836 | |
837 if (! error_state) | |
838 retval = map_fm_f (hypotf, x, y); | |
839 } | |
840 } | |
841 else | |
842 { | |
843 NDArray x; | |
844 if (arg_x.is_complex_type ()) | |
845 x = arg_x.complex_array_value ().abs (); | |
846 else | |
847 x = arg_x.array_value (); | |
848 | |
849 if (! error_state) | |
850 { | |
851 double y; | |
852 if (arg_y.is_complex_type ()) | |
853 y = abs (arg_y.complex_value ()); | |
854 else | |
855 y = arg_y.double_value (); | |
856 | |
857 if (! error_state) | |
858 retval = map_m_d (hypot, x, y); | |
859 } | |
860 } | |
861 } | |
862 else if (x_is_scalar) | |
863 { | |
864 if (is_float) | |
865 { | |
866 float x; | |
867 if (arg_x.is_complex_type ()) | |
868 x = abs (arg_x.float_complex_value ()); | |
869 else | |
870 x = arg_x.float_value (); | |
871 | |
872 if (! error_state) | |
873 { | |
874 FloatNDArray y; | |
875 if (arg_y.is_complex_type ()) | |
876 y = arg_y.float_complex_array_value ().abs (); | |
877 else | |
878 y = arg_y.float_array_value (); | |
879 | |
880 if (! error_state) | |
881 retval = map_f_fm (hypotf, x, y); | |
882 } | |
883 } | |
884 else | |
885 { | |
886 double x; | |
887 if (arg_x.is_complex_type ()) | |
888 x = abs (arg_x.complex_value ()); | |
889 else | |
890 x = arg_x.double_value (); | |
891 | |
892 if (! error_state) | |
893 { | |
894 NDArray y; | |
895 if (arg_y.is_complex_type ()) | |
896 y = arg_y.complex_array_value ().abs (); | |
897 else | |
898 y = arg_y.array_value (); | |
899 | |
900 if (! error_state) | |
901 retval = map_d_m (hypot, x, y); | |
902 } | |
903 } | |
904 } | |
905 else if (y_dims == x_dims) | |
906 { | |
907 if (arg_x.is_sparse_type () && arg_y.is_sparse_type ()) | |
908 { | |
909 SparseMatrix x; | |
910 if (arg_x.is_complex_type ()) | |
911 x = arg_x.sparse_complex_matrix_value ().abs (); | |
912 else | |
913 x = arg_x.sparse_matrix_value (); | |
914 | |
915 if (! error_state) | |
916 { | |
917 SparseMatrix y; | |
918 if (arg_y.is_complex_type ()) | |
919 y = arg_y.sparse_complex_matrix_value ().abs (); | |
920 else | |
921 y = arg_y.sparse_matrix_value (); | |
922 | |
923 if (! error_state) | |
924 retval = map_s_s (hypot, x, y); | |
925 } | |
926 } | |
927 else if (is_float) | |
928 { | |
929 FloatNDArray x; | |
930 if (arg_x.is_complex_type ()) | |
931 x = arg_x.float_complex_array_value ().abs (); | |
932 else | |
933 x = arg_x.float_array_value (); | |
934 | |
935 if (! error_state) | |
936 { | |
937 FloatNDArray y; | |
938 if (arg_y.is_complex_type ()) | |
939 y = arg_y.float_complex_array_value ().abs (); | |
940 else | |
941 y = arg_y.float_array_value (); | |
942 | |
943 if (! error_state) | |
944 retval = map_fm_fm (hypotf, x, y); | |
945 } | |
946 } | |
947 else | |
948 { | |
949 NDArray x; | |
950 if (arg_x.is_complex_type ()) | |
951 x = arg_x.complex_array_value ().abs (); | |
952 else | |
953 x = arg_x.array_value (); | |
954 | |
955 if (! error_state) | |
956 { | |
957 NDArray y; | |
958 if (arg_y.is_complex_type ()) | |
959 y = arg_y.complex_array_value ().abs (); | |
960 else | |
961 y = arg_y.array_value (); | |
962 | |
963 if (! error_state) | |
964 retval = map_m_m (hypot, x, y); | |
965 } | |
966 } | |
967 } | |
968 else | |
969 error ("hypot: nonconformant matrices"); | |
970 } | 344 } |
971 } | 345 } |
972 else | 346 else |
973 print_usage (); | 347 print_usage (); |
974 | 348 |
1110 { | 484 { |
1111 octave_value retval; | 485 octave_value retval; |
1112 | 486 |
1113 int nargin = args.length (); | 487 int nargin = args.length (); |
1114 | 488 |
1115 if (nargin == 2 && args(0).is_defined () && args(1).is_defined ()) | 489 if (nargin == 2) |
1116 { | 490 { |
1117 octave_value arg_x = args(0); | 491 if (! args(0).is_numeric_type ()) |
1118 octave_value arg_y = args(1); | 492 gripe_wrong_type_arg ("fmod", args(0)); |
1119 | 493 else if (! args(1).is_numeric_type ()) |
1120 dim_vector y_dims = arg_y.dims (); | 494 gripe_wrong_type_arg ("fmod", args(1)); |
1121 dim_vector x_dims = arg_x.dims (); | 495 else if (args(0).is_complex_type () || args(1).is_complex_type ()) |
1122 | 496 error ("fmod: not defined for complex numbers"); |
1123 bool y_is_scalar = y_dims.all_ones (); | 497 else if (args(0).is_single_type () || args(1).is_single_type ()) |
1124 bool x_is_scalar = x_dims.all_ones (); | 498 { |
1125 | 499 if (args(0).is_scalar_type () && args(1).is_scalar_type ()) |
1126 bool is_float = arg_y.is_single_type () || arg_x.is_single_type (); | 500 retval = fmodf (args(0).float_value (), args(1).float_value ()); |
1127 | 501 else |
1128 if (y_is_scalar && x_is_scalar) | |
1129 { | |
1130 if (is_float) | |
1131 { | 502 { |
1132 float y = arg_y.float_value (); | 503 FloatNDArray a0 = args(0).float_array_value (); |
1133 | 504 FloatNDArray a1 = args(1).float_array_value (); |
1134 if (! error_state) | 505 retval = binmap<float> (a0, a1, ::fmodf, "fmod"); |
1135 { | 506 } |
1136 float x = arg_x.float_value (); | 507 } |
1137 | 508 else |
1138 if (! error_state) | 509 { |
1139 retval = fmod (x, y); | 510 bool a0_scalar = args(0).is_scalar_type (); |
1140 } | 511 bool a1_scalar = args(1).is_scalar_type (); |
512 if (a0_scalar && a1_scalar) | |
513 retval = fmod (args(0).scalar_value (), args(1).scalar_value ()); | |
514 else if ((a0_scalar || args(0).is_sparse_type ()) | |
515 && (a1_scalar || args(1).is_sparse_type ())) | |
516 { | |
517 SparseMatrix m0 = args(0).sparse_matrix_value (); | |
518 SparseMatrix m1 = args(1).sparse_matrix_value (); | |
519 retval = binmap<double> (m0, m1, ::fmod, "fmod"); | |
1141 } | 520 } |
1142 else | 521 else |
1143 { | 522 { |
1144 double y = arg_y.double_value (); | 523 NDArray a0 = args(0).array_value (); |
1145 | 524 NDArray a1 = args(1).array_value (); |
1146 if (! error_state) | 525 retval = binmap<double> (a0, a1, ::fmod, "fmod"); |
1147 { | |
1148 double x = arg_x.double_value (); | |
1149 | |
1150 if (! error_state) | |
1151 retval = fmod (x, y); | |
1152 } | |
1153 } | 526 } |
1154 } | 527 } |
1155 else if (y_is_scalar) | |
1156 { | |
1157 if (is_float) | |
1158 { | |
1159 float y = arg_y.float_value (); | |
1160 | |
1161 if (! error_state) | |
1162 { | |
1163 FloatNDArray x = arg_x.float_array_value (); | |
1164 | |
1165 if (! error_state) | |
1166 retval = map_fm_f (fmodf, x, y); | |
1167 } | |
1168 } | |
1169 else | |
1170 { | |
1171 double y = arg_y.double_value (); | |
1172 | |
1173 if (! error_state) | |
1174 { | |
1175 if (arg_x.is_sparse_type ()) | |
1176 { | |
1177 SparseMatrix x = arg_x.sparse_matrix_value (); | |
1178 | |
1179 if (! error_state) | |
1180 retval = map_s_d (fmod, x, y); | |
1181 } | |
1182 else | |
1183 { | |
1184 NDArray x = arg_x.array_value (); | |
1185 | |
1186 if (! error_state) | |
1187 retval = map_m_d (fmod, x, y); | |
1188 } | |
1189 } | |
1190 } | |
1191 } | |
1192 else if (x_is_scalar) | |
1193 { | |
1194 if (arg_y.is_sparse_type ()) | |
1195 { | |
1196 SparseMatrix y = arg_y.sparse_matrix_value (); | |
1197 | |
1198 if (! error_state) | |
1199 { | |
1200 double x = arg_x.double_value (); | |
1201 | |
1202 if (! error_state) | |
1203 retval = map_d_s (fmod, x, y); | |
1204 } | |
1205 } | |
1206 else if (is_float) | |
1207 { | |
1208 FloatNDArray y = arg_y.float_array_value (); | |
1209 | |
1210 if (! error_state) | |
1211 { | |
1212 float x = arg_x.float_value (); | |
1213 | |
1214 if (! error_state) | |
1215 retval = map_f_fm (fmodf, x, y); | |
1216 } | |
1217 } | |
1218 else | |
1219 { | |
1220 NDArray y = arg_y.array_value (); | |
1221 | |
1222 if (! error_state) | |
1223 { | |
1224 double x = arg_x.double_value (); | |
1225 | |
1226 if (! error_state) | |
1227 retval = map_d_m (fmod, x, y); | |
1228 } | |
1229 } | |
1230 } | |
1231 else if (y_dims == x_dims) | |
1232 { | |
1233 if (arg_y.is_sparse_type () || arg_x.is_sparse_type ()) | |
1234 { | |
1235 SparseMatrix y = arg_y.sparse_matrix_value (); | |
1236 | |
1237 if (! error_state) | |
1238 { | |
1239 SparseMatrix x = arg_x.sparse_matrix_value (); | |
1240 | |
1241 if (! error_state) | |
1242 retval = map_s_s (fmod, x, y); | |
1243 } | |
1244 } | |
1245 else if (is_float) | |
1246 { | |
1247 FloatNDArray y = arg_y.float_array_value (); | |
1248 | |
1249 if (! error_state) | |
1250 { | |
1251 FloatNDArray x = arg_x.float_array_value (); | |
1252 | |
1253 if (! error_state) | |
1254 retval = map_fm_fm (fmodf, x, y); | |
1255 } | |
1256 } | |
1257 else | |
1258 { | |
1259 NDArray y = arg_y.array_value (); | |
1260 | |
1261 if (! error_state) | |
1262 { | |
1263 NDArray x = arg_x.array_value (); | |
1264 | |
1265 if (! error_state) | |
1266 retval = map_m_m (fmod, x, y); | |
1267 } | |
1268 } | |
1269 } | |
1270 else | |
1271 error ("fmod: nonconformant matrices"); | |
1272 } | 528 } |
1273 else | 529 else |
1274 print_usage (); | 530 print_usage (); |
1275 | 531 |
1276 return retval; | 532 return retval; |