kml_fft(f/h)_plan_many_r2r
建立howmany组数据序列n维R2R变换的plan。其中,单个FFT的数据序列不需要是连续的,可以以跨步的形式提供。
接口定义
C interface:
kml_fft_plan kml_fft_plan_many_r2r(int rank, const int *n, int howmany, double *in, const int *inembed, int istride, int idist, double *out, const int *onembed, int ostride, int odist, const kml_fft_r2r_kind *kind, unsigned flags);
kml_fftf_plan kml_fftf_plan_many_r2r(int rank, const int *n, int howmany, float *in, const int *inembed, int istride, int idist, float *out, const int *onembed, int ostride, int odist, const kml_fftf_r2r_kind *kind, unsigned flags);
kml_ffth_plan kml_ffth_plan_many_r2r(int rank, const int *n, int howmany, __fp16 *in, const int *inembed, int istride, int idist, __fp16 *out, const int *onembed, int ostride, int odist, const kml_ffth_r2r_kind *kind, unsigned flags);
Fortran interface:
RES = KML_FFT_PLAN_MANY_DFT_R2R(RANK, N, HOWMANY, IN, INEMBED, ISTRIDE, IDIST, OUT, ONEMBED, OSTRIDE, ODIST, KIND, FLAGS);
RES = KML_FFTF_PLAN_MANY_DFT_R2R(RANK, N, HOWMANY, IN, INEMBED, ISTRIDE, IDIST, OUT, ONEMBED, OSTRIDE, ODIST, KIND, FLAGS);
RES = KML_FFTH_PLAN_MANY_DFT_R2R(RANK, N, HOWMANY, IN, INEMBED, ISTRIDE, IDIST, OUT, ONEMBED, OSTRIDE, ODIST, KIND, FLAGS);
KML_FFT_REDFT11和KML_FFT_RODFT11只支持长度为4的整数倍大小的序列,其他变换类型只支持长度为2的整数倍大小的序列。
返回值
函数返回一个kml_fft(f/h)_plan类型的结构体指针。将该对象作为参数传入kml_fft(f/h)_execute函数中使用,将对当前提供的输入in和输出out执行FFT变换;另外,也可以通过将该对象作为参数传入kml_fft(f/h)_execute_r2r函数中以对新的输入in和输出out执行FFT变换。
如果函数返回非空指针,则表示plan执行成功,否则表示执行失败。
参数
| 参数名 | 数据类型 | 描述 | 输入/输出 | 
|---|---|---|---|
| rank | int | FFT变换的维度是rank,约束:1 ≤ rank ≤ 3。 | 输入 | 
| n | const int* | n是维度为rank的数组,包含FFT序列每一维度的大小,约束:n[i] ≥ 1, for i in 0 to rank - 1。 | 输入 | 
| howmany | int | howmany表示要多少个多维FFT变换。 | 输入 | 
| in | 
 | 输入待变换的数据。 | 输入 | 
| inembed | const int* | inembed是大小为rank的数组或者NULL,该数组表示FFT数据in存储的更大空间的每一维度的大小。 约束:inembed[i] ≥ n[i] for i in 0, rank-1 或者:如果inembed == NULL,则可以认为inembed与n相等。 | 输入 | 
| istride | int | FFT输入序列的相继元素之间的间隔。 | 输入 | 
| idist | int | idist表示输入的每个FFT序列的间隔。 | 输入 | 
| out | 
 | 输出快速傅里叶变换后的数据。 | 输出 | 
| onembed | const int* | onembed是大小为rank的数组或者NULL,该数组表示FFT数据out存储的更大空间的每一维度的大小。 约束:onembed[i] ≥ n[i] for i in 0, rank-1 或者:如果onembed == NULL,则可以认为onembed与n相等。 | 输入 | 
| ostride | int | FFT输出序列的相继元素之间的间隔。 | 输入 | 
| odist | int | odist表示输出的每个FFT序列的间隔。 | 输入 | 
| kind | 
 | kind是大小为rank的数组,包含FFT序列每一维度的R2R变换类型,kind[i] (for i in 0 to rank - 1)有以下可选值: 
 | 输入 | 
| flags | unsigned int | planning选项,描述ESTIMATE模式或PATIENT模式。 KML_FFT_ESTIMATE:ESTIMATE模式 KML_FFT_PATIENT:PATIENT模式 | 输入 | 
依赖
C: "kfft.h"
示例
C interface:
    int rank = 1; 
    int *n; 
    n = (int*)kml_fft_malloc(sizeof(int) * rank); 
    n[0] = 2; 
    int howmany = 3; 
    int istride = 1; 
    int ostride = 1; 
    int idist = 2; 
    int odist = 2; 
    double init[6] = {120, 1, 8, 8, 120, 1}; 
    double *in; 
    in = (double*)kml_fft_malloc(sizeof(double) * 6); 
    for (int i = 0; i < 6; i++) { 
        in[i] = init[i]; 
    } 
    double *out; 
    out = (double*)kml_fft_malloc(sizeof(double) * 6); 
    kml_fft_r2r_kind *kind; 
    kind = (kml_fft_r2r_kind*)kml_fft_malloc(sizeof(kml_fft_r2r_kind) * rank); 
    kind[0] = KML_FFT_DHT;     
    kml_fft_plan plan; 
    plan = kml_fft_plan_many_r2r(rank, n, howmany, in, NULL, istride, idist, out, NULL, ostride, odist, kind, KML_FFT_ESTIMATE); 
    kml_fft_execute_r2r(plan, in, out); 
 
    kml_fft_destroy_plan(plan); 
    kml_fft_free(n); 
    kml_fft_free(kind); 
    kml_fft_free(in); 
    kml_fft_free(out); 
 
    /* 
     * out = {1.210000e+02, 1.190000e+02, 1.600000e+01, 0.000000e+00, 
     *        1.210000e+02, 1.190000e+02} 
     */
Fortran interface:
    INTEGER(C_INT) :: RANK = 1 
    INTEGER(C_INT) :: N(1) 
    INTEGER(C_INT) :: KIND(1) 
    REAL(C_DOUBLE), DIMENSION(6) :: INIT 
    TYPE(C_DOUBLE), POINTER :: IN(:) 
    TYPE(KML_FFT_COMPLEX), POINTER :: OUT(:) 
    INTEGER, POINTER :: INEMBED(:), ONEMBED(:) 
    TYPE(C_PTR) :: PIN, POUT 
    INTEGER(C_INT) :: HOWMANY = 3 
    INTEGER(C_INT) :: ISTRIDE = 1 
    INTEGER(C_INT) :: OSTRIDE = 1 
    INTEGER(C_INT) :: IDIST = 2 
    INTEGER(C_INT) :: ODIST = 2 
    INTEGER(C_SIZE_T) :: SIZE 
    SIZE = 6 
    PIN = KML_FFT_MALLOC(8 * SIZE) 
    POUT = KML_FFT_MALLOC(8 * SIZE) 
    CALL C_F_POINTER(PIN, IN, SHAPE=[6]) 
    CALL C_F_POINTER(POUT, OUT, SHAPE=[6]) 
    CALL C_F_POINTER(C_NULL_PTR, INEMBED, SHAPE=[0]) 
    CALL C_F_POINTER(C_NULL_PTR, ONEMBED, SHAPE=[0]) 
    N(1) = 2 
    DATA INIT/120, 1, 8, 8, 120, 1/ 
    INTEGER :: I 
    DO WHILE(I <= 6) 
        IN(I) = INIT(I) 
    END DO 
    KIND(1) = KML_FFT_DHT 
 
    TYPE(C_PTR) :: PLAN 
    PLAN = KML_FFT_PLAN_MANY_DFT_R2R(RANK, N, HOWMANY, IN, INEMBED, ISTRIDE, IDIST, OUT, ONEMBED, OSTRIDE, ODIST, KIND, KML_FFT_ESTIMATE); 
    CALL KML_FFT_EXECUTE_DFT_R2R(PLAN, IN, OUT); 
 
    CALL KML_FFT_DESTROY_PLAN(PLAN) 
    CALL KML_FFT_FREE(PIN) 
    CALL KML_FFT_FREE(POUT) 
    ! 
    ! OUT = /1.210000E+02, 1.190000E+02, 1.600000E+01, 0.000000E+00, 
    !        1.210000E+02, 1.190000E+02/ 
    !