Skip to content

Commit

Permalink
New typemap for in-place arrays of arbitrary number of dimensions:
Browse files Browse the repository at this point in the history
  (DATA_TYPE* INPLACE_ARRAY_FLAT, DIM_TYPE DIM_FLAT)
Added unittests, updated documentation.
  • Loading branch information
Tom Krauss committed May 25, 2015
1 parent 0844302 commit a49ad15
Show file tree
Hide file tree
Showing 10 changed files with 406 additions and 26 deletions.
11 changes: 11 additions & 0 deletions doc/source/reference/swig.interface-file.rst
Expand Up @@ -320,6 +320,17 @@ signatures are
These typemaps now check to make sure that the ``INPLACE_ARRAY``
arguments use native byte ordering. If not, an exception is raised.

There is also a "flat" in-place array for situations in which
you would like to modify or process each element, regardless of the
number of dimensions. One example is a "quantization" function that
quantizes each element of an array in-place, be it 1D, 2D or whatever.
This form checks for continuity but allows either C or Fortran ordering.

ND:

* ``(DATA_TYPE* INPLACE_ARRAY_FLAT, DIM_TYPE DIM_FLAT)``


Argout Arrays
`````````````

Expand Down
1 change: 1 addition & 0 deletions doc/source/reference/swig.testing.rst
Expand Up @@ -57,6 +57,7 @@ Two-dimensional arrays are tested in exactly the same manner. The
above description applies, but with ``Matrix`` substituted for
``Vector``. For three-dimensional tests, substitute ``Tensor`` for
``Vector``. For four-dimensional tests, substitute ``SuperTensor``
for ``Vector``. For flat in-place array tests, substitute ``Flat``
for ``Vector``.
For the descriptions that follow, we will reference the
``Vector`` tests, but the same information applies to ``Matrix``,
Expand Down
52 changes: 31 additions & 21 deletions tools/swig/README
Expand Up @@ -37,6 +37,11 @@ The files related to testing are are in the test subdirectory::
SuperTensor.i
testSuperTensor.py

Flat.h
Flat.cxx
Flat.i
testFlat.py

The header files contain prototypes for functions that illustrate the
wrapping issues we wish to address. Right now, this consists of
functions with argument signatures of the following forms. Vector.h::
Expand Down Expand Up @@ -89,9 +94,12 @@ SuperTensor.h::

(type ARGOUT_ARRAY4[ANY][ANY][ANY][ANY])

Flat.h::
(type* INPLACE_ARRAY_FLAT, int DIM_FLAT)

These function signatures take a pointer to an array of type "type",
whose length is specified by the integer(s) DIM1 (and DIM2, and DIM3,
and DIM4).
and DIM4, or DIM_FLAT).

The objective for the IN_ARRAY signatures is for SWIG to generate
python wrappers that take a container that constitutes a valid
Expand All @@ -105,30 +113,32 @@ The objective for the INPLACE_ARRAY signatures is for SWIG to generate
python wrappers that accept a numpy array of any of the above-listed
types.

The source files Vector.cxx, Matrix.cxx Tensor.cxx and SuperTensor.cxx
contain the actual implementations of the functions described in
Vector.h, Matrix.h Tensor.h and SuperTensor.h. The python scripts
testVector.py, testMatrix.py testTensor.py and testSuperTensor.py
test the resulting python wrappers using the unittest module.

The SWIG interface files Vector.i, Matrix.i Tensor.i and SuperTensor.i
are used to generate the wrapper code. The SWIG_FILE_WITH_INIT macro
allows numpy.i to be used with multiple python modules. If it is
specified, then the %init block found in Vector.i, Matrix.i Tensor.i
and SuperTensor.i are required. The other things done in Vector.i,
Matrix.i Tensor.i and SuperTensor.i are the inclusion of the
appropriate header file and numpy.i file, and the "%apply" directives
to force the functions to use the typemaps.
The source files Vector.cxx, Matrix.cxx, Tensor.cxx, SuperTensor.cxx
and Flat.cxx contain the actual implementations of the functions
described in Vector.h, Matrix.h Tensor.h, SuperTensor.h and Flat.h.
The python scripts testVector.py, testMatrix.py testTensor.py,
testSuperTensor.py and testFlat.py test the resulting python wrappers
using the unittest module.

The SWIG interface files Vector.i, Matrix.i, Tensor.i, SuperTensor.i
and Flat.i are used to generate the wrapper code. The
SWIG_FILE_WITH_INIT macro allows numpy.i to be used with multiple
python modules. If it is specified, then the %init block found in
Vector.i, Matrix.i Tensor.i, SuperTensor.i and Flat.i are required.
The other things done in Vector.i, Matrix.i, Tensor.i, SuperTensor.i
and Flat.i are the inclusion of the appropriate header file and
numpy.i file, and the "%apply" directives to force the functions to
use the typemaps.

The setup.py script is a standard python distutils script. It defines
_Vector, _Matrix _Tensor and _SuperTensor extension modules and Vector
, Matrix, Tensor and SuperTensor python modules. The Makefile
_Vector, _Matrix, _Tensor, _SuperTensor, _Flat extension modules and Vector,
Matrix, Tensor, SuperTensor and Flat python modules. The Makefile
automates everything, setting up the dependencies, calling swig to
generate the wrappers, and calling setup.py to compile the wrapper
code and generate the shared objects.
Targets "all" (default), "test", "doc" and "clean" are supported. The
"doc" target creates HTML documentation (with make target "html"), and
PDF documentation (with make targets "tex" and "pdf").
code and generate the shared objects. Targets "all" (default), "test",
"doc" and "clean" are supported. The "doc" target creates HTML
documentation (with make target "html"), and PDF documentation
(with make targets "tex" and "pdf").

To build and run the test code, simply execute from the shell::

Expand Down
46 changes: 45 additions & 1 deletion tools/swig/numpy.i
Expand Up @@ -381,6 +381,22 @@
return contiguous;
}

/* Test whether a python object is (C_ or F_) contiguous. If array is
* contiguous, return 1. Otherwise, set the python error string and
* return 0.
*/
int require_c_or_f_contiguous(PyArrayObject* ary)
{
int contiguous = 1;
if (!(array_is_contiguous(ary) || array_is_fortran(ary)))
{
PyErr_SetString(PyExc_TypeError,
"Array must be contiguous (C_ or F_). A non-contiguous array was given");
contiguous = 0;
}
return contiguous;
}

/* Require that a numpy array is not byte-swapped. If the array is
* not byte-swapped, return 1. Otherwise, set the python error string
* and return 0.
Expand Down Expand Up @@ -543,7 +559,7 @@

/* %numpy_typemaps() macro
*
* This macro defines a family of 74 typemaps that allow C arguments
* This macro defines a family of 75 typemaps that allow C arguments
* of the form
*
* 1. (DATA_TYPE IN_ARRAY1[ANY])
Expand Down Expand Up @@ -640,6 +656,8 @@
* 73. (DATA_TYPE** ARGOUTVIEWM_FARRAY4, DIM_TYPE* DIM1, DIM_TYPE* DIM2, DIM_TYPE* DIM3, DIM_TYPE* DIM4)
* 74. (DIM_TYPE* DIM1, DIM_TYPE* DIM2, DIM_TYPE* DIM3, DIM_TYPE* DIM4, DATA_TYPE** ARGOUTVIEWM_FARRAY4)
*
* 75. (DATA_TYPE* INPLACE_ARRAY_FLAT, DIM_TYPE DIM_FLAT)
*
* where "DATA_TYPE" is any type supported by the NumPy module, and
* "DIM_TYPE" is any int-like type suitable for specifying dimensions.
* The difference between "ARRAY" typemaps and "FARRAY" typemaps is
Expand Down Expand Up @@ -3072,6 +3090,32 @@
$result = SWIG_Python_AppendOutput($result,obj);
}

/**************************************/
/* In-Place Array Typemap - flattened */
/**************************************/

/* Typemap suite for (DATA_TYPE* INPLACE_ARRAY_FLAT, DIM_TYPE DIM_FLAT)
*/
%typecheck(SWIG_TYPECHECK_DOUBLE_ARRAY,
fragment="NumPy_Macros")
(DATA_TYPE* INPLACE_ARRAY_FLAT, DIM_TYPE DIM_FLAT)
{
$1 = is_array($input) && PyArray_EquivTypenums(array_type($input),
DATA_TYPECODE);
}
%typemap(in,
fragment="NumPy_Fragments")
(DATA_TYPE* INPLACE_ARRAY_FLAT, DIM_TYPE DIM_FLAT)
(PyArrayObject* array=NULL, int i=1)
{
array = obj_to_array_no_conversion($input, DATA_TYPECODE);
if (!array || !require_c_or_f_contiguous(array)
|| !require_native(array)) SWIG_fail;
$1 = (DATA_TYPE*) array_data(array);
$2 = 1;
for (i=0; i < array_numdims(array); ++i) $2 *= array_size(array,i);
}

%enddef /* %numpy_typemaps() macro */
/* *************************************************************** */

Expand Down
36 changes: 36 additions & 0 deletions tools/swig/test/Flat.cxx
@@ -0,0 +1,36 @@
#include <stdlib.h>
#include <math.h>
#include <iostream>
#include "Flat.h"

// The following macro defines a family of functions that work with 1D
// arrays with the forms
//
// void SNAMEProcess(TYPE * array, int size);
//
// for any specified type TYPE (for example: short, unsigned int, long
// long, etc.) with given short name SNAME (for example: short, uint,
// longLong, etc.). The macro is then expanded for the given
// TYPE/SNAME pairs. The resulting functions are for testing numpy
// interfaces for:
//
// * in-place arrays (arbitrary number of dimensions) with a fixed number of elements
//
#define TEST_FUNCS(TYPE, SNAME) \
\
void SNAME ## Process(TYPE * array, int size) { \
for (int i=0; i<size; ++i) array[i] += 1; \
}

TEST_FUNCS(signed char , schar )
TEST_FUNCS(unsigned char , uchar )
TEST_FUNCS(short , short )
TEST_FUNCS(unsigned short , ushort )
TEST_FUNCS(int , int )
TEST_FUNCS(unsigned int , uint )
TEST_FUNCS(long , long )
TEST_FUNCS(unsigned long , ulong )
TEST_FUNCS(long long , longLong )
TEST_FUNCS(unsigned long long, ulongLong)
TEST_FUNCS(float , float )
TEST_FUNCS(double , double )
34 changes: 34 additions & 0 deletions tools/swig/test/Flat.h
@@ -0,0 +1,34 @@
#ifndef FLAT_H
#define FLAT_H

// The following macro defines the prototypes for a family of
// functions that work with arrays with the forms
//
// void SNAMEProcess(TYPE * array, int size);
//
// for any specified type TYPE (for example: short, unsigned int, long
// long, etc.) with given short name SNAME (for example: short, uint,
// longLong, etc.). The macro is then expanded for the given
// TYPE/SNAME pairs. The resulting functions are for testing numpy
// interfaces for:
//
// * in-place arrays (arbitrary number of dimensions) with a fixed number of elements
//
#define TEST_FUNC_PROTOS(TYPE, SNAME) \
\
void SNAME ## Process(TYPE * array, int size); \

TEST_FUNC_PROTOS(signed char , schar )
TEST_FUNC_PROTOS(unsigned char , uchar )
TEST_FUNC_PROTOS(short , short )
TEST_FUNC_PROTOS(unsigned short , ushort )
TEST_FUNC_PROTOS(int , int )
TEST_FUNC_PROTOS(unsigned int , uint )
TEST_FUNC_PROTOS(long , long )
TEST_FUNC_PROTOS(unsigned long , ulong )
TEST_FUNC_PROTOS(long long , longLong )
TEST_FUNC_PROTOS(unsigned long long, ulongLong)
TEST_FUNC_PROTOS(float , float )
TEST_FUNC_PROTOS(double , double )

#endif
36 changes: 36 additions & 0 deletions tools/swig/test/Flat.i
@@ -0,0 +1,36 @@
// -*- c++ -*-
%module Flat

%{
#define SWIG_FILE_WITH_INIT
#include "Flat.h"
%}

// Get the NumPy typemaps
%include "../numpy.i"

%init %{
import_array();
%}

%define %apply_numpy_typemaps(TYPE)

%apply (TYPE* INPLACE_ARRAY_FLAT, int DIM_FLAT) {(TYPE* array, int size)};

%enddef /* %apply_numpy_typemaps() macro */

%apply_numpy_typemaps(signed char )
%apply_numpy_typemaps(unsigned char )
%apply_numpy_typemaps(short )
%apply_numpy_typemaps(unsigned short )
%apply_numpy_typemaps(int )
%apply_numpy_typemaps(unsigned int )
%apply_numpy_typemaps(long )
%apply_numpy_typemaps(unsigned long )
%apply_numpy_typemaps(long long )
%apply_numpy_typemaps(unsigned long long)
%apply_numpy_typemaps(float )
%apply_numpy_typemaps(double )

// Include the header file to be wrapped
%include "Flat.h"
6 changes: 4 additions & 2 deletions tools/swig/test/Makefile
@@ -1,13 +1,14 @@
# SWIG
INTERFACES = Array.i Farray.i Vector.i Matrix.i Tensor.i Fortran.i
INTERFACES = Array.i Farray.i Vector.i Matrix.i Tensor.i Fortran.i Flat.i
WRAPPERS = $(INTERFACES:.i=_wrap.cxx)
PROXIES = $(INTERFACES:.i=.py )

# Default target: build the tests
.PHONY : all
all: $(WRAPPERS) Array1.cxx Array1.h Array2.cxx Array2.h ArrayZ.cxx ArrayZ.h \
Farray.cxx Farray.h Vector.cxx Vector.h \
Matrix.cxx Matrix.h Tensor.cxx Tensor.h Fortran.h Fortran.cxx
Matrix.cxx Matrix.h Tensor.cxx Tensor.h Fortran.h Fortran.cxx \
Flat.h Flat.cxx
./setup.py build_ext -i

# Test target: run the tests
Expand All @@ -19,6 +20,7 @@ test: all
python testArray.py
python testFarray.py
python testFortran.py
python testFlat.py

# Rule: %.i -> %_wrap.cxx
%_wrap.cxx: %.i %.h ../numpy.i
Expand Down
10 changes: 8 additions & 2 deletions tools/swig/test/setup.py
Expand Up @@ -54,12 +54,18 @@
include_dirs = [numpy_include],
)

_Flat = Extension("_Flat",
["Flat_wrap.cxx",
"Flat.cxx"],
include_dirs = [numpy_include],
)

# NumyTypemapTests setup
setup(name = "NumpyTypemapTests",
description = "Functions that work on arrays",
author = "Bill Spotz",
py_modules = ["Array", "Farray", "Vector", "Matrix", "Tensor",
"Fortran"],
"Fortran", "Flat"],
ext_modules = [_Array, _Farray, _Vector, _Matrix, _Tensor,
_Fortran]
_Fortran, _Flat]
)

0 comments on commit a49ad15

Please sign in to comment.