Skip to content

MAINT: use copy=False in a few astype() calls #5909

New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Merged
merged 3 commits into from
May 22, 2015
Merged
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
54 changes: 34 additions & 20 deletions numpy/linalg/linalg.py
Original file line number Diff line number Diff line change
Expand Up @@ -23,7 +23,7 @@
csingle, cdouble, inexact, complexfloating, newaxis, ravel, all, Inf, dot,
add, multiply, sqrt, maximum, fastCopyAndTranspose, sum, isfinite, size,
finfo, errstate, geterrobj, longdouble, rollaxis, amin, amax, product, abs,
broadcast, atleast_2d, intp, asanyarray
broadcast, atleast_2d, intp, asanyarray, isscalar
)
from numpy.lib import triu, asfarray
from numpy.linalg import lapack_lite, _umath_linalg
Expand Down Expand Up @@ -382,7 +382,7 @@ def solve(a, b):
extobj = get_linalg_error_extobj(_raise_linalgerror_singular)
r = gufunc(a, b, signature=signature, extobj=extobj)

return wrap(r.astype(result_t))
return wrap(r.astype(result_t, copy=False))


def tensorinv(a, ind=2):
Expand Down Expand Up @@ -522,7 +522,7 @@ def inv(a):
signature = 'D->D' if isComplexType(t) else 'd->d'
extobj = get_linalg_error_extobj(_raise_linalgerror_singular)
ainv = _umath_linalg.inv(a, signature=signature, extobj=extobj)
return wrap(ainv.astype(result_t))
return wrap(ainv.astype(result_t, copy=False))


# Cholesky decomposition
Expand Down Expand Up @@ -606,7 +606,8 @@ def cholesky(a):
_assertNdSquareness(a)
t, result_t = _commonType(a)
signature = 'D->D' if isComplexType(t) else 'd->d'
return wrap(gufunc(a, signature=signature, extobj=extobj).astype(result_t))
r = gufunc(a, signature=signature, extobj=extobj)
return wrap(r.astype(result_t, copy=False))

# QR decompostion

Expand Down Expand Up @@ -781,7 +782,7 @@ def qr(a, mode='reduced'):

if mode == 'economic':
if t != result_t :
a = a.astype(result_t)
a = a.astype(result_t, copy=False)
return wrap(a.T)

# generate q from a
Expand Down Expand Up @@ -908,7 +909,7 @@ def eigvals(a):
else:
result_t = _complexType(result_t)

return w.astype(result_t)
return w.astype(result_t, copy=False)

def eigvalsh(a, UPLO='L'):
"""
Expand Down Expand Up @@ -978,7 +979,7 @@ def eigvalsh(a, UPLO='L'):
t, result_t = _commonType(a)
signature = 'D->d' if isComplexType(t) else 'd->d'
w = gufunc(a, signature=signature, extobj=extobj)
return w.astype(_realType(result_t))
return w.astype(_realType(result_t), copy=False)

def _convertarray(a):
t, result_t = _commonType(a)
Expand Down Expand Up @@ -1124,8 +1125,8 @@ def eig(a):
else:
result_t = _complexType(result_t)

vt = vt.astype(result_t)
return w.astype(result_t), wrap(vt)
vt = vt.astype(result_t, copy=False)
return w.astype(result_t, copy=False), wrap(vt)


def eigh(a, UPLO='L'):
Expand Down Expand Up @@ -1232,8 +1233,8 @@ def eigh(a, UPLO='L'):

signature = 'D->dD' if isComplexType(t) else 'd->dd'
w, vt = gufunc(a, signature=signature, extobj=extobj)
w = w.astype(_realType(result_t))
vt = vt.astype(result_t)
w = w.astype(_realType(result_t), copy=False)
vt = vt.astype(result_t, copy=False)
return w, wrap(vt)


Expand Down Expand Up @@ -1344,9 +1345,9 @@ def svd(a, full_matrices=1, compute_uv=1):

signature = 'D->DdD' if isComplexType(t) else 'd->ddd'
u, s, vt = gufunc(a, signature=signature, extobj=extobj)
u = u.astype(result_t)
s = s.astype(_realType(result_t))
vt = vt.astype(result_t)
u = u.astype(result_t, copy=False)
s = s.astype(_realType(result_t), copy=False)
vt = vt.astype(result_t, copy=False)
return wrap(u), s, wrap(vt)
else:
if m < n:
Expand All @@ -1356,7 +1357,7 @@ def svd(a, full_matrices=1, compute_uv=1):

signature = 'D->d' if isComplexType(t) else 'd->d'
s = gufunc(a, signature=signature, extobj=extobj)
s = s.astype(_realType(result_t))
s = s.astype(_realType(result_t), copy=False)
return s

def cond(x, p=None):
Expand Down Expand Up @@ -1695,7 +1696,15 @@ def slogdet(a):
real_t = _realType(result_t)
signature = 'D->Dd' if isComplexType(t) else 'd->dd'
sign, logdet = _umath_linalg.slogdet(a, signature=signature)
return sign.astype(result_t), logdet.astype(real_t)
if isscalar(sign):
sign = sign.astype(result_t)
else:
sign = sign.astype(result_t, copy=False)
if isscalar(logdet):
logdet = logdet.astype(real_t)
else:
logdet = logdet.astype(real_t, copy=False)
return sign, logdet

def det(a):
"""
Expand Down Expand Up @@ -1749,7 +1758,12 @@ def det(a):
_assertNdSquareness(a)
t, result_t = _commonType(a)
signature = 'D->D' if isComplexType(t) else 'd->d'
return _umath_linalg.det(a, signature=signature).astype(result_t)
r = _umath_linalg.det(a, signature=signature)
if isscalar(r):
r = r.astype(result_t)
else:
r = r.astype(result_t, copy=False)
return r

# Linear Least Squares

Expand Down Expand Up @@ -1905,12 +1919,12 @@ def lstsq(a, b, rcond=-1):
if results['rank'] == n and m > n:
if isComplexType(t):
resids = sum(abs(transpose(bstar)[n:,:])**2, axis=0).astype(
result_real_t)
result_real_t, copy=False)
else:
resids = sum((transpose(bstar)[n:,:])**2, axis=0).astype(
result_real_t)
result_real_t, copy=False)

st = s[:min(n, m)].copy().astype(result_real_t)
st = s[:min(n, m)].astype(result_real_t, copy=True)
return wrap(x), wrap(resids), results['rank'], st


Expand Down