FastTransforms
Loading...
Searching...
No Matches
fasttransforms.h
Go to the documentation of this file.
1
5#ifndef FASTTRANSFORMS_H
6#define FASTTRANSFORMS_H
7
8#include <cblas.h>
9#include <fftw3.h>
10
11typedef double ft_complex[2];
12
13#ifdef _OPENMP
14 #include <omp.h>
15 #define FT_GET_THREAD_NUM() omp_get_thread_num()
16 #define FT_GET_NUM_THREADS() omp_get_num_threads()
17 #define FT_GET_MAX_THREADS() omp_get_max_threads()
18 #define FT_SET_NUM_THREADS(x) omp_set_num_threads(x)
19#else
20 #define FT_GET_THREAD_NUM() 0
21 #define FT_GET_NUM_THREADS() 1
22 #define FT_GET_MAX_THREADS() 1
23 #define FT_SET_NUM_THREADS(x)
24#endif
25
26#define FT_CONCAT(prefix, name, suffix) prefix ## name ## suffix
27
28void ft_horner(const int n, const double * c, const int incc, const int m, double * x, double * f);
29void ft_hornerf(const int n, const float * c, const int incc, const int m, float * x, float * f);
30
31void ft_clenshaw(const int n, const double * c, const int incc, const int m, double * x, double * f);
32void ft_clenshawf(const int n, const float * c, const int incc, const int m, float * x, float * f);
33
34void ft_orthogonal_polynomial_clenshaw(const int n, const double * c, const int incc, const double * A, const double * B, const double * C, const int m, double * x, double * phi0, double * f);
35void ft_orthogonal_polynomial_clenshawf(const int n, const float * c, const int incc, const float * A, const float * B, const float * C, const int m, float * x, float * phi0, float * f);
36
37void ft_eigen_eval(const int n, const double * c, const int incc, const double * A, const double * B, const double * C, const int m, double * x, const int sign, double * f);
38void ft_eigen_evalf(const int n, const float * c, const int incc, const float * A, const float * B, const float * C, const int m, float * x, const int sign, float * f);
39void ft_eigen_evall(const int n, const long double * c, const int incc, const long double * A, const long double * B, const long double * C, const int m, long double * x, const int sign, long double * f);
40#if defined(FT_QUADMATH)
41 #include <quadmath.h>
42 typedef __float128 quadruple;
43 void ft_eigen_evalq(const int n, const quadruple * c, const int incc, const quadruple * A, const quadruple * B, const quadruple * C, const int m, quadruple * x, const int sign, quadruple * f);
44#endif
45
46#define FT_SN (1U << 0)
47#define FT_CN (1U << 1)
48#define FT_DN (1U << 2)
49
50#define FT_MODIFIED_NMAX (1U << 30)
51
52#include "tdc.h"
53
62ft_tb_eigen_FMM * ft_plan_legendre_to_chebyshev(const int normleg, const int normcheb, const int n);
71ft_tb_eigen_FMM * ft_plan_chebyshev_to_legendre(const int normcheb, const int normleg, const int n);
80ft_tb_eigen_FMM * ft_plan_ultraspherical_to_ultraspherical(const int norm1, const int norm2, const int n, const double lambda, const double mu);
89ft_tb_eigen_FMM * ft_plan_jacobi_to_jacobi(const int norm1, const int norm2, const int n, const double alpha, const double beta, const double gamma, const double delta);
98ft_tb_eigen_FMM * ft_plan_laguerre_to_laguerre(const int norm1, const int norm2, const int n, const double alpha, const double beta);
107ft_tb_eigen_FMM * ft_plan_jacobi_to_ultraspherical(const int normjac, const int normultra, const int n, const double alpha, const double beta, const double lambda);
116ft_tb_eigen_FMM * ft_plan_ultraspherical_to_jacobi(const int normultra, const int normjac, const int n, const double lambda, const double alpha, const double beta);
125ft_tb_eigen_FMM * ft_plan_jacobi_to_chebyshev(const int normjac, const int normcheb, const int n, const double alpha, const double beta);
134ft_tb_eigen_FMM * ft_plan_chebyshev_to_jacobi(const int normcheb, const int normjac, const int n, const double alpha, const double beta);
143ft_tb_eigen_FMM * ft_plan_ultraspherical_to_chebyshev(const int normultra, const int normcheb, const int n, const double lambda);
152ft_tb_eigen_FMM * ft_plan_chebyshev_to_ultraspherical(const int normcheb, const int normultra, const int n, const double lambda);
153
154ft_btb_eigen_FMM * ft_plan_associated_jacobi_to_jacobi(const int norm1, const int norm2, const int n, const int c, const double alpha, const double beta, const double gamma, const double delta);
155ft_btb_eigen_FMM * ft_plan_associated_laguerre_to_laguerre(const int norm1, const int norm2, const int n, const int c, const double alpha, const double beta);
156ft_btb_eigen_FMM * ft_plan_associated_hermite_to_hermite(const int norm1, const int norm2, const int n, const int c);
157
165ft_modified_plan * ft_plan_modified_jacobi_to_jacobi(const int n, const double alpha, const double beta, const int nu, const double * u, const int nv, const double * v, const int verbose);
173ft_modified_plan * ft_plan_modified_laguerre_to_laguerre(const int n, const double alpha, const int nu, const double * u, const int nv, const double * v, const int verbose);
181ft_modified_plan * ft_plan_modified_hermite_to_hermite(const int n, const int nu, const double * u, const int nv, const double * v, const int verbose);
182
184ft_tb_eigen_FMMf * ft_plan_legendre_to_chebyshevf(const int normleg, const int normcheb, const int n);
186ft_tb_eigen_FMMf * ft_plan_chebyshev_to_legendref(const int normcheb, const int normleg, const int n);
188ft_tb_eigen_FMMf * ft_plan_ultraspherical_to_ultrasphericalf(const int norm1, const int norm2, const int n, const float lambda, const float mu);
190ft_tb_eigen_FMMf * ft_plan_jacobi_to_jacobif(const int norm1, const int norm2, const int n, const float alpha, const float beta, const float gamma, const float delta);
192ft_tb_eigen_FMMf * ft_plan_laguerre_to_laguerref(const int norm1, const int norm2, const int n, const float alpha, const float beta);
194ft_tb_eigen_FMMf * ft_plan_jacobi_to_ultrasphericalf(const int normjac, const int normultra, const int n, const float alpha, const float beta, const float lambda);
196ft_tb_eigen_FMMf * ft_plan_ultraspherical_to_jacobif(const int normultra, const int normjac, const int n, const float lambda, const float alpha, const float beta);
198ft_tb_eigen_FMMf * ft_plan_jacobi_to_chebyshevf(const int normjac, const int normcheb, const int n, const float alpha, const float beta);
200ft_tb_eigen_FMMf * ft_plan_chebyshev_to_jacobif(const int normcheb, const int normjac, const int n, const float alpha, const float beta);
202ft_tb_eigen_FMMf * ft_plan_ultraspherical_to_chebyshevf(const int normultra, const int normcheb, const int n, const float lambda);
204ft_tb_eigen_FMMf * ft_plan_chebyshev_to_ultrasphericalf(const int normcheb, const int normultra, const int n, const float lambda);
205
206ft_btb_eigen_FMMf * ft_plan_associated_jacobi_to_jacobif(const int norm1, const int norm2, const int n, const int c, const float alpha, const float beta, const float gamma, const float delta);
207ft_btb_eigen_FMMf * ft_plan_associated_laguerre_to_laguerref(const int norm1, const int norm2, const int n, const int c, const float alpha, const float beta);
208ft_btb_eigen_FMMf * ft_plan_associated_hermite_to_hermitef(const int norm1, const int norm2, const int n, const int c);
209
211ft_modified_planf * ft_plan_modified_jacobi_to_jacobif(const int n, const float alpha, const float beta, const int nu, const float * u, const int nv, const float * v, const int verbose);
213ft_modified_planf * ft_plan_modified_laguerre_to_laguerref(const int n, const float alpha, const int nu, const float * u, const int nv, const float * v, const int verbose);
215ft_modified_planf * ft_plan_modified_hermite_to_hermitef(const int n, const int nu, const float * u, const int nv, const float * v, const int verbose);
216
218ft_tb_eigen_FMMl * ft_plan_legendre_to_chebyshevl(const int normleg, const int normcheb, const int n);
220ft_tb_eigen_FMMl * ft_plan_chebyshev_to_legendrel(const int normcheb, const int normleg, const int n);
222ft_tb_eigen_FMMl * ft_plan_ultraspherical_to_ultrasphericall(const int norm1, const int norm2, const int n, const long double lambda, const long double mu);
224ft_tb_eigen_FMMl * ft_plan_jacobi_to_jacobil(const int norm1, const int norm2, const int n, const long double alpha, const long double beta, const long double gamma, const long double delta);
226ft_tb_eigen_FMMl * ft_plan_laguerre_to_laguerrel(const int norm1, const int norm2, const int n, const long double alpha, const long double beta);
228ft_tb_eigen_FMMl * ft_plan_jacobi_to_ultrasphericall(const int normjac, const int normultra, const int n, const long double alpha, const long double beta, const long double lambda);
230ft_tb_eigen_FMMl * ft_plan_ultraspherical_to_jacobil(const int normultra, const int normjac, const int n, const long double lambda, const long double alpha, const long double beta);
232ft_tb_eigen_FMMl * ft_plan_jacobi_to_chebyshevl(const int normjac, const int normcheb, const int n, const long double alpha, const long double beta);
234ft_tb_eigen_FMMl * ft_plan_chebyshev_to_jacobil(const int normcheb, const int normjac, const int n, const long double alpha, const long double beta);
236ft_tb_eigen_FMMl * ft_plan_ultraspherical_to_chebyshevl(const int normultra, const int normcheb, const int n, const long double lambda);
238ft_tb_eigen_FMMl * ft_plan_chebyshev_to_ultrasphericall(const int normcheb, const int normultra, const int n, const long double lambda);
239
240ft_btb_eigen_FMMl * ft_plan_associated_jacobi_to_jacobil(const int norm1, const int norm2, const int n, const int c, const long double alpha, const long double beta, const long double gamma, const long double delta);
241ft_btb_eigen_FMMl * ft_plan_associated_laguerre_to_laguerrel(const int norm1, const int norm2, const int n, const int c, const long double alpha, const long double beta);
242ft_btb_eigen_FMMl * ft_plan_associated_hermite_to_hermitel(const int norm1, const int norm2, const int n, const int c);
243
245ft_modified_planl * ft_plan_modified_jacobi_to_jacobil(const int n, const long double alpha, const long double beta, const int nu, const long double * u, const int nv, const long double * v, const int verbose);
247ft_modified_planl * ft_plan_modified_laguerre_to_laguerrel(const int n, const long double alpha, const int nu, const long double * u, const int nv, const long double * v, const int verbose);
249ft_modified_planl * ft_plan_modified_hermite_to_hermitel(const int n, const int nu, const long double * u, const int nv, const long double * v, const int verbose);
250
251#include <mpfr.h>
252
253typedef struct {
254 mpfr_t * data;
255 int n;
256 int b;
258
259void ft_mpfr_destroy_plan(mpfr_t * A, int n);
260void ft_mpfr_trmv(char TRANS, int n, mpfr_t * A, int LDA, mpfr_t * x, mpfr_rnd_t rnd);
261void ft_mpfr_trsv(char TRANS, int n, mpfr_t * A, int LDA, mpfr_t * x, mpfr_rnd_t rnd);
262void ft_mpfr_trmm(char TRANS, int n, mpfr_t * A, int LDA, mpfr_t * B, int LDB, int N, mpfr_rnd_t rnd);
263void ft_mpfr_trsm(char TRANS, int n, mpfr_t * A, int LDA, mpfr_t * B, int LDB, int N, mpfr_rnd_t rnd);
264// C -- Julia interoperability. Julia `BigFloat` does not have the same size as `mpfr_t`.
265// So we give all Julia-owned data its own address, and dereference to retrieve the number.
266void ft_mpfr_trmv_ptr(char TRANS, int n, mpfr_t * A, int LDA, mpfr_t ** x, mpfr_rnd_t rnd);
267void ft_mpfr_trsv_ptr(char TRANS, int n, mpfr_t * A, int LDA, mpfr_t ** x, mpfr_rnd_t rnd);
268void ft_mpfr_trmm_ptr(char TRANS, int n, mpfr_t * A, int LDA, mpfr_t ** B, int LDB, int N, mpfr_rnd_t rnd);
269void ft_mpfr_trsm_ptr(char TRANS, int n, mpfr_t * A, int LDA, mpfr_t ** B, int LDB, int N, mpfr_rnd_t rnd);
270
272mpfr_t * ft_mpfr_plan_legendre_to_chebyshev(const int normleg, const int normcheb, const int n, mpfr_prec_t prec, mpfr_rnd_t rnd);
274mpfr_t * ft_mpfr_plan_chebyshev_to_legendre(const int normcheb, const int normleg, const int n, mpfr_prec_t prec, mpfr_rnd_t rnd);
276mpfr_t * ft_mpfr_plan_ultraspherical_to_ultraspherical(const int norm1, const int norm2, const int n, mpfr_srcptr lambda, mpfr_srcptr mu, mpfr_prec_t prec, mpfr_rnd_t rnd);
278mpfr_t * ft_mpfr_plan_jacobi_to_jacobi(const int norm1, const int norm2, const int n, mpfr_srcptr alpha, mpfr_srcptr beta, mpfr_srcptr gamma, mpfr_srcptr delta, mpfr_prec_t prec, mpfr_rnd_t rnd);
280mpfr_t * ft_mpfr_plan_laguerre_to_laguerre(const int norm1, const int norm2, const int n, mpfr_srcptr alpha, mpfr_srcptr beta, mpfr_prec_t prec, mpfr_rnd_t rnd);
282mpfr_t * ft_mpfr_plan_jacobi_to_ultraspherical(const int normjac, const int normultra, const int n, mpfr_srcptr alpha, mpfr_srcptr beta, mpfr_srcptr lambda, mpfr_prec_t prec, mpfr_rnd_t rnd);
284mpfr_t * ft_mpfr_plan_ultraspherical_to_jacobi(const int normultra, const int normjac, const int n, mpfr_srcptr lambda, mpfr_srcptr alpha, mpfr_srcptr beta, mpfr_prec_t prec, mpfr_rnd_t rnd);
286mpfr_t * ft_mpfr_plan_jacobi_to_chebyshev(const int normjac, const int normcheb, const int n, mpfr_srcptr alpha, mpfr_srcptr beta, mpfr_prec_t prec, mpfr_rnd_t rnd);
288mpfr_t * ft_mpfr_plan_chebyshev_to_jacobi(const int normcheb, const int normjac, const int n, mpfr_srcptr alpha, mpfr_srcptr beta, mpfr_prec_t prec, mpfr_rnd_t rnd);
290mpfr_t * ft_mpfr_plan_ultraspherical_to_chebyshev(const int normultra, const int normcheb, const int n, mpfr_srcptr lambda, mpfr_prec_t prec, mpfr_rnd_t rnd);
292mpfr_t * ft_mpfr_plan_chebyshev_to_ultraspherical(const int normcheb, const int normultra, const int n, mpfr_srcptr lambda, mpfr_prec_t prec, mpfr_rnd_t rnd);
293
295void ft_set_num_threads(const int n);
296
298typedef struct ft_rotation_plan_s ft_rotation_plan;
299
302
303ft_rotation_plan * ft_plan_rotsphere(const int n);
304
306void ft_kernel_sph_hi2lo(const ft_rotation_plan * RP, const int m1, const int m2, double * A, const int S);
308void ft_kernel_sph_lo2hi(const ft_rotation_plan * RP, const int m1, const int m2, double * A, const int S);
309
310ft_rotation_plan * ft_plan_rottriangle(const int n, const double alpha, const double beta, const double gamma);
311
313void ft_kernel_tri_hi2lo(const ft_rotation_plan * RP, const int m1, const int m2, double * A, const int S);
315void ft_kernel_tri_lo2hi(const ft_rotation_plan * RP, const int m1, const int m2, double * A, const int S);
316
317ft_rotation_plan * ft_plan_rotdisk(const int n, const double alpha, const double beta);
318ft_rotation_plan * ft_plan_rotannulus(const int n, const double alpha, const double beta, const double gamma, const double rho);
319
321void ft_kernel_disk_hi2lo(const ft_rotation_plan * RP, const int m1, const int m2, double * A, const int S);
323void ft_kernel_disk_lo2hi(const ft_rotation_plan * RP, const int m1, const int m2, double * A, const int S);
324
325ft_rotation_plan * ft_plan_rotrectdisk(const int n, const double beta);
326
328void ft_kernel_rectdisk_hi2lo(const ft_rotation_plan * RP, const int m1, const int m2, double * A, const int S);
330void ft_kernel_rectdisk_lo2hi(const ft_rotation_plan * RP, const int m1, const int m2, double * A, const int S);
331
332void ft_kernel_tet_hi2lo(const ft_rotation_plan * RP, const int L, const int m, double * A);
333void ft_kernel_tet_lo2hi(const ft_rotation_plan * RP, const int L, const int m, double * A);
334
336typedef struct ft_spin_rotation_plan_s ft_spin_rotation_plan;
337
340
341ft_spin_rotation_plan * ft_plan_rotspinsphere(const int n, const int s);
342
344void ft_kernel_spinsph_hi2lo(const ft_spin_rotation_plan * SRP, const int m, ft_complex * A, const int S);
346void ft_kernel_spinsph_lo2hi(const ft_spin_rotation_plan * SRP, const int m, ft_complex * A, const int S);
347
348
349void ft_execute_sph_hi2lo(const ft_rotation_plan * RP, double * A, double * B, const int M);
350void ft_execute_sph_lo2hi(const ft_rotation_plan * RP, double * A, double * B, const int M);
351
352void ft_execute_sphv_hi2lo(const ft_rotation_plan * RP, double * A, double * B, const int M);
353void ft_execute_sphv_lo2hi(const ft_rotation_plan * RP, double * A, double * B, const int M);
354
355void ft_execute_tri_hi2lo(const ft_rotation_plan * RP, double * A, double * B, const int M);
356void ft_execute_tri_lo2hi(const ft_rotation_plan * RP, double * A, double * B, const int M);
357
358void ft_execute_disk_hi2lo(const ft_rotation_plan * RP, double * A, double * B, const int M);
359void ft_execute_disk_lo2hi(const ft_rotation_plan * RP, double * A, double * B, const int M);
360
361void ft_execute_rectdisk_hi2lo(const ft_rotation_plan * RP, double * A, double * B, const int M);
362void ft_execute_rectdisk_lo2hi(const ft_rotation_plan * RP, double * A, double * B, const int M);
363
364void ft_execute_tet_hi2lo(const ft_rotation_plan * RP1, const ft_rotation_plan * RP2, double * A, const int L, const int M);
365void ft_execute_tet_lo2hi(const ft_rotation_plan * RP1, const ft_rotation_plan * RP2, double * A, const int L, const int M);
366
367void ft_execute_spinsph_hi2lo(const ft_spin_rotation_plan * SRP, ft_complex * A, ft_complex * B, const int M);
368void ft_execute_spinsph_lo2hi(const ft_spin_rotation_plan * SRP, ft_complex * A, ft_complex * B, const int M);
369
371typedef struct {
372 ft_rotation_plan ** RP;
373 ft_modified_plan ** MP;
374 double * B;
375 double ** P;
376 double ** Pinv;
377 double alpha;
378 double beta;
379 double gamma;
380 double delta;
381 double rho;
382 int NRP;
383 int NMP;
384 int NP;
386
389
392
394void ft_execute_sph2fourier(const char TRANS, const ft_harmonic_plan * P, double * A, const int N, const int M);
396void ft_execute_fourier2sph(const char TRANS, const ft_harmonic_plan * P, double * A, const int N, const int M);
397
398void ft_execute_sphv2fourier(const char TRANS, const ft_harmonic_plan * P, double * A, const int N, const int M);
399void ft_execute_fourier2sphv(const char TRANS, const ft_harmonic_plan * P, double * A, const int N, const int M);
400
402ft_harmonic_plan * ft_plan_tri2cheb(const int n, const double alpha, const double beta, const double gamma);
403
405void ft_execute_tri2cheb(const char TRANS, const ft_harmonic_plan * P, double * A, const int N, const int M);
407void ft_execute_cheb2tri(const char TRANS, const ft_harmonic_plan * P, double * A, const int N, const int M);
408
410ft_harmonic_plan * ft_plan_disk2cxf(const int n, const double alpha, const double beta);
411
413void ft_execute_disk2cxf(const char TRANS, const ft_harmonic_plan * P, double * A, const int N, const int M);
415void ft_execute_cxf2disk(const char TRANS, const ft_harmonic_plan * P, double * A, const int N, const int M);
416
418ft_harmonic_plan * ft_plan_ann2cxf(const int n, const double alpha, const double beta, const double gamma, const double rho);
419
421void ft_execute_ann2cxf(const char TRANS, const ft_harmonic_plan * P, double * A, const int N, const int M);
423void ft_execute_cxf2ann(const char TRANS, const ft_harmonic_plan * P, double * A, const int N, const int M);
424
426ft_harmonic_plan * ft_plan_rectdisk2cheb(const int n, const double beta);
427
429void ft_execute_rectdisk2cheb(const char TRANS, const ft_harmonic_plan * P, double * A, const int N, const int M);
431void ft_execute_cheb2rectdisk(const char TRANS, const ft_harmonic_plan * P, double * A, const int N, const int M);
432
434ft_harmonic_plan * ft_plan_tet2cheb(const int n, const double alpha, const double beta, const double gamma, const double delta);
435
437void ft_execute_tet2cheb(const char TRANS, const ft_harmonic_plan * P, double * A, const int N, const int L, const int M);
439void ft_execute_cheb2tet(const char TRANS, const ft_harmonic_plan * P, double * A, const int N, const int L, const int M);
440
442typedef struct {
444 ft_complex * B;
445 ft_complex * P1;
446 ft_complex * P2;
447 ft_complex * P1inv;
448 ft_complex * P2inv;
449 int s;
451
454
457
459void ft_execute_spinsph2fourier(const char TRANS, const ft_spin_harmonic_plan * P, ft_complex * A, const int N, const int M);
461void ft_execute_fourier2spinsph(const char TRANS, const ft_spin_harmonic_plan * P, ft_complex * A, const int N, const int M);
462
463
464int ft_fftw_init_threads(void);
465void ft_fftw_plan_with_nthreads(const int n);
466
467typedef struct {
468 fftw_plan plantheta1;
469 fftw_plan plantheta2;
470 fftw_plan plantheta3;
471 fftw_plan plantheta4;
472 fftw_plan planphi;
473 double * Y;
475
478
479ft_sphere_fftw_plan * ft_plan_sph_with_kind(const int N, const int M, const fftw_r2r_kind kind[3][1], const unsigned flags);
481ft_sphere_fftw_plan * ft_plan_sph_synthesis(const int N, const int M, const unsigned flags);
483ft_sphere_fftw_plan * ft_plan_sph_analysis(const int N, const int M, const unsigned flags);
484ft_sphere_fftw_plan * ft_plan_sphv_synthesis(const int N, const int M, const unsigned flags);
485ft_sphere_fftw_plan * ft_plan_sphv_analysis(const int N, const int M, const unsigned flags);
486
488void ft_execute_sph_synthesis(const char TRANS, const ft_sphere_fftw_plan * P, double * X, const int N, const int M);
490void ft_execute_sph_analysis(const char TRANS, const ft_sphere_fftw_plan * P, double * X, const int N, const int M);
491
492void ft_execute_sphv_synthesis(const char TRANS, const ft_sphere_fftw_plan * P, double * X, const int N, const int M);
493void ft_execute_sphv_analysis(const char TRANS, const ft_sphere_fftw_plan * P, double * X, const int N, const int M);
494
495typedef struct {
496 fftw_plan planxy;
498
501
502ft_triangle_fftw_plan * ft_plan_tri_with_kind(const int N, const int M, const fftw_r2r_kind kind0, const fftw_r2r_kind kind1, const unsigned flags);
504ft_triangle_fftw_plan * ft_plan_tri_synthesis(const int N, const int M, const unsigned flags);
506ft_triangle_fftw_plan * ft_plan_tri_analysis(const int N, const int M, const unsigned flags);
507
509void ft_execute_tri_synthesis(const char TRANS, const ft_triangle_fftw_plan * P, double * X, const int N, const int M);
511void ft_execute_tri_analysis(const char TRANS, const ft_triangle_fftw_plan * P, double * X, const int N, const int M);
512
513typedef struct {
514 fftw_plan planxyz;
516
517void ft_destroy_tetrahedron_fftw_plan(ft_tetrahedron_fftw_plan * P);
518
519ft_tetrahedron_fftw_plan * ft_plan_tet_with_kind(const int N, const int L, const int M, const fftw_r2r_kind kind0, const fftw_r2r_kind kind1, const fftw_r2r_kind kind2, const unsigned flags);
520ft_tetrahedron_fftw_plan * ft_plan_tet_synthesis(const int N, const int L, const int M, const unsigned flags);
521ft_tetrahedron_fftw_plan * ft_plan_tet_analysis(const int N, const int L, const int M, const unsigned flags);
522
523void ft_execute_tet_synthesis(const char TRANS, const ft_tetrahedron_fftw_plan * P, double * X, const int N, const int L, const int M);
524void ft_execute_tet_analysis(const char TRANS, const ft_tetrahedron_fftw_plan * P, double * X, const int N, const int L, const int M);
525
526typedef struct {
527 fftw_plan planr1;
528 fftw_plan planr2;
529 fftw_plan planr3;
530 fftw_plan planr4;
531 fftw_plan plantheta;
532 double * Y;
534
537
538ft_disk_fftw_plan * ft_plan_disk_with_kind(const int N, const int M, const fftw_r2r_kind kind[3][1], const unsigned flags);
540ft_disk_fftw_plan * ft_plan_disk_synthesis(const int N, const int M, const unsigned flags);
542ft_disk_fftw_plan * ft_plan_disk_analysis(const int N, const int M, const unsigned flags);
543
545void ft_execute_disk_synthesis(const char TRANS, const ft_disk_fftw_plan * P, double * X, const int N, const int M);
547void ft_execute_disk_analysis(const char TRANS, const ft_disk_fftw_plan * P, double * X, const int N, const int M);
548
549typedef struct {
550 fftw_plan planr;
551 fftw_plan plantheta;
552 double rho;
553 double * w;
554 double * Y;
556
559
560double ft_get_rho_annulus_fftw_plan(ft_annulus_fftw_plan * P);
561
562ft_annulus_fftw_plan * ft_plan_annulus_with_kind(const int N, const int M, const double rho, const fftw_r2r_kind kind[3][1], const unsigned flags);
564ft_annulus_fftw_plan * ft_plan_annulus_synthesis(const int N, const int M, const double rho, const unsigned flags);
566ft_annulus_fftw_plan * ft_plan_annulus_analysis(const int N, const int M, const double rho, const unsigned flags);
567
569void ft_execute_annulus_synthesis(const char TRANS, const ft_annulus_fftw_plan * P, double * X, const int N, const int M);
571void ft_execute_annulus_analysis(const char TRANS, const ft_annulus_fftw_plan * P, double * X, const int N, const int M);
572
573typedef struct {
574 fftw_plan planx1;
575 fftw_plan planx2;
576 fftw_plan plany;
578
581
582ft_rectdisk_fftw_plan * ft_plan_rectdisk_with_kind(const int N, const int M, const fftw_r2r_kind kind[3][1], const unsigned flags);
584ft_rectdisk_fftw_plan * ft_plan_rectdisk_synthesis(const int N, const int M, const unsigned flags);
586ft_rectdisk_fftw_plan * ft_plan_rectdisk_analysis(const int N, const int M, const unsigned flags);
587
589void ft_execute_rectdisk_synthesis(const char TRANS, const ft_rectdisk_fftw_plan * P, double * X, const int N, const int M);
591void ft_execute_rectdisk_analysis(const char TRANS, const ft_rectdisk_fftw_plan * P, double * X, const int N, const int M);
592
593typedef struct {
594 fftw_plan plantheta1;
595 fftw_plan plantheta2;
596 fftw_plan plantheta3;
597 fftw_plan plantheta4;
598 fftw_plan planphi;
599 double * Y;
600 int S;
602
605
606int ft_get_spin_spinsphere_fftw_plan(const ft_spinsphere_fftw_plan * P);
607
608ft_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, const unsigned flags);
610ft_spinsphere_fftw_plan * ft_plan_spinsph_synthesis(const int N, const int M, const int S, const unsigned flags);
612ft_spinsphere_fftw_plan * ft_plan_spinsph_analysis(const int N, const int M, const int S, const unsigned flags);
613
615void ft_execute_spinsph_synthesis(const char TRANS, const ft_spinsphere_fftw_plan * P, ft_complex * X, const int N, const int M);
617void ft_execute_spinsph_analysis(const char TRANS, const ft_spinsphere_fftw_plan * P, ft_complex * X, const int N, const int M);
618
619typedef struct {
620 ft_banded ** B;
621 ft_triangular_banded ** T;
622 int n;
624
625void ft_destroy_gradient_plan(ft_gradient_plan * P);
626
627ft_gradient_plan * ft_plan_sph_gradient(const int n);
628
629void ft_execute_sph_gradient(ft_gradient_plan * P, double * U, double * Ut, double * Up, const int N, const int M);
630void ft_execute_sph_curl(ft_gradient_plan * P, double * U, double * Ut, double * Up, const int N, const int M);
631
632typedef struct {
633 ft_triangular_banded ** T;
634 ft_banded_qr ** F;
635 double * X;
636 int n;
638
639void ft_destroy_helmholtzhodge_plan(ft_helmholtzhodge_plan * P);
640
641ft_helmholtzhodge_plan * ft_plan_sph_helmholtzhodge(const int n);
642
643void ft_execute_sph_helmholtzhodge(ft_helmholtzhodge_plan * P, double * U1, double * U2, double * V1, double * V2, const int N, const int M);
644
649typedef struct {
650 double Q[9];
652
672typedef struct {
673 double s[3];
674 double c[3];
675 int sign;
676} ft_ZYZR;
677
681typedef struct {
682 double w[3];
684
685ft_ZYZR ft_create_ZYZR(ft_orthogonal_transformation Q);
686
687void ft_apply_ZYZR(ft_ZYZR Q, ft_orthogonal_transformation * U);
688
689void ft_apply_reflection(ft_reflection Q, ft_orthogonal_transformation * U);
690
691void ft_execute_sph_polar_rotation(double * A, const int N, const int M, double s, double c);
692
693void ft_execute_sph_polar_reflection(double * A, const int N, const int M);
694
695typedef struct {
696 ft_symmetric_tridiagonal_symmetric_eigen * F11;
697 ft_symmetric_tridiagonal_symmetric_eigen * F21;
698 ft_symmetric_tridiagonal_symmetric_eigen * F12;
699 ft_symmetric_tridiagonal_symmetric_eigen * F22;
700 int l;
702
703typedef struct {
705 int n;
707
708void ft_destroy_partial_sph_isometry_plan(ft_partial_sph_isometry_plan * F);
709
710void ft_destroy_sph_isometry_plan(ft_sph_isometry_plan * F);
711
712ft_partial_sph_isometry_plan * ft_plan_partial_sph_isometry(const int l);
713
714ft_sph_isometry_plan * ft_plan_sph_isometry(const int n);
715
716void ft_execute_sph_yz_axis_exchange(ft_sph_isometry_plan * J, double * A, const int N, const int M);
717
718void ft_execute_sph_isometry(ft_sph_isometry_plan * J, ft_ZYZR Q, double * A, const int N, const int M);
719
720void ft_execute_sph_rotation(ft_sph_isometry_plan * J, const double alpha, const double beta, const double gamma, double * A, const int N, const int M);
721
722void ft_execute_sph_reflection(ft_sph_isometry_plan * J, ft_reflection W, double * A, const int N, const int M);
723
724void ft_execute_sph_orthogonal_transformation(ft_sph_isometry_plan * J, ft_orthogonal_transformation Q, double * A, const int N, const int M);
725
726#endif // FASTTRANSFORMS_H
ft_disk_fftw_plan * ft_plan_disk_analysis(const int N, const int M, const unsigned flags)
Plan FFTW analysis on the disk.
ft_tb_eigen_FMMl * ft_plan_ultraspherical_to_ultrasphericall(const int norm1, const int norm2, const int n, const long double lambda, const long double mu)
A long double precision version of ft_plan_ultraspherical_to_ultraspherical.
void ft_execute_spinsph2fourier(const char TRANS, const ft_spin_harmonic_plan *P, ft_complex *A, const int N, const int M)
Transform a spin-weighted spherical harmonic expansion to a bivariate Fourier series.
ft_tb_eigen_FMMl * ft_plan_chebyshev_to_jacobil(const int normcheb, const int normjac, const int n, const long double alpha, const long double beta)
A long double precision version of ft_plan_chebyshev_to_jacobi.
ft_tb_eigen_FMMf * ft_plan_jacobi_to_jacobif(const int norm1, const int norm2, const int n, const float alpha, const float beta, const float gamma, const float delta)
A single precision version of ft_plan_jacobi_to_jacobi.
void ft_execute_cheb2rectdisk(const char TRANS, const ft_harmonic_plan *P, double *A, const int N, const int M)
Transform a Chebyshev series to a rectangularized disk harmonic expansion.
ft_tb_eigen_FMM * ft_plan_laguerre_to_laguerre(const int norm1, const int norm2, const int n, const double alpha, const double beta)
Pre-compute a factorization of the connection coefficients between Laguerre polynomials in double pre...
ft_sphere_fftw_plan * ft_plan_sph_synthesis(const int N, const int M, const unsigned flags)
Plan FFTW synthesis on the sphere.
ft_tb_eigen_FMM * ft_plan_jacobi_to_chebyshev(const int normjac, const int normcheb, const int n, const double alpha, const double beta)
Pre-compute a factorization of the connection coefficients between Jacobi and Chebyshev polynomials i...
void ft_execute_annulus_analysis(const char TRANS, const ft_annulus_fftw_plan *P, double *X, const int N, const int M)
Execute FFTW analysis on the annulus.
void ft_kernel_disk_hi2lo(const ft_rotation_plan *RP, const int m1, const int m2, double *A, const int S)
Convert a single vector of disk harmonic coefficients in A with stride S from order m2 down to order ...
ft_tb_eigen_FMMl * ft_plan_jacobi_to_chebyshevl(const int normjac, const int normcheb, const int n, const long double alpha, const long double beta)
A long double precision version of ft_plan_jacobi_to_chebyshev.
void ft_destroy_spinsphere_fftw_plan(ft_spinsphere_fftw_plan *P)
Destroy a ft_spinsphere_fftw_plan.
ft_modified_planf * ft_plan_modified_hermite_to_hermitef(const int n, const int nu, const float *u, const int nv, const float *v, const int verbose)
A single precision version of ft_plan_modified_hermite_to_hermite.
ft_modified_plan * ft_plan_modified_jacobi_to_jacobi(const int n, const double alpha, const double beta, const int nu, const double *u, const int nv, const double *v, const int verbose)
Pre-compute a factorization of the connection coefficients between modified Jacobi polynomials,...
void ft_execute_disk2cxf(const char TRANS, const ft_harmonic_plan *P, double *A, const int N, const int M)
Transform a disk harmonic expansion to a Chebyshev–Fourier series.
mpfr_t * ft_mpfr_plan_jacobi_to_chebyshev(const int normjac, const int normcheb, const int n, mpfr_srcptr alpha, mpfr_srcptr beta, mpfr_prec_t prec, mpfr_rnd_t rnd)
A multi-precision version of ft_plan_jacobi_to_chebyshev that returns a dense array of connection coe...
void ft_destroy_spin_rotation_plan(ft_spin_rotation_plan *SRP)
Destroy a ft_spin_rotation_plan.
ft_tb_eigen_FMMf * ft_plan_chebyshev_to_jacobif(const int normcheb, const int normjac, const int n, const float alpha, const float beta)
A single precision version of ft_plan_chebyshev_to_jacobi.
ft_rectdisk_fftw_plan * ft_plan_rectdisk_analysis(const int N, const int M, const unsigned flags)
Plan FFTW analysis on the rectangularized disk.
void ft_execute_rectdisk_analysis(const char TRANS, const ft_rectdisk_fftw_plan *P, double *X, const int N, const int M)
Execute FFTW analysis on the rectangularized disk.
void ft_execute_fourier2spinsph(const char TRANS, const ft_spin_harmonic_plan *P, ft_complex *A, const int N, const int M)
Transform a bivariate Fourier series to a spin-weighted spherical harmonic expansion.
ft_tb_eigen_FMMl * ft_plan_ultraspherical_to_jacobil(const int normultra, const int normjac, const int n, const long double lambda, const long double alpha, const long double beta)
A long double precision version of ft_plan_ultraspherical_to_jacobi.
ft_harmonic_plan * ft_plan_tri2cheb(const int n, const double alpha, const double beta, const double gamma)
Plan a triangular harmonic transform.
struct ft_spin_rotation_plan_s ft_spin_rotation_plan
Data structure to store sines and cosines of Givens rotations for spin-weighted harmonics.
Definition fasttransforms.h:336
void ft_execute_sph_analysis(const char TRANS, const ft_sphere_fftw_plan *P, double *X, const int N, const int M)
Execute FFTW analysis on the sphere.
ft_tb_eigen_FMMf * ft_plan_ultraspherical_to_ultrasphericalf(const int norm1, const int norm2, const int n, const float lambda, const float mu)
A single precision version of ft_plan_ultraspherical_to_ultraspherical.
ft_tb_eigen_FMMl * ft_plan_jacobi_to_ultrasphericall(const int normjac, const int normultra, const int n, const long double alpha, const long double beta, const long double lambda)
A long double precision version of ft_plan_jacobi_to_ultraspherical.
void ft_kernel_tri_lo2hi(const ft_rotation_plan *RP, const int m1, const int m2, double *A, const int S)
Convert a single vector of triangular harmonic coefficients in A with stride S from order m1 up to or...
ft_modified_planl * ft_plan_modified_jacobi_to_jacobil(const int n, const long double alpha, const long double beta, const int nu, const long double *u, const int nv, const long double *v, const int verbose)
A long double precision version of ft_plan_modified_jacobi_to_jacobi.
ft_modified_planf * ft_plan_modified_laguerre_to_laguerref(const int n, const float alpha, const int nu, const float *u, const int nv, const float *v, const int verbose)
A single precision version of ft_plan_modified_laguerre_to_laguerre.
void ft_execute_tri2cheb(const char TRANS, const ft_harmonic_plan *P, double *A, const int N, const int M)
Transform a triangular harmonic expansion to a bivariate Chebyshev series.
void ft_execute_spinsph_analysis(const char TRANS, const ft_spinsphere_fftw_plan *P, ft_complex *X, const int N, const int M)
Execute FFTW analysis on the sphere with spin.
ft_modified_plan * ft_plan_modified_laguerre_to_laguerre(const int n, const double alpha, const int nu, const double *u, const int nv, const double *v, const int verbose)
Pre-compute a factorization of the connection coefficients between modified Laguerre polynomials,...
ft_harmonic_plan * ft_plan_tet2cheb(const int n, const double alpha, const double beta, const double gamma, const double delta)
Plan a tetrahedral harmonic transform.
void ft_execute_tri_synthesis(const char TRANS, const ft_triangle_fftw_plan *P, double *X, const int N, const int M)
Execute FFTW synthesis on the triangle.
void ft_kernel_spinsph_hi2lo(const ft_spin_rotation_plan *SRP, const int m, ft_complex *A, const int S)
Convert a single vector of spin-weighted spherical harmonic coefficients in A with stride S from orde...
ft_tb_eigen_FMMf * ft_plan_chebyshev_to_legendref(const int normcheb, const int normleg, const int n)
A single precision version of ft_plan_chebyshev_to_legendre.
void ft_kernel_tri_hi2lo(const ft_rotation_plan *RP, const int m1, const int m2, double *A, const int S)
Convert a single vector of triangular harmonic coefficients in A with stride S from order m2 down to ...
void ft_execute_sph2fourier(const char TRANS, const ft_harmonic_plan *P, double *A, const int N, const int M)
Transform a spherical harmonic expansion to a bivariate Fourier series.
ft_tb_eigen_FMMl * ft_plan_chebyshev_to_ultrasphericall(const int normcheb, const int normultra, const int n, const long double lambda)
A long double precision version of ft_plan_chebyshev_to_ultraspherical.
ft_modified_planl * ft_plan_modified_laguerre_to_laguerrel(const int n, const long double alpha, const int nu, const long double *u, const int nv, const long double *v, const int verbose)
A long double precision version of ft_plan_modified_laguerre_to_laguerre.
ft_tb_eigen_FMM * ft_plan_jacobi_to_jacobi(const int norm1, const int norm2, const int n, const double alpha, const double beta, const double gamma, const double delta)
Pre-compute a factorization of the connection coefficients between Jacobi polynomials in double preci...
ft_spinsphere_fftw_plan * ft_plan_spinsph_analysis(const int N, const int M, const int S, const unsigned flags)
Plan FFTW analysis on the sphere with spin.
ft_spin_harmonic_plan * ft_plan_spinsph2fourier(const int n, const int s)
Plan a spin-weighted spherical harmonic transform.
ft_tb_eigen_FMMf * ft_plan_laguerre_to_laguerref(const int norm1, const int norm2, const int n, const float alpha, const float beta)
A single precision version of ft_plan_laguerre_to_laguerre.
mpfr_t * ft_mpfr_plan_ultraspherical_to_chebyshev(const int normultra, const int normcheb, const int n, mpfr_srcptr lambda, mpfr_prec_t prec, mpfr_rnd_t rnd)
A multi-precision version of ft_plan_ultraspherical_to_chebyshev that returns a dense array of connec...
ft_harmonic_plan * ft_plan_sph2fourier(const int n)
Plan a spherical harmonic transform.
ft_tb_eigen_FMMl * ft_plan_ultraspherical_to_chebyshevl(const int normultra, const int normcheb, const int n, const long double lambda)
A long double precision version of ft_plan_ultraspherical_to_chebyshev.
ft_triangle_fftw_plan * ft_plan_tri_analysis(const int N, const int M, const unsigned flags)
Plan FFTW analysis on the triangle.
void ft_destroy_annulus_fftw_plan(ft_annulus_fftw_plan *P)
Destroy a ft_annulus_fftw_plan.
mpfr_t * ft_mpfr_plan_chebyshev_to_legendre(const int normcheb, const int normleg, const int n, mpfr_prec_t prec, mpfr_rnd_t rnd)
A multi-precision version of ft_plan_chebyshev_to_legendre that returns a dense array of connection c...
void ft_kernel_disk_lo2hi(const ft_rotation_plan *RP, const int m1, const int m2, double *A, const int S)
Convert a single vector of disk harmonic coefficients in A with stride S from order m1 up to order m2...
void ft_execute_cheb2tet(const char TRANS, const ft_harmonic_plan *P, double *A, const int N, const int L, const int M)
Transform a trivariate Chebyshev series to a tetrahedral harmonic expansion.
void ft_destroy_harmonic_plan(ft_harmonic_plan *P)
Destroy a ft_harmonic_plan.
ft_tb_eigen_FMM * ft_plan_chebyshev_to_ultraspherical(const int normcheb, const int normultra, const int n, const double lambda)
Pre-compute a factorization of the connection coefficients between Chebyshev and ultraspherical polyn...
void ft_execute_annulus_synthesis(const char TRANS, const ft_annulus_fftw_plan *P, double *X, const int N, const int M)
Execute FFTW synthesis on the annulus.
void ft_kernel_rectdisk_lo2hi(const ft_rotation_plan *RP, const int m1, const int m2, double *A, const int S)
Convert a single vector of rectangularized disk harmonic coefficients in A with stride S from order m...
mpfr_t * ft_mpfr_plan_legendre_to_chebyshev(const int normleg, const int normcheb, const int n, mpfr_prec_t prec, mpfr_rnd_t rnd)
A multi-precision version of ft_plan_legendre_to_chebyshev that returns a dense array of connection c...
ft_tb_eigen_FMM * ft_plan_ultraspherical_to_chebyshev(const int normultra, const int normcheb, const int n, const double lambda)
Pre-compute a factorization of the connection coefficients between ultraspherical and Chebyshev polyn...
void ft_execute_disk_analysis(const char TRANS, const ft_disk_fftw_plan *P, double *X, const int N, const int M)
Execute FFTW analysis on the disk.
void ft_kernel_sph_lo2hi(const ft_rotation_plan *RP, const int m1, const int m2, double *A, const int S)
Convert a single vector of spherical harmonic coefficients in A with stride S from order m1 up to ord...
ft_sphere_fftw_plan * ft_plan_sph_analysis(const int N, const int M, const unsigned flags)
Plan FFTW analysis on the sphere.
ft_tb_eigen_FMM * ft_plan_jacobi_to_ultraspherical(const int normjac, const int normultra, const int n, const double alpha, const double beta, const double lambda)
Pre-compute a factorization of the connection coefficients between Jacobi and ultraspherical polynomi...
ft_tb_eigen_FMM * ft_plan_chebyshev_to_jacobi(const int normcheb, const int normjac, const int n, const double alpha, const double beta)
Pre-compute a factorization of the connection coefficients between Chebyshev and Jacobi polynomials i...
ft_harmonic_plan * ft_plan_ann2cxf(const int n, const double alpha, const double beta, const double gamma, const double rho)
Plan an annulus harmonic transform.
mpfr_t * ft_mpfr_plan_ultraspherical_to_jacobi(const int normultra, const int normjac, const int n, mpfr_srcptr lambda, mpfr_srcptr alpha, mpfr_srcptr beta, mpfr_prec_t prec, mpfr_rnd_t rnd)
A multi-precision version of ft_plan_ultraspherical_to_jacobi that returns a dense array of connectio...
void ft_destroy_disk_fftw_plan(ft_disk_fftw_plan *P)
Destroy a ft_disk_fftw_plan.
ft_modified_planf * ft_plan_modified_jacobi_to_jacobif(const int n, const float alpha, const float beta, const int nu, const float *u, const int nv, const float *v, const int verbose)
A single precision version of ft_plan_modified_jacobi_to_jacobi.
void ft_destroy_rectdisk_fftw_plan(ft_rectdisk_fftw_plan *P)
Destroy a ft_rectdisk_fftw_plan.
ft_tb_eigen_FMMf * ft_plan_chebyshev_to_ultrasphericalf(const int normcheb, const int normultra, const int n, const float lambda)
A single precision version of ft_plan_chebyshev_to_ultraspherical.
mpfr_t * ft_mpfr_plan_jacobi_to_jacobi(const int norm1, const int norm2, const int n, mpfr_srcptr alpha, mpfr_srcptr beta, mpfr_srcptr gamma, mpfr_srcptr delta, mpfr_prec_t prec, mpfr_rnd_t rnd)
A multi-precision version of ft_plan_jacobi_to_jacobi that returns a dense array of connection coeffi...
mpfr_t * ft_mpfr_plan_chebyshev_to_jacobi(const int normcheb, const int normjac, const int n, mpfr_srcptr alpha, mpfr_srcptr beta, mpfr_prec_t prec, mpfr_rnd_t rnd)
A multi-precision version of ft_plan_chebyshev_to_jacobi that returns a dense array of connection coe...
ft_tb_eigen_FMMf * ft_plan_legendre_to_chebyshevf(const int normleg, const int normcheb, const int n)
A single precision version of ft_plan_legendre_to_chebyshev.
ft_annulus_fftw_plan * ft_plan_annulus_synthesis(const int N, const int M, const double rho, const unsigned flags)
Plan FFTW synthesis on the annulus.
ft_harmonic_plan * ft_plan_disk2cxf(const int n, const double alpha, const double beta)
Plan a disk harmonic transform.
void ft_destroy_triangle_fftw_plan(ft_triangle_fftw_plan *P)
Destroy a ft_triangle_fftw_plan.
void ft_kernel_rectdisk_hi2lo(const ft_rotation_plan *RP, const int m1, const int m2, double *A, const int S)
Convert a single vector of rectangularized disk harmonic coefficients in A with stride S from order m...
struct ft_rotation_plan_s ft_rotation_plan
Data structure to store sines and cosines of Givens rotations.
Definition fasttransforms.h:298
void ft_destroy_spin_harmonic_plan(ft_spin_harmonic_plan *P)
Destroy a ft_spin_harmonic_plan.
ft_spinsphere_fftw_plan * ft_plan_spinsph_synthesis(const int N, const int M, const int S, const unsigned flags)
Plan FFTW synthesis on the sphere with spin.
void ft_execute_cxf2disk(const char TRANS, const ft_harmonic_plan *P, double *A, const int N, const int M)
Transform a Chebyshev–Fourier series to a disk harmonic expansion.
void ft_execute_fourier2sph(const char TRANS, const ft_harmonic_plan *P, double *A, const int N, const int M)
Transform a bivariate Fourier series to a spherical harmonic expansion.
ft_tb_eigen_FMMl * ft_plan_laguerre_to_laguerrel(const int norm1, const int norm2, const int n, const long double alpha, const long double beta)
A long double precision version of ft_plan_laguerre_to_laguerre.
void ft_execute_ann2cxf(const char TRANS, const ft_harmonic_plan *P, double *A, const int N, const int M)
Transform an annulus harmonic expansion to a Chebyshev–Fourier series.
ft_tb_eigen_FMMl * ft_plan_chebyshev_to_legendrel(const int normcheb, const int normleg, const int n)
A long double precision version of ft_plan_chebyshev_to_legendre.
ft_tb_eigen_FMMl * ft_plan_legendre_to_chebyshevl(const int normleg, const int normcheb, const int n)
A long double precision version of ft_plan_legendre_to_chebyshev.
ft_triangle_fftw_plan * ft_plan_tri_synthesis(const int N, const int M, const unsigned flags)
Plan FFTW synthesis on the triangle.
ft_annulus_fftw_plan * ft_plan_annulus_analysis(const int N, const int M, const double rho, const unsigned flags)
Plan FFTW analysis on the annulus.
mpfr_t * ft_mpfr_plan_ultraspherical_to_ultraspherical(const int norm1, const int norm2, const int n, mpfr_srcptr lambda, mpfr_srcptr mu, mpfr_prec_t prec, mpfr_rnd_t rnd)
A multi-precision version of ft_plan_ultraspherical_to_ultraspherical that returns a dense array of c...
void ft_kernel_spinsph_lo2hi(const ft_spin_rotation_plan *SRP, const int m, ft_complex *A, const int S)
Convert a single vector of spin-weighted spherical harmonic coefficients in A with stride S from orde...
ft_disk_fftw_plan * ft_plan_disk_synthesis(const int N, const int M, const unsigned flags)
Plan FFTW synthesis on the disk.
ft_tb_eigen_FMM * ft_plan_ultraspherical_to_ultraspherical(const int norm1, const int norm2, const int n, const double lambda, const double mu)
Pre-compute a factorization of the connection coefficients between ultraspherical polynomials in doub...
void ft_kernel_sph_hi2lo(const ft_rotation_plan *RP, const int m1, const int m2, double *A, const int S)
Convert a single vector of spherical harmonic coefficients in A with stride S from order m2 down to o...
ft_tb_eigen_FMM * ft_plan_legendre_to_chebyshev(const int normleg, const int normcheb, const int n)
Pre-compute a factorization of the connection coefficients between Legendre and Chebyshev polynomials...
ft_harmonic_plan * ft_plan_rectdisk2cheb(const int n, const double beta)
Plan a rectangularized disk harmonic transform.
ft_rectdisk_fftw_plan * ft_plan_rectdisk_synthesis(const int N, const int M, const unsigned flags)
Plan FFTW synthesis on the rectangularized disk.
mpfr_t * ft_mpfr_plan_laguerre_to_laguerre(const int norm1, const int norm2, const int n, mpfr_srcptr alpha, mpfr_srcptr beta, mpfr_prec_t prec, mpfr_rnd_t rnd)
A multi-precision version of ft_plan_laguerre_to_laguerre that returns a dense array of connection co...
void ft_destroy_sphere_fftw_plan(ft_sphere_fftw_plan *P)
Destroy a ft_sphere_fftw_plan.
void ft_execute_rectdisk2cheb(const char TRANS, const ft_harmonic_plan *P, double *A, const int N, const int M)
Transform a rectangularized disk harmonic expansion to a Chebyshev series.
mpfr_t * ft_mpfr_plan_jacobi_to_ultraspherical(const int normjac, const int normultra, const int n, mpfr_srcptr alpha, mpfr_srcptr beta, mpfr_srcptr lambda, mpfr_prec_t prec, mpfr_rnd_t rnd)
A multi-precision version of ft_plan_jacobi_to_ultraspherical that returns a dense array of connectio...
ft_tb_eigen_FMM * ft_plan_ultraspherical_to_jacobi(const int normultra, const int normjac, const int n, const double lambda, const double alpha, const double beta)
Pre-compute a factorization of the connection coefficients between ultraspherical and Jacobi polynomi...
void ft_execute_tet2cheb(const char TRANS, const ft_harmonic_plan *P, double *A, const int N, const int L, const int M)
Transform a tetrahedral harmonic expansion to a trivariate Chebyshev series.
mpfr_t * ft_mpfr_plan_chebyshev_to_ultraspherical(const int normcheb, const int normultra, const int n, mpfr_srcptr lambda, mpfr_prec_t prec, mpfr_rnd_t rnd)
A multi-precision version of ft_plan_chebyshev_to_ultraspherical that returns a dense array of connec...
void ft_execute_cheb2tri(const char TRANS, const ft_harmonic_plan *P, double *A, const int N, const int M)
Transform a bivariate Chebyshev series to a triangular harmonic expansion.
ft_modified_plan * ft_plan_modified_hermite_to_hermite(const int n, const int nu, const double *u, const int nv, const double *v, const int verbose)
Pre-compute a factorization of the connection coefficients between modified Hermite polynomials,...
void ft_execute_disk_synthesis(const char TRANS, const ft_disk_fftw_plan *P, double *X, const int N, const int M)
Execute FFTW synthesis on the disk.
void ft_execute_cxf2ann(const char TRANS, const ft_harmonic_plan *P, double *A, const int N, const int M)
Transform a Chebyshev–Fourier series to an annulus harmonic expansion.
void ft_execute_tri_analysis(const char TRANS, const ft_triangle_fftw_plan *P, double *X, const int N, const int M)
Execute FFTW analysis on the triangle.
void ft_execute_sph_synthesis(const char TRANS, const ft_sphere_fftw_plan *P, double *X, const int N, const int M)
Execute FFTW synthesis on the sphere.
void ft_destroy_rotation_plan(ft_rotation_plan *RP)
Destroy a ft_rotation_plan.
ft_tb_eigen_FMMf * ft_plan_ultraspherical_to_chebyshevf(const int normultra, const int normcheb, const int n, const float lambda)
A single precision version of ft_plan_ultraspherical_to_chebyshev.
ft_tb_eigen_FMMf * ft_plan_jacobi_to_chebyshevf(const int normjac, const int normcheb, const int n, const float alpha, const float beta)
A single precision version of ft_plan_jacobi_to_chebyshev.
void ft_execute_spinsph_synthesis(const char TRANS, const ft_spinsphere_fftw_plan *P, ft_complex *X, const int N, const int M)
Execute FFTW synthesis on the sphere with spin.
ft_tb_eigen_FMMf * ft_plan_jacobi_to_ultrasphericalf(const int normjac, const int normultra, const int n, const float alpha, const float beta, const float lambda)
A single precision version of ft_plan_jacobi_to_ultraspherical.
void ft_set_num_threads(const int n)
Set the number of OpenMP threads.
ft_tb_eigen_FMMl * ft_plan_jacobi_to_jacobil(const int norm1, const int norm2, const int n, const long double alpha, const long double beta, const long double gamma, const long double delta)
A long double precision version of ft_plan_jacobi_to_jacobi.
void ft_execute_rectdisk_synthesis(const char TRANS, const ft_rectdisk_fftw_plan *P, double *X, const int N, const int M)
Execute FFTW synthesis on the rectangularized disk.
ft_tb_eigen_FMM * ft_plan_chebyshev_to_legendre(const int normcheb, const int normleg, const int n)
Pre-compute a factorization of the connection coefficients between Chebyshev and Legendre polynomials...
ft_tb_eigen_FMMf * ft_plan_ultraspherical_to_jacobif(const int normultra, const int normjac, const int n, const float lambda, const float alpha, const float beta)
A single precision version of ft_plan_ultraspherical_to_jacobi.
ft_modified_planl * ft_plan_modified_hermite_to_hermitel(const int n, const int nu, const long double *u, const int nv, const long double *v, const int verbose)
A long double precision version of ft_plan_modified_hermite_to_hermite.
Every orthogonal matrix can be decomposed as a product of Euler angles and, if necessary,...
Definition fasttransforms.h:672
Definition fasttransforms.h:549
Definition fasttransforms.h:526
Definition fasttransforms.h:619
Data structure to store ft_rotation_plans, arrays to represent 1D orthogonal polynomial transforms an...
Definition fasttransforms.h:371
Definition fasttransforms.h:632
Definition fasttransforms.h:253
A static struct to store an orthogonal matrix , such that . has column-major storage.
Definition fasttransforms.h:649
Definition fasttransforms.h:695
Definition fasttransforms.h:573
A static struct to store a reflection about the plane in .
Definition fasttransforms.h:681
Definition fasttransforms.h:703
Definition fasttransforms.h:467
Data structure to store a ft_spin_rotation_plan, and various arrays to represent 1D orthogonal polyno...
Definition fasttransforms.h:442
Definition fasttransforms.h:593
Definition fasttransforms.h:513
Definition fasttransforms.h:495