Mercurial > octave-antonio
annotate scripts/general/curl.m @ 14846:460a3c6d8bf1
maint: Use Octave coding convention for cuddled parenthis in function calls with empty argument lists.
Example: func() => func ()
* dynamic.txi, func.txi, oop.txi, var.txi, embedded.cc, fortdemo.cc,
funcdemo.cc, paramdemo.cc, stringdemo.cc, unwinddemo.cc, Array.cc, Array.h,
CColVector.cc, CDiagMatrix.h, CMatrix.cc, CNDArray.cc, CRowVector.cc,
CSparse.cc, CmplxGEPBAL.cc, EIG.cc, MSparse.cc, MatrixType.cc,
Sparse-op-defs.h, Sparse-perm-op-defs.h, Sparse.cc, Sparse.h,
SparseCmplxCHOL.cc, SparseCmplxCHOL.h, SparseCmplxLU.cc, SparseCmplxQR.cc,
SparseCmplxQR.h, SparseQR.cc, SparseQR.h, SparsedbleCHOL.cc, SparsedbleCHOL.h,
SparsedbleLU.cc, SparsedbleLU.h, base-lu.cc, cmd-hist.cc, dColVector.cc,
dDiagMatrix.h, dMatrix.cc, dNDArray.cc, dRowVector.cc, dSparse.cc, dbleCHOL.cc,
dbleGEPBAL.cc, dim-vector.cc, eigs-base.cc, f2c-main.c, fCColVector.cc,
fCDiagMatrix.h, fCMatrix.cc, fCNDArray.cc, fCRowVector.cc, fCmplxGEPBAL.cc,
fColVector.cc, fDiagMatrix.h, fEIG.cc, fMatrix.cc, fNDArray.cc, fRowVector.cc,
file-ops.cc, file-stat.cc, floatCHOL.cc, floatGEPBAL.cc, idx-vector.h,
lo-specfun.cc, lo-sysdep.cc, mx-inlines.cc, oct-binmap.h, oct-convn.cc,
oct-md5.cc, oct-mem.h, oct-rand.cc, oct-syscalls.cc, randgamma.c, randmtzig.c,
sparse-base-chol.cc, sparse-base-chol.h, sparse-base-lu.cc, sparse-dmsolve.cc,
tempname.c, curl.m, divergence.m, randi.m, dlmwrite.m, edit.m, getappdata.m,
what.m, getarchdir.m, install.m, installed_packages.m, repackage.m,
unload_packages.m, colorbar.m, figure.m, isosurface.m, legend.m, loglog.m,
plot.m, plot3.m, plotyy.m, polar.m, __errplot__.m, __ghostscript__.m,
__marching_cube__.m, __plt__.m, __scatter__.m, semilogx.m, semilogy.m,
trimesh.m, trisurf.m, demo.m, test.m, datetick.m, __delaunayn__.cc,
__dsearchn__.cc, __fltk_uigetfile__.cc, __glpk__.cc, __init_fltk__.cc,
__lin_interpn__.cc, __magick_read__.cc, __pchip_deriv__.cc, balance.cc,
bsxfun.cc, ccolamd.cc, cellfun.cc, chol.cc, daspk.cc, dasrt.cc, dassl.cc,
dmperm.cc, eig.cc, eigs.cc, fftw.cc, filter.cc, find.cc, kron.cc, lookup.cc,
lsode.cc, matrix_type.cc, md5sum.cc, mgorth.cc, qr.cc, quad.cc, rand.cc,
regexp.cc, symbfact.cc, tril.cc, urlwrite.cc, op-bm-bm.cc, op-cdm-cdm.cc,
op-cell.cc, op-chm.cc, op-cm-cm.cc, op-cm-scm.cc, op-cm-sm.cc, op-cs-scm.cc,
op-cs-sm.cc, op-dm-dm.cc, op-dm-scm.cc, op-dm-sm.cc, op-fcdm-fcdm.cc,
op-fcm-fcm.cc, op-fdm-fdm.cc, op-fm-fm.cc, op-int.h, op-m-m.cc, op-m-scm.cc,
op-m-sm.cc, op-pm-pm.cc, op-pm-scm.cc, op-pm-sm.cc, op-range.cc, op-s-scm.cc,
op-s-sm.cc, op-sbm-sbm.cc, op-scm-cm.cc, op-scm-cs.cc, op-scm-m.cc,
op-scm-s.cc, op-scm-scm.cc, op-scm-sm.cc, op-sm-cm.cc, op-sm-cs.cc, op-sm-m.cc,
op-sm-s.cc, op-sm-scm.cc, op-sm-sm.cc, op-str-str.cc, op-struct.cc, bitfcns.cc,
data.cc, debug.cc, dynamic-ld.cc, error.cc, gl-render.cc, graphics.cc,
graphics.in.h, load-path.cc, ls-hdf5.cc, ls-mat5.cc, ls-mat5.h,
ls-oct-ascii.cc, ls-oct-ascii.h, mex.cc, mk-errno-list, oct-map.cc, oct-obj.h,
oct-parse.yy, octave-config.in.cc, ov-base-int.cc, ov-base-mat.cc, ov-base.cc,
ov-bool-mat.cc, ov-bool-sparse.cc, ov-bool.cc, ov-cell.cc, ov-class.cc,
ov-class.h, ov-cx-mat.cc, ov-cx-sparse.cc, ov-fcn-handle.cc, ov-flt-cx-mat.cc,
ov-flt-re-mat.cc, ov-intx.h, ov-range.h, ov-re-mat.cc, ov-re-sparse.cc,
ov-str-mat.cc, ov-struct.cc, ov-usr-fcn.h, ov.h, pr-output.cc, pt-id.cc,
pt-id.h, pt-mat.cc, pt-select.cc, sparse.cc, symtab.cc, symtab.h, syscalls.cc,
toplev.cc, txt-eng-ft.cc, variables.cc, zfstream.cc, zfstream.h, Dork.m,
getStash.m, myStash.m, Gork.m, Pork.m, myStash.m, getStash.m, myStash.m,
getStash.m, myStash.m, fntests.m: Use Octave coding convention for
cuddled parenthis in function calls with empty argument lists.
author | Rik <octave@nomad.inbox5.com> |
---|---|
date | Sun, 08 Jul 2012 11:28:50 -0700 |
parents | f3d52523cde1 |
children | d63878346099 |
rev | line source |
---|---|
14138
72c96de7a403
maint: update copyright notices for 2012
John W. Eaton <jwe@octave.org>
parents:
12520
diff
changeset
|
1 ## Copyright (C) 2009-2012 Kai Habel |
11428 | 2 ## |
3 ## This file is part of Octave. | |
4 ## | |
5 ## Octave is free software; you can redistribute it and/or modify it | |
6 ## under the terms of the GNU General Public License as published by | |
7 ## the Free Software Foundation; either version 3 of the License, or (at | |
8 ## your option) any later version. | |
9 ## | |
10 ## Octave is distributed in the hope that it will be useful, but | |
11 ## WITHOUT ANY WARRANTY; without even the implied warranty of | |
12 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU | |
13 ## General Public License for more details. | |
14 ## | |
15 ## You should have received a copy of the GNU General Public License | |
16 ## along with Octave; see the file COPYING. If not, see | |
17 ## <http://www.gnu.org/licenses/>. | |
18 | |
19 ## -*- texinfo -*- | |
11563
3c6e8aaa9555
Grammarcheck m-files before 3.4 release.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
20 ## @deftypefn {Function File} {[@var{cx}, @var{cy}, @var{cz}, @var{v}] =} curl (@var{x}, @var{y}, @var{z}, @var{fx}, @var{fy}, @var{fz}) |
11428 | 21 ## @deftypefnx {Function File} {[@var{cz}, @var{v}] =} curl (@var{x}, @var{y}, @var{fx}, @var{fy}) |
22 ## @deftypefnx {Function File} {[@dots{}] =} curl (@var{fx}, @var{fy}, @var{fz}) | |
23 ## @deftypefnx {Function File} {[@dots{}] =} curl (@var{fx}, @var{fy}) | |
24 ## @deftypefnx {Function File} {@var{v} =} curl (@dots{}) | |
11563
3c6e8aaa9555
Grammarcheck m-files before 3.4 release.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
25 ## Calculate curl of vector field given by the arrays @var{fx}, @var{fy}, and |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11575
diff
changeset
|
26 ## @var{fz} or @var{fx}, @var{fy} respectively. |
11428 | 27 ## @tex |
28 ## $$ curl F(x,y,z) = \left( {\partial{d} \over \partial{y}} F_z - {\partial{d} \over \partial{z}} F_y, {\partial{d} \over \partial{z}} F_x - {\partial{d} \over \partial{x}} F_z, {\partial{d} \over \partial{x}} F_y - {\partial{d} \over \partial{y}} F_x \right)$$ | |
29 ## @end tex | |
30 ## @ifnottex | |
11563
3c6e8aaa9555
Grammarcheck m-files before 3.4 release.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
31 ## |
11428 | 32 ## @example |
33 ## @group | |
34 ## / d d d d d d \ | |
35 ## curl F(x,y,z) = | -- Fz - -- Fy, -- Fx - -- Fz, -- Fy - -- Fx | | |
36 ## \ dy dz dz dx dx dy / | |
37 ## @end group | |
38 ## @end example | |
11563
3c6e8aaa9555
Grammarcheck m-files before 3.4 release.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
39 ## |
11428 | 40 ## @end ifnottex |
11563
3c6e8aaa9555
Grammarcheck m-files before 3.4 release.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
41 ## The coordinates of the vector field can be given by the arguments @var{x}, |
3c6e8aaa9555
Grammarcheck m-files before 3.4 release.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
42 ## @var{y}, @var{z} or @var{x}, @var{y} respectively. @var{v} calculates the |
3c6e8aaa9555
Grammarcheck m-files before 3.4 release.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
43 ## scalar component of the angular velocity vector in direction of the z-axis |
3c6e8aaa9555
Grammarcheck m-files before 3.4 release.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
44 ## for two-dimensional input. For three-dimensional input the scalar |
11428 | 45 ## rotation is calculated at each grid point in direction of the vector field |
46 ## at that point. | |
12520
ad05e1547398
Add function chop to documentation.
Rik <octave@nomad.inbox5.com>
parents:
12176
diff
changeset
|
47 ## @seealso{divergence, gradient, del2, cross} |
11428 | 48 ## @end deftypefn |
49 | |
50 ## Author: Kai Habel <kai.habel@gmx.de> | |
51 | |
52 function varargout = curl (varargin) | |
53 | |
54 fidx = 1; | |
55 if (nargin == 2) | |
56 sz = size (varargin{fidx}); | |
57 dx = (1:sz(2))(:); | |
58 dy = (1:sz(1))(:); | |
59 elseif (nargin == 3) | |
60 sz = size (varargin{fidx}); | |
61 dx = (1:sz(2))(:); | |
62 dy = (1:sz(1))(:); | |
63 dz = (1:sz(3))(:); | |
64 elseif (nargin == 4) | |
65 fidx = 3; | |
66 dx = varargin{1}(1,:); | |
67 dy = varargin{2}(:,1); | |
68 elseif (nargin == 6) | |
69 fidx = 4; | |
70 dx = varargin{1}(1,:,1)(:); | |
71 dy = varargin{2}(:,1,1)(:); | |
72 dz = varargin{3}(1,1,:)(:); | |
73 else | |
14846
460a3c6d8bf1
maint: Use Octave coding convention for cuddled parenthis in function calls with empty argument lists.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
74 print_usage (); |
11428 | 75 endif |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11575
diff
changeset
|
76 |
11428 | 77 if ((nargin == 4) || (nargin == 2)) |
78 if (!size_equal (varargin{fidx}, varargin{fidx + 1})) | |
11588
d5bd2766c640
style fixes for warning and error messages in script files
John W. Eaton <jwe@octave.org>
parents:
11587
diff
changeset
|
79 error ("curl: size of X and Y must match"); |
11428 | 80 elseif (ndims (varargin{fidx}) != 2) |
11588
d5bd2766c640
style fixes for warning and error messages in script files
John W. Eaton <jwe@octave.org>
parents:
11587
diff
changeset
|
81 error ("curl: expected two-dimensional matrices X and Y"); |
11428 | 82 elseif ((length (dx) != columns (varargin{fidx})) |
83 || (length (dy) != rows (varargin{fidx}))) | |
11472
1740012184f9
Use uppercase for variable names in error() strings to match Info documentation. Only m-files done.
Rik <octave@nomad.inbox5.com>
parents:
11428
diff
changeset
|
84 error ("curl: size of dx and dy must match the respective dimension of X and Y"); |
11428 | 85 endif |
86 | |
87 dFx_dy = gradient (varargin{fidx}.', dy, dx).'; | |
88 dFy_dx = gradient (varargin{fidx + 1}, dx, dy); | |
89 rot_z = dFy_dx - dFx_dy; | |
90 av = rot_z / 2; | |
12173
dd2af7b8dafe
curl: simplify processing of output values
John W. Eaton <jwe@octave.org>
parents:
11589
diff
changeset
|
91 if (nargout == 0 || nargout == 1) |
11428 | 92 varargout{1} = av; |
12176
515446c8fe23
curl.m: fix thinko in previous change
John W. Eaton <jwe@octave.org>
parents:
12173
diff
changeset
|
93 else |
11428 | 94 varargout{1} = rot_z; |
95 varargout{2} = av; | |
96 endif | |
97 | |
98 elseif ((nargin == 6) || (nargin == 3)) | |
99 if (!size_equal (varargin{fidx}, varargin{fidx + 1}, varargin{fidx + 2})) | |
11589
b0084095098e
missing semicolons in script files
John W. Eaton <jwe@octave.org>
parents:
11588
diff
changeset
|
100 error ("curl: size of X, Y, and Z must match"); |
11428 | 101 elseif (ndims (varargin{fidx}) != 3) |
11588
d5bd2766c640
style fixes for warning and error messages in script files
John W. Eaton <jwe@octave.org>
parents:
11587
diff
changeset
|
102 error ("curl: expected two-dimensional matrices X, Y, and Z"); |
11428 | 103 elseif ((length (dx) != size (varargin{fidx}, 2)) |
104 || (length (dy) != size (varargin{fidx}, 1)) | |
105 || (length (dz) != size (varargin{fidx}, 3))) | |
11588
d5bd2766c640
style fixes for warning and error messages in script files
John W. Eaton <jwe@octave.org>
parents:
11587
diff
changeset
|
106 error ("curl: size of dx, dy, and dz must match the respective dimesion of X, Y, and Z"); |
11428 | 107 endif |
108 | |
109 [~, dFx_dy, dFx_dz] = gradient (varargin{fidx}, dx, dy, dz); | |
110 [dFy_dx, ~, dFy_dz] = gradient (varargin{fidx + 1}, dx, dy, dz); | |
111 [dFz_dx, dFz_dy] = gradient (varargin{fidx + 2}, dx, dy, dz); | |
112 rot_x = dFz_dy - dFy_dz; | |
113 rot_y = dFx_dz - dFz_dx; | |
114 rot_z = dFy_dx - dFx_dy; | |
115 l = sqrt(varargin{fidx}.^2 + varargin{fidx + 1}.^2 + varargin{fidx + 2}.^2); | |
116 av = (rot_x .* varargin{fidx} + | |
117 rot_y .* varargin{fidx + 1} + | |
118 rot_z .* varargin{fidx + 2}) ./ (2 * l); | |
119 | |
12173
dd2af7b8dafe
curl: simplify processing of output values
John W. Eaton <jwe@octave.org>
parents:
11589
diff
changeset
|
120 if (nargout == 0 || nargout == 1) |
11428 | 121 varargout{1} = av; |
12173
dd2af7b8dafe
curl: simplify processing of output values
John W. Eaton <jwe@octave.org>
parents:
11589
diff
changeset
|
122 else |
11428 | 123 varargout{1} = rot_x; |
124 varargout{2} = rot_y; | |
125 varargout{3} = rot_z; | |
126 varargout{4} = av; | |
127 endif | |
128 endif | |
129 | |
130 endfunction | |
131 | |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
132 |
11428 | 133 %!test |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
134 %! [X,Y] = meshgrid (-20:20,-22:22); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
135 %! av = curl (2*(X-Y), Y); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
136 %! assert (all (av(:) == 1)); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
137 %! [cz,av] = curl (2*(X-Y), Y); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
138 %! assert (all (cz(:) == 2)); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
139 %! assert (all (av(:) == 1)); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
140 %! [cz,av] = curl (X/2, Y/2, 2*(X-Y), Y); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
141 %! assert (all (cz(:) == 4)); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
142 %! assert (all (av(:) == 2)); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
143 %! assert (size_equal (X,Y,cz,av)); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
144 |