changeset 11662:fbd81057dab0 octave-forge

image: new spatial transformation functions * maketform: accept input/output control points as second and third argument * imtransform: new function * findsbounds: new function * private/istform: add comment on what function does * NEWS/INDEX: update list of new functions
author carandraug
date Sun, 28 Apr 2013 02:31:13 +0000
parents 1b4e81051b66
children 3863d5c25ae3
files main/image/INDEX main/image/NEWS main/image/inst/findbounds.m main/image/inst/imtransform.m main/image/inst/maketform.m main/image/inst/private/istform.m
diffstat 6 files changed, 482 insertions(+), 26 deletions(-) [+]
line wrap: on
line diff
--- a/main/image/INDEX	Sat Apr 27 19:17:33 2013 +0000
+++ b/main/image/INDEX	Sun Apr 28 02:31:13 2013 +0000
@@ -63,8 +63,10 @@
  ordfiltn
  stretchlim
 Filtering and Transforms
+ findbounds
  fspecial
  imfilter
+ imtransform
  iradon
  nonmax_supress
  phantom
--- a/main/image/NEWS	Sat Apr 27 19:17:33 2013 +0000
+++ b/main/image/NEWS	Sun Apr 28 02:31:13 2013 +0000
@@ -5,6 +5,8 @@
 
       checkerboard
       cp2tform
+      findbounds
+      imtransform
       maketform
       strel
       tformfwd
--- /dev/null	Thu Jan 01 00:00:00 1970 +0000
+++ b/main/image/inst/findbounds.m	Sun Apr 28 02:31:13 2013 +0000
@@ -0,0 +1,65 @@
+## Copyright (C) 2012 Pantxo Diribarne
+##
+## This program 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.
+##
+## This program 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} {@var{outbnd} =} findbounds (@var{T}, @var{inbnd})
+## Estimate bounds for spatial transformation.
+##
+## Given a transformation structure @var{T} (see e.g. maketform) 
+## and bounds @var{inbnd} (2-by-ndims_in) in an input space, returns
+## an estimation of the bounds in the output space @var{outbnd} 
+## (2-by-ndims_out). For instance two dimensionnal bounds could 
+## be represented as : [xmin ymin; xmax ymax]. If @var{T} does not
+## define a forward trasform (i.e. for 'polynomial'), the output
+## bounds are infered using fsolve and the inverse transform.
+##
+## @seealso{maketform, cp2tform, tformfwd}
+## @end deftypefn
+
+## Author: Pantxo Diribarne <pantxo@dibona>
+
+function [outbnd] = findbounds (T, inbnd)
+  if (nargin != 2)
+    print_usage ();
+  elseif (! istform (T))
+    error ("imtransform: T must be a transformation structure (see `maketform')");
+  elseif (! all (size (inbnd) == [2 T.ndims_in]))
+    error ("imtransform: INBDN must have same size as T.ndims_in");
+  endif
+
+  ## Control points grid
+  if (columns (inbnd) == 2)
+    [xcp, ycp] = meshgrid (linspace (inbnd(1,1), inbnd(2,1), 3),
+                           linspace (inbnd(1,2), inbnd(2,2), 3));
+    xcp = reshape (xcp, numel (xcp), 1);
+    ycp = reshape (ycp, numel (ycp), 1);
+    xycp = [xcp, ycp];
+  else
+    error ("findbounds: support only 2D inputbounds.");
+  endif
+
+  ## Output bounds
+  if (!is_function_handle (T.forward_fcn))
+    outbnd = zeros (size (xycp));
+    for ii = 1:rows (xycp)
+      fun = @(x) tforminv (T, x) - xycp(ii,:);
+      outbnd(ii,:) = fsolve (fun, xycp(ii,:));
+    endfor
+  else
+    outbnd = tformfwd (T, xycp);
+  endif
+  outbnd = [min(outbnd); max(outbnd)];
+endfunction
--- /dev/null	Thu Jan 01 00:00:00 1970 +0000
+++ b/main/image/inst/imtransform.m	Sun Apr 28 02:31:13 2013 +0000
@@ -0,0 +1,350 @@
+## Copyright (C) 2012 Pantxo Diribarne
+##
+## This program 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.
+##
+## This program 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} {@var{B} =} imtransform (@var{A}, @var{T})
+## @deftypefnx {Function File} {@var{B} =} imtransform (@var{A}, @var{T}, @var{interp})
+## @deftypefnx {Function File} {@var{B} =} imtransform (@dots{}, @var{prop}, @var{val})
+## Transform image.
+##
+## Given an image @var{A} in one space, returns an image @var{B}
+## resulting from the forward transform defined in the transformation
+## structure @var{T}.  An additionnal input argument @var{interp},
+## 'bicubic', 'bilinear' (default) or 'nearest',
+## specifies the interpolation method to be used.  Finally, the
+## transformation can be tuned using  @var{prop}/@var{val} pairs.  The
+## following properties are supported:
+## 
+## @table @asis
+## @item "udata"
+## Specifies the input space horizontal limits. Value must be a two
+## elements vector [minval maxval]. Default: [1 columns(@var{A})]
+##
+## @item "vdata"
+## Specifies the input space vertical limits. Value must be a two
+## elements vector [minval maxval]. Default: [1 rows(@var{A})]
+##
+## @item "xdata"
+## Specifies the requiered output space horizontal limits. Value must
+## be a two elements vector [minval maxval]. Default: estimated using
+## udata, vdata and findbounds function.
+##
+## @item "ydata"
+## Specifies the requiered output space vertical limits. Value must
+## be a two elements vector [minval maxval]. Default: estimated using
+## udata, vdata and findbounds function.
+##
+## @item "xyscale"
+## Specifies the output scale in outputspace_units/pixel. If a scalar
+## is provided, both vertical and horizontal dimensions are scaled the
+## same way. If @var{val} is a two element vector, it must indicate
+## consecutively horizontal and vertical scales. Default value is
+## computed using the input space scale provided
+## that the number of pixel of any dimension of the output image does
+## not exceed 20000.
+##
+## @item "size"
+## Size of the output image (1-by-2 vector). Overrides the effect of
+## "xyscale" property. 
+##
+## @item "fillvalues"
+## Color of the areas where no interpolation where possible, e.g. when
+## coorfinates of points in the output space are out of the limits of
+## the input space. @var{val} must be coherent with the input image
+## format: for grayscale and indexed images (2D) @var{val} must be
+## scalar, for RGB (n-by-m-by-3) @var{val} must be a 3 element vector.
+## e.g.:
+##
+## @end table
+##
+## @seealso{maketform, cp2tform, tforminv, tformfwd, findbounds}
+## @end deftypefn
+
+## Author: Pantxo Diribarne <pantxo@dibona>
+
+function [varargout] = imtransform (im, T, varargin)
+
+  if (nargin < 2)
+    print_usage ();
+  elseif (! istform (T))
+    error ("imtransform: T must be a transformation structure (see `maketform')");
+  endif
+
+  ## Parameters
+  interp = "linear";
+  imdepth = size (im, 3);
+  maximsize = [20000 20000];
+
+  udata = [1; columns(im)];
+  vdata = [1; rows(im)];
+  xdata = ydata = [];
+  xyscale = [];
+  imsize = [];
+  fillvalues = ones (1, imdepth) * NaN;
+
+  if (isempty (varargin))
+    xydata = findbounds (T, [udata vdata]);
+    xdata = xydata(:,1);
+    ydata = xydata(:,2);
+  else
+    ## interp
+    if (floor (numel (varargin)/2) != numel (varargin)/2)
+      allowed = {"bicubic", "bilinear", "nearest"};
+      tst = strcmp (varargin{1}, allowed);
+      if (!any (tst))
+        error ("imtransform: expect one of %s as interp method", disp (allowed));
+      else
+        interp = {"pchip", "linear", "nearest"}{find (tst)};
+      endif
+      varargin = varargin(2:end);
+    endif
+
+    ## options
+    allowed = {"udata", "vdata", "xdata", "ydata", ...
+               "xyscale", "size", "fillvalues"};
+    props = varargin(1:2:end);
+    vals = varargin(2:2:end);
+    np = numel (props);
+    if (!all (cellfun (@ischar, props)))
+      error ("imtransform: expect property/value pairs.");
+    endif
+
+    props = tolower (props);
+    tst = cellfun (@(x) any (strcmp (x, allowed)), props);
+    if (!all (tst))
+      error ("imtransform: unknown property %s", disp (props{!tst}));
+    endif
+
+    ## u(vxy)data
+    iolims = allowed(1:4);
+    for ii = 1:numel (iolims)
+      tst = cellfun (@(x) any (strcmp (x, iolims{ii})), props);
+      if (any (tst))
+        prop = props{find (tst)(1)};
+        val = vals{find (tst)(1)};
+        if (isnumeric (val) && numel (val) == 2)
+          if (isrow (val))
+            val = val';
+          endif
+          eval (sprintf ("%s = val;", prop),
+                "error (\"imtransform: %s\n\", lasterr ());");
+        else
+          error ("imtransform: expect 2 elements real vector for %s", prop)
+        endif
+      endif
+    endfor
+    if (isempty (xdata) && isempty (ydata))
+      xydata = findbounds (T, [udata vdata]);
+      xdata = xydata(:,1);
+      ydata = xydata(:,2);
+    elseif (isempty (xdata))
+      xydata = findbounds (T, [udata vdata]);
+      xdata = xydata(:,1);
+    elseif (isempty (ydata))
+      xydata = findbounds (T, [udata vdata]);
+      ydata = xydata(:,2);
+    endif
+
+    ## size and xyscale
+    tst = strcmp ("size", props);
+    if (any (tst))
+      val = vals{find (tst)(1)};
+      if (isnumeric (val) && numel (val) == 2 &&
+          all (val > 0))
+        imsize = val;
+      else
+        error ("imtransform: expect 2 elements real vector for size");
+      endif
+    elseif (any (tst = strcmp ("xyscale", props)))
+      val = vals{find (tst)(1)};
+      if (isnumeric (val) && all (val > 0))
+        if (numel (val) == 1)
+          xyscale(1:2) = val;
+        elseif (numel (val) == 2)
+          xyscale = val;
+        else
+          error ("imtransform: expect 1 or 2 element(s) real vector for xyscale");
+        endif
+      else
+        error ("imtransform: expect 1 or 2 elements real vector for xyscale");
+      endif
+    else
+      xyscale = [(diff (udata) / columns (im)) (diff (vdata) / rows (im))];
+    endif
+
+    ## Fillvalues
+    tst = strcmp ("fillvalues", props);
+    if (any (tst))
+      val = vals{find (tst)(1)};
+      if (isnumeric (val) && numel (val) == 1)
+        fillvalues(1:end) = val;
+      elseif (isnumeric (val) && numel (val) == 3)
+        fillvalues = val;
+      else
+        error ("imtransform: expect 1 or 3 elements real vector for `fillvalues'");
+      endif
+    endif
+  endif
+
+  ## Ouput/Input pixels
+  if (isempty (imsize))
+    if (isempty (xyscale))
+      xyscale = [(diff (udata) / columns (im)) (diff (vdata) / rows (im))];
+    endif
+    xscale = xyscale(1);
+    yscale = xyscale(2);
+    xsize = floor (diff (xdata) / xscale);
+    ysize = floor (diff (ydata) / yscale);
+    if (xsize > maximsize(2) || ysize > maximsize(1))
+      if (xsize >= ysize)
+        scalefactor = (diff (xdata) / maximsize(2)) / xscale;
+      else
+        scalefactor = (diff (ydata) / maximsize(1)) / yscale;
+      endif
+      xscale *= scalefactor
+      yscale *= scalefactor
+      xsize = floor (diff (xdata) / xscale);
+      ysize = floor (diff (ydata) / yscale);
+      warning ("imtransform: output image two large, adjusting the largest dimension to %d", maximsize);
+    endif
+    imsize = [ysize xsize];
+  endif
+  [xx yy] = meshgrid (linspace (xdata(1), xdata(2), imsize(2)),
+                      linspace (ydata(1), ydata(2), imsize(1)));
+
+  [uu vv] = meshgrid (linspace (udata(1), udata(2), size(im)(2)),
+                      linspace (vdata(1), vdata(2), size(im)(1)));
+
+  ## Input coordinates
+  [uui, vvi] = tforminv (T, reshape (xx, numel (xx), 1),
+                         reshape (yy, numel (yy), 1));
+  uui = reshape (uui, size (xx));
+  vvi = reshape (vvi, size (yy));
+  ## Interpolation
+  for layer = 1:imdepth
+    imout(:,:,layer) = interp2 (uu, vv, im(:,:,layer), ...
+                                uui, vvi, interp, fillvalues(layer));
+  endfor
+  if (nargout == 1)
+    varargout{1} = imout;
+  else
+    varargout = {imout, xdata, ydata};
+  endif
+endfunction
+
+%!demo
+%! ## Various linear transforms
+%! figure (); 
+%! im = [checkerboard(20, 2, 4); checkerboard(40, 1, 2)];
+%! %input space corners
+%! incp = [1 1; 160 1; 160 160; 1 160];
+%! udata = [min(incp(:,1)) max(incp(:,1))];
+%! vdata = [min(incp(:,2)) max(incp(:,2))];
+%! subplot (2,3,1); 
+%! imshow (im)
+%! hold on
+%! plot (incp(:,1), incp(:,2), 'ob')
+%! axis on
+%! xlabel ('Original')
+%! 
+%! % Translation and scaling
+%! outcp = incp * 2;
+%! outcp(:,1) += 200; 
+%! outcp(:,2) += 500;
+%! T = maketform ('affine', incp(1:3,:), outcp(1:3,:));
+%! subplot (2,3,2); 
+%! [im2 xdata ydata] = imtransform (im, T, 'udata', udata,
+%!                                  'vdata', vdata, 'fillvalues', 1);
+%! imh = imshow (im2); set (imh, 'xdata', xdata, 'ydata', ydata) 
+%! set (gca, 'xlim', xdata, 'ylim', ydata)
+%! axis on, hold on, xlabel ('Translation / Scaling');
+%! plot (outcp(:,1), outcp(:,2), 'or')
+%! 
+%! % Shear
+%! outcp = [1 1; 160 1; 140 160; -19 160]; % affine only needs 3 control points
+%! T = maketform ('affine', incp(1:3,:), outcp(1:3,:));
+%! subplot (2,3,3); 
+%! [im2 xdata ydata] = imtransform (im, T, 'udata', udata,
+%!                                  'vdata', vdata, 'fillvalues', 1);
+%! imh = imshow (im2); set (imh, 'xdata', xdata, 'ydata', ydata) 
+%! set (gca, 'xlim', xdata, 'ylim', ydata)
+%! axis on, hold on, xlabel ('Shear');
+%! plot (outcp(:,1), outcp(:,2), 'or')
+%! 
+%! % Rotation 
+%! theta = pi/4;
+%! T = maketform ('affine', [cos(theta) -sin(theta); ...
+%!                           sin(theta) cos(theta); 0 0]);
+%! outcp = tformfwd (T, incp);
+%! subplot (2,3,4); 
+%! [im2 xdata ydata] = imtransform (im, T, 'udata', udata,
+%!                                  'vdata', vdata, 'fillvalues', 1 );
+%! imh = imshow (im2); set (imh, 'xdata', xdata, 'ydata', ydata) 
+%! set (gca, 'xlim', xdata, 'ylim', ydata)
+%! axis on, hold on, xlabel ('Rotation');
+%! plot (outcp(:,1), outcp(:,2), 'or')
+%!
+%! % Reflection around x axis
+%! outcp = incp;
+%! outcp(:,2) *= -1;
+%! T = cp2tform (incp, outcp, 'similarity');
+%! subplot (2,3,5); 
+%! [im2 xdata ydata] = imtransform (im, T, 'udata', udata,
+%!                                  'vdata', vdata, 'fillvalues', 1 );
+%! imh = imshow (im2); set (imh, 'xdata', xdata, 'ydata', ydata) 
+%! set (gca, 'xlim', xdata, 'ylim', ydata)
+%! axis on, hold on, xlabel ('Reflection');
+%! plot (outcp(:,1), outcp(:,2), 'or')
+%!
+%! % Projection
+%! outcp = [1 1; 160 -40; 220 220; 12 140];
+%! T = maketform ('projective', incp, outcp);
+%! subplot (2,3,6); 
+%! [im2 xdata ydata] = imtransform (im, T, 'udata', udata,
+%!                                  'vdata', vdata, 'fillvalues', 1 );
+%! imh = imshow (im2); set (imh, 'xdata', xdata, 'ydata', ydata) 
+%! set (gca, 'xlim', xdata, 'ylim', ydata)
+%! axis on, hold on, xlabel ('Projection');
+%! plot (outcp(:,1), outcp(:,2), 'or')
+
+%!demo
+%! ## Streched image
+%! rad = 2; % minimum value: 4/pi
+%! [uu vv] = meshgrid ((-2:2)/rad, (-2:2)/rad);
+%! rescfactor = sin ((uu.^2 + vv.^2).^.5);
+%! inpts = [(reshape (uu, numel (uu), 1)), (reshape (vv, numel (uu), 1))]; 
+%! xx = rescfactor .* sign(uu);
+%! yy = rescfactor .* sign(vv);
+%! outpts = [reshape(xx, numel (xx), 1) reshape(yy, numel (yy), 1)];
+%! 
+%! T = cp2tform (inpts, outpts, "polynomial", 4);
+%! figure;
+%! subplot (1,2,1)
+%! im = zeros (800, 800, 3);
+%! im(:,:,1) = checkerboard (100) > 0.2;
+%! im(:,:,3) = checkerboard (100) < 0.2;
+%! [im2 xdata ydata] = imtransform (im, T, 'udata', [-2 2],
+%!                                  'vdata', [-2 2], 'fillvalues',
+%!                                  [0 1 0]);
+%! imh = imshow (im2);
+%! set (imh, 'xdata', xdata, 'ydata', ydata)
+%! set (gca, 'xlim', xdata, 'ylim', ydata)
+%! [im cmap] = imread ('default.img');
+%! subplot (1,2,2)
+%! [im2 xdata ydata] = imtransform (im, T, 'udata', [-1 1],
+%!                                  'vdata', [-1 1], 'fillvalues',
+%!                                  round (length (cmap) / 2));
+%! imh = imshow (im2, cmap);
--- a/main/image/inst/maketform.m	Sat Apr 27 19:17:33 2013 +0000
+++ b/main/image/inst/maketform.m	Sun Apr 28 02:31:13 2013 +0000
@@ -1,46 +1,68 @@
 ## Copyright (C) 2012 Pantxo Diribarne
-## 
+##
 ## This program 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.
-## 
+##
 ## This program 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} {@var{T} =} maketform (@var{ttype}, @var{tmat})
+## @deftypefnx {Function File} {@var{T} =} maketform (@var{ttype}, @var{inc}, @var{outc})
 ## @deftypefnx {Function File} {@var{T} =} maketform ("custom", @var{ndims_in}, @var{ndims_out}, @var{forward_fcn}, @var{inverse_fcn}, @var{tdata})
+## Create structure for spatial transformations.
+##
 ## Returns a transform structure containing fields @var{ndims_in},
-## @var{ndims_out}, @var{forward_fcn}, @var{inverse_fcn} and @var{tdata}. The content
-## of each field depends on the requested transform type @var{ttype}:
+## @var{ndims_out}, @var{forward_fcn}, @var{inverse_fcn} and @var{tdata}.  The
+## content of each field depends on the requested transform type @var{ttype}:
+##
 ## @table @asis
 ## @item "projective"
-## A ndims_in = N -> ndims_out = N projective transformation structure
+## A ndims_in = N -> @var{ndims_out} = N projective transformation structure
 ## is returned.
 ## The second input argument @var{tmat} must be a (N+1)-by-(N+1)
-## transformation matrix. The
-## (N+1)th column must contain projection coefficients. As an example a two
+## transformation matrix.  The
+## (N+1)th column must contain projection coefficients.  As an example a two
 ## dimentionnal transform from [x y] coordinates to [u v] coordinates
 ## is represented by a transformation matrix defined so that:
+##
 ## @example
 ## [xx yy zz] = [u v 1] * [a d g;
 ##                         b e h;
 ##                         c f i]
 ## [x y] =  [xx./zz yy./zz];
 ## @end example
+## 
+## Alternatively the transform can be specified using a quadilateral
+## coordinates (typically the 4 corners of the
+## image) in the input space (@var{inc}, 4-by-ndims_in matrix) and in
+## the output space (@var{outc}, 4-by-ndims_out matrix).  This is
+## equivalent to building the transform using
+## @code{T = cp2tform (@var{inc}, @var{outc}, "projective")}.
+##
 ## @item "affine"
-## Affine is a subset of projective transform (see above). A ndims_in = N ->
-## ndims_out = N  affine transformation structure is returned.
+## Affine is a subset of projective transform (see above).  A
+## @var{ndims_in} = N -> @var{ndims_out} = N affine transformation structure is
+## returned.
 ## The second input argument @var{tmat} must be a (N+1)-by-(N+1) or
 ## (N+1)-by-(N) transformation matrix. If present, the (N+1)th column  must
 ## contain [zeros(N,1); 1] so that projection is suppressed.
+##
+## Alternatively the transform can be specified a using a triangle
+## coordinates (typically the 3 corners of the
+## image)  in the input space (@var{inc}, 3-by-ndims_in matrix) and in
+## the  output space (@var{outc}, 3-by-ndims_out matrix). This is
+## equivalent to building the transform using "T = cp2tform (@var{inc}, @var{outc},
+## 'affine')".
+## 
 ## @item "custom"
 ## For user defined transforms every field of the transform structure
 ## must be supplied. The prototype of the transform functions,
@@ -50,17 +72,18 @@
 ## The argument T is the transformation structure which will contain
 ## the user supplied transformation matrix @var{tdata}. 
 ## @end table
+##
 ## @seealso{tformfwd, tforminv, cp2tform}
 ## @end deftypefn
 
 ## Author: Pantxo Diribarne <pantxo@dibona>
-## Created: 2012-09-05
 
 function T = maketform (ttype, varargin)
 
-  if (nargin < 2 || ! any (strcmp (ttype, {"affine", "projective", "custom"})))
+  if (nargin < 2 || ! any (strcmpi (ttype, {"affine", "projective", "custom"})))
     print_usage ();
   endif
+
   if (numel (varargin) == 1)
     tmat = varargin {1};
     ndin = rows (tmat) - 1;
@@ -94,6 +117,26 @@
     T.inverse_fcn = inverse_fcn;
     T.tdata.T = tmat;
     T.tdata.Tinv = inv (tmat);
+
+  elseif (numel (varargin) == 2)
+    inc = varargin{1};
+    outc = varargin{2};
+    if (strcmp (ttype, "affine"))
+      if (all (size (inc) == size (outc)) &&
+          all (size (inc) == [3 2]))
+        T = cp2tform (inc, outc, ttype);
+      else
+        error ("maketform: expect INC and OUTC to be 3-by-2 vectors.");
+      endif
+    elseif (strcmp (ttype, "projective"))
+      if (all (size (inc) == size (outc)) &&
+          all (size (inc) == [4 2]))
+        T = cp2tform (inc, outc, ttype);
+      else
+        error ("maketform: expect INC and OUTC to be 4-by-2 vectors.");
+      endif
+    endif
+
   elseif (numel (varargin) == 5 && strcmpi (ttype, "custom"))
     if (isscalar (varargin{1}) && isscalar (varargin{2})
         && varargin{1} > 0 && varargin{2} > 0)
@@ -114,6 +157,7 @@
     endif
     
     T.tdata = varargin{5};
+
   else
     print_usage ();
   endif
--- a/main/image/inst/private/istform.m	Sat Apr 27 19:17:33 2013 +0000
+++ b/main/image/inst/private/istform.m	Sun Apr 28 02:31:13 2013 +0000
@@ -1,29 +1,22 @@
 ## Copyright (C) 2012 Pantxo Diribarne
-## 
+##
 ## This program 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.
-## 
+##
 ## This program 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} {[@var{UV}] =} tformfwd (@var{T}, @var{XY})
-## @deftypefnx {Function File} {[@var{U} @var{V}] =} tformfwd (@var{T}, @var{X},@var{Y})
-## 
-## Given to dimensionnal coordinates from one space, returns two 
-## dimensionnal coordinates in the other space, as defined in 
-## the transform structure @var{T}. Input and output coordinates 
-## may be gigen either as a n-by-2 arrays, or as two n-by-1 vectors.
-## @seealso{maketform, cp2tform, tforminv}
-## @end deftypefn
+## Private internal function to check if a argument is a transformation
+## structure (as created by maketform and to be use by findbounds,
+## imtransform, and the like)
 
 ## Author: Pantxo Diribarne <pantxo@dibona>
 
@@ -35,7 +28,7 @@
     required = {"ndims_in";"ndims_out"; ...
                 "forward_fcn"; "inverse_fcn"; ...
                 "tdata"};
-  
+
     fields = fieldnames (T);
     tst = cellfun (@(x) any (strcmp (fields, x)), required);
     if (! all (tst))