пятница, 10 марта 2017 г.

GFortran: вызываем функции LAPACK из DLL динамически

Как известно, существует два способа связывания DLL-библиотек с остальной частью программы:

  1. Статическое, или неявное связывание. Все заботы берет на себя компилятор (точнее, линковщик), вы просто вызываете функцию, содержащуюся в DLL. Для этого ему нужно иметь соответствующую (сопоставленную) нашей DLL библиотеку *.a (или *.lib), которая содержит необходимую для линковщика информацию. Все связанные так DLL подключаются сразу после запуска программа.

  2. Динамическое или явное связывание. Может быть выполнено в любое время работы программы и выполняется программистом «руками». Библиотеки *.a (или *.lib) – не требуются.

Недостаток первого метода для нас заключается в том, что все DLL должны иметься на момент запуска программы. Если какой-то DLL нет, то программа не запустится. А нам иногда очень хочется подключить ту или иную DLL в зависимости от настроек. Кроме того, первый способ не позволяет указать путь к нужным DLL, из-за чего может получиться, что подключатся какие-то не те версии.

Например, как мы уже отмечали, библиотека линейной алгебры LAPACK имеется как минимум в трех реализациях: классической, ObenBLAS и MKL (вообще-то, их куда больше). Для экспериментов нам хотелось бы иметь возможность переключатся между ними «на лету». Пусть эти библиотеки имеются в отдельных DLL. (На самом деле, в каждом случае нужна не одна, а несколько связанных библиотек). Напишем модуль, реализующий их динамическое подключение (на примере трех функций LAPACK: dgesv, dposv и dsysv). На Фортране он получается чуть-чуть мудреным, но в принципе, все не так страшно:

!> ******   Модуль с интерфейсом подключения процедур LAPACK из DLL   ******
module LapackFromDll
    use, intrinsic :: iso_c_binding, only:  c_f_procpointer, &
         c_funptr, c_intptr_t, c_null_char, c_char, c_associated, c_bool

    implicit none

    !> Интерфейсы функции WinAPI, работающих с DLL
    interface
        ! Интерфейс загрузки DLL (из Kernel32.lib)
        function LoadLibrary(lpFileName) bind(c, name='LoadLibraryA')
            use, intrinsic :: iso_c_binding, only: c_intptr_t, c_char
            implicit none 
            character(kind=c_char)          :: lpFileName(*) 
            !GCC$ attributes stdcall        :: LoadLibrary 
            integer(c_intptr_t)             :: LoadLibrary 
        endfunction 
        ! Интерфейс получения адреса функции DLL (из Kernel32.dll)
        function GetProcAddress(hModule, lpProcName) bind(c, name='GetProcAddress')
            use, intrinsic :: iso_c_binding, only: c_funptr, c_intptr_t, c_char
            implicit none
            !GCC$ attributes stdcall        :: GetProcAddress
            type(c_funptr)                  :: GetProcAddress
            integer(c_intptr_t), value      :: hModule
            character(kind=c_char)          :: lpProcName(*)
        endfunction 
        ! Интерфейс освобождения DLL (из Kernel32.dll)
        function FreeLibrary(hModule) bind(c, name='FreeLibrary')
            use, intrinsic :: iso_c_binding, only: c_funptr, c_intptr_t, c_bool
            implicit none
            !GCC$ attributes stdcall        :: FreeLibrary
            logical(c_bool)                 :: FreeLibrary
            integer(c_intptr_t), value      :: hModule
        endfunction  
    endinterface 

    !> Интерфейсы функций LAPACK
    abstract interface
        subroutine dgesv_(N, NRHS, A, LDA, IPIV, B, LDB, INFO) 
            integer            INFO, LDA, LDB, N, NRHS
            integer            IPIV(*)
            double precision   A(lda,*), B(ldb,*)
        endsubroutine
        subroutine dposv_(UPLO, N, NRHS, A, LDA, B, LDB, INFO) 
            character          UPLO
            integer            INFO, LDA, LDB, N, NRHS
            double precision   A(LDA,*), B(LDB,*)
        endsubroutine   
        subroutine dsysv_(UPLO, N, NRHS, A, LDA, IPIV, B, LDB, WORK, LWORK, INFO)
            character          UPLO
            integer            INFO, LDA, LDB, LWORK, N, NRHS
            integer            IPIV(*)
            double precision   A(LDA,*), B(LDB,*), WORK(*)
        endsubroutine  
    endinterface

    !> Процедуры-указатели на функции LAPACK
    procedure(dgesv_), pointer  :: dgesv        
    procedure(dposv_), pointer  :: dposv
    procedure(dsysv_), pointer  :: dsysv        

 contains  ! ============== Обертки функций для работы с DLL =================

    !> Загружает dll и связывает с ним Dll_handle
    subroutine LoadDll(Dll_handle, Dll_name)
        integer(c_intptr_t)     :: Dll_handle
        character(*, c_char)    :: Dll_name
        !------------------------------------------------------
        Dll_handle = LoadLibrary(Dll_name // c_null_char)
        if (Dll_handle == 0) stop 'Не удалось загрузить ' // Dll_name
    endsubroutine

    !> Ассоциирует процедуру-указатель Proc с функцией из Dll
    subroutine DllProcAssoc(Dll_handl, Proc, Proc_name)
        integer(c_intptr_t)     :: Dll_handl
        character(*, c_char)    :: Proc_name
        procedure(), pointer    :: Proc 
        type(c_funptr)          :: proc_address
        !------------------------------------------------------
        proc_address = GetProcAddress(Dll_handl, Proc_name // c_null_char)
        if (.not. c_associated(proc_address)) stop 'Не могу найти адрес процедуры '&
                                                    // Proc_name // " в DLL" 
        call c_f_procpointer(proc_address, Proc)
    endsubroutine

    !> Освобождает связанные с DLL ресурсы
    subroutine FreeDll(Dll_handle)
        integer(c_intptr_t)       :: Dll_handle
        logical(c_bool)           :: FrDll
        !------------------------------------------------------
        FrDll = FreeLibrary( Dll_handle )
        if (.not. FrDll) stop 'Не удалось выгрузить DLL'
    endsubroutine    

endmodule

Интерфейсы функций LAPACK можно копировать прямо из исходников с сайта www.netlib.org/lapack/.

Проверяем. Вызовем функцию dgesv «на лету» последовательно из трех реализаций LAPACK – классической, ObenBLAS и MKL:

program LinTest
    use, intrinsic :: iso_c_binding, only: c_intptr_t
    use  LapackFromDll
    implicit none
    integer(c_intptr_t)      :: LAPACK_dll ! Хендл LAPACK_dll
    integer, parameter       :: n = 10     ! Размерность матриц
    real(8), dimension(n,n)  :: A, A0      ! Матрица коэффициентов
    real(8), dimension(n,1)  :: B, X       ! Векторы правых частей / рез.
    integer                  :: i,j        ! Вспомогательный индекс
    integer, dimension(n)    :: ipiv       ! Вспомогательные для LAPACK
    integer                  :: info       ! Инфо об успешности решения
    character(100)           :: FRM        ! Формат для вывода матриц
    !---------------------------------------------------------------------
    print *, "Проверка динамической загрузки процедур LAPACK из DLL"
    forall(i=1:n, j=1:n) A0(i,j) = 1-(dble(i-j)/n)**2  ! Формируем матрицу
    forall(i=1:n, j=1:n, i==j) A0(i,j) = 10.
    forall(i = 1:n) B(i,1) = 100*i        ! Заполняем вектор B
    write(FRM, "(a,i10,a)") "(",n,"f8.2)" ! Формат для вывода матриц
    print *, "Матрица коэффициентов A:"
    print FRM, transpose(A0)
    print *, "Вектор правых частей B:"
    print FRM, B

    call LoadDll(LAPACK_dll, "liblapack.dll")
    call DllProcAssoc(LAPACK_dll, dgesv, "dgesv_")
    A = A0 ; X = B 
    call dgesv(n, 1, A, n, ipiv, X, n, info)  
    print "(' *** Решение с помощью liblapack.dll (info=',i1,'):')", info
    print FRM, X
    call FreeDll(LAPACK_dll)

    call LoadDll(LAPACK_dll, "libopenblas.dll")
    call DllProcAssoc(LAPACK_dll, dgesv, "dgesv_")
    A = A0 ; X = B 
    call dgesv(n, 1, A, n, ipiv, X, n, info)  
    print "(' *** Решение с помощью libopenblas.dll (info=',i1,'):')", info
    print FRM, X
    call FreeDll(LAPACK_dll)

    call LoadDll(LAPACK_dll, "mkl_rt.dll")
    call DllProcAssoc(LAPACK_dll, dgesv, "dgesv_")
    A = A0 ; X = B 
    call dgesv(n, 1, A, n, ipiv, X, n, info)  
    print "(' *** Решение с помощью mkl_rt.dll (info=',i1,'):')", info
    print FRM, X
    call FreeDll(LAPACK_dll)    
endprogram 

Компилируем:

gfortran -c LapackFromDll.F90
gfortran -c LinTest.F90
gfortran *.o -o LinTest.exe

Получаем:

 Проверка динамической загрузки процедур LAPACK из DLL
 Матрица коэффициентов A:
  10.00    0.99    0.96    0.91    0.84    0.75    0.64    0.51    0.36    0.19
   0.99   10.00    0.99    0.96    0.91    0.84    0.75    0.64    0.51    0.36
   0.96    0.99   10.00    0.99    0.96    0.91    0.84    0.75    0.64    0.51
   0.91    0.96    0.99   10.00    0.99    0.96    0.91    0.84    0.75    0.64
   0.84    0.91    0.96    0.99   10.00    0.99    0.96    0.91    0.84    0.75
   0.75    0.84    0.91    0.96    0.99   10.00    0.99    0.96    0.91    0.84
   0.64    0.75    0.84    0.91    0.96    0.99   10.00    0.99    0.96    0.91
   0.51    0.64    0.75    0.84    0.91    0.96    0.99   10.00    0.99    0.96
   0.36    0.51    0.64    0.75    0.84    0.91    0.96    0.99   10.00    0.99
   0.19    0.36    0.51    0.64    0.75    0.84    0.91    0.96    0.99   10.00
 Вектор правых частей B:
 100.00  200.00  300.00  400.00  500.00  600.00  700.00  800.00  900.00 1000.00
 *** Решение с помощью liblapack.dll (info=0):
  -6.20    0.36    7.63   15.60   24.29   33.68   43.77   54.58   66.09   78.30
 *** Решение с помощью libopenblas.dll (info=0):
  -6.20    0.36    7.63   15.60   24.29   33.68   43.77   54.58   66.09   78.30
 *** Решение с помощью mkl_rt.dll (info=0):
  -6.20    0.36    7.63   15.60   24.29   33.68   43.77   54.58   66.09   78.30

1 комментарий:

  1. As claimed by Stanford Medical, It is in fact the ONLY reason this country's women live 10 years longer and weigh on average 19 KG lighter than us.

    (And by the way, it has totally NOTHING to do with genetics or some secret diet and absolutely EVERYTHING to related to "how" they are eating.)

    BTW, What I said is "HOW", and not "WHAT"...

    Tap on this link to discover if this easy test can help you discover your real weight loss possibility

    ОтветитьУдалить