annotate libinterp/corefcn/bsxfun.cc @ 19697:4197fc428c7d

maint: Update copyright notices for 2015.
author John W. Eaton <jwe@octave.org>
date Wed, 11 Feb 2015 14:19:08 -0500
parents b560bac0fca2
children 3fa35defe495
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
1 /*
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
2
19697
4197fc428c7d maint: Update copyright notices for 2015.
John W. Eaton <jwe@octave.org>
parents: 18112
diff changeset
3 Copyright (C) 2007-2015 David Bateman
9783
119d97db51f0 avoid repeated table init in bsxfun
Jaroslav Hajek <highegg@gmail.com>
parents: 9743
diff changeset
4 Copyright (C) 2009 VZLU Prague
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
5
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
6 This file is part of Octave.
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
7
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
8 Octave is free software; you can redistribute it and/or modify it
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
9 under the terms of the GNU General Public License as published by the
7016
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 7001
diff changeset
10 Free Software Foundation; either version 3 of the License, or (at your
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 7001
diff changeset
11 option) any later version.
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
12
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
13 Octave is distributed in the hope that it will be useful, but WITHOUT
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
14 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
15 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
16 for more details.
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
17
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
18 You should have received a copy of the GNU General Public License
7016
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 7001
diff changeset
19 along with Octave; see the file COPYING. If not, see
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 7001
diff changeset
20 <http://www.gnu.org/licenses/>.
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
21
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
22 */
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
23
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
24 #ifdef HAVE_CONFIG_H
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
25 #include <config.h>
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
26 #endif
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
27
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
28 #include <string>
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
29 #include <vector>
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
30 #include <list>
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
31
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
32 #include "lo-mappers.h"
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
33
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
34 #include "oct-map.h"
15039
e753177cde93 maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents: 14854
diff changeset
35 #include "defun.h"
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
36 #include "parse.h"
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
37 #include "variables.h"
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
38 #include "ov-colon.h"
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
39 #include "unwind-prot.h"
9743
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
40 #include "ov-fcn-handle.h"
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
41
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
42 // Optimized bsxfun operations
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
43 enum bsxfun_builtin_op
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
44 {
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
45 bsxfun_builtin_plus = 0,
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
46 bsxfun_builtin_minus,
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
47 bsxfun_builtin_times,
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
48 bsxfun_builtin_divide,
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
49 bsxfun_builtin_max,
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
50 bsxfun_builtin_min,
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
51 bsxfun_builtin_eq,
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
52 bsxfun_builtin_ne,
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
53 bsxfun_builtin_lt,
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
54 bsxfun_builtin_le,
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
55 bsxfun_builtin_gt,
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
56 bsxfun_builtin_ge,
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
57 bsxfun_builtin_and,
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
58 bsxfun_builtin_or,
9827
c15a5ed0da58 optimize bsxfun (@power, ...)
Jaroslav Hajek <highegg@gmail.com>
parents: 9783
diff changeset
59 bsxfun_builtin_power,
9743
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
60 bsxfun_builtin_unknown,
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
61 bsxfun_num_builtin_ops = bsxfun_builtin_unknown
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
62 };
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
63
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
64 const char *bsxfun_builtin_names[] =
9743
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
65 {
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
66 "plus",
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
67 "minus",
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
68 "times",
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
69 "rdivide",
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
70 "max",
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
71 "min",
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
72 "eq",
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
73 "ne",
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
74 "lt",
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
75 "le",
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
76 "gt",
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
77 "ge",
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
78 "and",
9827
c15a5ed0da58 optimize bsxfun (@power, ...)
Jaroslav Hajek <highegg@gmail.com>
parents: 9783
diff changeset
79 "or",
c15a5ed0da58 optimize bsxfun (@power, ...)
Jaroslav Hajek <highegg@gmail.com>
parents: 9783
diff changeset
80 "power"
9743
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
81 };
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
82
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
83 static bsxfun_builtin_op
9743
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
84 bsxfun_builtin_lookup (const std::string& name)
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
85 {
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
86 for (int i = 0; i < bsxfun_num_builtin_ops; i++)
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
87 if (name == bsxfun_builtin_names[i])
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
88 return static_cast<bsxfun_builtin_op> (i);
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
89 return bsxfun_builtin_unknown;
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
90 }
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
91
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
92 typedef octave_value (*bsxfun_handler) (const octave_value&,
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
93 const octave_value&);
9743
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
94
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
95 // Static table of handlers.
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
96 bsxfun_handler bsxfun_handler_table[bsxfun_num_builtin_ops][btyp_num_types];
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
97
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
98 template <class NDA, NDA (bsxfun_op) (const NDA&, const NDA&)>
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
99 static octave_value
9743
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
100 bsxfun_forward_op (const octave_value& x, const octave_value& y)
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
101 {
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
102 NDA xa = octave_value_extract<NDA> (x);
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
103 NDA ya = octave_value_extract<NDA> (y);
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
104 return octave_value (bsxfun_op (xa, ya));
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
105 }
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
106
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
107 template <class NDA, boolNDArray (bsxfun_rel) (const NDA&, const NDA&)>
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
108 static octave_value
9743
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
109 bsxfun_forward_rel (const octave_value& x, const octave_value& y)
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
110 {
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
111 NDA xa = octave_value_extract<NDA> (x);
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
112 NDA ya = octave_value_extract<NDA> (y);
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
113 return octave_value (bsxfun_rel (xa, ya));
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
114 }
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
115
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
116 // pow() needs a special handler for reals
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
117 // because of the potentially complex result.
9827
c15a5ed0da58 optimize bsxfun (@power, ...)
Jaroslav Hajek <highegg@gmail.com>
parents: 9783
diff changeset
118 template <class NDA, class CNDA>
c15a5ed0da58 optimize bsxfun (@power, ...)
Jaroslav Hajek <highegg@gmail.com>
parents: 9783
diff changeset
119 static octave_value
c15a5ed0da58 optimize bsxfun (@power, ...)
Jaroslav Hajek <highegg@gmail.com>
parents: 9783
diff changeset
120 do_bsxfun_real_pow (const octave_value& x, const octave_value& y)
c15a5ed0da58 optimize bsxfun (@power, ...)
Jaroslav Hajek <highegg@gmail.com>
parents: 9783
diff changeset
121 {
c15a5ed0da58 optimize bsxfun (@power, ...)
Jaroslav Hajek <highegg@gmail.com>
parents: 9783
diff changeset
122 NDA xa = octave_value_extract<NDA> (x);
c15a5ed0da58 optimize bsxfun (@power, ...)
Jaroslav Hajek <highegg@gmail.com>
parents: 9783
diff changeset
123 NDA ya = octave_value_extract<NDA> (y);
c15a5ed0da58 optimize bsxfun (@power, ...)
Jaroslav Hajek <highegg@gmail.com>
parents: 9783
diff changeset
124 if (! ya.all_integers () && xa.any_element_is_negative ())
c15a5ed0da58 optimize bsxfun (@power, ...)
Jaroslav Hajek <highegg@gmail.com>
parents: 9783
diff changeset
125 return octave_value (bsxfun_pow (CNDA (xa), ya));
c15a5ed0da58 optimize bsxfun (@power, ...)
Jaroslav Hajek <highegg@gmail.com>
parents: 9783
diff changeset
126 else
c15a5ed0da58 optimize bsxfun (@power, ...)
Jaroslav Hajek <highegg@gmail.com>
parents: 9783
diff changeset
127 return octave_value (bsxfun_pow (xa, ya));
c15a5ed0da58 optimize bsxfun (@power, ...)
Jaroslav Hajek <highegg@gmail.com>
parents: 9783
diff changeset
128 }
c15a5ed0da58 optimize bsxfun (@power, ...)
Jaroslav Hajek <highegg@gmail.com>
parents: 9783
diff changeset
129
9743
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
130 static void maybe_fill_table (void)
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
131 {
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
132 static bool filled = false;
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
133 if (filled)
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
134 return;
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
135
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
136 #define REGISTER_OP_HANDLER(OP, BTYP, NDA, FUNOP) \
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
137 bsxfun_handler_table[OP][BTYP] = bsxfun_forward_op<NDA, FUNOP>
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
138 #define REGISTER_REL_HANDLER(REL, BTYP, NDA, FUNREL) \
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
139 bsxfun_handler_table[REL][BTYP] = bsxfun_forward_rel<NDA, FUNREL>
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
140 #define REGISTER_STD_HANDLERS(BTYP, NDA) \
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
141 REGISTER_OP_HANDLER (bsxfun_builtin_plus, BTYP, NDA, bsxfun_add); \
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
142 REGISTER_OP_HANDLER (bsxfun_builtin_minus, BTYP, NDA, bsxfun_sub); \
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
143 REGISTER_OP_HANDLER (bsxfun_builtin_times, BTYP, NDA, bsxfun_mul); \
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
144 REGISTER_OP_HANDLER (bsxfun_builtin_divide, BTYP, NDA, bsxfun_div); \
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
145 REGISTER_OP_HANDLER (bsxfun_builtin_max, BTYP, NDA, bsxfun_max); \
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
146 REGISTER_OP_HANDLER (bsxfun_builtin_min, BTYP, NDA, bsxfun_min); \
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
147 REGISTER_REL_HANDLER (bsxfun_builtin_eq, BTYP, NDA, bsxfun_eq); \
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
148 REGISTER_REL_HANDLER (bsxfun_builtin_ne, BTYP, NDA, bsxfun_ne); \
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
149 REGISTER_REL_HANDLER (bsxfun_builtin_lt, BTYP, NDA, bsxfun_lt); \
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
150 REGISTER_REL_HANDLER (bsxfun_builtin_le, BTYP, NDA, bsxfun_le); \
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
151 REGISTER_REL_HANDLER (bsxfun_builtin_gt, BTYP, NDA, bsxfun_gt); \
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
152 REGISTER_REL_HANDLER (bsxfun_builtin_ge, BTYP, NDA, bsxfun_ge)
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
153
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
154 REGISTER_STD_HANDLERS (btyp_double, NDArray);
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
155 REGISTER_STD_HANDLERS (btyp_float, FloatNDArray);
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
156 REGISTER_STD_HANDLERS (btyp_complex, ComplexNDArray);
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
157 REGISTER_STD_HANDLERS (btyp_float_complex, FloatComplexNDArray);
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
158 REGISTER_STD_HANDLERS (btyp_int8, int8NDArray);
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
159 REGISTER_STD_HANDLERS (btyp_int16, int16NDArray);
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
160 REGISTER_STD_HANDLERS (btyp_int32, int32NDArray);
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
161 REGISTER_STD_HANDLERS (btyp_int64, int64NDArray);
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
162 REGISTER_STD_HANDLERS (btyp_uint8, uint8NDArray);
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
163 REGISTER_STD_HANDLERS (btyp_uint16, uint16NDArray);
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
164 REGISTER_STD_HANDLERS (btyp_uint32, uint32NDArray);
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
165 REGISTER_STD_HANDLERS (btyp_uint64, uint64NDArray);
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
166
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
167 // For bools, we register and/or.
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
168 REGISTER_OP_HANDLER (bsxfun_builtin_and, btyp_bool, boolNDArray, bsxfun_and);
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
169 REGISTER_OP_HANDLER (bsxfun_builtin_or, btyp_bool, boolNDArray, bsxfun_or);
9783
119d97db51f0 avoid repeated table init in bsxfun
Jaroslav Hajek <highegg@gmail.com>
parents: 9743
diff changeset
170
9827
c15a5ed0da58 optimize bsxfun (@power, ...)
Jaroslav Hajek <highegg@gmail.com>
parents: 9783
diff changeset
171 // Register power handlers.
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
172 bsxfun_handler_table[bsxfun_builtin_power][btyp_double] =
9827
c15a5ed0da58 optimize bsxfun (@power, ...)
Jaroslav Hajek <highegg@gmail.com>
parents: 9783
diff changeset
173 do_bsxfun_real_pow<NDArray, ComplexNDArray>;
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
174 bsxfun_handler_table[bsxfun_builtin_power][btyp_float] =
9827
c15a5ed0da58 optimize bsxfun (@power, ...)
Jaroslav Hajek <highegg@gmail.com>
parents: 9783
diff changeset
175 do_bsxfun_real_pow<FloatNDArray, FloatComplexNDArray>;
c15a5ed0da58 optimize bsxfun (@power, ...)
Jaroslav Hajek <highegg@gmail.com>
parents: 9783
diff changeset
176
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
177 REGISTER_OP_HANDLER (bsxfun_builtin_power, btyp_complex, ComplexNDArray,
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
178 bsxfun_pow);
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
179 REGISTER_OP_HANDLER (bsxfun_builtin_power, btyp_float_complex,
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
180 FloatComplexNDArray, bsxfun_pow);
9827
c15a5ed0da58 optimize bsxfun (@power, ...)
Jaroslav Hajek <highegg@gmail.com>
parents: 9783
diff changeset
181
10107
fd262afea1d1 optimize bsxfun for chars
Jaroslav Hajek <highegg@gmail.com>
parents: 9827
diff changeset
182 // For chars, we want just relational handlers.
fd262afea1d1 optimize bsxfun for chars
Jaroslav Hajek <highegg@gmail.com>
parents: 9827
diff changeset
183 REGISTER_REL_HANDLER (bsxfun_builtin_eq, btyp_char, charNDArray, bsxfun_eq);
fd262afea1d1 optimize bsxfun for chars
Jaroslav Hajek <highegg@gmail.com>
parents: 9827
diff changeset
184 REGISTER_REL_HANDLER (bsxfun_builtin_ne, btyp_char, charNDArray, bsxfun_ne);
fd262afea1d1 optimize bsxfun for chars
Jaroslav Hajek <highegg@gmail.com>
parents: 9827
diff changeset
185 REGISTER_REL_HANDLER (bsxfun_builtin_lt, btyp_char, charNDArray, bsxfun_lt);
fd262afea1d1 optimize bsxfun for chars
Jaroslav Hajek <highegg@gmail.com>
parents: 9827
diff changeset
186 REGISTER_REL_HANDLER (bsxfun_builtin_le, btyp_char, charNDArray, bsxfun_le);
fd262afea1d1 optimize bsxfun for chars
Jaroslav Hajek <highegg@gmail.com>
parents: 9827
diff changeset
187 REGISTER_REL_HANDLER (bsxfun_builtin_gt, btyp_char, charNDArray, bsxfun_gt);
fd262afea1d1 optimize bsxfun for chars
Jaroslav Hajek <highegg@gmail.com>
parents: 9827
diff changeset
188 REGISTER_REL_HANDLER (bsxfun_builtin_ge, btyp_char, charNDArray, bsxfun_ge);
fd262afea1d1 optimize bsxfun for chars
Jaroslav Hajek <highegg@gmail.com>
parents: 9827
diff changeset
189
9783
119d97db51f0 avoid repeated table init in bsxfun
Jaroslav Hajek <highegg@gmail.com>
parents: 9743
diff changeset
190 filled = true;
9743
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
191 }
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
192
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
193 static octave_value
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
194 maybe_optimized_builtin (const std::string& name,
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
195 const octave_value& a, const octave_value& b)
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
196 {
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
197 octave_value retval;
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
198
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
199 maybe_fill_table ();
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
200
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
201 bsxfun_builtin_op op = bsxfun_builtin_lookup (name);
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
202 if (op != bsxfun_builtin_unknown)
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
203 {
18100
6a71e5030df5 Follow coding convention of defining and initializing only 1 variable per line in liboctinterp.
Rik <rik@octave.org>
parents: 17787
diff changeset
204 builtin_type_t btyp_a = a.builtin_type ();
6a71e5030df5 Follow coding convention of defining and initializing only 1 variable per line in liboctinterp.
Rik <rik@octave.org>
parents: 17787
diff changeset
205 builtin_type_t btyp_b = b.builtin_type ();
9743
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
206
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
207 // Simplify single/double combinations.
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
208 if (btyp_a == btyp_float && btyp_b == btyp_double)
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
209 btyp_b = btyp_float;
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
210 else if (btyp_a == btyp_double && btyp_b == btyp_float)
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
211 btyp_a = btyp_float;
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
212 else if (btyp_a == btyp_float_complex && btyp_b == btyp_complex)
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
213 btyp_b = btyp_float_complex;
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
214 else if (btyp_a == btyp_complex && btyp_b == btyp_float_complex)
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
215 btyp_a = btyp_float_complex;
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
216
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
217 if (btyp_a == btyp_b && btyp_a != btyp_unknown)
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
218 {
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
219 bsxfun_handler handler = bsxfun_handler_table[op][btyp_a];
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
220 if (handler)
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
221 retval = handler (a, b);
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
222 }
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
223 }
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
224
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
225 return retval;
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
226 }
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
227
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
228 static bool
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
229 maybe_update_column (octave_value& Ac, const octave_value& A,
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
230 const dim_vector& dva, const dim_vector& dvc,
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
231 octave_idx_type i, octave_value_list &idx)
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
232 {
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
233 octave_idx_type nd = dva.length ();
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
234
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
235 if (i == 0)
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
236 {
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
237 idx(0) = octave_value (':');
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
238 for (octave_idx_type j = 1; j < nd; j++)
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
239 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
240 if (dva (j) == 1)
14854
5ae9f0f77635 maint: Use Octave coding conventions for coddling parenthis is DLD-FUNCTIONS directory
Rik <octave@nomad.inbox5.com>
parents: 14846
diff changeset
241 idx(j) = octave_value (1);
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
242 else
14854
5ae9f0f77635 maint: Use Octave coding conventions for coddling parenthis is DLD-FUNCTIONS directory
Rik <octave@nomad.inbox5.com>
parents: 14846
diff changeset
243 idx(j) = octave_value ((i % dvc(j)) + 1);
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
244
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
245 i = i / dvc (j);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
246 }
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
247
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
248 Ac = A;
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
249 Ac = Ac.single_subsref ("(", idx);
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
250 return true;
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
251 }
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
252 else
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
253 {
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
254 bool is_changed = false;
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
255 octave_idx_type k = i;
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
256 octave_idx_type k1 = i - 1;
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
257 for (octave_idx_type j = 1; j < nd; j++)
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
258 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
259 if (dva(j) != 1 && k % dvc (j) != k1 % dvc (j))
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
260 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
261 idx (j) = octave_value ((k % dvc(j)) + 1);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
262 is_changed = true;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
263 }
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
264
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
265 k = k / dvc (j);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
266 k1 = k1 / dvc (j);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
267 }
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
268
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
269 if (is_changed)
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
270 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
271 Ac = A;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
272 Ac = Ac.single_subsref ("(", idx);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
273 return true;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
274 }
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
275 else
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
276 return false;
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
277 }
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
278 }
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
279
6959
47f4f4e88166 [project @ 2007-10-04 20:43:32 by jwe]
jwe
parents: 6881
diff changeset
280 #if 0
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
281 // FIXME: this function is not used; is it OK to delete it?
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
282 static void
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
283 update_index (octave_value_list& idx, const dim_vector& dv, octave_idx_type i)
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
284 {
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
285 octave_idx_type nd = dv.length ();
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
286
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
287 if (i == 0)
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
288 {
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
289 for (octave_idx_type j = nd - 1; j > 0; j--)
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
290 idx(j) = octave_value (static_cast<double>(1));
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
291 idx(0) = octave_value (':');
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
292 }
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
293 else
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
294 {
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
295 for (octave_idx_type j = 1; j < nd; j++)
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
296 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
297 idx (j) = octave_value (i % dv (j) + 1);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
298 i = i / dv (j);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
299 }
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
300 }
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
301 }
6959
47f4f4e88166 [project @ 2007-10-04 20:43:32 by jwe]
jwe
parents: 6881
diff changeset
302 #endif
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
303
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
304 static void
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
305 update_index (Array<int>& idx, const dim_vector& dv, octave_idx_type i)
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
306 {
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
307 octave_idx_type nd = dv.length ();
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
308
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
309 idx(0) = 0;
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
310 for (octave_idx_type j = 1; j < nd; j++)
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
311 {
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
312 idx (j) = i % dv (j);
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
313 i = i / dv (j);
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
314 }
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
315 }
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
316
15039
e753177cde93 maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents: 14854
diff changeset
317 DEFUN (bsxfun, args, ,
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
318 "-*- texinfo -*-\n\
15039
e753177cde93 maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory
Rik <rik@octave.org>
parents: 14854
diff changeset
319 @deftypefn {Built-in Function} {} bsxfun (@var{f}, @var{A}, @var{B})\n\
14116
951eacaf9381 Initial documentation for broadcasting and general vectorization guidelines
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 13929
diff changeset
320 The binary singleton expansion function applier performs broadcasting,\n\
951eacaf9381 Initial documentation for broadcasting and general vectorization guidelines
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 13929
diff changeset
321 that is, applies a binary function @var{f} element-by-element to two\n\
12972
e4f82a337d66 bsxfun.cc: Expand cryptic bsxfun acronym in docstring and explain it a bit more
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12639
diff changeset
322 array arguments @var{A} and @var{B}, and expands as necessary\n\
13929
9cae456085c2 Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents: 13013
diff changeset
323 singleton dimensions in either input argument. @var{f} is a function\n\
12972
e4f82a337d66 bsxfun.cc: Expand cryptic bsxfun acronym in docstring and explain it a bit more
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12639
diff changeset
324 handle, inline function, or string containing the name of the function\n\
13929
9cae456085c2 Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents: 13013
diff changeset
325 to evaluate. The function @var{f} must be capable of accepting two\n\
12972
e4f82a337d66 bsxfun.cc: Expand cryptic bsxfun acronym in docstring and explain it a bit more
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12639
diff changeset
326 column-vector arguments of equal length, or one column vector argument\n\
e4f82a337d66 bsxfun.cc: Expand cryptic bsxfun acronym in docstring and explain it a bit more
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12639
diff changeset
327 and a scalar.\n\
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
328 \n\
13929
9cae456085c2 Grammarcheck of documentation before 3.6.0 release.
Rik <octave@nomad.inbox5.com>
parents: 13013
diff changeset
329 The dimensions of @var{A} and @var{B} must be equal or singleton. The\n\
12972
e4f82a337d66 bsxfun.cc: Expand cryptic bsxfun acronym in docstring and explain it a bit more
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12639
diff changeset
330 singleton dimensions of the arrays will be expanded to the same\n\
e4f82a337d66 bsxfun.cc: Expand cryptic bsxfun acronym in docstring and explain it a bit more
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12639
diff changeset
331 dimensionality as the other array.\n\
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
332 @seealso{arrayfun, cellfun}\n\
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
333 @end deftypefn")
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
334 {
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
335 int nargin = args.length ();
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
336 octave_value_list retval;
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
337
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
338 if (nargin != 3)
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
339 print_usage ();
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
340 else
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
341 {
9450
cf714e75c656 implement overloaded function handles
Jaroslav Hajek <highegg@gmail.com>
parents: 9041
diff changeset
342 octave_value func = args(0);
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
343
9450
cf714e75c656 implement overloaded function handles
Jaroslav Hajek <highegg@gmail.com>
parents: 9041
diff changeset
344 if (func.is_string ())
cf714e75c656 implement overloaded function handles
Jaroslav Hajek <highegg@gmail.com>
parents: 9041
diff changeset
345 {
cf714e75c656 implement overloaded function handles
Jaroslav Hajek <highegg@gmail.com>
parents: 9041
diff changeset
346 std::string name = func.string_value ();
cf714e75c656 implement overloaded function handles
Jaroslav Hajek <highegg@gmail.com>
parents: 9041
diff changeset
347 func = symbol_table::find_function (name);
cf714e75c656 implement overloaded function handles
Jaroslav Hajek <highegg@gmail.com>
parents: 9041
diff changeset
348 if (func.is_undefined ())
cf714e75c656 implement overloaded function handles
Jaroslav Hajek <highegg@gmail.com>
parents: 9041
diff changeset
349 error ("bsxfun: invalid function name: %s", name.c_str ());
cf714e75c656 implement overloaded function handles
Jaroslav Hajek <highegg@gmail.com>
parents: 9041
diff changeset
350 }
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
351 else if (! (args(0).is_function_handle ()
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
352 || args(0).is_inline_function ()))
11553
01f703952eff Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
353 error ("bsxfun: F must be a string or function handle");
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
354
18112
b560bac0fca2 maint: Don't use space between 'args' and '(' when doing indexing.
Rik <rik@octave.org>
parents: 18100
diff changeset
355 const octave_value A = args(1);
b560bac0fca2 maint: Don't use space between 'args' and '(' when doing indexing.
Rik <rik@octave.org>
parents: 18100
diff changeset
356 const octave_value B = args(2);
9743
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
357
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
358 if (func.is_builtin_function ()
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
359 || (func.is_function_handle ()
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
360 && ! A.is_object () && ! B.is_object ()))
9743
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
361 {
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
362 // This may break if the default behavior is overriden. But if you
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
363 // override arithmetic operators for builtin classes, you should
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
364 // expect mayhem anyway (constant folding etc). Querying
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
365 // is_overloaded() may not be exactly what we need here.
9743
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
366 octave_function *fcn_val = func.function_value ();
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
367 if (fcn_val)
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
368 {
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
369 octave_value tmp = maybe_optimized_builtin (fcn_val->name (),
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
370 A, B);
9743
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
371 if (tmp.is_defined ())
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
372 retval(0) = tmp;
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
373 }
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
374 }
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
375
26abff55f6fe optimize bsxfun for common built-in operations
Jaroslav Hajek <highegg@gmail.com>
parents: 9450
diff changeset
376 if (! error_state && retval.empty ())
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
377 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
378 dim_vector dva = A.dims ();
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
379 octave_idx_type nda = dva.length ();
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
380 dim_vector dvb = B.dims ();
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
381 octave_idx_type ndb = dvb.length ();
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
382 octave_idx_type nd = nda;
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
383
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
384 if (nda > ndb)
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
385 dvb.resize (nda, 1);
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
386 else if (nda < ndb)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
387 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
388 dva.resize (ndb, 1);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
389 nd = ndb;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
390 }
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
391
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
392 for (octave_idx_type i = 0; i < nd; i++)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
393 if (dva (i) != dvb (i) && dva (i) != 1 && dvb (i) != 1)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
394 {
11553
01f703952eff Improve docstrings for functions in DLD-FUNCTIONS directory.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
395 error ("bsxfun: dimensions of A and B must match");
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
396 break;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
397 }
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
398
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
399 if (!error_state)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
400 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
401 // Find the size of the output
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
402 dim_vector dvc;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
403 dvc.resize (nd);
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
404
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
405 for (octave_idx_type i = 0; i < nd; i++)
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
406 dvc (i) = (dva (i) < 1 ? dva (i)
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
407 : (dvb (i) < 1 ? dvb (i)
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
408 : (dva (i) > dvb (i)
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
409 ? dva (i) : dvb (i))));
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
410
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
411 if (dva == dvb || dva.numel () == 1 || dvb.numel () == 1)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
412 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
413 octave_value_list inputs;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
414 inputs (0) = A;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
415 inputs (1) = B;
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
416 retval = func.do_multi_index_op (1, inputs);
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
417 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
418 else if (dvc.numel () < 1)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
419 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
420 octave_value_list inputs;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
421 inputs (0) = A.resize (dvc);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
422 inputs (1) = B.resize (dvc);
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
423 retval = func.do_multi_index_op (1, inputs);
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
424 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
425 else
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
426 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
427 octave_idx_type ncount = 1;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
428 for (octave_idx_type i = 1; i < nd; i++)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
429 ncount *= dvc (i);
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
430
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
431 #define BSXDEF(T) \
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
432 T result_ ## T; \
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
433 bool have_ ## T = false;
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
434
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
435 BSXDEF(NDArray);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
436 BSXDEF(ComplexNDArray);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
437 BSXDEF(FloatNDArray);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
438 BSXDEF(FloatComplexNDArray);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
439 BSXDEF(boolNDArray);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
440 BSXDEF(int8NDArray);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
441 BSXDEF(int16NDArray);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
442 BSXDEF(int32NDArray);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
443 BSXDEF(int64NDArray);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
444 BSXDEF(uint8NDArray);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
445 BSXDEF(uint16NDArray);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
446 BSXDEF(uint32NDArray);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
447 BSXDEF(uint64NDArray);
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
448
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
449 octave_value Ac ;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
450 octave_value_list idxA;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
451 octave_value Bc;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
452 octave_value_list idxB;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
453 octave_value C;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
454 octave_value_list inputs;
14846
460a3c6d8bf1 maint: Use Octave coding convention for cuddled parenthis in function calls with empty argument lists.
Rik <octave@nomad.inbox5.com>
parents: 14844
diff changeset
455 Array<int> ra_idx (dim_vector (dvc.length (), 1), 0);
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
456
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
457
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
458 for (octave_idx_type i = 0; i < ncount; i++)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
459 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
460 if (maybe_update_column (Ac, A, dva, dvc, i, idxA))
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
461 inputs (0) = Ac;
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
462
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
463 if (maybe_update_column (Bc, B, dvb, dvc, i, idxB))
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
464 inputs (1) = Bc;
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
465
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
466 octave_value_list tmp = func.do_multi_index_op (1,
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
467 inputs);
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
468
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
469 if (error_state)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
470 break;
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
471
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
472 #define BSXINIT(T, CLS, EXTRACTOR) \
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
473 (result_type == CLS) \
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
474 { \
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
475 have_ ## T = true; \
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
476 result_ ## T = \
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
477 tmp (0). EXTRACTOR ## _array_value (); \
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
478 result_ ## T .resize (dvc); \
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
479 }
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
480
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
481 if (i == 0)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
482 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
483 if (! tmp(0).is_sparse_type ())
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
484 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
485 std::string result_type = tmp(0).class_name ();
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
486 if (result_type == "double")
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
487 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
488 if (tmp(0).is_real_type ())
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
489 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
490 have_NDArray = true;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
491 result_NDArray = tmp(0).array_value ();
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
492 result_NDArray.resize (dvc);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
493 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
494 else
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
495 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
496 have_ComplexNDArray = true;
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
497 result_ComplexNDArray =
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
498 tmp(0).complex_array_value ();
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
499 result_ComplexNDArray.resize (dvc);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
500 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
501 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
502 else if (result_type == "single")
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
503 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
504 if (tmp(0).is_real_type ())
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
505 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
506 have_FloatNDArray = true;
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
507 result_FloatNDArray
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
508 = tmp(0).float_array_value ();
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
509 result_FloatNDArray.resize (dvc);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
510 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
511 else
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
512 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
513 have_ComplexNDArray = true;
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
514 result_ComplexNDArray =
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
515 tmp(0).complex_array_value ();
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
516 result_ComplexNDArray.resize (dvc);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
517 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
518 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
519 else if BSXINIT(boolNDArray, "logical", bool)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
520 else if BSXINIT(int8NDArray, "int8", int8)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
521 else if BSXINIT(int16NDArray, "int16", int16)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
522 else if BSXINIT(int32NDArray, "int32", int32)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
523 else if BSXINIT(int64NDArray, "int64", int64)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
524 else if BSXINIT(uint8NDArray, "uint8", uint8)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
525 else if BSXINIT(uint16NDArray, "uint16", uint16)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
526 else if BSXINIT(uint32NDArray, "uint32", uint32)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
527 else if BSXINIT(uint64NDArray, "uint64", uint64)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
528 else
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
529 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
530 C = tmp (0);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
531 C = C.resize (dvc);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
532 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
533 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
534 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
535 else
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
536 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
537 update_index (ra_idx, dvc, i);
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
538
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
539 if (have_FloatNDArray ||
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
540 have_FloatComplexNDArray)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
541 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
542 if (! tmp(0).is_float_type ())
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
543 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
544 if (have_FloatNDArray)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
545 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
546 have_FloatNDArray = false;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
547 C = result_FloatNDArray;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
548 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
549 else
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
550 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
551 have_FloatComplexNDArray = false;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
552 C = result_FloatComplexNDArray;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
553 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
554 C = do_cat_op (C, tmp(0), ra_idx);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
555 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
556 else if (tmp(0).is_double_type ())
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
557 {
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
558 if (tmp(0).is_complex_type () &&
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
559 have_FloatNDArray)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
560 {
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
561 result_ComplexNDArray =
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
562 ComplexNDArray (result_FloatNDArray);
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
563 result_ComplexNDArray.insert
14846
460a3c6d8bf1 maint: Use Octave coding convention for cuddled parenthis in function calls with empty argument lists.
Rik <octave@nomad.inbox5.com>
parents: 14844
diff changeset
564 (tmp(0).complex_array_value (), ra_idx);
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
565 have_FloatComplexNDArray = false;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
566 have_ComplexNDArray = true;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
567 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
568 else
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
569 {
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
570 result_NDArray =
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
571 NDArray (result_FloatNDArray);
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
572 result_NDArray.insert
14846
460a3c6d8bf1 maint: Use Octave coding convention for cuddled parenthis in function calls with empty argument lists.
Rik <octave@nomad.inbox5.com>
parents: 14844
diff changeset
573 (tmp(0).array_value (), ra_idx);
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
574 have_FloatNDArray = false;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
575 have_NDArray = true;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
576 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
577 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
578 else if (tmp(0).is_real_type ())
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
579 result_FloatNDArray.insert
14846
460a3c6d8bf1 maint: Use Octave coding convention for cuddled parenthis in function calls with empty argument lists.
Rik <octave@nomad.inbox5.com>
parents: 14844
diff changeset
580 (tmp(0).float_array_value (), ra_idx);
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
581 else
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
582 {
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
583 result_FloatComplexNDArray =
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
584 FloatComplexNDArray (result_FloatNDArray);
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
585 result_FloatComplexNDArray.insert
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
586 (tmp(0).float_complex_array_value (),
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
587 ra_idx);
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
588 have_FloatNDArray = false;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
589 have_FloatComplexNDArray = true;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
590 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
591 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
592 else if (have_NDArray)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
593 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
594 if (! tmp(0).is_float_type ())
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
595 {
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
596 have_NDArray = false;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
597 C = result_NDArray;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
598 C = do_cat_op (C, tmp(0), ra_idx);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
599 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
600 else if (tmp(0).is_real_type ())
14846
460a3c6d8bf1 maint: Use Octave coding convention for cuddled parenthis in function calls with empty argument lists.
Rik <octave@nomad.inbox5.com>
parents: 14844
diff changeset
601 result_NDArray.insert (tmp(0).array_value (),
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
602 ra_idx);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
603 else
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
604 {
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
605 result_ComplexNDArray =
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
606 ComplexNDArray (result_NDArray);
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
607 result_ComplexNDArray.insert
14846
460a3c6d8bf1 maint: Use Octave coding convention for cuddled parenthis in function calls with empty argument lists.
Rik <octave@nomad.inbox5.com>
parents: 14844
diff changeset
608 (tmp(0).complex_array_value (), ra_idx);
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
609 have_NDArray = false;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
610 have_ComplexNDArray = true;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
611 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
612 }
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
613
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
614 #define BSXLOOP(T, CLS, EXTRACTOR) \
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
615 (have_ ## T) \
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
616 { \
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
617 if (tmp (0).class_name () != CLS) \
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
618 { \
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
619 have_ ## T = false; \
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
620 C = result_ ## T; \
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
621 C = do_cat_op (C, tmp (0), ra_idx); \
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
622 } \
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
623 else \
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
624 result_ ## T .insert \
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
625 (tmp(0). EXTRACTOR ## _array_value (), \
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
626 ra_idx); \
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
627 }
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
628
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
629 else if BSXLOOP(ComplexNDArray, "double", complex)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
630 else if BSXLOOP(boolNDArray, "logical", bool)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
631 else if BSXLOOP(int8NDArray, "int8", int8)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
632 else if BSXLOOP(int16NDArray, "int16", int16)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
633 else if BSXLOOP(int32NDArray, "int32", int32)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
634 else if BSXLOOP(int64NDArray, "int64", int64)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
635 else if BSXLOOP(uint8NDArray, "uint8", uint8)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
636 else if BSXLOOP(uint16NDArray, "uint16", uint16)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
637 else if BSXLOOP(uint32NDArray, "uint32", uint32)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
638 else if BSXLOOP(uint64NDArray, "uint64", uint64)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
639 else
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
640 C = do_cat_op (C, tmp(0), ra_idx);
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
641 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
642 }
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
643
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
644 #define BSXEND(T) \
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
645 (have_ ## T) \
14844
5bc9b9cb4362 maint: Use Octave coding conventions for cuddled parenthesis in retval assignments.
Rik <octave@nomad.inbox5.com>
parents: 14501
diff changeset
646 retval(0) = result_ ## T;
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
647
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
648 if BSXEND(NDArray)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
649 else if BSXEND(ComplexNDArray)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
650 else if BSXEND(FloatNDArray)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
651 else if BSXEND(FloatComplexNDArray)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
652 else if BSXEND(boolNDArray)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
653 else if BSXEND(int8NDArray)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
654 else if BSXEND(int16NDArray)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
655 else if BSXEND(int32NDArray)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
656 else if BSXEND(int64NDArray)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
657 else if BSXEND(uint8NDArray)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
658 else if BSXEND(uint16NDArray)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
659 else if BSXEND(uint32NDArray)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
660 else if BSXEND(uint64NDArray)
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
661 else
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
662 retval(0) = C;
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
663 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
664 }
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 10141
diff changeset
665 }
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11570
diff changeset
666 }
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
667
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
668 return retval;
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
669 }
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
670
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
671 /*
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
672
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
673 %!shared a, b, c, f
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
674 %! a = randn (4, 4);
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
675 %! b = mean (a, 1);
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
676 %! c = mean (a, 2);
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
677 %! f = @minus;
14501
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
678 %!error (bsxfun (f))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
679 %!error (bsxfun (f, a))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
680 %!error (bsxfun (a, b))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
681 %!error (bsxfun (a, b, c))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
682 %!error (bsxfun (f, a, b, c))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
683 %!error (bsxfun (f, ones (4, 0), ones (4, 4)))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
684 %!assert (bsxfun (f, ones (4, 0), ones (4, 1)), zeros (4, 0))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
685 %!assert (bsxfun (f, ones (1, 4), ones (4, 1)), zeros (4, 4))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
686 %!assert (bsxfun (f, a, b), a - repmat (b, 4, 1))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
687 %!assert (bsxfun (f, a, c), a - repmat (c, 1, 4))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
688 %!assert (bsxfun ("minus", ones (1, 4), ones (4, 1)), zeros (4, 4))
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
689
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
690 %!shared a, b, c, f
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
691 %! a = randn (4, 4);
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
692 %! a(1) *= 1i;
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
693 %! b = mean (a, 1);
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
694 %! c = mean (a, 2);
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
695 %! f = @minus;
14501
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
696 %!error (bsxfun (f))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
697 %!error (bsxfun (f, a))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
698 %!error (bsxfun (a, b))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
699 %!error (bsxfun (a, b, c))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
700 %!error (bsxfun (f, a, b, c))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
701 %!error (bsxfun (f, ones (4, 0), ones (4, 4)))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
702 %!assert (bsxfun (f, ones (4, 0), ones (4, 1)), zeros (4, 0))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
703 %!assert (bsxfun (f, ones (1, 4), ones (4, 1)), zeros (4, 4))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
704 %!assert (bsxfun (f, a, b), a - repmat (b, 4, 1))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
705 %!assert (bsxfun (f, a, c), a - repmat (c, 1, 4))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
706 %!assert (bsxfun ("minus", ones (1, 4), ones (4, 1)), zeros (4, 4))
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
707
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
708 %!shared a, b, c, f
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
709 %! a = randn (4, 4);
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
710 %! a(end) *= 1i;
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
711 %! b = mean (a, 1);
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
712 %! c = mean (a, 2);
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
713 %! f = @minus;
14501
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
714 %!error (bsxfun (f))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
715 %!error (bsxfun (f, a))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
716 %!error (bsxfun (a, b))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
717 %!error (bsxfun (a, b, c))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
718 %!error (bsxfun (f, a, b, c))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
719 %!error (bsxfun (f, ones (4, 0), ones (4, 4)))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
720 %!assert (bsxfun (f, ones (4, 0), ones (4, 1)), zeros (4, 0))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
721 %!assert (bsxfun (f, ones (1, 4), ones (4, 1)), zeros (4, 4))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
722 %!assert (bsxfun (f, a, b), a - repmat (b, 4, 1))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
723 %!assert (bsxfun (f, a, c), a - repmat (c, 1, 4))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
724 %!assert (bsxfun ("minus", ones (1, 4), ones (4, 1)), zeros (4, 4))
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
725
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
726 %!shared a, b, c, f
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
727 %! a = randn (4, 4);
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
728 %! b = a (1, :);
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
729 %! c = a (:, 1);
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
730 %! f = @(x, y) x == y;
14501
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
731 %!error (bsxfun (f))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
732 %!error (bsxfun (f, a))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
733 %!error (bsxfun (a, b))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
734 %!error (bsxfun (a, b, c))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
735 %!error (bsxfun (f, a, b, c))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
736 %!error (bsxfun (f, ones (4, 0), ones (4, 4)))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
737 %!assert (bsxfun (f, ones (4, 0), ones (4, 1)), zeros (4, 0, "logical"))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
738 %!assert (bsxfun (f, ones (1, 4), ones (4, 1)), ones (4, 4, "logical"))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
739 %!assert (bsxfun (f, a, b), a == repmat (b, 4, 1))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
740 %!assert (bsxfun (f, a, c), a == repmat (c, 1, 4))
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
741
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
742 %!shared a, b, c, d, f
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
743 %! a = randn (4, 4, 4);
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
744 %! b = mean (a, 1);
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
745 %! c = mean (a, 2);
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
746 %! d = mean (a, 3);
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
747 %! f = @minus;
14501
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
748 %!error (bsxfun (f, ones ([4, 0, 4]), ones ([4, 4, 4])))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
749 %!assert (bsxfun (f, ones ([4, 0, 4]), ones ([4, 1, 4])), zeros ([4, 0, 4]))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
750 %!assert (bsxfun (f, ones ([4, 4, 0]), ones ([4, 1, 1])), zeros ([4, 4, 0]))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
751 %!assert (bsxfun (f, ones ([1, 4, 4]), ones ([4, 1, 4])), zeros ([4, 4, 4]))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
752 %!assert (bsxfun (f, ones ([4, 4, 1]), ones ([4, 1, 4])), zeros ([4, 4, 4]))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
753 %!assert (bsxfun (f, ones ([4, 1, 4]), ones ([1, 4, 4])), zeros ([4, 4, 4]))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
754 %!assert (bsxfun (f, ones ([4, 1, 4]), ones ([1, 4, 1])), zeros ([4, 4, 4]))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
755 %!assert (bsxfun (f, a, b), a - repmat (b, [4, 1, 1]))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
756 %!assert (bsxfun (f, a, c), a - repmat (c, [1, 4, 1]))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
757 %!assert (bsxfun (f, a, d), a - repmat (d, [1, 1, 4]))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
758 %!assert (bsxfun ("minus", ones ([4, 0, 4]), ones ([4, 1, 4])), zeros ([4, 0, 4]))
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
759
14501
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
760 %% The test below is a very hard case to treat
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
761 %!assert (bsxfun (f, ones ([4, 1, 4, 1]), ones ([1, 4, 1, 4])), zeros ([4, 4, 4, 4]));
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
762
10141
e409546ac0a8 add bsxfun tests
Jaroslav Hajek <highegg@gmail.com>
parents: 10107
diff changeset
763 %!shared a, b, aa, bb
e409546ac0a8 add bsxfun tests
Jaroslav Hajek <highegg@gmail.com>
parents: 10107
diff changeset
764 %! a = randn (3, 1, 3);
e409546ac0a8 add bsxfun tests
Jaroslav Hajek <highegg@gmail.com>
parents: 10107
diff changeset
765 %! aa = a(:, ones (1, 3), :, ones (1, 3));
e409546ac0a8 add bsxfun tests
Jaroslav Hajek <highegg@gmail.com>
parents: 10107
diff changeset
766 %! b = randn (1, 3, 3, 3);
e409546ac0a8 add bsxfun tests
Jaroslav Hajek <highegg@gmail.com>
parents: 10107
diff changeset
767 %! bb = b(ones (1, 3), :, :, :);
14501
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
768 %!assert (bsxfun (@plus, a, b), aa + bb)
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
769 %!assert (bsxfun (@minus, a, b), aa - bb)
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
770 %!assert (bsxfun (@times, a, b), aa .* bb)
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
771 %!assert (bsxfun (@rdivide, a, b), aa ./ bb)
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
772 %!assert (bsxfun (@ldivide, a, b), aa .\ bb)
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
773 %!assert (bsxfun (@power, a, b), aa .^ bb)
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
774 %!assert (bsxfun (@power, abs (a), b), abs (aa) .^ bb)
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
775 %!assert (bsxfun (@eq, round (a), round (b)), round (aa) == round (bb))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
776 %!assert (bsxfun (@ne, round (a), round (b)), round (aa) != round (bb))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
777 %!assert (bsxfun (@lt, a, b), aa < bb)
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
778 %!assert (bsxfun (@le, a, b), aa <= bb)
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
779 %!assert (bsxfun (@gt, a, b), aa > bb)
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
780 %!assert (bsxfun (@ge, a, b), aa >= bb)
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
781 %!assert (bsxfun (@min, a, b), min (aa, bb))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
782 %!assert (bsxfun (@max, a, b), max (aa, bb))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
783 %!assert (bsxfun (@and, a > 0, b > 0), (aa > 0) & (bb > 0))
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
784 %!assert (bsxfun (@or, a > 0, b > 0), (aa > 0) | (bb > 0))
13013
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
785
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
786 %% Test automatic bsxfun
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
787 %
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
788 %!test
17245
7babcdb9bc13 Use ... instead of \ for line continuation marker.
Stefan Mahr <dac922@gmx.de>
parents: 15195
diff changeset
789 %! funs = {@plus, @minus, @times, @rdivide, @ldivide, @power, @max, @min, ...
7babcdb9bc13 Use ... instead of \ for line continuation marker.
Stefan Mahr <dac922@gmx.de>
parents: 15195
diff changeset
790 %! @rem, @mod, @atan2, @hypot, @eq, @ne, @lt, @le, @gt, @ge, ...
13013
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
791 %! @and, @or, @xor };
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
792 %!
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
793 %! float_types = {@single, @double};
17245
7babcdb9bc13 Use ... instead of \ for line continuation marker.
Stefan Mahr <dac922@gmx.de>
parents: 15195
diff changeset
794 %! int_types = {@int8, @int16, @int32, @int64, ...
14501
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
795 %! @uint8, @uint16, @uint32, @uint64};
13013
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
796 %!
14501
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
797 %! x = rand (3) * 10-5;
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
798 %! y = rand (3,1) * 10-5;
13013
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
799 %!
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
800 %! for i=1:length (funs)
14854
5ae9f0f77635 maint: Use Octave coding conventions for coddling parenthis is DLD-FUNCTIONS directory
Rik <octave@nomad.inbox5.com>
parents: 14846
diff changeset
801 %! for j = 1:length (float_types)
5ae9f0f77635 maint: Use Octave coding conventions for coddling parenthis is DLD-FUNCTIONS directory
Rik <octave@nomad.inbox5.com>
parents: 14846
diff changeset
802 %! for k = 1:length (int_types)
13013
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
803 %!
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
804 %! fun = funs{i};
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
805 %! f_type = float_types{j};
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
806 %! i_type = int_types{k};
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
807 %!
17245
7babcdb9bc13 Use ... instead of \ for line continuation marker.
Stefan Mahr <dac922@gmx.de>
parents: 15195
diff changeset
808 %! assert (bsxfun (fun, f_type (x), i_type (y)), ...
13013
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
809 %! fun (f_type(x), i_type (y)));
17245
7babcdb9bc13 Use ... instead of \ for line continuation marker.
Stefan Mahr <dac922@gmx.de>
parents: 15195
diff changeset
810 %! assert (bsxfun (fun, f_type (y), i_type (x)), ...
13013
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
811 %! fun (f_type(y), i_type (x)));
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
812 %!
17245
7babcdb9bc13 Use ... instead of \ for line continuation marker.
Stefan Mahr <dac922@gmx.de>
parents: 15195
diff changeset
813 %! assert (bsxfun (fun, i_type (x), i_type (y)), ...
13013
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
814 %! fun (i_type (x), i_type (y)));
17245
7babcdb9bc13 Use ... instead of \ for line continuation marker.
Stefan Mahr <dac922@gmx.de>
parents: 15195
diff changeset
815 %! assert (bsxfun (fun, i_type (y), i_type (x)), ...
13013
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
816 %! fun (i_type (y), i_type (x)));
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
817 %!
17245
7babcdb9bc13 Use ... instead of \ for line continuation marker.
Stefan Mahr <dac922@gmx.de>
parents: 15195
diff changeset
818 %! assert (bsxfun (fun, f_type (x), f_type (y)), ...
13013
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
819 %! fun (f_type (x), f_type (y)));
17245
7babcdb9bc13 Use ... instead of \ for line continuation marker.
Stefan Mahr <dac922@gmx.de>
parents: 15195
diff changeset
820 %! assert (bsxfun (fun, f_type(y), f_type(x)), ...
13013
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
821 %! fun (f_type (y), f_type (x)));
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
822 %! endfor
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
823 %! endfor
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
824 %! endfor
5307bf69cd88 Add tests for automatic bsxfun
Jordi Gutiérrez Hermoso <jordigh@gmail.com>
parents: 12972
diff changeset
825 %!
6869
f9c893831e68 [project @ 2007-09-06 16:38:44 by dbateman]
dbateman
parents:
diff changeset
826 */