1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211
|
#!/usr/bin/env python
#
# Created by: Pearu Peterson, April 2002
#
__usage__ = """
Build linalg:
python setup_linalg.py build
Run tests if scipy is installed:
python -c 'import scipy;scipy.linalg.test(<level>)'
Run tests if linalg is not installed:
python tests/test_blas.py [<level>]
"""
from numpy import arange, add, array
import math
import sys
from numpy.testing import *
set_package_path()
from linalg import fblas
from linalg import cblas
restore_path()
class test_cblas1_simple(ScipyTestCase):
def check_axpy(self):
for p in 'sd':
f = getattr(cblas,p+'axpy',None)
if f is None: continue
assert_array_almost_equal(f(5,[1,2,3],[2,-1,3]),[7,9,18])
for p in 'cz':
f = getattr(cblas,p+'axpy',None)
if f is None: continue
assert_array_almost_equal(f(5,[1,2j,3],[2,-1,3]),[7,10j-1,18])
class test_fblas1_simple(ScipyTestCase):
def check_axpy(self):
for p in 'sd':
f = getattr(fblas,p+'axpy',None)
if f is None: continue
assert_array_almost_equal(f([1,2,3],[2,-1,3],a=5),[7,9,18])
for p in 'cz':
f = getattr(fblas,p+'axpy',None)
if f is None: continue
assert_array_almost_equal(f([1,2j,3],[2,-1,3],a=5),[7,10j-1,18])
def check_copy(self):
for p in 'sd':
f = getattr(fblas,p+'copy',None)
if f is None: continue
assert_array_almost_equal(f([3,4,5],[8]*3),[3,4,5])
for p in 'cz':
f = getattr(fblas,p+'copy',None)
if f is None: continue
assert_array_almost_equal(f([3,4j,5+3j],[8]*3),[3,4j,5+3j])
def check_asum(self):
for p in 'sd':
f = getattr(fblas,p+'asum',None)
if f is None: continue
assert_almost_equal(f([3,-4,5]),12)
for p in ['sc','dz']:
f = getattr(fblas,p+'asum',None)
if f is None: continue
assert_almost_equal(f([3j,-4,3-4j]),14)
def check_dot(self):
for p in 'sd':
f = getattr(fblas,p+'dot',None)
if f is None: continue
assert_almost_equal(f([3,-4,5],[2,5,1]),-9)
for p in 'cz':
f = getattr(fblas,p+'dotu',None)
if f is None: continue
assert_almost_equal(f([3j,-4,3-4j],[2,3,1]),-9+2j)
f = getattr(fblas,p+'dotc')
assert_almost_equal(f([3j,-4,3-4j],[2,3j,1]),3-14j)
def check_nrm2(self):
for p in 'sd':
f = getattr(fblas,p+'nrm2',None)
if f is None: continue
assert_almost_equal(f([3,-4,5]),math.sqrt(50))
for p in ['sc','dz']:
f = getattr(fblas,p+'nrm2',None)
if f is None: continue
assert_almost_equal(f([3j,-4,3-4j]),math.sqrt(50))
def check_scal(self):
for p in 'sd':
f = getattr(fblas,p+'scal',None)
if f is None: continue
assert_array_almost_equal(f(2,[3,-4,5]),[6,-8,10])
for p in 'cz':
f = getattr(fblas,p+'scal',None)
if f is None: continue
assert_array_almost_equal(f(3j,[3j,-4,3-4j]),[-9,-12j,12+9j])
for p in ['cs','zd']:
f = getattr(fblas,p+'scal',None)
if f is None: continue
assert_array_almost_equal(f(3,[3j,-4,3-4j]),[9j,-12,9-12j])
def check_swap(self):
for p in 'sd':
f = getattr(fblas,p+'swap',None)
if f is None: continue
x,y = [2,3,1],[-2,3,7]
x1,y1 = f(x,y)
assert_array_almost_equal(x1,y)
assert_array_almost_equal(y1,x)
for p in 'cz':
f = getattr(fblas,p+'swap',None)
if f is None: continue
x,y = [2,3j,1],[-2,3,7-3j]
x1,y1 = f(x,y)
assert_array_almost_equal(x1,y)
assert_array_almost_equal(y1,x)
def check_amax(self):
for p in 'sd':
f = getattr(fblas,'i'+p+'amax')
assert_equal(f([-2,4,3]),1)
for p in 'cz':
f = getattr(fblas,'i'+p+'amax')
assert_equal(f([-5,4+3j,6]),1)
#XXX: need tests for rot,rotm,rotg,rotmg
class test_fblas2_simple(ScipyTestCase):
def check_gemv(self):
for p in 'sd':
f = getattr(fblas,p+'gemv',None)
if f is None: continue
assert_array_almost_equal(f(3,[[3]],[-4]),[-36])
assert_array_almost_equal(f(3,[[3]],[-4],3,[5]),[-21])
for p in 'cz':
f = getattr(fblas,p+'gemv',None)
if f is None: continue
assert_array_almost_equal(f(3j,[[3-4j]],[-4]),[-48-36j])
assert_array_almost_equal(f(3j,[[3-4j]],[-4],3,[5j]),[-48-21j])
def check_ger(self):
for p in 'sd':
f = getattr(fblas,p+'ger',None)
if f is None: continue
assert_array_almost_equal(f(1,[1,
2],[3,4]),[[3,4],[6,8]])
assert_array_almost_equal(f(2,[1,
2,
3],[3,4]),[[6,8],[12,16],[18,24]])
assert_array_almost_equal(f(1,[1,
2],[3,4],
a=[[1,2],[3,4]]
),[[4,6],[9,12]])
for p in 'cz':
f = getattr(fblas,p+'geru',None)
if f is None: continue
assert_array_almost_equal(f(1,[1j,
2],[3,4]),[[3j,4j],[6,8]])
assert_array_almost_equal(f(-2,[1j,
2j,
3j],[3j,4j]),[[6,8],[12,16],[18,24]])
for p in 'cz':
f = getattr(fblas,p+'gerc',None)
if f is None: continue
assert_array_almost_equal(f(1,[1j,
2],[3,4]),[[3j,4j],[6,8]])
assert_array_almost_equal(f(2,[1j,
2j,
3j],[3j,4j]),[[6,8],[12,16],[18,24]])
class test_fblas3_simple(ScipyTestCase):
def check_gemm(self):
for p in 'sd':
f = getattr(fblas,p+'gemm',None)
if f is None: continue
assert_array_almost_equal(f(3,[3],[-4]),[[-36]])
assert_array_almost_equal(f(3,[3],[-4],3,[5]),[-21])
for p in 'cz':
f = getattr(fblas,p+'gemm',None)
if f is None: continue
assert_array_almost_equal(f(3j,[3-4j],[-4]),[[-48-36j]])
assert_array_almost_equal(f(3j,[3-4j],[-4],3,[5j]),[-48-21j])
class test_blas(ScipyTestCase):
def check_fblas(self):
if hasattr(fblas,'empty_module'):
print """
****************************************************************
WARNING: fblas module is empty.
-----------
See scipy/INSTALL.txt for troubleshooting.
****************************************************************
"""
def check_cblas(self):
if hasattr(cblas,'empty_module'):
print """
****************************************************************
WARNING: cblas module is empty
-----------
See scipy/INSTALL.txt for troubleshooting.
Notes:
* If atlas library is not found by numpy/distutils/system_info.py,
then scipy uses fblas instead of cblas.
****************************************************************
"""
if __name__ == "__main__":
ScipyTest().run()
|