summaryrefslogtreecommitdiff
diff options
context:
space:
mode:
-rw-r--r--doc/source/reference/c-api.array.rst8
-rw-r--r--doc/source/reference/routines.array-manipulation.rst7
-rw-r--r--numpy/add_newdocs.py53
-rw-r--r--numpy/core/code_generators/numpy_api.py5
-rw-r--r--numpy/core/numeric.py3
-rw-r--r--numpy/core/src/multiarray/ctors.c275
-rw-r--r--numpy/core/src/multiarray/item_selection.c5
-rw-r--r--numpy/core/src/multiarray/multiarraymodule.c77
-rw-r--r--numpy/core/src/umath/ufunc_type_resolution.c4
-rw-r--r--numpy/core/tests/test_api.py36
10 files changed, 468 insertions, 5 deletions
diff --git a/doc/source/reference/c-api.array.rst b/doc/source/reference/c-api.array.rst
index 02aa03ff0..51802c436 100644
--- a/doc/source/reference/c-api.array.rst
+++ b/doc/source/reference/c-api.array.rst
@@ -475,6 +475,10 @@ From other objects
indicates that the array should be aligned in the sense that the
strides are multiples of the element size.
+ In versions 1.6 and earlier of NumPy, the following flags
+ did not have the _ARRAY_ macro namespace in them. That form
+ of the constant names is deprecated in 1.7.
+
.. cvar:: NPY_ARRAY_NOTSWAPPED
Make sure the returned array has a data-type descriptor that is in
@@ -1210,6 +1214,10 @@ might not be writeable. It might be in Fortan-contiguous order. The
array flags are used to indicate what can be said about data
associated with an array.
+In versions 1.6 and earlier of NumPy, the following flags
+did not have the _ARRAY_ macro namespace in them. That form
+of the constant names is deprecated in 1.7.
+
.. cvar:: NPY_ARRAY_C_CONTIGUOUS
The data area is in C-style contiguous order (last index varies the
diff --git a/doc/source/reference/routines.array-manipulation.rst b/doc/source/reference/routines.array-manipulation.rst
index 2c1a5b200..ca97bb479 100644
--- a/doc/source/reference/routines.array-manipulation.rst
+++ b/doc/source/reference/routines.array-manipulation.rst
@@ -3,6 +3,13 @@ Array manipulation routines
.. currentmodule:: numpy
+Basic operations
+================
+.. autosummary::
+ :toctree: generated/
+
+ copyto
+
Changing array shape
====================
.. autosummary::
diff --git a/numpy/add_newdocs.py b/numpy/add_newdocs.py
index 2f9071378..c73845ee7 100644
--- a/numpy/add_newdocs.py
+++ b/numpy/add_newdocs.py
@@ -3184,6 +3184,10 @@ add_newdoc('numpy.core.multiarray', 'ndarray', ('copy',
order. If order is 'A' ('Any'), then the result has the same order
as the input.
+ See also
+ --------
+ numpy.copyto
+
Examples
--------
>>> x = np.array([[1,2,3],[4,5,6]], order='F')
@@ -3690,11 +3694,50 @@ add_newdoc('numpy.core.multiarray', 'ndarray', ('put',
"""))
+add_newdoc('numpy.core.multiarray', 'copyto',
+ """
+ copyto(dst, src, casting='same_kind', where=None)
+
+ Copies values from `src` into `dst`, broadcasting as necessary.
+ Raises a TypeError if the casting rule is violated, and if
+ `where` is provided, it selects which elements to copy.
+
+ .. versionadded:: 1.7.0
+
+ Parameters
+ ----------
+ dst : ndarray
+ The array into which values are copied.
+ src : array_like
+ The array from which values are copied.
+ casting : {'no', 'equiv', 'safe', 'same_kind', 'unsafe'}, optional
+ Controls what kind of data casting may occur when copying.
+
+ * 'no' means the data types should not be cast at all.
+ * 'equiv' means only byte-order changes are allowed.
+ * 'safe' means only casts which can preserve values are allowed.
+ * 'same_kind' means only safe casts or casts within a kind,
+ like float64 to float32, are allowed.
+ * 'unsafe' means any data conversions may be done.
+ where : array_like of bool
+ A boolean array which is broadcasted to match the dimensions
+ of `dst`, and selects elements to copy from `src` to `dst`
+ wherever it contains the value True.
+
+ Returns
+ -------
+ out : ndarray
+ Returns the array `dst`.
+
+ """)
add_newdoc('numpy.core.multiarray', 'putmask',
"""
putmask(a, mask, values)
+ This function is deprecated as of NumPy 1.7. Use the function
+ ``np.copy(a, values, where=mask)`` to achieve this functionality.
+
Changes elements of an array based on conditional and input values.
Sets ``a.flat[n] = values[n]`` for each n where ``mask.flat[n]==True``.
@@ -3714,7 +3757,7 @@ add_newdoc('numpy.core.multiarray', 'putmask',
See Also
--------
- place, put, take
+ place, put, take, copyto
Examples
--------
@@ -5959,6 +6002,8 @@ add_newdoc('numpy.core.multiarray', 'busdaycalendar',
information defining business days for the business
day-related functions.
+ .. versionadded:: 1.7.0
+
Parameters
----------
weekmask : str or array_like of bool, optional
@@ -6018,6 +6063,8 @@ add_newdoc('numpy.core.multiarray', 'is_busday',
Calculates which of the given dates are valid business days, and
which are not.
+ .. versionadded:: 1.7.0
+
Parameters
----------
dates : array_like of datetime64[D]
@@ -6070,6 +6117,8 @@ add_newdoc('numpy.core.multiarray', 'busday_offset',
the ``roll`` rule, then applies offsets to the given dates
counted in business days.
+ .. versionadded:: 1.7.0
+
Parameters
----------
dates : array_like of datetime64[D]
@@ -6158,6 +6207,8 @@ add_newdoc('numpy.core.multiarray', 'busday_count',
Counts the number of business days between `begindates` and
`enddates`, not including the day of `enddates`.
+ .. versionadded:: 1.7.0
+
Parameters
----------
begindates : array_like of datetime64[D]
diff --git a/numpy/core/code_generators/numpy_api.py b/numpy/core/code_generators/numpy_api.py
index b8cde7b3b..d83371319 100644
--- a/numpy/core/code_generators/numpy_api.py
+++ b/numpy/core/code_generators/numpy_api.py
@@ -256,7 +256,7 @@ multiarray_funcs_api = {
'PyArray_TimedeltaToTimedeltaStruct': 221,
'PyArray_DatetimeStructToDatetime': 222,
'PyArray_TimedeltaStructToTimedelta': 223,
- # New Iterator API
+ # NDIter API
'NpyIter_New': 224,
'NpyIter_MultiNew': 225,
'NpyIter_AdvancedNew': 226,
@@ -315,6 +315,9 @@ multiarray_funcs_api = {
'PyArray_GetArrayParamsFromObject': 278,
'PyArray_ConvertClipmodeSequence': 279,
'PyArray_MatrixProduct2': 280,
+ # End 1.6 API
+ 'PyArray_MaskedCopyInto': 281,
+ 'PyArray_MaskedMoveInto': 282,
}
ufunc_types_api = {
diff --git a/numpy/core/numeric.py b/numpy/core/numeric.py
index 67235e4d1..a28b6f3c4 100644
--- a/numpy/core/numeric.py
+++ b/numpy/core/numeric.py
@@ -1,7 +1,7 @@
__all__ = ['newaxis', 'ndarray', 'flatiter', 'nditer', 'nested_iters', 'ufunc',
'arange', 'array', 'zeros', 'count_nonzero', 'empty', 'broadcast',
'dtype', 'fromstring', 'fromfile', 'frombuffer',
- 'int_asbuffer', 'where', 'argwhere',
+ 'int_asbuffer', 'where', 'argwhere', 'copyto',
'concatenate', 'fastCopyAndTranspose', 'lexsort', 'set_numeric_ops',
'can_cast', 'promote_types', 'min_scalar_type', 'result_type',
'asarray', 'asanyarray', 'ascontiguousarray', 'asfortranarray',
@@ -58,6 +58,7 @@ nditer = multiarray.nditer
nested_iters = multiarray.nested_iters
broadcast = multiarray.broadcast
dtype = multiarray.dtype
+copyto = multiarray.copyto
ufunc = type(sin)
diff --git a/numpy/core/src/multiarray/ctors.c b/numpy/core/src/multiarray/ctors.c
index 1d9b83208..1d8270b3e 100644
--- a/numpy/core/src/multiarray/ctors.c
+++ b/numpy/core/src/multiarray/ctors.c
@@ -21,6 +21,7 @@
#include "buffer.h"
#include "numpymemoryview.h"
#include "lowlevel_strided_loops.h"
+#include "methods.h"
#include "_datetime.h"
#include "datetime_strings.h"
@@ -505,6 +506,91 @@ PyArray_MoveInto(PyArrayObject *dst, PyArrayObject *src)
}
}
+/*NUMPY_API
+ * Copy the memory of one array into another, allowing for overlapping data
+ * and selecting which elements to move based on a mask.
+ *
+ * Precisely handling the overlapping data is in general a difficult
+ * problem to solve efficiently, because strides can be negative.
+ * Consider "a = np.arange(3); a[::-1] = a", which previously produced
+ * the incorrect [0, 1, 0].
+ *
+ * Instead of trying to be fancy, we simply check for overlap and make
+ * a temporary copy when one exists.
+ *
+ * Returns 0 on success, negative on failure.
+ */
+NPY_NO_EXPORT int
+PyArray_MaskedMoveInto(PyArrayObject *dst, PyArrayObject *src,
+ PyArrayObject *mask, NPY_CASTING casting)
+{
+ /*
+ * Performance fix for expresions like "a[1000:6000] += x". In this
+ * case, first an in-place add is done, followed by an assignment,
+ * equivalently expressed like this:
+ *
+ * tmp = a[1000:6000] # Calls array_subscript_nice in mapping.c
+ * np.add(tmp, x, tmp)
+ * a[1000:6000] = tmp # Calls array_ass_sub in mapping.c
+ *
+ * In the assignment the underlying data type, shape, strides, and
+ * data pointers are identical, but src != dst because they are separately
+ * generated slices. By detecting this and skipping the redundant
+ * copy of values to themselves, we potentially give a big speed boost.
+ *
+ * Note that we don't call EquivTypes, because usually the exact same
+ * dtype object will appear, and we don't want to slow things down
+ * with a complicated comparison. The comparisons are ordered to
+ * try and reject this with as little work as possible.
+ */
+ if (PyArray_DATA(src) == PyArray_DATA(dst) &&
+ PyArray_DESCR(src) == PyArray_DESCR(dst) &&
+ PyArray_NDIM(src) == PyArray_NDIM(dst) &&
+ PyArray_CompareLists(PyArray_DIMS(src),
+ PyArray_DIMS(dst),
+ PyArray_NDIM(src)) &&
+ PyArray_CompareLists(PyArray_STRIDES(src),
+ PyArray_STRIDES(dst),
+ PyArray_NDIM(src))) {
+ /*printf("Redundant copy operation detected\n");*/
+ return 0;
+ }
+
+ /*
+ * A special case is when there is just one dimension with positive
+ * strides, and we pass that to CopyInto, which correctly handles
+ * it for most cases. It may still incorrectly handle copying of
+ * partially-overlapping data elements, where the data pointer was offset
+ * by a fraction of the element size.
+ */
+ if ((PyArray_NDIM(dst) == 1 &&
+ PyArray_NDIM(src) == 1 &&
+ PyArray_STRIDE(dst, 0) > 0 &&
+ PyArray_STRIDE(src, 0) > 0) ||
+ !_arrays_overlap(dst, src)) {
+ return PyArray_MaskedCopyInto(dst, src, mask, casting);
+ }
+ else {
+ PyArrayObject *tmp;
+ int ret;
+
+ /*
+ * Allocate a temporary copy array.
+ */
+ tmp = (PyArrayObject *)PyArray_NewLikeArray(dst,
+ NPY_KEEPORDER, NULL, 0);
+ if (tmp == NULL) {
+ return -1;
+ }
+ ret = PyArray_CopyInto(tmp, src);
+ if (ret == 0) {
+ ret = PyArray_MaskedCopyInto(dst, tmp, mask, casting);
+ }
+ Py_DECREF(tmp);
+ return ret;
+ }
+}
+
/* adapted from Numarray */
@@ -2798,6 +2884,195 @@ PyArray_CopyInto(PyArrayObject *dst, PyArrayObject *src)
}
}
+/*NUMPY_API
+ * Copy an Array into another array, wherever the mask specifies.
+ * The memory of src and dst must not overlap.
+ *
+ * Broadcast to the destination shape if necessary.
+ *
+ * Returns 0 on success, -1 on failure.
+ */
+NPY_NO_EXPORT int
+PyArray_MaskedCopyInto(PyArrayObject *dst, PyArrayObject *src,
+ PyArrayObject *mask, NPY_CASTING casting)
+{
+ PyArray_MaskedStridedTransferFn *stransfer = NULL;
+ NpyAuxData *transferdata = NULL;
+ NPY_BEGIN_THREADS_DEF;
+
+ if (!PyArray_ISWRITEABLE(dst)) {
+ PyErr_SetString(PyExc_RuntimeError,
+ "cannot write to array");
+ return -1;
+ }
+
+ if (!PyArray_CanCastArrayTo(src, PyArray_DESCR(dst), casting)) {
+ PyObject *errmsg;
+ errmsg = PyUString_FromString("Cannot cast array data from ");
+ PyUString_ConcatAndDel(&errmsg,
+ PyObject_Repr((PyObject *)PyArray_DESCR(src)));
+ PyUString_ConcatAndDel(&errmsg,
+ PyUString_FromString(" to "));
+ PyUString_ConcatAndDel(&errmsg,
+ PyObject_Repr((PyObject *)PyArray_DESCR(dst)));
+ PyUString_ConcatAndDel(&errmsg,
+ PyUString_FromFormat(" according to the rule %s",
+ npy_casting_to_string(casting)));
+ PyErr_SetObject(PyExc_TypeError, errmsg);
+ return -1;
+ }
+
+
+ if (PyArray_NDIM(dst) >= PyArray_NDIM(src) &&
+ PyArray_NDIM(dst) >= PyArray_NDIM(mask) &&
+ PyArray_TRIVIALLY_ITERABLE_TRIPLE(dst, src, mask)) {
+ char *dst_data, *src_data, *mask_data;
+ npy_intp count, dst_stride, src_stride, src_itemsize, mask_stride;
+
+ int needs_api = 0;
+
+ PyArray_PREPARE_TRIVIAL_TRIPLE_ITERATION(dst, src, mask, count,
+ dst_data, src_data, mask_data,
+ dst_stride, src_stride, mask_stride);
+
+ /*
+ * Check for overlap with positive strides, and if found,
+ * possibly reverse the order
+ */
+ if (dst_data > src_data && src_stride > 0 && dst_stride > 0 &&
+ (dst_data < src_data+src_stride*count) &&
+ (src_data < dst_data+dst_stride*count)) {
+ dst_data += dst_stride*(count-1);
+ src_data += src_stride*(count-1);
+ mask_data += mask_stride*(count-1);
+ dst_stride = -dst_stride;
+ src_stride = -src_stride;
+ mask_stride = -mask_stride;
+ }
+
+ if (PyArray_GetMaskedDTypeTransferFunction(
+ PyArray_ISALIGNED(src) && PyArray_ISALIGNED(dst),
+ src_stride, dst_stride, mask_stride,
+ PyArray_DESCR(src),
+ PyArray_DESCR(dst),
+ PyArray_DESCR(mask),
+ 0,
+ &stransfer, &transferdata,
+ &needs_api) != NPY_SUCCEED) {
+ return -1;
+ }
+
+ src_itemsize = PyArray_DESCR(src)->elsize;
+
+ if (!needs_api) {
+ NPY_BEGIN_THREADS;
+ }
+
+ stransfer(dst_data, dst_stride, src_data, src_stride,
+ (npy_uint8 *)mask_data, mask_stride,
+ count, src_itemsize, transferdata);
+
+ if (!needs_api) {
+ NPY_END_THREADS;
+ }
+
+ NPY_AUXDATA_FREE(transferdata);
+
+ return PyErr_Occurred() ? -1 : 0;
+ }
+ else {
+ PyArrayObject *op[3];
+ npy_uint32 op_flags[3];
+ NpyIter *iter;
+
+ NpyIter_IterNextFunc *iternext;
+ char **dataptr;
+ npy_intp *stride;
+ npy_intp *countptr;
+ npy_intp src_itemsize;
+ int needs_api;
+
+ op[0] = dst;
+ op[1] = src;
+ op[2] = mask;
+ /*
+ * TODO: In NumPy 2.0, renable NPY_ITER_NO_BROADCAST. This
+ * was removed during NumPy 1.6 testing for compatibility
+ * with NumPy 1.5, as per Travis's -10 veto power.
+ */
+ /*op_flags[0] = NPY_ITER_WRITEONLY|NPY_ITER_NO_BROADCAST;*/
+ op_flags[0] = NPY_ITER_WRITEONLY;
+ op_flags[1] = NPY_ITER_READONLY;
+ op_flags[2] = NPY_ITER_READONLY;
+
+ iter = NpyIter_MultiNew(3, op,
+ NPY_ITER_EXTERNAL_LOOP|
+ NPY_ITER_REFS_OK|
+ NPY_ITER_ZEROSIZE_OK,
+ NPY_KEEPORDER,
+ NPY_NO_CASTING,
+ op_flags,
+ NULL);
+ if (iter == NULL) {
+ return -1;
+ }
+
+ iternext = NpyIter_GetIterNext(iter, NULL);
+ if (iternext == NULL) {
+ NpyIter_Deallocate(iter);
+ return -1;
+ }
+ dataptr = NpyIter_GetDataPtrArray(iter);
+ stride = NpyIter_GetInnerStrideArray(iter);
+ countptr = NpyIter_GetInnerLoopSizePtr(iter);
+ src_itemsize = PyArray_DESCR(src)->elsize;
+
+ needs_api = NpyIter_IterationNeedsAPI(iter);
+
+ /*
+ * Because buffering is disabled in the iterator, the inner loop
+ * strides will be the same throughout the iteration loop. Thus,
+ * we can pass them to this function to take advantage of
+ * contiguous strides, etc.
+ */
+ if (PyArray_GetMaskedDTypeTransferFunction(
+ PyArray_ISALIGNED(src) && PyArray_ISALIGNED(dst),
+ stride[1], stride[0], stride[2],
+ PyArray_DESCR(src),
+ PyArray_DESCR(dst),
+ PyArray_DESCR(mask),
+ 0,
+ &stransfer, &transferdata,
+ &needs_api) != NPY_SUCCEED) {
+ NpyIter_Deallocate(iter);
+ return -1;
+ }
+
+
+ if (NpyIter_GetIterSize(iter) != 0) {
+ if (!needs_api) {
+ NPY_BEGIN_THREADS;
+ }
+
+ do {
+ stransfer(dataptr[0], stride[0],
+ dataptr[1], stride[1],
+ (npy_uint8 *)dataptr[2], stride[2],
+ *countptr, src_itemsize, transferdata);
+ } while(iternext(iter));
+
+ if (!needs_api) {
+ NPY_END_THREADS;
+ }
+ }
+
+ NPY_AUXDATA_FREE(transferdata);
+ NpyIter_Deallocate(iter);
+
+ return PyErr_Occurred() ? -1 : 0;
+ }
+}
+
/*NUMPY_API
* PyArray_CheckAxis
diff --git a/numpy/core/src/multiarray/item_selection.c b/numpy/core/src/multiarray/item_selection.c
index b4aef4ad6..74203774c 100644
--- a/numpy/core/src/multiarray/item_selection.c
+++ b/numpy/core/src/multiarray/item_selection.c
@@ -393,6 +393,11 @@ PyArray_PutMask(PyArrayObject *self, PyObject* values0, PyObject* mask0)
char *src, *dest;
int copied = 0;
+ if (DEPRECATE("putmask has been deprecated. Use copyto with 'where' as "
+ "the mask instead") < 0) {
+ return NULL;
+ }
+
mask = NULL;
values = NULL;
if (!PyArray_Check(self)) {
diff --git a/numpy/core/src/multiarray/multiarraymodule.c b/numpy/core/src/multiarray/multiarraymodule.c
index 3b2b8f583..452ddec20 100644
--- a/numpy/core/src/multiarray/multiarraymodule.c
+++ b/numpy/core/src/multiarray/multiarraymodule.c
@@ -45,6 +45,7 @@ NPY_NO_EXPORT int NPY_NUMUSERTYPES = 0;
#include "numpymemoryview.h"
#include "convert_datatype.h"
#include "nditer_pywrap.h"
+#include "methods.h"
#include "_datetime.h"
#include "datetime_strings.h"
#include "datetime_busday.h"
@@ -1649,6 +1650,79 @@ clean_type:
}
static PyObject *
+array_copyto(PyObject *NPY_UNUSED(ignored), PyObject *args, PyObject *kwds)
+{
+
+ static char *kwlist[] = {"dst","src","casting","where",NULL};
+ PyObject *wheremask_in = NULL;
+ PyArrayObject *dst = NULL, *src = NULL, *wheremask = NULL;
+ NPY_CASTING casting = NPY_SAME_KIND_CASTING;
+
+ if (!PyArg_ParseTupleAndKeywords(args, kwds, "O!O&|O&O", kwlist,
+ &PyArray_Type, &dst,
+ &PyArray_Converter, &src,
+ &PyArray_CastingConverter, &casting,
+ &wheremask_in)) {
+ goto fail;
+ }
+
+ if (wheremask_in != NULL) {
+ /* Get the boolean where mask */
+ PyArray_Descr *dtype = PyArray_DescrFromType(NPY_BOOL);
+ if (dtype == NULL) {
+ goto fail;
+ }
+ wheremask = (PyArrayObject *)PyArray_FromAny(wheremask_in,
+ dtype, 0, 0, 0, NULL);
+ if (wheremask == NULL) {
+ goto fail;
+ }
+
+ /* Use the 'move' function which handles overlapping */
+ if (PyArray_MaskedMoveInto(dst, src, wheremask, casting) < 0) {
+ goto fail;
+ }
+ }
+ else {
+ /*
+ * MoveInto doesn't accept a casting rule, must check it
+ * ourselves.
+ */
+ if (!PyArray_CanCastArrayTo(src, PyArray_DESCR(dst), casting)) {
+ PyObject *errmsg;
+ errmsg = PyUString_FromString("Cannot cast array data from ");
+ PyUString_ConcatAndDel(&errmsg,
+ PyObject_Repr((PyObject *)PyArray_DESCR(src)));
+ PyUString_ConcatAndDel(&errmsg,
+ PyUString_FromString(" to "));
+ PyUString_ConcatAndDel(&errmsg,
+ PyObject_Repr((PyObject *)PyArray_DESCR(dst)));
+ PyUString_ConcatAndDel(&errmsg,
+ PyUString_FromFormat(" according to the rule %s",
+ npy_casting_to_string(casting)));
+ PyErr_SetObject(PyExc_TypeError, errmsg);
+ goto fail;
+ }
+
+ /* Use the 'move' function which handles overlapping */
+ if (PyArray_MoveInto(dst, src) < 0) {
+ goto fail;
+ }
+ }
+
+ Py_XDECREF(src);
+ Py_XDECREF(wheremask);
+
+ Py_INCREF(Py_None);
+ return Py_None;
+
+ fail:
+ Py_XDECREF(src);
+ Py_XDECREF(wheremask);
+ return NULL;
+}
+
+static PyObject *
array_empty(PyObject *NPY_UNUSED(ignored), PyObject *args, PyObject *kwds)
{
@@ -3363,6 +3437,9 @@ static struct PyMethodDef array_module_methods[] = {
{"array",
(PyCFunction)_array_fromobject,
METH_VARARGS|METH_KEYWORDS, NULL},
+ {"copyto",
+ (PyCFunction)array_copyto,
+ METH_VARARGS|METH_KEYWORDS, NULL},
{"nested_iters",
(PyCFunction)NpyIter_NestedIters,
METH_VARARGS|METH_KEYWORDS, NULL},
diff --git a/numpy/core/src/umath/ufunc_type_resolution.c b/numpy/core/src/umath/ufunc_type_resolution.c
index 5fb055b28..8fb17a441 100644
--- a/numpy/core/src/umath/ufunc_type_resolution.c
+++ b/numpy/core/src/umath/ufunc_type_resolution.c
@@ -1407,7 +1407,7 @@ unmasked_ufunc_loop_as_masked(
do {
/* Skip masked values */
subloopsize = 0;
- while (subloopsize < loopsize && (*(npy_uint8 *)mask)&0x01 == 0) {
+ while (subloopsize < loopsize && ((*(npy_uint8 *)mask)&0x01) == 0) {
++subloopsize;
mask += mask_stride;
}
@@ -1420,7 +1420,7 @@ unmasked_ufunc_loop_as_masked(
* mess with the 'args' pointer values)
*/
subloopsize = 0;
- while (subloopsize < loopsize && (*(npy_uint8 *)mask)&0x01 != 0) {
+ while (subloopsize < loopsize && ((*(npy_uint8 *)mask)&0x01) != 0) {
++subloopsize;
mask += mask_stride;
}
diff --git a/numpy/core/tests/test_api.py b/numpy/core/tests/test_api.py
index 2c9b7d4d0..7ebcb932b 100644
--- a/numpy/core/tests/test_api.py
+++ b/numpy/core/tests/test_api.py
@@ -84,5 +84,41 @@ def test_array_astype():
assert_(not (a is b))
assert_(type(b) != np.matrix)
+def test_copyto():
+ a = np.arange(6, dtype='i4').reshape(2,3)
+
+ # Simple copy
+ np.copyto(a, [[3,1,5], [6,2,1]])
+ assert_equal(a, [[3,1,5], [6,2,1]])
+
+ # Overlapping copy should work
+ np.copyto(a[:,:2], a[::-1, 1::-1])
+ assert_equal(a, [[2,6,5], [1,3,1]])
+
+ # Defaults to 'same_kind' casting
+ assert_raises(TypeError, np.copyto, a, 1.5)
+
+ # Force a copy with 'unsafe' casting, truncating 1.5 to 1
+ np.copyto(a, 1.5, casting='unsafe')
+ assert_equal(a, 1)
+
+ # Copying with a mask
+ np.copyto(a, 3, where=[True,False,True])
+ assert_equal(a, [[3,1,3],[3,1,3]])
+
+ # Casting rule still applies with a mask
+ assert_raises(TypeError, np.copyto, a, 3.5, where=[True,False,True])
+
+ # Lists of integer 0's and 1's is ok too
+ np.copyto(a, 4, where=[[0,1,1], [1,0,0]])
+ assert_equal(a, [[3,4,4], [4,1,3]])
+
+ # Overlapping copy with mask should work
+ np.copyto(a[:,:2], a[::-1, 1::-1], where=[[0,1],[1,1]])
+ assert_equal(a, [[3,4,4], [4,3,3]])
+
+ # 'dst' must be an array
+ assert_raises(TypeError, np.copyto, [1,2,3], [2,3,4])
+
if __name__ == "__main__":
run_module_suite()