Skip to content

Commit

Permalink
Merge pull request #2694 from seberg/cflags2
Browse files Browse the repository at this point in the history
ENH: Make 1-dimensional axes not matter for contiguous flags
  • Loading branch information
njsmith committed Oct 25, 2012
2 parents c5ccca9 + 02ebf8b commit a890a85
Show file tree
Hide file tree
Showing 13 changed files with 154 additions and 177 deletions.
7 changes: 4 additions & 3 deletions numpy/core/include/numpy/ndarraytypes.h
Expand Up @@ -756,8 +756,9 @@ typedef int (PyArray_FinalizeFunc)(PyArrayObject *, PyObject *);
#define NPY_ARRAY_F_CONTIGUOUS 0x0002 #define NPY_ARRAY_F_CONTIGUOUS 0x0002


/* /*
* Note: all 0-d arrays are C_CONTIGUOUS and F_CONTIGUOUS. If a * Note: all 0-d arrays are C_CONTIGUOUS and F_CONTIGUOUS. An N-d
* 1-d array is C_CONTIGUOUS it is also F_CONTIGUOUS * array that is C_CONTIGUOUS is also F_CONTIGUOUS if only
* one axis has a dimension different from one (ie. a 1x3x1 array).
*/ */


/* /*
Expand Down Expand Up @@ -1370,7 +1371,7 @@ PyArrayNeighborhoodIter_Next2D(PyArrayNeighborhoodIterObject* iter);
PyArray_CHKFLAGS(m, NPY_ARRAY_F_CONTIGUOUS)) PyArray_CHKFLAGS(m, NPY_ARRAY_F_CONTIGUOUS))


#define PyArray_ISFORTRAN(m) (PyArray_CHKFLAGS(m, NPY_ARRAY_F_CONTIGUOUS) && \ #define PyArray_ISFORTRAN(m) (PyArray_CHKFLAGS(m, NPY_ARRAY_F_CONTIGUOUS) && \
(PyArray_NDIM(m) > 1)) (!PyArray_CHKFLAGS(m, NPY_ARRAY_C_CONTIGUOUS)))


#define PyArray_FORTRAN_IF(m) ((PyArray_CHKFLAGS(m, NPY_ARRAY_F_CONTIGUOUS) ? \ #define PyArray_FORTRAN_IF(m) ((PyArray_CHKFLAGS(m, NPY_ARRAY_F_CONTIGUOUS) ? \
NPY_ARRAY_F_CONTIGUOUS : 0)) NPY_ARRAY_F_CONTIGUOUS : 0))
Expand Down
4 changes: 2 additions & 2 deletions numpy/core/numeric.py
Expand Up @@ -538,7 +538,7 @@ def require(a, dtype=None, requirements=None):
def isfortran(a): def isfortran(a):
""" """
Returns True if array is arranged in Fortran-order in memory Returns True if array is arranged in Fortran-order in memory
and dimension > 1. and not C-order.
Parameters Parameters
---------- ----------
Expand Down Expand Up @@ -584,7 +584,7 @@ def isfortran(a):
>>> np.isfortran(b) >>> np.isfortran(b)
True True
1-D arrays always evaluate as False. C-ordered arrays evaluate as False even if they are also FORTRAN-ordered.
>>> np.isfortran(np.array([1, 2], order='FORTRAN')) >>> np.isfortran(np.array([1, 2], order='FORTRAN'))
False False
Expand Down
4 changes: 2 additions & 2 deletions numpy/core/src/multiarray/convert.c
Expand Up @@ -265,8 +265,8 @@ PyArray_ToString(PyArrayObject *self, NPY_ORDER order)
*/ */


numbytes = PyArray_NBYTES(self); numbytes = PyArray_NBYTES(self);
if ((PyArray_ISCONTIGUOUS(self) && (order == NPY_CORDER)) if ((PyArray_IS_C_CONTIGUOUS(self) && (order == NPY_CORDER))
|| (PyArray_ISFORTRAN(self) && (order == NPY_FORTRANORDER))) { || (PyArray_IS_F_CONTIGUOUS(self) && (order == NPY_FORTRANORDER))) {
ret = PyBytes_FromStringAndSize(PyArray_DATA(self), (Py_ssize_t) numbytes); ret = PyBytes_FromStringAndSize(PyArray_DATA(self), (Py_ssize_t) numbytes);
} }
else { else {
Expand Down
41 changes: 30 additions & 11 deletions numpy/core/src/multiarray/ctors.c
Expand Up @@ -1112,7 +1112,6 @@ PyArray_NewLikeArray(PyArrayObject *prototype, NPY_ORDER order,
int idim; int idim;


PyArray_CreateSortedStridePerm(PyArray_NDIM(prototype), PyArray_CreateSortedStridePerm(PyArray_NDIM(prototype),
PyArray_SHAPE(prototype),
PyArray_STRIDES(prototype), PyArray_STRIDES(prototype),
strideperm); strideperm);


Expand Down Expand Up @@ -1825,9 +1824,6 @@ PyArray_FromArray(PyArrayObject *arr, PyArray_Descr *newtype, int flags)
} }


arrflags = PyArray_FLAGS(arr); arrflags = PyArray_FLAGS(arr);
if (PyArray_NDIM(arr) <= 1 && (flags & NPY_ARRAY_F_CONTIGUOUS)) {
flags |= NPY_ARRAY_C_CONTIGUOUS;
}
/* If a guaranteed copy was requested */ /* If a guaranteed copy was requested */
copy = (flags & NPY_ARRAY_ENSURECOPY) || copy = (flags & NPY_ARRAY_ENSURECOPY) ||
/* If C contiguous was requested, and arr is not */ /* If C contiguous was requested, and arr is not */
Expand All @@ -1837,9 +1833,8 @@ PyArray_FromArray(PyArrayObject *arr, PyArray_Descr *newtype, int flags)
((flags & NPY_ARRAY_ALIGNED) && ((flags & NPY_ARRAY_ALIGNED) &&
(!(arrflags & NPY_ARRAY_ALIGNED))) || (!(arrflags & NPY_ARRAY_ALIGNED))) ||
/* If a Fortran contiguous array was requested, and arr is not */ /* If a Fortran contiguous array was requested, and arr is not */
(PyArray_NDIM(arr) > 1 && ((flags & NPY_ARRAY_F_CONTIGUOUS) &&
((flags & NPY_ARRAY_F_CONTIGUOUS) && (!(arrflags & NPY_ARRAY_F_CONTIGUOUS))) ||
(!(arrflags & NPY_ARRAY_F_CONTIGUOUS)))) ||
/* If a writeable array was requested, and arr is not */ /* If a writeable array was requested, and arr is not */
((flags & NPY_ARRAY_WRITEABLE) && ((flags & NPY_ARRAY_WRITEABLE) &&
(!(arrflags & NPY_ARRAY_WRITEABLE))) || (!(arrflags & NPY_ARRAY_WRITEABLE))) ||
Expand Down Expand Up @@ -3570,14 +3565,33 @@ _array_fill_strides(npy_intp *strides, npy_intp *dims, int nd, size_t itemsize,
int inflag, int *objflags) int inflag, int *objflags)
{ {
int i; int i;
npy_bool not_cf_contig = 0;
npy_bool nod = 0; /* A dim != 1 was found */

/* Check if new array is both F- and C-contiguous */
for (i = 0; i < nd; i++) {
if (dims[i] != 1) {
if (nod) {
not_cf_contig = 1;
break;
}
nod = 1;
}
}

/* Only make Fortran strides if not contiguous as well */ /* Only make Fortran strides if not contiguous as well */
if ((inflag & (NPY_ARRAY_F_CONTIGUOUS|NPY_ARRAY_C_CONTIGUOUS)) == if ((inflag & (NPY_ARRAY_F_CONTIGUOUS|NPY_ARRAY_C_CONTIGUOUS)) ==
NPY_ARRAY_F_CONTIGUOUS) { NPY_ARRAY_F_CONTIGUOUS) {
for (i = 0; i < nd; i++) { for (i = 0; i < nd; i++) {
strides[i] = itemsize; strides[i] = itemsize;
itemsize *= dims[i] ? dims[i] : 1; if (dims[i]) {
itemsize *= dims[i];
}
else {
not_cf_contig = 0;
}
} }
if (nd > 1) { if (not_cf_contig) {
*objflags = ((*objflags)|NPY_ARRAY_F_CONTIGUOUS) & *objflags = ((*objflags)|NPY_ARRAY_F_CONTIGUOUS) &
~NPY_ARRAY_C_CONTIGUOUS; ~NPY_ARRAY_C_CONTIGUOUS;
} }
Expand All @@ -3588,9 +3602,14 @@ _array_fill_strides(npy_intp *strides, npy_intp *dims, int nd, size_t itemsize,
else { else {
for (i = nd - 1; i >= 0; i--) { for (i = nd - 1; i >= 0; i--) {
strides[i] = itemsize; strides[i] = itemsize;
itemsize *= dims[i] ? dims[i] : 1; if (dims[i]) {
itemsize *= dims[i];
}
else {
not_cf_contig = 0;
}
} }
if (nd > 1) { if (not_cf_contig) {
*objflags = ((*objflags)|NPY_ARRAY_C_CONTIGUOUS) & *objflags = ((*objflags)|NPY_ARRAY_C_CONTIGUOUS) &
~NPY_ARRAY_F_CONTIGUOUS; ~NPY_ARRAY_F_CONTIGUOUS;
} }
Expand Down
6 changes: 3 additions & 3 deletions numpy/core/src/multiarray/dtype_transfer.c
Expand Up @@ -3921,7 +3921,7 @@ PyArray_PrepareOneRawArrayIter(int ndim, npy_intp *shape,
} }


/* Sort the axes based on the destination strides */ /* Sort the axes based on the destination strides */
PyArray_CreateSortedStridePerm(ndim, shape, strides, strideperm); PyArray_CreateSortedStridePerm(ndim, strides, strideperm);
for (i = 0; i < ndim; ++i) { for (i = 0; i < ndim; ++i) {
int iperm = strideperm[ndim - i - 1].perm; int iperm = strideperm[ndim - i - 1].perm;
out_shape[i] = shape[iperm]; out_shape[i] = shape[iperm];
Expand Down Expand Up @@ -4051,7 +4051,7 @@ PyArray_PrepareTwoRawArrayIter(int ndim, npy_intp *shape,
} }


/* Sort the axes based on the destination strides */ /* Sort the axes based on the destination strides */
PyArray_CreateSortedStridePerm(ndim, shape, stridesA, strideperm); PyArray_CreateSortedStridePerm(ndim, stridesA, strideperm);
for (i = 0; i < ndim; ++i) { for (i = 0; i < ndim; ++i) {
int iperm = strideperm[ndim - i - 1].perm; int iperm = strideperm[ndim - i - 1].perm;
out_shape[i] = shape[iperm]; out_shape[i] = shape[iperm];
Expand Down Expand Up @@ -4185,7 +4185,7 @@ PyArray_PrepareThreeRawArrayIter(int ndim, npy_intp *shape,
} }


/* Sort the axes based on the destination strides */ /* Sort the axes based on the destination strides */
PyArray_CreateSortedStridePerm(ndim, shape, stridesA, strideperm); PyArray_CreateSortedStridePerm(ndim, stridesA, strideperm);
for (i = 0; i < ndim; ++i) { for (i = 0; i < ndim; ++i) {
int iperm = strideperm[ndim - i - 1].perm; int iperm = strideperm[ndim - i - 1].perm;
out_shape[i] = shape[iperm]; out_shape[i] = shape[iperm];
Expand Down
100 changes: 34 additions & 66 deletions numpy/core/src/multiarray/flagsobject.c
Expand Up @@ -15,11 +15,8 @@


#include "common.h" #include "common.h"


static int static void
_IsContiguous(PyArrayObject *ap); _UpdateContiguousFlags(PyArrayObject *ap);

static int
_IsFortranContiguous(PyArrayObject *ap);


/*NUMPY_API /*NUMPY_API
* *
Expand Down Expand Up @@ -62,28 +59,9 @@ PyArray_NewFlagsObject(PyObject *obj)
NPY_NO_EXPORT void NPY_NO_EXPORT void
PyArray_UpdateFlags(PyArrayObject *ret, int flagmask) PyArray_UpdateFlags(PyArrayObject *ret, int flagmask)
{ {

/* Always update both, as its not trivial to guess one from the other */
if (flagmask & NPY_ARRAY_F_CONTIGUOUS) { if (flagmask & (NPY_ARRAY_F_CONTIGUOUS | NPY_ARRAY_C_CONTIGUOUS)) {
if (_IsFortranContiguous(ret)) { _UpdateContiguousFlags(ret);
PyArray_ENABLEFLAGS(ret, NPY_ARRAY_F_CONTIGUOUS);
if (PyArray_NDIM(ret) > 1) {
PyArray_CLEARFLAGS(ret, NPY_ARRAY_C_CONTIGUOUS);
}
}
else {
PyArray_CLEARFLAGS(ret, NPY_ARRAY_F_CONTIGUOUS);
}
}
if (flagmask & NPY_ARRAY_C_CONTIGUOUS) {
if (_IsContiguous(ret)) {
PyArray_ENABLEFLAGS(ret, NPY_ARRAY_C_CONTIGUOUS);
if (PyArray_NDIM(ret) > 1) {
PyArray_CLEARFLAGS(ret, NPY_ARRAY_F_CONTIGUOUS);
}
}
else {
PyArray_CLEARFLAGS(ret, NPY_ARRAY_C_CONTIGUOUS);
}
} }
if (flagmask & NPY_ARRAY_ALIGNED) { if (flagmask & NPY_ARRAY_ALIGNED) {
if (_IsAligned(ret)) { if (_IsAligned(ret)) {
Expand All @@ -110,66 +88,56 @@ PyArray_UpdateFlags(PyArrayObject *ret, int flagmask)


/* /*
* Check whether the given array is stored contiguously * Check whether the given array is stored contiguously
* (row-wise) in memory. * in memory. And update the passed in ap flags apropriately.
* *
* 0-strided arrays are not contiguous (even if dimension == 1) * A dimension == 1 stride is ignored for contiguous flags and a 0-sized array
* is always both C- and F-Contiguous. 0-strided arrays are not contiguous.
*/ */
static int static void
_IsContiguous(PyArrayObject *ap) _UpdateContiguousFlags(PyArrayObject *ap)
{ {
npy_intp sd; npy_intp sd;
npy_intp dim; npy_intp dim;
int i; int i;
npy_bool is_c_contig = 1;


if (PyArray_NDIM(ap) == 0) {
return 1;
}
sd = PyArray_DESCR(ap)->elsize; sd = PyArray_DESCR(ap)->elsize;
if (PyArray_NDIM(ap) == 1) {
return PyArray_DIMS(ap)[0] == 1 || sd == PyArray_STRIDES(ap)[0];
}
for (i = PyArray_NDIM(ap) - 1; i >= 0; --i) { for (i = PyArray_NDIM(ap) - 1; i >= 0; --i) {
dim = PyArray_DIMS(ap)[i]; dim = PyArray_DIMS(ap)[i];
/* contiguous by definition */ /* contiguous by definition */
if (dim == 0) { if (dim == 0) {
return 1; PyArray_ENABLEFLAGS(ap, NPY_ARRAY_C_CONTIGUOUS);
PyArray_ENABLEFLAGS(ap, NPY_ARRAY_F_CONTIGUOUS);
return;
} }
if (PyArray_STRIDES(ap)[i] != sd) { if (dim != 1) {
return 0; if (PyArray_STRIDES(ap)[i] != sd) {
is_c_contig = 0;
}
sd *= dim;
} }
sd *= dim;
} }
return 1; if (is_c_contig) {
} PyArray_ENABLEFLAGS(ap, NPY_ARRAY_C_CONTIGUOUS);


/* 0-strided arrays are not contiguous (even if dimension == 1) */
static int
_IsFortranContiguous(PyArrayObject *ap)
{
npy_intp sd;
npy_intp dim;
int i;

if (PyArray_NDIM(ap) == 0) {
return 1;
} }
sd = PyArray_DESCR(ap)->elsize; else {
if (PyArray_NDIM(ap) == 1) { PyArray_CLEARFLAGS(ap, NPY_ARRAY_C_CONTIGUOUS);
return PyArray_DIMS(ap)[0] == 1 || sd == PyArray_STRIDES(ap)[0];
} }

/* check if fortran contiguous */
sd = PyArray_DESCR(ap)->elsize;
for (i = 0; i < PyArray_NDIM(ap); ++i) { for (i = 0; i < PyArray_NDIM(ap); ++i) {
dim = PyArray_DIMS(ap)[i]; dim = PyArray_DIMS(ap)[i];
/* fortran contiguous by definition */ if (dim != 1) {
if (dim == 0) { if (PyArray_STRIDES(ap)[i] != sd) {
return 1; PyArray_CLEARFLAGS(ap, NPY_ARRAY_F_CONTIGUOUS);
} return;
if (PyArray_STRIDES(ap)[i] != sd) { }
return 0; sd *= dim;
} }
sd *= dim;
} }
return 1; PyArray_ENABLEFLAGS(ap, NPY_ARRAY_F_CONTIGUOUS);
return;
} }


static void static void
Expand Down
18 changes: 5 additions & 13 deletions numpy/core/src/multiarray/multiarraymodule.c
Expand Up @@ -1517,26 +1517,18 @@ PyArray_EquivTypenums(int typenum1, int typenum2)
/*** END C-API FUNCTIONS **/ /*** END C-API FUNCTIONS **/


static PyObject * static PyObject *
_prepend_ones(PyArrayObject *arr, int nd, int ndmin, NPY_ORDER order) _prepend_ones(PyArrayObject *arr, int nd, int ndmin)
{ {
npy_intp newdims[NPY_MAXDIMS]; npy_intp newdims[NPY_MAXDIMS];
npy_intp newstrides[NPY_MAXDIMS]; npy_intp newstrides[NPY_MAXDIMS];
npy_intp newstride;
int i, k, num; int i, k, num;
PyArrayObject *ret; PyArrayObject *ret;
PyArray_Descr *dtype; PyArray_Descr *dtype;


if (order == NPY_FORTRANORDER || PyArray_ISFORTRAN(arr) || PyArray_NDIM(arr) == 0) {
newstride = PyArray_DESCR(arr)->elsize;
}
else {
newstride = PyArray_STRIDES(arr)[0] * PyArray_DIMS(arr)[0];
}

num = ndmin - nd; num = ndmin - nd;
for (i = 0; i < num; i++) { for (i = 0; i < num; i++) {
newdims[i] = 1; newdims[i] = 1;
newstrides[i] = newstride; newstrides[i] = PyArray_DESCR(arr)->elsize;
} }
for (i = num; i < ndmin; i++) { for (i = num; i < ndmin; i++) {
k = i - num; k = i - num;
Expand Down Expand Up @@ -1566,8 +1558,8 @@ _prepend_ones(PyArrayObject *arr, int nd, int ndmin, NPY_ORDER order)
#define STRIDING_OK(op, order) \ #define STRIDING_OK(op, order) \
((order) == NPY_ANYORDER || \ ((order) == NPY_ANYORDER || \
(order) == NPY_KEEPORDER || \ (order) == NPY_KEEPORDER || \
((order) == NPY_CORDER && PyArray_ISCONTIGUOUS(op)) || \ ((order) == NPY_CORDER && PyArray_IS_C_CONTIGUOUS(op)) || \
((order) == NPY_FORTRANORDER && PyArray_ISFORTRAN(op))) ((order) == NPY_FORTRANORDER && PyArray_IS_F_CONTIGUOUS(op)))


static PyObject * static PyObject *
_array_fromobject(PyObject *NPY_UNUSED(ignored), PyObject *args, PyObject *kws) _array_fromobject(PyObject *NPY_UNUSED(ignored), PyObject *args, PyObject *kws)
Expand Down Expand Up @@ -1677,7 +1669,7 @@ _array_fromobject(PyObject *NPY_UNUSED(ignored), PyObject *args, PyObject *kws)
* create a new array from the same data with ones in the shape * create a new array from the same data with ones in the shape
* steals a reference to ret * steals a reference to ret
*/ */
return _prepend_ones(ret, nd, ndmin, order); return _prepend_ones(ret, nd, ndmin);


clean_type: clean_type:
Py_XDECREF(type); Py_XDECREF(type);
Expand Down

0 comments on commit a890a85

Please sign in to comment.