Skip to content

Commit ecababf

Browse files
add const char TRANS to remaining fftw transforms
1 parent 081adf2 commit ecababf

File tree

6 files changed

+383
-191
lines changed

6 files changed

+383
-191
lines changed

examples/calculus.c

Lines changed: 3 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -84,7 +84,7 @@ int main(void) {
8484
printmat("F", FMT, F, N, M);
8585
printf("\n");
8686

87-
ft_execute_tri_analysis(PA, F, N, M);
87+
ft_execute_tri_analysis('N', PA, F, N, M);
8888
ft_execute_cheb2tri('N', P, F, N, M);
8989

9090
printf("Its Proriol-"MAGENTA("(α,β,γ)")" coefficients are:\n\n");
@@ -114,7 +114,7 @@ int main(void) {
114114
printf("to "MAGENTA("y")", the analogous formulae result in a Proriol-"MAGENTA("(α,β+1,γ+1)")" series.\n");
115115
printf("These expressions are adapted from "CYAN("https://doi.org/10.1137/19M1245888")"\n");
116116

117-
ft_execute_tri_analysis(PA, Fx, N, M);
117+
ft_execute_tri_analysis('N', PA, Fx, N, M);
118118
ft_execute_cheb2tri('N', Px, Fx, N, M);
119119

120120
printmat("Ux from sampling", FMT, Fx, N, M);
@@ -137,7 +137,7 @@ int main(void) {
137137
printmat("Ux from U", FMT, Gx, N, M);
138138
printf("\n");
139139

140-
ft_execute_tri_analysis(PA, Fy, N, M);
140+
ft_execute_tri_analysis('N', PA, Fy, N, M);
141141
ft_execute_cheb2tri('N', Py, Fy, N, M);
142142

143143
printmat("Uy from sampling", FMT, Fy, N, M);

examples/holomorphic.c

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -88,7 +88,7 @@ int main(void) {
8888
ft_harmonic_plan * P = ft_plan_disk2cxf(N, alpha, beta);
8989
ft_disk_fftw_plan * PA = ft_plan_disk_analysis(N, M);
9090

91-
ft_execute_disk_analysis(PA, F, N, M);
91+
ft_execute_disk_analysis('N', PA, F, N, M);
9292
ft_execute_cxf2disk('N', P, F, N, M);
9393

9494
printf("Its Zernike coefficients are:\n\n");
@@ -144,7 +144,7 @@ int main(void) {
144144
P = ft_plan_rectdisk2cheb(N, beta);
145145
ft_rectdisk_fftw_plan * QA = ft_plan_rectdisk_analysis(N, M);
146146

147-
ft_execute_rectdisk_analysis(QA, F, N, M);
147+
ft_execute_rectdisk_analysis('N', QA, F, N, M);
148148
ft_execute_cheb2rectdisk('N', P, F, N, M);
149149

150150
printf("Its Dunkl-Xu coefficients are:\n\n");

examples/spinweighted.c

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -70,7 +70,7 @@ int main(void) {
7070
ft_spin_harmonic_plan * P = ft_plan_spinsph2fourier(N, 0);
7171
ft_spinsphere_fftw_plan * PA = ft_plan_spinsph_analysis(N, M, 0);
7272

73-
ft_execute_spinsph_analysis(PA, F, N, M);
73+
ft_execute_spinsph_analysis('N', PA, F, N, M);
7474
ft_execute_fourier2spinsph('N', P, F, N, M);
7575

7676
printf("Its spin-0 spherical harmonic coefficients are:\n\n");
@@ -108,7 +108,7 @@ int main(void) {
108108
P = ft_plan_spinsph2fourier(N, 1);
109109
PA = ft_plan_spinsph_analysis(N, M, 1);
110110

111-
ft_execute_spinsph_analysis(PA, F, N, M);
111+
ft_execute_spinsph_analysis('N', PA, F, N, M);
112112
ft_execute_fourier2spinsph('N', P, F, N, M);
113113

114114
printmat("U¹sampling", FMT, (double *) F, 2*N, M);

src/fasttransforms.h

Lines changed: 12 additions & 10 deletions
Original file line numberDiff line numberDiff line change
@@ -453,9 +453,9 @@ ft_triangle_fftw_plan * ft_plan_tri_synthesis(const int N, const int M);
453453
ft_triangle_fftw_plan * ft_plan_tri_analysis(const int N, const int M);
454454

455455
/// Execute FFTW synthesis on the triangle.
456-
void ft_execute_tri_synthesis(const ft_triangle_fftw_plan * P, double * X, const int N, const int M);
456+
void ft_execute_tri_synthesis(const char TRANS, const ft_triangle_fftw_plan * P, double * X, const int N, const int M);
457457
/// Execute FFTW analysis on the triangle.
458-
void ft_execute_tri_analysis(const ft_triangle_fftw_plan * P, double * X, const int N, const int M);
458+
void ft_execute_tri_analysis(const char TRANS, const ft_triangle_fftw_plan * P, double * X, const int N, const int M);
459459

460460
typedef struct {
461461
fftw_plan planxyz;
@@ -467,8 +467,8 @@ ft_tetrahedron_fftw_plan * ft_plan_tet_with_kind(const int N, const int L, const
467467
ft_tetrahedron_fftw_plan * ft_plan_tet_synthesis(const int N, const int L, const int M);
468468
ft_tetrahedron_fftw_plan * ft_plan_tet_analysis(const int N, const int L, const int M);
469469

470-
void ft_execute_tet_synthesis(const ft_tetrahedron_fftw_plan * P, double * X, const int N, const int L, const int M);
471-
void ft_execute_tet_analysis(const ft_tetrahedron_fftw_plan * P, double * X, const int N, const int L, const int M);
470+
void ft_execute_tet_synthesis(const char TRANS, const ft_tetrahedron_fftw_plan * P, double * X, const int N, const int L, const int M);
471+
void ft_execute_tet_analysis(const char TRANS, const ft_tetrahedron_fftw_plan * P, double * X, const int N, const int L, const int M);
472472

473473
typedef struct {
474474
fftw_plan planr1;
@@ -489,9 +489,9 @@ ft_disk_fftw_plan * ft_plan_disk_synthesis(const int N, const int M);
489489
ft_disk_fftw_plan * ft_plan_disk_analysis(const int N, const int M);
490490

491491
/// Execute FFTW synthesis on the disk.
492-
void ft_execute_disk_synthesis(const ft_disk_fftw_plan * P, double * X, const int N, const int M);
492+
void ft_execute_disk_synthesis(const char TRANS, const ft_disk_fftw_plan * P, double * X, const int N, const int M);
493493
/// Execute FFTW analysis on the disk.
494-
void ft_execute_disk_analysis(const ft_disk_fftw_plan * P, double * X, const int N, const int M);
494+
void ft_execute_disk_analysis(const char TRANS, const ft_disk_fftw_plan * P, double * X, const int N, const int M);
495495

496496
typedef struct {
497497
fftw_plan planx1;
@@ -509,9 +509,9 @@ ft_rectdisk_fftw_plan * ft_plan_rectdisk_synthesis(const int N, const int M);
509509
ft_rectdisk_fftw_plan * ft_plan_rectdisk_analysis(const int N, const int M);
510510

511511
/// Execute FFTW synthesis on the rectangularized disk.
512-
void ft_execute_rectdisk_synthesis(const ft_rectdisk_fftw_plan * P, double * X, const int N, const int M);
512+
void ft_execute_rectdisk_synthesis(const char TRANS, const ft_rectdisk_fftw_plan * P, double * X, const int N, const int M);
513513
/// Execute FFTW analysis on the rectangularized disk.
514-
void ft_execute_rectdisk_analysis(const ft_rectdisk_fftw_plan * P, double * X, const int N, const int M);
514+
void ft_execute_rectdisk_analysis(const char TRANS, const ft_rectdisk_fftw_plan * P, double * X, const int N, const int M);
515515

516516
typedef struct {
517517
fftw_plan plantheta1;
@@ -526,16 +526,18 @@ typedef struct {
526526
/// Destroy a \ref ft_spinsphere_fftw_plan.
527527
void ft_destroy_spinsphere_fftw_plan(ft_spinsphere_fftw_plan * P);
528528

529+
int ft_get_spin_spinsphere_fftw_plan(const ft_spinsphere_fftw_plan * P);
530+
529531
ft_spinsphere_fftw_plan * ft_plan_spinsph_with_kind(const int N, const int M, const int S, const fftw_r2r_kind kind[2][1], const int sign);
530532
/// Plan FFTW synthesis on the sphere with spin.
531533
ft_spinsphere_fftw_plan * ft_plan_spinsph_synthesis(const int N, const int M, const int S);
532534
/// Plan FFTW analysis on the sphere with spin.
533535
ft_spinsphere_fftw_plan * ft_plan_spinsph_analysis(const int N, const int M, const int S);
534536

535537
/// Execute FFTW synthesis on the sphere with spin.
536-
void ft_execute_spinsph_synthesis(const ft_spinsphere_fftw_plan * P, ft_complex * X, const int N, const int M);
538+
void ft_execute_spinsph_synthesis(const char TRANS, const ft_spinsphere_fftw_plan * P, ft_complex * X, const int N, const int M);
537539
/// Execute FFTW analysis on the sphere with spin.
538-
void ft_execute_spinsph_analysis(const ft_spinsphere_fftw_plan * P, ft_complex * X, const int N, const int M);
540+
void ft_execute_spinsph_analysis(const char TRANS, const ft_spinsphere_fftw_plan * P, ft_complex * X, const int N, const int M);
539541

540542
typedef struct {
541543
ft_banded ** B;

0 commit comments

Comments
 (0)