Mercurial > hg > octave-avbm
view scripts/general/cplxpair.m @ 14868:5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
* lin2mu.m, loadaudio.m, wavread.m, accumarray.m, bicubic.m, celldisp.m,
colon.m, cplxpair.m, dblquad.m, divergence.m, genvarname.m, gradient.m,
int2str.m, interp1.m, interp1q.m, interp2.m, interpn.m, loadobj.m, nthargout.m,
__isequal__.m, __splinen__.m, quadgk.m, quadl.m, quadv.m, rat.m, rot90.m,
rotdim.m, saveobj.m, subsindex.m, triplequad.m, delaunay3.m, griddata.m,
inpolygon.m, tsearchn.m, voronoi.m, get_first_help_sentence.m, which.m,
gray2ind.m, pink.m, dlmwrite.m, strread.m, textread.m, textscan.m, housh.m,
ishermitian.m, issymmetric.m, krylov.m, logm.m, null.m, rref.m,
compare_versions.m, copyfile.m, dump_prefs.m, edit.m, fileparts.m,
getappdata.m, isappdata.m, movefile.m, orderfields.m, parseparams.m,
__xzip__.m, rmappdata.m, setappdata.m, swapbytes.m, unpack.m, ver.m, fminbnd.m,
fminunc.m, fsolve.m, glpk.m, lsqnonneg.m, qp.m, sqp.m, configure_make.m,
copy_files.m, describe.m, get_description.m, get_forge_pkg.m, install.m,
installed_packages.m, is_architecture_dependent.m, load_package_dirs.m,
print_package_description.m, rebuild.m, repackage.m, save_order.m, shell.m,
allchild.m, ancestor.m, area.m, axes.m, axis.m, clabel.m, close.m, colorbar.m,
comet.m, comet3.m, contour.m, cylinder.m, ezmesh.m, ezsurf.m, findobj.m,
fplot.m, hist.m, isocolors.m, isonormals.m, isosurface.m, isprop.m, legend.m,
mesh.m, meshz.m, pareto.m, pcolor.m, peaks.m, plot3.m, plotmatrix.m, plotyy.m,
polar.m, print.m, __add_datasource__.m, __add_default_menu__.m,
__axes_limits__.m, __bar__.m, __clabel__.m, __contour__.m, __errcomm__.m,
__errplot__.m, __ezplot__.m, __file_filter__.m, __fltk_print__.m,
__ghostscript__.m, __gnuplot_print__.m, __go_draw_axes__.m,
__go_draw_figure__.m, __interp_cube__.m, __marching_cube__.m, __patch__.m,
__pie__.m, __plt__.m, __print_parse_opts__.m, __quiver__.m, __scatter__.m,
__stem__.m, __tight_eps_bbox__.m, __uigetdir_fltk__.m, __uigetfile_fltk__.m,
__uiputfile_fltk__.m, quiver.m, quiver3.m, rectangle.m, refreshdata.m,
ribbon.m, scatter.m, semilogy.m, shading.m, slice.m, subplot.m, surface.m,
surfl.m, surfnorm.m, text.m, uigetfile.m, uiputfile.m, whitebg.m, deconv.m,
mkpp.m, pchip.m, polyaffine.m, polyder.m, polygcd.m, polyout.m, polyval.m,
ppint.m, ppjumps.m, ppval.m, residue.m, roots.m, spline.m, splinefit.m,
addpref.m, getpref.m, setpref.m, ismember.m, setxor.m, arch_fit.m, arch_rnd.m,
arch_test.m, autoreg_matrix.m, diffpara.m, fftconv.m, filter2.m, hanning.m,
hurst.m, periodogram.m, triangle_sw.m, sinc.m, spectral_xdf.m, spencer.m,
stft.m, synthesis.m, unwrap.m, yulewalker.m, bicgstab.m, gmres.m, pcg.m, pcr.m,
__sprand_impl__.m, speye.m, spfun.m, sprandn.m, spstats.m, svds.m,
treelayout.m, treeplot.m, bessel.m, factor.m, legendre.m, perms.m, primes.m,
magic.m, toeplitz.m, corr.m, cov.m, mean.m, median.m, mode.m, qqplot.m,
quantile.m, ranks.m, zscore.m, logistic_regression_likelihood.m,
bartlett_test.m, chisquare_test_homogeneity.m, chisquare_test_independence.m,
kolmogorov_smirnov_test.m, run_test.m, u_test.m, wilcoxon_test.m, z_test.m,
z_test_2.m, bin2dec.m, dec2base.m, mat2str.m, strcat.m, strchr.m, strjust.m,
strtok.m, substr.m, untabify.m, assert.m, demo.m, example.m, fail.m, speed.m,
test.m, now.m: Use Octave coding conventions for cuddling parentheses in
scripts directory.
author | Rik <octave@nomad.inbox5.com> |
---|---|
date | Tue, 17 Jul 2012 07:08:39 -0700 |
parents | 86854d032a37 |
children |
line wrap: on
line source
## Copyright (C) 2000-2012 Paul Kienzle ## ## This file is part of Octave. ## ## Octave is free software; you can redistribute it and/or modify it ## under the terms of the GNU General Public License as published by ## the Free Software Foundation; either version 3 of the License, or (at ## your option) any later version. ## ## Octave is distributed in the hope that it will be useful, but ## WITHOUT ANY WARRANTY; without even the implied warranty of ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU ## General Public License for more details. ## ## You should have received a copy of the GNU General Public License ## along with Octave; see the file COPYING. If not, see ## <http://www.gnu.org/licenses/>. ## -*- texinfo -*- ## @deftypefn {Function File} {} cplxpair (@var{z}) ## @deftypefnx {Function File} {} cplxpair (@var{z}, @var{tol}) ## @deftypefnx {Function File} {} cplxpair (@var{z}, @var{tol}, @var{dim}) ## Sort the numbers @var{z} into complex conjugate pairs ordered by ## increasing real part. Place the negative imaginary complex number ## first within each pair. Place all the real numbers (those with ## @code{abs (imag (@var{z}) / @var{z}) < @var{tol})}) after the ## complex pairs. ## ## If @var{tol} is unspecified the default value is 100*@code{eps}. ## ## By default the complex pairs are sorted along the first non-singleton ## dimension of @var{z}. If @var{dim} is specified, then the complex ## pairs are sorted along this dimension. ## ## Signal an error if some complex numbers could not be paired. Signal an ## error if all complex numbers are not exact conjugates (to within ## @var{tol}). Note that there is no defined order for pairs with identical ## real parts but differing imaginary parts. ## @c Set example in small font to prevent overfull line ## ## @smallexample ## cplxpair (exp(2i*pi*[0:4]'/5)) == exp(2i*pi*[3; 2; 4; 1; 0]/5) ## @end smallexample ## @end deftypefn ## FIXME: subsort returned pairs by imaginary magnitude ## FIXME: Why doesn't exp (2i*pi*[0:4]'/5) produce exact conjugates. Does ## FIXME: it in Matlab? The reason is that complex pairs are supposed ## FIXME: to be exact conjugates, and not rely on a tolerance test. ## 2006-05-12 David Bateman - Modified for NDArrays function y = cplxpair (z, tol, dim) if (nargin < 1 || nargin > 3) print_usage (); endif if (length (z) == 0) y = zeros (size (z)); return; endif if (nargin < 2 || isempty (tol)) if (isa (z, "single")) tol = 100 * eps("single"); else tol = 100*eps; endif endif nd = ndims (z); orig_dims = size (z); if (nargin < 3) ## Find the first singleton dimension. dim = 0; while (dim < nd && orig_dims(dim+1) == 1) dim++; endwhile dim++; if (dim > nd) dim = 1; endif else dim = floor (dim); if (dim < 1 || dim > nd) error ("cplxpair: invalid dimension along which to sort"); endif endif ## Move dimension to treat first, and convert to a 2-D matrix. perm = [dim:nd, 1:dim-1]; z = permute (z, perm); sz = size (z); n = sz (1); m = prod (sz) / n; z = reshape (z, n, m); ## Sort the sequence in terms of increasing real values. [q, idx] = sort (real (z), 1); z = z(idx + n * ones (n, 1) * [0:m-1]); ## Put the purely real values at the end of the returned list. cls = "double"; if (isa (z, "single")) cls = "single"; endif [idxi, idxj] = find (abs (imag (z)) ./ (abs (z) + realmin (cls)) < tol); q = sparse (idxi, idxj, 1, n, m); nr = sum (q, 1); [q, idx] = sort (q, 1); z = z(idx); y = z; ## For each remaining z, place the value and its conjugate at the ## start of the returned list, and remove them from further ## consideration. for j = 1:m p = n - nr(j); for i = 1:2:p if (i+1 > p) error ("cplxpair: could not pair all complex numbers"); endif [v, idx] = min (abs (z(i+1:p) - conj (z(i)))); if (v > tol) error ("cplxpair: could not pair all complex numbers"); endif if (imag (z(i)) < 0) y([i, i+1]) = z([i, idx+i]); else y([i, i+1]) = z([idx+i, i]); endif z(idx+i) = z(i+1); endfor endfor ## Reshape the output matrix. y = ipermute (reshape (y, sz), perm); endfunction %!demo %! [ cplxpair(exp(2i*pi*[0:4]'/5)), exp(2i*pi*[3; 2; 4; 1; 0]/5) ] %!assert (isempty (cplxpair ([]))) %!assert (cplxpair (1), 1) %!assert (cplxpair ([1+1i, 1-1i]), [1-1i, 1+1i]) %!assert (cplxpair ([1+1i, 1+1i, 1, 1-1i, 1-1i, 2]), ... %! [1-1i, 1+1i, 1-1i, 1+1i, 1, 2]) %!assert (cplxpair ([1+1i; 1+1i; 1; 1-1i; 1-1i; 2]), ... %! [1-1i; 1+1i; 1-1i; 1+1i; 1; 2]) %!assert (cplxpair ([0, 1, 2]), [0, 1, 2]) %!shared z %! z = exp (2i*pi*[4; 3; 5; 2; 6; 1; 0]/7); %!assert (cplxpair (z(randperm (7))), z) %!assert (cplxpair (z(randperm (7))), z) %!assert (cplxpair (z(randperm (7))), z) %!assert (cplxpair ([z(randperm(7)),z(randperm(7))]), [z,z]) %!assert (cplxpair ([z(randperm(7)),z(randperm(7))],[],1), [z,z]) %!assert (cplxpair ([z(randperm(7)).';z(randperm(7)).'],[],2), [z.';z.']) %!## tolerance test %!assert (cplxpair ([1i, -1i, 1+(1i*eps)],2*eps), [-1i, 1i, 1+(1i*eps)])