summaryrefslogtreecommitdiff
path: root/numpy/linalg
diff options
context:
space:
mode:
Diffstat (limited to 'numpy/linalg')
-rw-r--r--numpy/linalg/lapack_lite/clapack_scrub.py2
-rwxr-xr-xnumpy/linalg/lapack_lite/make_lite.py6
-rw-r--r--numpy/linalg/linalg.py32
-rw-r--r--numpy/linalg/setup.py2
-rw-r--r--numpy/linalg/umath_linalg.cpp81
5 files changed, 77 insertions, 46 deletions
diff --git a/numpy/linalg/lapack_lite/clapack_scrub.py b/numpy/linalg/lapack_lite/clapack_scrub.py
index fffd70910..d4da9070d 100644
--- a/numpy/linalg/lapack_lite/clapack_scrub.py
+++ b/numpy/linalg/lapack_lite/clapack_scrub.py
@@ -297,7 +297,7 @@ def scrubSource(source, nsteps=None, verbose=False):
if __name__ == '__main__':
filename = sys.argv[1]
outfilename = os.path.join(sys.argv[2], os.path.basename(filename))
- with open(filename, 'r') as fo:
+ with open(filename) as fo:
source = fo.read()
if len(sys.argv) > 3:
diff --git a/numpy/linalg/lapack_lite/make_lite.py b/numpy/linalg/lapack_lite/make_lite.py
index ca8d4c62c..20b212792 100755
--- a/numpy/linalg/lapack_lite/make_lite.py
+++ b/numpy/linalg/lapack_lite/make_lite.py
@@ -253,7 +253,7 @@ def dumpRoutineNames(library, output_dir):
def concatenateRoutines(routines, output_file):
with open(output_file, 'w') as output_fo:
for r in routines:
- with open(r.filename, 'r') as fo:
+ with open(r.filename) as fo:
source = fo.read()
output_fo.write(source)
@@ -296,7 +296,7 @@ def create_name_header(output_dir):
if not fn.endswith('.f'):
continue
- with open(fn, 'r') as f:
+ with open(fn) as f:
for line in f:
m = routine_re.match(line)
if m:
@@ -304,7 +304,7 @@ def create_name_header(output_dir):
# f2c symbols
f2c_symbols = set()
- with open('f2c.h', 'r') as f:
+ with open('f2c.h') as f:
for line in f:
m = extern_re.match(line)
if m:
diff --git a/numpy/linalg/linalg.py b/numpy/linalg/linalg.py
index 0f06c8520..b838b9397 100644
--- a/numpy/linalg/linalg.py
+++ b/numpy/linalg/linalg.py
@@ -24,7 +24,7 @@ from numpy.core import (
array, asarray, zeros, empty, empty_like, intc, single, double,
csingle, cdouble, inexact, complexfloating, newaxis, all, Inf, dot,
add, multiply, sqrt, sum, isfinite,
- finfo, errstate, geterrobj, moveaxis, amin, amax, product, abs,
+ finfo, errstate, geterrobj, moveaxis, amin, amax, prod, abs,
atleast_2d, intp, asanyarray, object_, matmul,
swapaxes, divide, count_nonzero, isnan, sign, argsort, sort,
reciprocal
@@ -65,11 +65,11 @@ fortran_int = intc
@set_module('numpy.linalg')
-class LinAlgError(Exception):
+class LinAlgError(ValueError):
"""
Generic Python-exception-derived object raised by linalg functions.
- General purpose exception class, derived from Python's exception.Exception
+ General purpose exception class, derived from Python's ValueError
class, programmatically raised in linalg functions when a Linear
Algebra-related condition would prevent further correct execution of the
function.
@@ -161,24 +161,24 @@ def _commonType(*arrays):
result_type = single
is_complex = False
for a in arrays:
- if issubclass(a.dtype.type, inexact):
- if isComplexType(a.dtype.type):
+ type_ = a.dtype.type
+ if issubclass(type_, inexact):
+ if isComplexType(type_):
is_complex = True
- rt = _realType(a.dtype.type, default=None)
- if rt is None:
+ rt = _realType(type_, default=None)
+ if rt is double:
+ result_type = double
+ elif rt is None:
# unsupported inexact scalar
raise TypeError("array type %s is unsupported in linalg" %
(a.dtype.name,))
else:
- rt = double
- if rt is double:
result_type = double
if is_complex:
- t = cdouble
result_type = _complex_types_map[result_type]
+ return cdouble, result_type
else:
- t = double
- return t, result_type
+ return double, result_type
def _to_native_byte_order(*arrays):
@@ -219,7 +219,7 @@ def _assert_finite(*arrays):
def _is_empty_2d(arr):
# check size first for efficiency
- return arr.size == 0 and product(arr.shape[-2:]) == 0
+ return arr.size == 0 and prod(arr.shape[-2:]) == 0
def transpose(a):
@@ -924,12 +924,12 @@ def qr(a, mode='reduced'):
msg = "".join((
"The 'full' option is deprecated in favor of 'reduced'.\n",
"For backward compatibility let mode default."))
- warnings.warn(msg, DeprecationWarning, stacklevel=3)
+ warnings.warn(msg, DeprecationWarning, stacklevel=2)
mode = 'reduced'
elif mode in ('e', 'economic'):
# 2013-04-01, 1.8
msg = "The 'economic' option is deprecated."
- warnings.warn(msg, DeprecationWarning, stacklevel=3)
+ warnings.warn(msg, DeprecationWarning, stacklevel=2)
mode = 'economic'
else:
raise ValueError(f"Unrecognized mode '{mode}'")
@@ -2308,7 +2308,7 @@ def lstsq(a, b, rcond="warn"):
"To use the future default and silence this warning "
"we advise to pass `rcond=None`, to keep using the old, "
"explicitly pass `rcond=-1`.",
- FutureWarning, stacklevel=3)
+ FutureWarning, stacklevel=2)
rcond = -1
if rcond is None:
rcond = finfo(t).eps * max(n, m)
diff --git a/numpy/linalg/setup.py b/numpy/linalg/setup.py
index 1c4e1295e..6f72635ab 100644
--- a/numpy/linalg/setup.py
+++ b/numpy/linalg/setup.py
@@ -4,7 +4,6 @@ import sysconfig
def configuration(parent_package='', top_path=None):
from numpy.distutils.misc_util import Configuration
- from numpy.distutils.ccompiler_opt import NPY_CXX_FLAGS
from numpy.distutils.system_info import get_info, system_info
config = Configuration('linalg', parent_package, top_path)
@@ -81,7 +80,6 @@ def configuration(parent_package='', top_path=None):
sources=['umath_linalg.cpp', get_lapack_lite_sources],
depends=['lapack_lite/f2c.h'],
extra_info=lapack_info,
- extra_cxx_compile_args=NPY_CXX_FLAGS,
libraries=['npymath'],
)
config.add_data_files('*.pyi')
diff --git a/numpy/linalg/umath_linalg.cpp b/numpy/linalg/umath_linalg.cpp
index bbb4bb896..68db2b2f1 100644
--- a/numpy/linalg/umath_linalg.cpp
+++ b/numpy/linalg/umath_linalg.cpp
@@ -22,6 +22,7 @@
#include <cstdio>
#include <cassert>
#include <cmath>
+#include <type_traits>
#include <utility>
@@ -1148,7 +1149,7 @@ slogdet(char **args,
void *NPY_UNUSED(func))
{
fortran_int m;
- npy_uint8 *tmp_buff = NULL;
+ char *tmp_buff = NULL;
size_t matrix_size;
size_t pivot_size;
size_t safe_m;
@@ -1162,10 +1163,11 @@ slogdet(char **args,
*/
INIT_OUTER_LOOP_3
m = (fortran_int) dimensions[0];
- safe_m = m;
+ /* avoid empty malloc (buffers likely unused) and ensure m is `size_t` */
+ safe_m = m != 0 ? m : 1;
matrix_size = safe_m * safe_m * sizeof(typ);
pivot_size = safe_m * sizeof(fortran_int);
- tmp_buff = (npy_uint8 *)malloc(matrix_size + pivot_size);
+ tmp_buff = (char *)malloc(matrix_size + pivot_size);
if (tmp_buff) {
LINEARIZE_DATA_t lin_data;
@@ -1182,6 +1184,13 @@ slogdet(char **args,
free(tmp_buff);
}
+ else {
+ /* TODO: Requires use of new ufunc API to indicate error return */
+ NPY_ALLOW_C_API_DEF
+ NPY_ALLOW_C_API;
+ PyErr_NoMemory();
+ NPY_DISABLE_C_API;
+ }
}
template<typename typ, typename basetyp>
@@ -1192,7 +1201,7 @@ det(char **args,
void *NPY_UNUSED(func))
{
fortran_int m;
- npy_uint8 *tmp_buff;
+ char *tmp_buff;
size_t matrix_size;
size_t pivot_size;
size_t safe_m;
@@ -1206,10 +1215,11 @@ det(char **args,
*/
INIT_OUTER_LOOP_2
m = (fortran_int) dimensions[0];
- safe_m = m;
+ /* avoid empty malloc (buffers likely unused) and ensure m is `size_t` */
+ safe_m = m != 0 ? m : 1;
matrix_size = safe_m * safe_m * sizeof(typ);
pivot_size = safe_m * sizeof(fortran_int);
- tmp_buff = (npy_uint8 *)malloc(matrix_size + pivot_size);
+ tmp_buff = (char *)malloc(matrix_size + pivot_size);
if (tmp_buff) {
LINEARIZE_DATA_t lin_data;
@@ -1230,6 +1240,13 @@ det(char **args,
free(tmp_buff);
}
+ else {
+ /* TODO: Requires use of new ufunc API to indicate error return */
+ NPY_ALLOW_C_API_DEF
+ NPY_ALLOW_C_API;
+ PyErr_NoMemory();
+ NPY_DISABLE_C_API;
+ }
}
@@ -3737,16 +3754,16 @@ scalar_trait)
fortran_int lda = fortran_int_max(1, m);
fortran_int ldb = fortran_int_max(1, fortran_int_max(m,n));
- mem_buff = (npy_uint8 *)malloc(a_size + b_size + s_size);
-
- if (!mem_buff)
- goto error;
+ size_t msize = a_size + b_size + s_size;
+ mem_buff = (npy_uint8 *)malloc(msize != 0 ? msize : 1);
+ if (!mem_buff) {
+ goto no_memory;
+ }
a = mem_buff;
b = a + a_size;
s = b + b_size;
-
params->M = m;
params->N = n;
params->NRHS = nrhs;
@@ -3766,9 +3783,9 @@ scalar_trait)
params->RWORK = NULL;
params->LWORK = -1;
- if (call_gelsd(params) != 0)
+ if (call_gelsd(params) != 0) {
goto error;
-
+ }
work_count = (fortran_int)work_size_query;
work_size = (size_t) work_size_query * sizeof(ftyp);
@@ -3776,9 +3793,9 @@ scalar_trait)
}
mem_buff2 = (npy_uint8 *)malloc(work_size + iwork_size);
- if (!mem_buff2)
- goto error;
-
+ if (!mem_buff2) {
+ goto no_memory;
+ }
work = mem_buff2;
iwork = work + work_size;
@@ -3788,12 +3805,18 @@ scalar_trait)
params->LWORK = work_count;
return 1;
+
+ no_memory:
+ NPY_ALLOW_C_API_DEF
+ NPY_ALLOW_C_API;
+ PyErr_NoMemory();
+ NPY_DISABLE_C_API;
+
error:
TRACE_TXT("%s failed init\n", __FUNCTION__);
free(mem_buff);
free(mem_buff2);
memset(params, 0, sizeof(*params));
-
return 0;
}
@@ -3857,16 +3880,17 @@ using frealtyp = basetype_t<ftyp>;
fortran_int lda = fortran_int_max(1, m);
fortran_int ldb = fortran_int_max(1, fortran_int_max(m,n));
- mem_buff = (npy_uint8 *)malloc(a_size + b_size + s_size);
+ size_t msize = a_size + b_size + s_size;
+ mem_buff = (npy_uint8 *)malloc(msize != 0 ? msize : 1);
- if (!mem_buff)
- goto error;
+ if (!mem_buff) {
+ goto no_memory;
+ }
a = mem_buff;
b = a + a_size;
s = b + b_size;
-
params->M = m;
params->N = n;
params->NRHS = nrhs;
@@ -3887,8 +3911,9 @@ using frealtyp = basetype_t<ftyp>;
params->RWORK = &rwork_size_query;
params->LWORK = -1;
- if (call_gelsd(params) != 0)
+ if (call_gelsd(params) != 0) {
goto error;
+ }
work_count = (fortran_int)work_size_query.r;
@@ -3898,8 +3923,9 @@ using frealtyp = basetype_t<ftyp>;
}
mem_buff2 = (npy_uint8 *)malloc(work_size + rwork_size + iwork_size);
- if (!mem_buff2)
- goto error;
+ if (!mem_buff2) {
+ goto no_memory;
+ }
work = mem_buff2;
rwork = work + work_size;
@@ -3911,6 +3937,13 @@ using frealtyp = basetype_t<ftyp>;
params->LWORK = work_count;
return 1;
+
+ no_memory:
+ NPY_ALLOW_C_API_DEF
+ NPY_ALLOW_C_API;
+ PyErr_NoMemory();
+ NPY_DISABLE_C_API;
+
error:
TRACE_TXT("%s failed init\n", __FUNCTION__);
free(mem_buff);