view libinterp/corefcn/__lin_interpn__.cc @ 17289:bc924baa2c4e

doc: Add new @qcode macro for code samples which are quoted. Macro handles options ("on") or properties ("position") more elegantly than @code{"text"}. * doc/interpreter/macros.texi: Add new @qcode macro. * doc/interpreter/tips.txi: Add documentation about @qcode macro. * doc/interpreter/basics.txi, doc/interpreter/container.txi, doc/interpreter/emacs.txi, doc/interpreter/errors.txi, doc/interpreter/eval.txi, doc/interpreter/expr.txi, doc/interpreter/external.txi, doc/interpreter/func.txi, doc/interpreter/grammar.txi, doc/interpreter/image.txi, doc/interpreter/install.txi, doc/interpreter/interp.txi, doc/interpreter/io.txi, doc/interpreter/matrix.txi, doc/interpreter/numbers.txi, doc/interpreter/oop.txi, doc/interpreter/package.txi, doc/interpreter/plot.txi, doc/interpreter/quad.txi, doc/interpreter/sparse.txi, doc/interpreter/strings.txi, doc/interpreter/system.txi, doc/interpreter/vectorize.txi, libinterp/corefcn/balance.cc, libinterp/corefcn/bitfcns.cc, libinterp/corefcn/cellfun.cc, libinterp/corefcn/conv2.cc, libinterp/corefcn/data.cc, libinterp/corefcn/debug.cc, libinterp/corefcn/defaults.cc, libinterp/corefcn/dirfns.cc, libinterp/corefcn/dlmread.cc, libinterp/corefcn/error.cc, libinterp/corefcn/file-io.cc, libinterp/corefcn/find.cc, libinterp/corefcn/gammainc.cc, libinterp/corefcn/graphics.cc, libinterp/corefcn/help.cc, libinterp/corefcn/hex2num.cc, libinterp/corefcn/input.cc, libinterp/corefcn/load-path.cc, libinterp/corefcn/load-save.cc, libinterp/corefcn/ls-oct-ascii.cc, libinterp/corefcn/lu.cc, libinterp/corefcn/luinc.cc, libinterp/corefcn/matrix_type.cc, libinterp/corefcn/oct-hist.cc, libinterp/corefcn/pager.cc, libinterp/corefcn/pr-output.cc, libinterp/corefcn/pt-jit.cc, libinterp/corefcn/qz.cc, libinterp/corefcn/rand.cc, libinterp/corefcn/regexp.cc, libinterp/corefcn/schur.cc, libinterp/corefcn/sighandlers.cc, libinterp/corefcn/sparse.cc, libinterp/corefcn/spparms.cc, libinterp/corefcn/str2double.cc, libinterp/corefcn/svd.cc, libinterp/corefcn/symtab.cc, libinterp/corefcn/syscalls.cc, libinterp/corefcn/toplev.cc, libinterp/corefcn/tril.cc, libinterp/corefcn/typecast.cc, libinterp/corefcn/utils.cc, libinterp/corefcn/variables.cc, libinterp/dldfcn/__init_fltk__.cc, libinterp/dldfcn/chol.cc, libinterp/dldfcn/colamd.cc, libinterp/dldfcn/fftw.cc, libinterp/dldfcn/qr.cc, libinterp/dldfcn/symbfact.cc, libinterp/octave-value/ov-base.cc, libinterp/octave-value/ov-fcn-handle.cc, libinterp/octave-value/ov-fcn-inline.cc, libinterp/octave-value/ov-java.cc, libinterp/octave-value/ov-range.cc, libinterp/octave-value/ov-struct.cc, libinterp/octave-value/ov-usr-fcn.cc, libinterp/parse-tree/oct-parse.in.yy, libinterp/parse-tree/pt-binop.cc, libinterp/parse-tree/pt-eval.cc, libinterp/parse-tree/pt-mat.cc, scripts/@ftp/ftp.m, scripts/deprecated/java_convert_matrix.m, scripts/deprecated/java_debug.m, scripts/deprecated/java_unsigned_conversion.m, scripts/deprecated/shell_cmd.m, scripts/general/dblquad.m, scripts/general/display.m, scripts/general/genvarname.m, scripts/general/idivide.m, scripts/general/interp1.m, scripts/general/interp2.m, scripts/general/interp3.m, scripts/general/interpn.m, scripts/general/isa.m, scripts/general/profexplore.m, scripts/general/profile.m, scripts/general/quadgk.m, scripts/general/randi.m, scripts/general/structfun.m, scripts/general/subsindex.m, scripts/general/triplequad.m, scripts/geometry/griddata.m, scripts/geometry/griddata3.m, scripts/geometry/griddatan.m, scripts/geometry/voronoi.m, scripts/help/help.m, scripts/help/lookfor.m, scripts/image/cmpermute.m, scripts/image/colormap.m, scripts/image/image.m, scripts/image/imagesc.m, scripts/image/imfinfo.m, scripts/image/imformats.m, scripts/image/imread.m, scripts/image/imshow.m, scripts/image/imwrite.m, scripts/image/ind2gray.m, scripts/image/lines.m, scripts/image/rgb2ind.m, scripts/image/spinmap.m, scripts/io/dlmwrite.m, scripts/io/strread.m, scripts/io/textread.m, scripts/io/textscan.m, scripts/java/javaclasspath.m, scripts/java/usejava.m, scripts/miscellaneous/bzip2.m, scripts/miscellaneous/computer.m, scripts/miscellaneous/copyfile.m, scripts/miscellaneous/debug.m, scripts/miscellaneous/dos.m, scripts/miscellaneous/edit.m, scripts/miscellaneous/gzip.m, scripts/miscellaneous/license.m, scripts/miscellaneous/mkoctfile.m, scripts/miscellaneous/movefile.m, scripts/miscellaneous/parseparams.m, scripts/miscellaneous/unix.m, scripts/optimization/fminbnd.m, scripts/optimization/fminsearch.m, scripts/optimization/fminunc.m, scripts/optimization/fsolve.m, scripts/optimization/fzero.m, scripts/optimization/glpk.m, scripts/optimization/lsqnonneg.m, scripts/optimization/optimset.m, scripts/optimization/pqpnonneg.m, scripts/pkg/pkg.m, scripts/plot/allchild.m, scripts/plot/ancestor.m, scripts/plot/area.m, scripts/plot/axis.m, scripts/plot/bar.m, scripts/plot/barh.m, scripts/plot/box.m, scripts/plot/caxis.m, scripts/plot/cla.m, scripts/plot/clabel.m, scripts/plot/clf.m, scripts/plot/close.m, scripts/plot/colorbar.m, scripts/plot/daspect.m, scripts/plot/ezmesh.m, scripts/plot/ezmeshc.m, scripts/plot/ezsurf.m, scripts/plot/ezsurfc.m, scripts/plot/findall.m, scripts/plot/findobj.m, scripts/plot/gcbo.m, scripts/plot/gcf.m, scripts/plot/gco.m, scripts/plot/grid.m, scripts/plot/guihandles.m, scripts/plot/hdl2struct.m, scripts/plot/hidden.m, scripts/plot/hold.m, scripts/plot/isonormals.m, scripts/plot/isosurface.m, scripts/plot/legend.m, scripts/plot/mesh.m, scripts/plot/meshc.m, scripts/plot/meshz.m, scripts/plot/newplot.m, scripts/plot/orient.m, scripts/plot/pareto.m, scripts/plot/patch.m, scripts/plot/pbaspect.m, scripts/plot/pcolor.m, scripts/plot/plot.m, scripts/plot/print.m, scripts/plot/private/__add_default_menu__.m, scripts/plot/quiver.m, scripts/plot/quiver3.m, scripts/plot/refreshdata.m, scripts/plot/saveas.m, scripts/plot/scatter.m, scripts/plot/scatter3.m, scripts/plot/shading.m, scripts/plot/shrinkfaces.m, scripts/plot/slice.m, scripts/plot/stem.m, scripts/plot/stem3.m, scripts/plot/struct2hdl.m, scripts/plot/subplot.m, scripts/plot/surf.m, scripts/plot/surfc.m, scripts/plot/surfl.m, scripts/plot/tetramesh.m, scripts/plot/uigetfile.m, scripts/plot/uimenu.m, scripts/plot/uiputfile.m, scripts/plot/waterfall.m, scripts/plot/whitebg.m, scripts/plot/xlim.m, scripts/plot/ylim.m, scripts/plot/zlim.m, scripts/polynomial/conv.m, scripts/polynomial/polyout.m, scripts/polynomial/splinefit.m, scripts/set/ismember.m, scripts/set/powerset.m, scripts/set/setdiff.m, scripts/set/union.m, scripts/set/unique.m, scripts/signal/detrend.m, scripts/signal/filter2.m, scripts/signal/freqz.m, scripts/signal/periodogram.m, scripts/signal/spectral_adf.m, scripts/signal/spectral_xdf.m, scripts/sparse/eigs.m, scripts/sparse/svds.m, scripts/specfun/legendre.m, scripts/special-matrix/gallery.m, scripts/statistics/base/mean.m, scripts/statistics/base/moment.m, scripts/statistics/tests/cor_test.m, scripts/statistics/tests/kolmogorov_smirnov_test.m, scripts/statistics/tests/kolmogorov_smirnov_test_2.m, scripts/statistics/tests/kruskal_wallis_test.m, scripts/statistics/tests/prop_test_2.m, scripts/statistics/tests/sign_test.m, scripts/statistics/tests/t_test.m, scripts/statistics/tests/t_test_2.m, scripts/statistics/tests/t_test_regression.m, scripts/statistics/tests/u_test.m, scripts/statistics/tests/var_test.m, scripts/statistics/tests/welch_test.m, scripts/statistics/tests/wilcoxon_test.m, scripts/statistics/tests/z_test.m, scripts/statistics/tests/z_test_2.m, scripts/strings/base2dec.m, scripts/strings/index.m, scripts/strings/isstrprop.m, scripts/strings/mat2str.m, scripts/strings/regexptranslate.m, scripts/strings/rindex.m, scripts/strings/str2num.m, scripts/strings/strcat.m, scripts/strings/strjust.m, scripts/strings/strmatch.m, scripts/strings/validatestring.m, scripts/testfun/demo.m, scripts/testfun/example.m, scripts/testfun/test.m, scripts/time/addtodate.m, scripts/time/asctime.m, scripts/time/datestr.m, scripts/time/datetick.m, scripts/time/weekday.m, scripts/ui/errordlg.m, scripts/ui/helpdlg.m, scripts/ui/inputdlg.m, scripts/ui/listdlg.m, scripts/ui/msgbox.m, scripts/ui/questdlg.m, scripts/ui/warndlg.m: Use new @qcode macro.
author Rik <rik@octave.org>
date Mon, 19 Aug 2013 20:46:38 -0700
parents 2fc554ffbc28
children
line wrap: on
line source

/*

Copyright (C) 2007-2012 Alexander Barth

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/>.

*/

#ifdef HAVE_CONFIG_H
#include <config.h>
#endif

#include "lo-ieee.h"
#include "dNDArray.h"
#include "oct-locbuf.h"

#include "defun.h"
#include "error.h"
#include "oct-obj.h"

// equivalent to isvector.m

template <class T>
bool
isvector (const T& array)
{
  const dim_vector dv = array.dims ();
  return dv.length () == 2 && (dv(0) == 1 || dv(1) == 1);
}

// lookup a value in a sorted table (lookup.m)
template <class T>
octave_idx_type
lookup (const T *x, octave_idx_type n, T y)
{
  octave_idx_type j;

  if (x[0] < x[n-1])
    {
      // increasing x

      if (y > x[n-1] || y < x[0])
        return -1;

#ifdef EXHAUSTIF
      for (j = 0; j < n - 1; j++)
        {
          if (x[j] <= y && y <= x[j+1])
            return j;
        }
#else
      octave_idx_type j0 = 0;
      octave_idx_type j1 = n - 1;

      while (true)
        {
          j = (j0+j1)/2;

          if (y <= x[j+1])
            {
              if (x[j] <= y)
                return j;

              j1 = j;
            }

          if (x[j] <= y)
            j0 = j;
        }
#endif
    }
  else
    {
      // decreasing x
      // previous code with x -> -x and y -> -y

      if (y > x[0] || y < x[n-1])
        return -1;

#ifdef EXHAUSTIF
      for (j = 0; j < n - 1; j++)
        {
          if (x[j+1] <= y && y <= x[j])
            return j;
        }
#else
      octave_idx_type j0 = 0;
      octave_idx_type j1 = n - 1;

      while (true)
        {
          j = (j0+j1)/2;

          if (y >= x[j+1])
            {
              if (x[j] >= y)
                return j;

              j1 = j;
            }

          if (x[j] >= y)
            j0 = j;
        }
#endif
    }
}

// n-dimensional linear interpolation

template <class T>
void
lin_interpn (int n, const octave_idx_type *size, const octave_idx_type *scale,
             octave_idx_type Ni, T extrapval, const T **x,
             const T *v, const T **y, T *vi)
{
  bool out = false;
  int bit;

  OCTAVE_LOCAL_BUFFER (T, coef, 2*n);
  OCTAVE_LOCAL_BUFFER (octave_idx_type, index, n);

  // loop over all points
  for (octave_idx_type m = 0; m < Ni; m++)
    {
      // loop over all dimensions
      for (int i = 0; i < n; i++)
        {
          index[i] = lookup (x[i], size[i], y[i][m]);
          out = index[i] == -1;

          if (out)
            break;
          else
            {
              octave_idx_type j = index[i];
              coef[2*i+1] = (y[i][m] - x[i][j])/(x[i][j+1] - x[i][j]);
              coef[2*i] = 1 - coef[2*i+1];
            }
        }


      if (out)
        vi[m] = extrapval;
      else
        {
          vi[m] = 0;

          // loop over all corners of hypercube (1<<n = 2^n)
          for (int i = 0; i < (1 << n); i++)
            {
              T c = 1;
              octave_idx_type l = 0;

              // loop over all dimensions
              for (int j = 0; j < n; j++)
                {
                  // test if the jth bit in i is set
                  bit = i >> j & 1;
                  l += scale[j] * (index[j] + bit);
                  c *= coef[2*j+bit];
                }

              vi[m] += c * v[l];
            }
        }
    }
}

template <class T, class M>
octave_value
lin_interpn (int n, M *X, const M V, M *Y)
{
  octave_value retval;

  M Vi = M (Y[0].dims ());

  OCTAVE_LOCAL_BUFFER (const T *, y, n);
  OCTAVE_LOCAL_BUFFER (octave_idx_type, size, n);

  for (int i = 0; i < n; i++)
    {
      y[i] = Y[i].data ();
      size[i] =  V.dims ()(i);
    }

  OCTAVE_LOCAL_BUFFER (const T *, x, n);
  OCTAVE_LOCAL_BUFFER (octave_idx_type, scale, n);

  const T *v = V.data ();
  T *vi = Vi.fortran_vec ();
  octave_idx_type Ni = Vi.numel ();

  T extrapval = octave_NA;

  // offset in memory of each dimension

  scale[0] = 1;

  for (int i = 1; i < n; i++)
    scale[i] = scale[i-1] * size[i-1];

  // tests if X[0] is a vector, if yes, assume that all elements of X are
  // in the ndgrid format.

  if (! isvector (X[0]))
    {
      for (int i = 0; i < n; i++)
        {
          if (X[i].dims () != V.dims ())
            {
              error ("interpn: incompatible size of argument number %d", i+1);
              return retval;
            }
          else
            {
              M tmp = M (dim_vector (size[i], 1));

              for (octave_idx_type j = 0; j < size[i]; j++)
                tmp(j) =  X[i](scale[i]*j);

              X[i] = tmp;
            }
        }
    }

  for (int i = 0; i < n; i++)
    {
      if (! isvector (X[i]) && X[i].numel () != size[i])
        {
          error ("interpn: incompatible size of argument number %d", i+1);
          return retval;
        }
      else
        x[i] = X[i].data ();
    }

  lin_interpn (n, size, scale, Ni, extrapval, x, v, y, vi);

  retval = Vi;

  return retval;
}

// Perform @var{n}-dimensional interpolation.  Each element of then
// @var{n}-dimensional array @var{v} represents a value at a location
// given by the parameters @var{x1}, @var{x2},...,@var{xn}. The parameters
// @var{x1}, @var{x2}, @dots{}, @var{xn} are either @var{n}-dimensional
// arrays of the same size as the array @var{v} in the \"ndgrid\" format
// or vectors.  The parameters @var{y1}, @var{y2}, @dots{}, @var{yn} are
// all @var{n}-dimensional arrays of the same size and represent the
// points at which the array @var{vi} is interpolated.
//
//This function only performs linear interpolation.

DEFUN (__lin_interpn__, args, ,
  "-*- texinfo -*-\n\
@deftypefn {Built-in Function} {@var{vi} =} __lin_interpn__ (@var{x1}, @var{x2}, @dots{}, @var{xn}, @var{v}, @var{y1}, @var{y2}, @dots{}, @var{yn})\n\
Undocumented internal function.\n\
@end deftypefn")
{
  octave_value retval;

  int nargin = args.length ();

  if (nargin < 2 ||  nargin % 2 == 0)
    {
      print_usage ();
      return retval;
    }

  // dimension of the problem
  int n = (nargin-1)/2;

  if (args(n).is_single_type ())
    {
      OCTAVE_LOCAL_BUFFER (FloatNDArray, X, n);
      OCTAVE_LOCAL_BUFFER (FloatNDArray, Y, n);

      const FloatNDArray V = args(n).float_array_value ();

      if (error_state)
        {
          print_usage ();
          return retval;
        }

      for (int i = 0; i < n; i++)
        {
          X[i] = args(i).float_array_value ();
          Y[i] = args(n+i+1).float_array_value ();

          if (error_state)
            {
              print_usage ();
              return retval;
            }

          if (Y[0].dims () != Y[i].dims ())
            {
              error ("interpn: incompatible size of argument number %d", n+i+2);
              return retval;
            }
        }

      retval = lin_interpn<float, FloatNDArray> (n, X, V, Y);
    }
  else
    {
      OCTAVE_LOCAL_BUFFER (NDArray, X, n);
      OCTAVE_LOCAL_BUFFER (NDArray, Y, n);

      const NDArray V = args(n).array_value ();

      if (error_state)
        {
          print_usage ();
          return retval;
        }

      for (int i = 0; i < n; i++)
        {
          X[i] = args(i).array_value ();
          Y[i] = args(n+i+1).array_value ();

          if (error_state)
            {
              print_usage ();
              return retval;
            }

          if (Y[0].dims () != Y[i].dims ())
            {
              error ("interpn: incompatible size of argument number %d", n+i+2);
              return retval;
            }
        }

      retval = lin_interpn<double, NDArray> (n, X, V, Y);
    }

  return retval;
}

/*
## No test needed for internal helper function.
%!assert (1)
*/