Отладка проблем, связанных с линейной алгеброй#
Отчеты об ошибках, связанных с линейной алгеброй, являются одними из самых сложных для диагностики и исправления. Это связано не только с тем, что линейная алгебра может быть сложной математически/алгоритмически (это верно для многих частей SciPy), но и потому, что библиотеки BLAS/LAPACK представляют собой сложную зависимость как на этапе сборки, так и во время выполнения — и мы поддерживаем значительное количество библиотек BLAS/LAPACK.
Этот документ призван предоставить руководство по отладке проблем линейной алгебры.
Если есть реальная ошибка, она может быть в одном из трёх мест:
В используемой библиотеке BLAS,
В привязках SciPy к BLAS или LAPACK (сгенерированных
numpy.f2pyи/или Cython),В алгоритмическом коде SciPy.
Ключевым первым шагом является определение, находится ли ошибка в SciPy или в библиотеке BLAS. Для этого наиболее эффективный способ разграничить их — настроить ваше окружение таким образом, чтобы можно было переключаться между различными библиотеками BLAS во время выполнения (что мы не поддерживаем из коробки и невозможно с колесами SciPy из PyPI).
Авторы библиотек BLAS предпочитают получать чистые воспроизводимые примеры (так же, как и мы), что означает: без участия Python. Поэтому это руководство также охватывает, как создавать воспроизводимые примеры на C или Fortran.
Поиск используемой библиотеки BLAS#
SciPy имеет одну функцию, show_config, для интроспекции конфигурации сборки
установленного пакета. Это позволяет запрашивать детали BLAS
и LAPACK. Например:
>>> 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/miniforge/envs/scipy-dev/include
lib directory: /home/user/miniforge/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/miniforge/envs/scipy-dev/lib/pkgconfig
Этот метод будет корректным для SciPy wheels и для локальных dev-сборок. Он может
может быть корректным для других установок, однако учтите, что дистрибутивы, такие как
conda-forge и Debian, могут собираться с использованием заглушек библиотек (обычно
blas.so/lapack.so) и затем установка другой библиотеки для пользователя -
в таких случаях, простой blas и lapack будет сообщено, даже если, например,
OpenBLAS или MKL установлены в среде. Для таких установок,
threadpoolctl обычно сможет
сообщить, какая именно библиотека BLAS используется (за исключением того, что не сообщает
о простом Netlib BLAS, см.
threadpoolctl#159):
$ python -m threadpoolctl -i scipy.linalg
[
{
"user_api": "blas",
"internal_api": "openblas",
"prefix": "libopenblas",
"filepath": "/home/user/miniforge/envs/dev/lib/libopenblasp-r0.3.21.so",
"version": "0.3.21",
"threading_layer": "pthreads",
"architecture": "SkylakeX",
"num_threads": 24
}
]
Другие способы интроспекции, которые могут быть полезны в локальных средах разработки, включают:
Проверка зависимостей общих библиотек:
$ ldd build/scipy/linalg/_fblas.cpython-*.so
...
libopenblas.so.0 => /home/user/miniforge/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)
Проверка, содержит ли связанная библиотека данный символ. Например, conda-forge устанавливает
libblas.soкоторый может быть любой поддерживаемой библиотекой:
$ nm -gD ~/miniforge/envs/scipy-dev/lib/libblas.so | rg openblas_set_num_threads
0000000000362990 T openblas_set_num_threads
% nm ~/miniforge/envs/scipy-dev/lib/libblas.3.dylib | rg openblas_set_num_threads
000000000015b6b0 T _openblas_set_num_threads
Настройка окружения для переключения библиотек BLAS#
Мы рассмотрим несколько способов переключения между различными библиотеками BLAS, потому что самый простой метод может зависеть от вашей ОС/дистрибутива и от того, хотите ли вы релизную версию SciPy или сборку для разработки.
Conda-forge#
Возможно, самый простой способ — использовать возможности переключения во время выполнения, предоставляемые дистрибутивами. Например, создать окружение conda с установленным последним выпуском SciPy, а затем переключаться между OpenBLAS, Netlib BLAS/LAPACK и MKL так же просто, как:
$ 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",
Это можно сделать и для сборок разработки при сборке через dev.py
точно так же, как в Рецепт сборки SciPy для conda-forge
(выводы опущены для краткости, они похожи на приведённые выше):
$ 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#
Ряд дистрибутивов Linux используют update-alternatives механизм для переключения между различными библиотеками BLAS через системный менеджер пакетов. Обратите внимание,
что это общий механизм для управления ситуациями с 'несколькими реализациями
одной и той же библиотеки или инструмента', а не что-то специфичное для
BLAS/LAPACK. Он похож на метод conda-forge выше, так как работает с scipy пакеты, а также для сборок разработки с использованием эталонного libblas/liblapack интерфейсы.
Интерфейс выглядит так:
$ update-alternatives --config libblas.so.3
$ update-alternatives --config liblapack.so.3
что откроет меню в терминале со всеми доступными библиотеками для выбора. Поскольку интерфейс и доступные опции, вероятно, различаются в разных дистрибутивах, мы ссылаемся здесь на документация Debian по переключению BLAS/LAPACK и избегать подробного документирования работы на других дистрибутивах.
Обратите внимание, что Fedora является исключением; это единственный дистрибутив, который поставляет FlexiBLAS (см. следующий раздел для подробностей) и позволяет устанавливать несколько библиотек BLAS параллельно, что делает возможным переключение во время выполнения без вызова менеджера пакетов системы. См. документация Fedora о выборе на уровне системы и пользователя для получения дополнительной информации.
FlexiBLAS#
FlexiBLAS предоставляет поддержку переключения во время выполнения
(среди прочего) для всех установленных библиотек BLAS, которые он
может обнаружить. На момент написания (март 2024) есть несколько ограничений,
в основном: нет поддержки Windows, нет поддержки macOS Accelerate (обновлённой
версии, с NEWLAPACK символы). Если эти ограничения для вас не важны, FlexiBLAS может быть весьма полезным инструментом для эффективной отладки, включая версии OpenBLAS и других библиотек BLAS, которые вам приходится собирать из исходного кода.
После того как все настроено, опыт разработки следующий:
$ 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
Вы также можете указать путь к собранной библиотеке BLAS (например,
FLEXIBLAS="libbhlas_atlas.so") - см. документация по использованию в его README
для получения дополнительной информации.
Если вы не используете Fedora, вам, скорее всего, придётся собирать FlexiBLAS из исходного кода, что требует некоторых усилий. Хорошая новость в том, что это должно работать независимо от того, используете ли вы Linux или macOS, и Python через virtualenv, conda-окружения или каким-либо другим способом. Мы рассмотрим, как собрать OpenBLAS и FlexiBLAS из исходного кода, чтобы проверить, отличается ли что-то в последней версии OpenBLAS от Netlib BLAS/LAPACK (или MKL) или нет.
Нижеследующее должно работать в любой среде, где можно собрать сам SciPy; единственный дополнительный инструмент, который нам нужен, — это CMake (установите, например, с помощью pip install
cmake).
Клонировать каждый репозиторий:
$ 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
Сборка 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 ../..
Сборка 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 ../..
Теперь мы готовы собрать SciPy с использованием 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
Теперь мы можем запустить тесты. Обратите внимание, что NETLIB опция собрана без необходимости её указания; это значение по умолчанию в FlexiBLAS, и исходники включены в его репозиторий:
$ 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)
Это переключение бэкенда также можно выполнить внутри интерпретатора Python с помощью
threadpoolctl (см. его README
подробности).
Другие библиотеки, доступные в системе, можно проверить с помощью:
$ ./flexiblas-setup/built-libs/bin/flexiblas list
Примечание
Использование локальных сборок эталонных BLAS/LAPACK или BLIS сложнее,
потому что FlexiBLAS требует единой общей библиотеки, содержащей все
необходимые символы. Это может быть выполнимо
использовать отдельный libblas и liblapack как «системная библиотека»,
но это оказалось более хрупким и сложным для сборки (так что это
зависит от ситуации). Если вы всё же хотите попробовать:
Сборка эталонных BLAS и 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
Затем добавьте следующие две строки в cmake .. команда настройки для
FlexiBLAS:
-DSYS_BLAS_LIBRARY=$PWD/../../built-libs/lib/libblas.a \
-DSYS_LAPACK_LIBRARY=$PWD/../../built-libs/lib/liblapack.a \
Создание воспроизводимых примеров на C или Fortran#
Наш опыт говорит, что подавляющее большинство ошибок находятся внутри SciPy, а не в OpenBLAS или другой библиотеке BLAS. Однако, если тестирование с разными библиотеками BLAS показывает, что проблема специфична для одной библиотеки BLAS (возможно, даже для одной версии этой библиотеки с регрессией), следующий шаг — создать воспроизводимый пример на C или Fortran; это необходимо для сообщения об ошибке вышестоящим разработчикам и значительно упрощает им решение проблемы.
Чтобы перейти от воспроизводимого примера на Python, который использует scipy функция с массивами NumPy
в качестве входных данных для воспроизводителя на C/Fortran, необходимо найти путь
кода, пройденный в SciPy, и определить, какая именно функция BLAS или LAPACK
вызывается, и с какими входными данными (примечание: ответ может содержаться в
.pyf.src файлы сигнатур f2py; анализ сгенерированного
_flapackmodule.c в каталоге сборки также может быть полезно). Это затем
можно воспроизвести на C/Fortran, определив некоторые целочисленные/вещественные переменные и массивы
(обычно достаточно небольших массивов с жёстко заданными числами).
Списки аргументов функций BLAS и LAPACK можно найти, например, в документация Netlib LAPACK или репозиторий Reference-LAPACK/lapack.
Ниже показан воспроизводимый пример для проблемы в эталонной LAPACK, которая была
сообщена как проблема SciPy в scipy#11577. Мы назовем файл
ggev_repro_gh_11577.c|f90:
#include
#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
Теперь нам нужно скомпилировать этот воспроизводитель локально и запустить его. Если мы вызываем компилятор напрямую, нам нужно вручную добавить необходимые флаги компиляции и линковки. Путь включения будет зависеть от вашей локальной установки, а флаги линковки будут зависеть от того, какую библиотеку вы тестируете. Например, для тестирования локальной сборки 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
Для эталонной BLAS/LAPACK, -lopenblas должен быть заменён на
-lblas -llapack.
Обратите внимание, что явные пути нужны только для библиотек в нестандартных
расположениях (как те, что мы собрали в этом руководстве). Для сборки против
библиотеки, установленной менеджером пакетов, для которой разделяемая библиотека и заголовки
находятся в обычном пути поиска компилятора (например, в /usr/lib и /usr/includeили внутри conda env при использовании компиляторов из того же env), их можно опустить:
$ 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
Альтернативно (и, вероятно, более надёжный способ) — использовать небольшой meson.build
файл для автоматизации этого и избежания ручных путей:
project('repro_gh_11577', 'c')
openblas_dep = dependency('openblas')
executable('repro_c',
'ggev_repro_gh_11577.c',
dependencies: openblas_dep
)
Чтобы затем собрать тест и запустить его:
$ 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
)
Чтобы затем собрать тест и запустить его:
$ 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
Предупреждение
Когда у вас есть несколько версий/сборок одной и той же библиотеки BLAS на вашем
компьютере, легко случайно выбрать неправильную во время сборки
(помните: -lopenblas только говорит "линковать против" некоторые
libopenblas.so). Если вы не уверены, используйте ldd в тестовом исполняемом файле, который вы собрали, чтобы проверить, с какой общей библиотекой он связан.
Отладка проблем с linalg с gdb#
При отладке проблем с linalg иногда полезно пройтись как по Python, так и по C коду. Вы можете использовать pdb для первого, и gdb для последнего.
Сначала подготовьте небольшой воспроизводимый пример на Python с точкой останова. Например:
$ 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)
Затем вам нужно будет запустить его под gdb и добавьте точку останова на уровне C в функции LAPACK. Таким образом, ваше выполнение остановится дважды: сначала на точке останова Python, а затем на точке останова C, и вы сможете пройти пошагово и проверить значения как переменных Python, так и C.
Чтобы узнать имя LAPACK, прочитайте исходный код функции SciPy на Python и
используйте nm на .so библиотеки, чтобы узнать точное название.
Для факторизации Холецкого выше, функция LAPACK - ?potrf, и имя C в Ubuntu Linux - dpotrf_ (он может писаться с конечным подчеркиванием или без него,
в верхнем или нижнем регистре, в зависимости от системы).
Вот пример gdb сессия:
$ 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)()
-> 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/miniforge/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=, capi_args=,
capi_keywds=, f2py_func=0x7ffff4c48820 )
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=, args=, nargs=1,
keywords=('lower', 'overwrite_a', 'clean'))
at /usr/local/src/conda/python-3.11.8/Objects/call.c:214
...
В зависимости от вашей системы, вам может потребоваться собрать SciPy с типом сборки debug, и сбросить переменные окружения CFLAGS/CXXFLAGS. См. Руководство по отладке NumPy для более подробной информации.