Debugging linear algebra related issues#
Linear algebra related bug reports are among the most challenging issues to diagnose and address. This is not only because linear algebra can be challenging mathematically/algorithmically (that is true for many parts of SciPy), but because BLAS/LAPACK libraries are a complex build-time as well as runtime dependency - and we support a significant number of BLAS/LAPACK libraries.
This document aims to provide guidance about how to go about debugging linear algebra issues.
If there is a real bug, it can be in one of three places:
In the BLAS library being used,
In SciPy’s bindings for BLAS or LAPACK (generated by
numpy.f2py
and/or Cython),In SciPy’s algorithmic code.
A key first step is to determine whether the bug is in SciPy or in the BLAS library. To do so, the most efficient way to disambiguate the two is to set up your environment in such a way that you can achieve runtime switching between different BLAS libraries (something we don’t support out of the box, and isn’t possible with SciPy’s wheels from PyPI).
Upstream BLAS library authors strongly prefer to get clean reproducers (just like we do), which means: no Python involved. So this guide will also cover how to create reproducers in C or Fortran.
Finding the BLAS library being used#
SciPy has one function, show_config
, to introspect the build
configuration of an installed package. This allows querying for details of BLAS
and LAPACK. E.g.:
>>> blas_dep = scipy.show_config(mode='dicts')['Build Dependencies']['blas']
>>> for key in blas_dep:
... print(f"{key}: {blas_dep[key]}")
...
name: openblas
found: True
version: 0.3.23
detection method: pkgconfig
include directory: /home/user/mambaforge/envs/scipy-dev/include
lib directory: /home/user/mambaforge/envs/scipy-dev/lib
openblas configuration: USE_64BITINT=0 DYNAMIC_ARCH=1 DYNAMIC_OLDER= NO_CBLAS= NO_LAPACK=0 NO_LAPACKE= NO_AFFINITY=1 USE_OPENMP=0 PRESCOTT MAX_THREADS=128
pc file directory: /home/user/mambaforge/envs/scipy-dev/lib/pkgconfig
This method will be correct for SciPy wheels and for local dev builds. It may
be correct for other installs, however keep in mind that distros like
conda-forge and Debian may be building against stub libraries (typically
blas.so
/lapack.so
) and then installing another library for the user -
in such cases, plain blas
and lapack
will be reported even if for
example OpenBLAS or MKL is installed in the environment. For such installs,
threadpoolctl will usually be
able to report what the actual BLAS library in use is (except it doesn’t report
on plain Netlib BLAS, see
threadpoolctl#159):
$ python -m threadpoolctl -i scipy.linalg
[
{
"user_api": "blas",
"internal_api": "openblas",
"prefix": "libopenblas",
"filepath": "/home/user/mambaforge/envs/dev/lib/libopenblasp-r0.3.21.so",
"version": "0.3.21",
"threading_layer": "pthreads",
"architecture": "SkylakeX",
"num_threads": 24
}
]
Other ways of introspecting that can be helpful in local dev environments include:
Checking dependencies of shared libraries:
$ ldd build/scipy/linalg/_fblas.cpython-*.so
...
libopenblas.so.0 => /home/user/mambaforge/envs/scipy-dev/lib/libopenblas.so.0 (0x0000780d6d000000)
% otool -L build/scipy/linalg/_fblas.cpython-310-darwin.so
build/scipy/linalg/_fblas.*.so:
/usr/lib/libSystem.B.dylib (compatibility version 1.0.0, current version 1336.61.1)
@rpath/libopenblas.0.dylib (compatibility version 0.0.0, current version 0.0.0)
Checking whether the linked library contains a given symbol. E.g., conda-forge installs a
libblas.so
that may be any supported library:
$ nm -gD ~/mambaforge/envs/scipy-dev/lib/libblas.so | rg openblas_set_num_threads
0000000000362990 T openblas_set_num_threads
% nm ~/mambaforge/envs/scipy-dev/lib/libblas.3.dylib | rg openblas_set_num_threads
000000000015b6b0 T _openblas_set_num_threads
Setting up your environment for switching BLAS libraries#
We’ll cover several ways of switching between different BLAS libraries, because the easiest method may depend on your OS/distro and on whether you want a released version of SciPy or a development build.
Conda-forge#
Perhaps the easiest way is to use the runtime switching abilities provided by distros. For example, to create a conda environment with the latest SciPy release installed and then switching between OpenBLAS, Netlib BLAS/LAPACK, and MKL is as simple as:
$ mamba create -n blas-switch scipy threadpoolctl
$ mamba activate blas-switch
$ python -m threadpoolctl -i scipy.linalg
...
"user_api": "blas",
"internal_api": "openblas",
$ mamba install "libblas=*=*netlib"
...
libblas 3.9.0-21_linux64_openblas --> 3.9.0-5_h92ddd45_netlib
libcblas 3.9.0-21_linux64_openblas --> 3.9.0-5_h92ddd45_netlib
liblapack 3.9.0-21_linux64_openblas --> 3.9.0-5_h92ddd45_netlib
$ mamba install "libblas=*=*mkl"
...
libblas 3.9.0-5_h92ddd45_netlib --> 3.9.0-21_linux64_mkl
libcblas 3.9.0-5_h92ddd45_netlib --> 3.9.0-21_linux64_mkl
liblapack 3.9.0-5_h92ddd45_netlib --> 3.9.0-21_linux64_mkl
$ python -m threadpoolctl -i scipy.linalg
...
"user_api": "blas",
"internal_api": "mkl",
This can be done for development builds as well, when building via dev.py
in the exact same way as in SciPy’s conda-forge build recipe
(outputs omitted for brevity, they’re similar to the ones above):
$ mamba env create -f environment.yml
$ mamba activate scipy-dev
$ mamba install "libblas=*=*netlib" # necessary, we need to build against blas/lapack
$ python dev.py build -C-Dblas=blas -C-Dlapack=lapack -C-Duse-g77-abi=true
$ python dev.py test -s linalg # run tests to verify
$ mamba install "libblas=*=*mkl"
$ python dev.py test -s linalg
$ mamba install "libblas=*=*openblas"
Linux distro package managers#
A number of Linux distros use the update-alternatives
mechanism to allow
switching between different BLAS libraries via the system package manager. Note
that this is a generic mechanism to manage “multiple implementations of the
same library or tool” situations, rather than something specific to
BLAS/LAPACK. It’s similar to the conda-forge method above, in that it works for
distro-provided scipy
packages as well as for development builds against
the reference libblas
/liblapack
interfaces.
The interface looks like:
$ update-alternatives --config libblas.so.3
$ update-alternatives --config liblapack.so.3
which will open a menu in the terminal with all available libraries to choose from. Because the interface and available options are likely to vary across distros, we link here to the Debian documentation for BLAS/LAPACK switching and avoid documenting in more detail how this works on other distros.
Note that Fedora is an exception; it is the only distro that ships FlexiBLAS (see the next section for more on that) and allows installing multiple BLAS libraries in parallel so true runtime switching without having to invoke the system package manager becomes possible. See the Fedora docs on system-level and user-level selection for more details.
FlexiBLAS#
FlexiBLAS provides runtime
switching support (among other things) for all installed BLAS libraries that it
can detect. There are a few limitations at the time of writing (March 2024),
primarily: no support for Windows, no support for macOS Accelerate (the updated
version, with NEWLAPACK
symbols). If those limitations don’t matter for
you, FlexiBLAS can be a quite useful tool for efficient debugging, including
for versions of OpenBLAS and other BLAS libraries that you have to build from
source.
Once you have everything set up, the development experience is:
$ python dev.py build -C-Dblas=flexiblas -C-Dlapack=flexiblas
$ FLEXIBLAS=NETLIB python dev.py test -s linalg
$ FLEXIBLAS=OpenBLAS python dev.py test -s linalg
# Or export the environment variable to make the selection stick:
$ export FLEXIBLAS=OpenBLAS
You can also provide a path to a built BLAS library (e.g.,
FLEXIBLAS="libbhlas_atlas.so"
) - see the usage docs in its README
for more details.
Unless you’re on Fedora, you will likely have to build FlexiBLAS from source, which is a bit of work. The good news is that this should work no matter if you’re on Linux or macOS, and use Python via virtualenvs, conda environments, or in some other way. We’ll go through how to build OpenBLAS and FlexiBLAS from source, to allow debugging whether something in the latest OpenBLAS version is different from Netlib BLAS/LAPACK (or MKL) or not.
The below should work in any environment where you can build SciPy itself; the
only additional tool we need is CMake (install with, for example, pip install
cmake
).
Clone each repository:
$ cd .. # starting from the root of the local `scipy` repo
$ mkdir flexiblas-setup && cd flexiblas-setup
$ git clone https://github.com/OpenMathLib/OpenBLAS.git openblas
$ git clone https://github.com/mpimd-csc/flexiblas.git
$ mkdir built-libs # our local prefix to install everything to
Build OpenBLAS:
$ cd openblas
$ mkdir build && cd build
$ cmake .. -DBUILD_SHARED_LIBS=ON -DCMAKE_INSTALL_PREFIX=$PWD/../../built-libs
$ cmake --build . -j
$ cmake --install . --prefix $PWD/../../built-libs
$ cd ../..
Build FlexiBLAS:
$ cd flexiblas
$ mkdir build && cd build
$ # Note: this will also pick up the libraries in your system/env libdir
$ cmake .. -DEXTRA="OpenBLAS" -DLAPACK_API_VERSION=3.9.0 \
-DOpenBLAS_LIBRARY=$PWD/../../built-libs/lib/libopenblas.so \
-DCMAKE_INSTALL_PREFIX=$PWD/../../built-libs
$ cmake --build . -j
$ cmake --install . --prefix $PWD/../../built-libs
$ cd ../..
We’re now ready to build SciPy against FlexiBLAS:
$ export PKG_CONFIG_PATH=$PWD/flexiblas-setup/built-libs/lib/pkgconfig/
$ cd scipy
$ python dev.py build -C-Dblas=flexiblas -C-Dlapack=flexiblas
...
Run-time dependency flexiblas found: YES 3.4.2
Now we can run the tests. Note that the NETLIB
option is built without
having to specify it; it’s the default in FlexiBLAS and sources are included in
its repository:
$ FLEXIBLAS=OpenBLAS python dev.py test -s linalg
$ FLEXIBLAS=NETLIB python dev.py test -s linalg
$ python dev.py test -s linalg # uses the default (NETLIB)
This backend switching can also be done inside a Python interpreter with
threadpoolctl
(see its README
for details).
Other libraries available on the system can be inspected with:
$ ./flexiblas-setup/built-libs/bin/flexiblas list
Note
Using local builds of reference BLAS/LAPACK or BLIS is more difficult,
because FlexiBLAS requires a single shared library which contains all
needed symbols. It may be feasible
to use a separate libblas
and liblapack
as the “system library”,
but this has proven to be more fragile and difficult to build (so this is
YMMV). In case you do want to try:
Build reference BLAS and LAPACK:
$ git clone Reference-LAPACK/lapack.git $ cd lapack $ mkdir build && cd build $ cmake -DCBLAS=ON -DBUILD_SHARED_LIBS=OFF .. $ cmake –build . -j $ cmake –install . –prefix $PWD/../../built-libs
Then add the following two lines to the cmake ..
configure command for
FlexiBLAS:
-DSYS_BLAS_LIBRARY=$PWD/../../built-libs/lib/libblas.a \
-DSYS_LAPACK_LIBRARY=$PWD/../../built-libs/lib/liblapack.a \
Creating reproducers in C or Fortran#
Our experience tells us that a large majority of bugs are inside SciPy rather than in OpenBLAS or another BLAS library. If the testing with different BLAS libraries tells us though that the problem is specific to a single BLAS library (maybe even a single version of that library with a regression), the next step is to produce a reproducer in C or Fortran; doing so is necessary for reporting the bug upstream, and makes it much easier for the BLAS library developers to address the problem.
To get from a Python reproducer which uses a scipy
function with NumPy
arrays as input to a C/Fortran reproducer, it is necessary to find the code
path taken in SciPy and determine which exact BLAS or LAPACK function is
called, and with what inputs (note: the answer may be contained in the
.pyf.src
f2py signature files; looking into the generated
_flapackmodule.c
in the build directory may be useful too). This can then
be reproduced in C/Fortran by defining some integer/float variables and arrays
(typically small arrays with hardcoded numbers are enough).
Argument lists of BLAS and LAPACK functions can be looked up in for example the Netlib LAPACK docs or the Reference-LAPACK/lapack repository.
Below a reproducer is shown for an issue in reference LAPACK, which was
reported as a SciPy issue in scipy#11577. We’ll name the file
ggev_repro_gh_11577.c|f90
:
#include <stdio.h>
#include "lapacke.h"
#define n 4
int main()
{
int lda=n, ldb=n, ldvr=n, ldvl=n, info;
char jobvl='V', jobvr='V';
double alphar[n], alphai[n], beta[n];
double vl[n*n], vr[n*n];
// int lwork = 156;
// double work[156]; /* cheat: 156 is the optimal lwork from the actual lapack call*/
double a[n*n] = {12.0, 28.0, 76.0, 220.0,
16.0, 32.0, 80.0, 224.0,
24.0, 40.0, 88.0, 232.0,
40.0, 56.0, 104.0, 248.0};
double b[n*n] = {2.0, 4.0, 10.0, 28.0,
3.0, 5.0, 11.0, 29.0,
5.0, 7.0, 13.0, 31.0,
9.0, 11.0, 17.0, 35.0};
info = LAPACKE_dggev(LAPACK_ROW_MAJOR, jobvl, jobvr, n, a, lda, b,
ldb, alphar, alphai, beta, vl, ldvl, vr, ldvr); //, work, lwork, info);
printf("info = %d\n", info);
printf("Re(eigv) = ");
for(int i=0; i < n; i++){
printf("%f , ", alphar[i] / beta[i] );
}
printf("\nIm(eigv = ");
for(int i=0; i < n; i++){
printf("%f , ", alphai[i] / beta[i] );
}
printf("\n");
}
!A = numpy.array([[12.0, 28.0, 76.0, 220.0], [16.0, 32.0, 80.0, 224.0], [24.0, 40.0, 88.0, 232.0], [40.0, 56.0, 104.0, 248.0]], dtype='float64')
! B = numpy.array([[2.0, 4.0, 10.0, 28.0], [3.0, 5.0, 11.0, 29.0], [5.0, 7.0, 13.0, 31.0], [9.0, 11.0, 17.0, 35.0]], dtype='float64')
! D, V = scipy.linalg.eig(A, B); D
implicit none
integer, parameter :: n = 4
integer :: lda, ldb, ldvr, ldvl, lwork, info
character*1 :: jobvl, jobvr
real*8 :: alphar(n)
real*8 :: alphai(n)
real*8 :: beta(n)
real*8 :: vl(n, n)
real*8 :: vr(n, n)
real*8, allocatable :: work(:)
real*8 :: a(n, n)
real*8 :: b(n, n)
a(1, :) = (/12.0, 28.0, 76.0, 220.0/)
a(2, :) = (/16.0, 32.0, 80.0, 224.0/)
a(3, :) = (/24.0, 40.0, 88.0, 232.0/)
a(4, :) = (/40.0, 56.0, 104.0, 248.0/)
b(1, :) = (/2.0, 4.0, 10.0, 28.0/)
b(2, :) = (/3.0, 5.0, 11.0, 29.0/)
b(3, :) = (/5.0, 7.0, 13.0, 31.0/)
b(4, :) = (/9.0, 11.0, 17.0, 35.0/)
lda = n
ldb = n
ldvr = n
ldvl = n
jobvr = 'V'
jobvl = 'V'
! workspace query
allocate(work(1:8*n)) ! min value
lwork = -1
print*, 'workspace query: lwork = ', lwork
call dggev(jobvl, jobvr, n, a, lda, b, ldb, alphar, alphai, beta, vl, ldvl, vr, ldvr, work, lwork, info)
print*, 'info = ', info
lwork = int(work(1))
print*, 'opt lwork =', lwork
! do the work
deallocate(work)
allocate(work(1:lwork))
call dggev(jobvl, jobvr, n, a, lda, b, ldb, alphar, alphai, beta, vl, ldvl, vr, ldvr, work, lwork, info)
print*
print*, 'info = ', info
print*, 'alphar = ', alphar
print*, 'alphai = ', alphai
print*, 'beta = ', beta
print*
print*, 'Re(eigv) = ', alphar / beta
print*, 'Im(eigv) = ', alphai / beta
deallocate(work)
end
Now we need to compile this reproducer locally and run it. If we’re invoking a compiler directly, we need to add the needed compile and link flags by hand. The include path will depend on your local install, and the link flags will depend on which library you’re testing. For example, to test against a local build of OpenBLAS:
$ gcc ggev_repro_gh_11577.c \
-I$PWD/../flexiblas-setup/built-libs/include/ \
-L$PWD/../flexiblas-setup/built-libs/lib -lopenblas
$ ./a.out # to run the reproducer
$ gfortran ggev_repro_gh_11577.f90 \
-I/$PWD/../flexiblas-setup/built-libs/include/ \
-L$PWD/../flexiblas-setup/built-libs/lib -lopenblas
$ ./a.out # to run the reproducer
For reference BLAS/LAPACK, the -lopenblas
should be replaced with
-lblas -llapack
.
Note that the explicit paths are only needed for libraries in non-standard
locations (like the ones we built in this guide). For building against a
package manager-installed library for which the shared library and headers are
on the normal compiler search path (e.g., in /usr/lib
and /usr/include
,
or inside a conda env when using compilers from the same env), they can be left out:
$ gcc ggev_repro_gh_11577.c -lopenblas
$ ./a.out # to run the reproducer
$ gfortran ggev_repro_gh_11577.f90 -lopenblas
$ ./a.out # to run the reproducer
Alternatively (and probably a more robust way), use a small meson.build
file to automate this and avoid the manual paths:
project('repro_gh_11577', 'c')
openblas_dep = dependency('openblas')
executable('repro_c',
'ggev_repro_gh_11577.c',
dependencies: openblas_dep
)
To then build the test and run it:
$ export PKG_CONFIG_PATH=~/code/tmp/flexiblas-setup/built-libs/lib/pkgconfig/
$ meson setup build
$ ninja -C build
$ ./build/repro_c # output may vary
info = 0
Re(eigv) = 4.000000 , 8.000000 , inf , -inf ,
Im(eigv = 0.000000 , 0.000000 , -nan , -nan ,
project('repro_gh_11577', 'fortran')
openblas_dep = dependency('openblas')
executable('repro_f90',
'ggev_repro_gh_11577.f90',
dependencies: openblas_dep
)
To then build the test and run it:
$ export PKG_CONFIG_PATH=~/code/tmp/flexiblas-setup/built-libs/lib/pkgconfig/
$ meson setup build
$ ninja -C build
$ ./build/repro_f90 # output may vary
workspace query: lwork = -1
info = 0
opt lwork = 156
info = 0
alphar = 1.0204501477442456 11.707793036240817 3.7423579363517347E-014 -1.1492523608519701E-014
alphai = 0.0000000000000000 0.0000000000000000 0.0000000000000000 0.0000000000000000
beta = 0.25511253693606051 1.4634741295300704 0.0000000000000000 0.0000000000000000
Re(eigv) = 4.0000000000000142 8.0000000000001741 Infinity -Infinity
Im(eigv) = 0.0000000000000000 0.0000000000000000 NaN NaN
Warning
When you have multiple versions/builds of the same BLAS library on your
machine, it’s easy to accidentally pick up the wrong one during the build
(remember: -lopenblas
only says “link against some
libopenblas.so
). If you’re not sure, use ldd
on the test executable
you built to inspect which shared library it’s linked against.
Debugging linalg issues with gdb
#
When debugging linalg issues, it is sometimes useful to step through both Python
and C code. You can use pdb
for the former, and gdb
for the latter.
First, prepare a small python reproducer, with a breakpoint. For example:
$ cat chol.py
import numpy as np
from scipy import linalg
n = 40
rng = np.random.default_rng(1234)
x = rng.uniform(size=n)
a = x[:, None] @ x[None, :] + np.identity(n)
breakpoint() # note a breakpoint
linalg.cholesky(a)
Then, you will need to run it under the gdb
and add a C-level breakpoint at
the LAPACK function. This way, your execution will stop twice: first on the
Python breakpoint and then on the C breakpoint, and you will be able to step
through and inspect values of both Python and C variables.
To find out the LAPACK name, read the python source of the SciPy function and
use nm
on the .so
library to find out the exact name.
For the Cholesky factorization above, the LAPACK function is ?potrf
, and the
C name on Ubuntu linux is dpotrf_
(it may be spelled with or without the
trailing underscore, in uppper case or lower case, depending on the system).
Here is an example gdb
session:
$ gdb --args python chol.py
...
(gdb) b dpotrf_ # this adds a C breakpoint (type "y" below)
Function "dpotrf_" not defined.
Make breakpoint pending on future shared library load? (y or [n]) y
Breakpoint 1 (dpotrf_) pending.
(gdb) run # run the python script
...
> /home/br/temp/chol/chol.py(12)<module>()
-> linalg.cholesky(a) # execution stopped at the python breakpoint
(Pdb) p a.shape # ... inspect the python state here
(40, 40)
(Pdb) c # continue execution until the C breakpoint
Thread 1 "python" hit Breakpoint 1, 0x00007ffff4c48820 in dpotrf_ ()
from /home/br/mambaforge/envs/scipy-dev/lib/python3.10/site-packages/numpy/core/../../../../libcblas.so.3
(gdb) s # step through the C function
Single stepping until exit from function dpotrf_,
which has no line number information.
f2py_rout__flapack_dpotrf (capi_self=<optimized out>, capi_args=<optimized out>,
capi_keywds=<optimized out>, f2py_func=0x7ffff4c48820 <dpotrf_>)
at scipy/linalg/_flapackmodule.c:63281
....
(gdb) p lda # inspect values of C variables
$1 = 40
# print out the C backtrace
(gdb) bt
#0 0x00007ffff3056b1e in f2py_rout.flapack_dpotrf ()
from /path/to/site-packages/scipy/linalg/_flapack.cpython-311-x86_64-linux-gnu.so
#1 0x0000555555734323 in _PyObject_MakeTpCall (tstate=0x555555ad0558 <_PyRuntime+166328>,
callable=<fortran at remote 0x7ffff40ffc00>, args=<optimized out>, nargs=1,
keywords=('lower', 'overwrite_a', 'clean'))
at /usr/local/src/conda/python-3.11.8/Objects/call.c:214
...
Depending on your system, you may need to build SciPy with debug build type, and unset CFLAGS/CXXFLAGS environment variables. See the NumPy debugging guide for more details.