樣板程式碼縮減與模板#

使用 FYPP 綁定通用介面#

f2py 目前不支援綁定介面區塊。然而,目前有一些替代方案。其中最廣為人知的是使用 tempita 來使用 .pyf.src 檔案,如同 scipy 的綁定 中所做的那樣。tempita 支援已被移除,且在任何情況下都不再建議使用。

注意

介面無法在 f2py 本身中被支援的原因,是因為它們不對應到編譯後程式庫中輸出的符號。

 nm gen.o
 0000000000000078 T __add_mod_MOD_add_complex
 0000000000000000 T __add_mod_MOD_add_complex_dp
 0000000000000150 T __add_mod_MOD_add_integer
 0000000000000124 T __add_mod_MOD_add_real
 00000000000000ee T __add_mod_MOD_add_real_dp

在此我們將討論一些技術,以結合 f2pyfypp 來模擬通用介面,並簡化多個(相似)函式的綁定。

基本範例:加法模組#

讓我們以使用者指南中的範例(F2PY 範例)為基礎,該範例是一個子程序,它接收兩個陣列並傳回它們的總和。

C
      SUBROUTINE ZADD(A,B,C,N)
C
      DOUBLE COMPLEX A(*)
      DOUBLE COMPLEX B(*)
      DOUBLE COMPLEX C(*)
      INTEGER N
      DO 20 J = 1, N
         C(J) = A(J)+B(J)
 20   CONTINUE
      END

我們將把這個範例改寫成現代 Fortran

module adder
    implicit none
contains

    subroutine zadd(a, b, c, n)
        integer, intent(in) :: n
        double complex, intent(in) :: a(n), b(n)
        double complex, intent(out) :: c(n)
        integer :: j
        do j = 1, n
            c(j) = a(j) + b(j)
        end do
    end subroutine zadd

end module adder

我們可以像原始範例一樣繼續下去,手動添加 intents 等等,但是在生產環境中,通常還有其他考量。其中之一是,我們可以透過 FYPP 模板化相似函式的建構。

module adder
    implicit none
contains

#:def add_subroutine(dtype_prefix, dtype)
    subroutine ${dtype_prefix}$add(a, b, c, n)
        integer, intent(in) :: n
        ${dtype}$, intent(in) :: a(n), b(n)
        ${dtype}$ :: c(n)
        integer :: j
        do j = 1, n
            c(j) = a(j) + b(j)
        end do
    end subroutine ${dtype_prefix}$add

#:enddef

#:for dtype_prefix, dtype in [('i', 'integer'), ('s', 'real'), ('d', 'real(kind=8)'), ('c', 'complex'), ('z', 'double complex')]
    @:add_subroutine(${dtype_prefix}$, ${dtype}$)
#:endfor

end module adder

這可以預先處理以產生完整的 Fortran 程式碼

 fypp gen_adder.f90.fypp > adder.f90

不出所料,這可以隨後被 f2py 包裝。

現在我們將考慮在單獨的檔案中維護綁定。請注意以下基本的 .pyf,它可以透過 f2py -m adder adder_base.f90 -h adder.pyf 為單個子程序生成

!    -*- f90 -*-
! Note: the context of this file is case sensitive.

python module adder ! in
    interface  ! in :adder
        module adder ! in :adder:adder_base.f90
            subroutine zadd(a,b,c,n) ! in :adder:adder_base.f90:adder
                double complex dimension(n),intent(in) :: a
                double complex dimension(n),intent(in),depend(n) :: b
                double complex dimension(n),intent(out),depend(n) :: c
                integer, optional,intent(in),check(shape(a, 0) == n),depend(a) :: n=shape(a, 0)
            end subroutine zadd
        end module adder
    end interface
end python module adder

! This file was auto-generated with f2py (version:2.0.0.dev0+git20240101.bab7280).
! See:
! https://web.archive.org/web/20140822061353/http://cens.ioc.ee/projects/f2py2e

帶有 docstring

c = zadd(a,b,[n])

Wrapper for ``zadd``.

Parameters
----------
a : input rank-1 array('D') with bounds (n)
b : input rank-1 array('D') with bounds (n)

Other Parameters
----------------
n : input int, optional
    Default: shape(a, 0)

Returns
-------
c : rank-1 array('D') with bounds (n)

這已經相當不錯了。然而,n 根本不應該被傳遞,所以我們將做一些小的調整。

!    -*- f90 -*-
! Note: the context of this file is case sensitive.

python module adder ! in
    interface  ! in :adder
        module adder ! in :adder:adder_base.f90
            subroutine zadd(a,b,c,n) ! in :adder:adder_base.f90:adder
                integer intent(hide),depend(a) :: n=len(a)
                double complex dimension(n),intent(in) :: a
                double complex dimension(n),intent(in),depend(n) :: b
                double complex dimension(n),intent(out),depend(n) :: c
            end subroutine zadd
        end module adder
    end interface
end python module adder

! This file was auto-generated with f2py (version:2.0.0.dev0+git20240101.bab7280).
! See:
! https://numpy.dev.org.tw/doc/stable/f2py/

這對應到

In [3]: ?adder.adder.zadd
Call signature: adder.adder.zadd(*args, **kwargs)
Type:           fortran
String form:    <fortran function zadd>
Docstring:
c = zadd(a,b)

Wrapper for ``zadd``.

Parameters
----------
a : input rank-1 array('D') with bounds (n)
b : input rank-1 array('D') with bounds (n)

Returns
-------
c : rank-1 array('D') with bounds (n)

最後,我們可以以類似的方式對其進行模板化,以達到最初的目標,即擁有利用 f2py 指令並具有最少多餘重複的綁定。

!    -*- f90 -*-
! Note: the context of this file is case sensitive.

python module adder ! in
    interface  ! in :adder
        module adder ! in :adder:adder_base.f90
#:def add_subroutine(dtype_prefix, dtype)
            subroutine ${dtype_prefix}$add(a,b,c,n) ! in :adder:adder_base.f90:adder
                integer intent(hide),depend(a) :: n=len(a)
                ${dtype}$ dimension(n),intent(in) :: a
                ${dtype}$ dimension(n),intent(in),depend(n) :: b
                ${dtype}$ dimension(n),intent(out),depend(n) :: c
            end subroutine ${dtype_prefix}$add

#:enddef

#:for dtype_prefix, dtype in [('i', 'integer'), ('s', 'real'), ('d', 'real(kind=8)'), ('c', 'complex'), ('z', 'complex(kind=8)')]
    @:add_subroutine(${dtype_prefix}$, ${dtype}$)
#:endfor
        end module adder
    end interface
end python module adder

! This file was auto-generated with f2py (version:2.0.0.dev0+git20240101.bab7280).
! See:
! https://numpy.dev.org.tw/doc/stable/f2py/

使用方式歸結為

fypp gen_adder.f90.fypp > adder.f90
fypp adder.pyf.fypp > adder.pyf
f2py -m adder -c adder.pyf adder.f90 --backend meson