FastTransforms
Loading...
Searching...
No Matches
Data Structures | Macros | Typedefs | Functions
fasttransforms.h File Reference

fasttransforms.h is the main header file for FastTransforms. More...

#include <cblas.h>
#include <fftw3.h>
#include "tdc.h"
#include <mpfr.h>

Go to the source code of this file.

Data Structures

struct  ft_mpfr_triangular_banded
 
struct  ft_harmonic_plan
 Data structure to store ft_rotation_plans, arrays to represent 1D orthogonal polynomial transforms and their inverses, and Greek parameters. More...
 
struct  ft_spin_harmonic_plan
 Data structure to store a ft_spin_rotation_plan, and various arrays to represent 1D orthogonal polynomial transforms. More...
 
struct  ft_sphere_fftw_plan
 
struct  ft_triangle_fftw_plan
 
struct  ft_tetrahedron_fftw_plan
 
struct  ft_disk_fftw_plan
 
struct  ft_annulus_fftw_plan
 
struct  ft_rectdisk_fftw_plan
 
struct  ft_spinsphere_fftw_plan
 
struct  ft_gradient_plan
 
struct  ft_helmholtzhodge_plan
 
struct  ft_orthogonal_transformation
 A static struct to store an orthogonal matrix \(Q \in \mathbb{R}^{3\times3}\), such that \(Q^\top Q = I\). \(Q\) has column-major storage. More...
 
struct  ft_ZYZR
 Every orthogonal matrix \(Q \in \mathbb{R}^{3\times3}\) can be decomposed as a product of \(ZYZ\) Euler angles and, if necessary, a reflection \(R\) about the \(xy\)-plane. More...
 
struct  ft_reflection
 A static struct to store a reflection about the plane \(w\cdot x = 0\) in \(\mathbb{R}^3\). More...
 
struct  ft_partial_sph_isometry_plan
 
struct  ft_sph_isometry_plan
 

Macros

#define FT_GET_THREAD_NUM()   0
 
#define FT_GET_NUM_THREADS()   1
 
#define FT_GET_MAX_THREADS()   1
 
#define FT_SET_NUM_THREADS(x)
 
#define FT_CONCAT(prefix, name, suffix)   prefix ## name ## suffix
 
#define FT_SN   (1U << 0)
 
#define FT_CN   (1U << 1)
 
#define FT_DN   (1U << 2)
 
#define FT_MODIFIED_NMAX   (1U << 30)
 

Typedefs

typedef double ft_complex[2]
 
typedef struct ft_rotation_plan_s ft_rotation_plan
 Data structure to store sines and cosines of Givens rotations.
 
typedef struct ft_spin_rotation_plan_s ft_spin_rotation_plan
 Data structure to store sines and cosines of Givens rotations for spin-weighted harmonics.
 

Functions

void ft_horner (const int n, const double *c, const int incc, const int m, double *x, double *f)
 
void ft_hornerf (const int n, const float *c, const int incc, const int m, float *x, float *f)
 
void ft_clenshaw (const int n, const double *c, const int incc, const int m, double *x, double *f)
 
void ft_clenshawf (const int n, const float *c, const int incc, const int m, float *x, float *f)
 
void 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)
 
void 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)
 
void 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)
 
void 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)
 
void 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)
 
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 in double precision so that ft_bfmv converts between expansions: More...
 
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 in double precision so that ft_bfmv converts between expansions: More...
 
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 double precision so that ft_bfmv converts between expansions: More...
 
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 precision so that ft_bfmv converts between expansions: More...
 
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 precision so that ft_bfmv converts between expansions: More...
 
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 polynomials in double precision so that ft_bfmv converts between expansions: More...
 
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 polynomials in double precision so that ft_bfmv converts between expansions: More...
 
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 in double precision so that ft_bfmv converts between expansions: More...
 
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 in double precision so that ft_bfmv converts between expansions: More...
 
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 polynomials in double precision so that ft_bfmv converts between expansions: More...
 
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 polynomials in double precision so that ft_bfmv converts between expansions: More...
 
ft_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)
 
ft_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)
 
ft_btb_eigen_FMM * ft_plan_associated_hermite_to_hermite (const int norm1, const int norm2, const int n, const int c)
 
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, orthonormal with respect to a rationally modified weight and orthonormal Jacobi polynomials. The rational function is expressed as a ratio of orthonormal Jacobi polynomial expansions: More...
 
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, orthonormal with respect to a rationally modified weight and orthonormal Laguerre polynomials. The rational function is expressed as a ratio of orthonormal Laguerre polynomial expansions: More...
 
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, orthonormal with respect to a rationally modified weight and orthonormal Hermite polynomials. The rational function is expressed as a ratio of orthonormal Hermite polynomial expansions: More...
 
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_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.
 
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_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.
 
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.
 
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.
 
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_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.
 
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_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_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.
 
ft_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)
 
ft_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)
 
ft_btb_eigen_FMMf * ft_plan_associated_hermite_to_hermitef (const int norm1, const int norm2, const int n, const int c)
 
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.
 
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.
 
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_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_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_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.
 
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.
 
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.
 
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.
 
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_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.
 
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_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_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_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)
 
ft_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)
 
ft_btb_eigen_FMMl * ft_plan_associated_hermite_to_hermitel (const int norm1, const int norm2, const int n, const int c)
 
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_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_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.
 
void ft_mpfr_destroy_plan (mpfr_t *A, int n)
 
void ft_mpfr_trmv (char TRANS, int n, mpfr_t *A, int LDA, mpfr_t *x, mpfr_rnd_t rnd)
 
void ft_mpfr_trsv (char TRANS, int n, mpfr_t *A, int LDA, mpfr_t *x, mpfr_rnd_t rnd)
 
void ft_mpfr_trmm (char TRANS, int n, mpfr_t *A, int LDA, mpfr_t *B, int LDB, int N, mpfr_rnd_t rnd)
 
void ft_mpfr_trsm (char TRANS, int n, mpfr_t *A, int LDA, mpfr_t *B, int LDB, int N, mpfr_rnd_t rnd)
 
void ft_mpfr_trmv_ptr (char TRANS, int n, mpfr_t *A, int LDA, mpfr_t **x, mpfr_rnd_t rnd)
 
void ft_mpfr_trsv_ptr (char TRANS, int n, mpfr_t *A, int LDA, mpfr_t **x, mpfr_rnd_t rnd)
 
void ft_mpfr_trmm_ptr (char TRANS, int n, mpfr_t *A, int LDA, mpfr_t **B, int LDB, int N, mpfr_rnd_t rnd)
 
void ft_mpfr_trsm_ptr (char TRANS, int n, mpfr_t *A, int LDA, mpfr_t **B, int LDB, int N, mpfr_rnd_t rnd)
 
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 coefficients.
 
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 coefficients.
 
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 connection coefficients.
 
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 coefficients.
 
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 coefficients.
 
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 connection coefficients.
 
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 connection coefficients.
 
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 coefficients.
 
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 coefficients.
 
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 connection coefficients.
 
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 connection coefficients.
 
void ft_set_num_threads (const int n)
 Set the number of OpenMP threads.
 
void ft_destroy_rotation_plan (ft_rotation_plan *RP)
 Destroy a ft_rotation_plan.
 
ft_rotation_planft_plan_rotsphere (const int n)
 
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 order m1.
 
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 order m2.
 
ft_rotation_planft_plan_rottriangle (const int n, const double alpha, const double beta, const double gamma)
 
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 order m1.
 
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 order m2.
 
ft_rotation_planft_plan_rotdisk (const int n, const double alpha, const double beta)
 
ft_rotation_planft_plan_rotannulus (const int n, const double alpha, const double beta, const double gamma, const double rho)
 
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 m1.
 
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.
 
ft_rotation_planft_plan_rotrectdisk (const int n, const double beta)
 
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 m2 down to order m1.
 
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 m1 up to order m2.
 
void ft_kernel_tet_hi2lo (const ft_rotation_plan *RP, const int L, const int m, double *A)
 
void ft_kernel_tet_lo2hi (const ft_rotation_plan *RP, const int L, const int m, double *A)
 
void ft_destroy_spin_rotation_plan (ft_spin_rotation_plan *SRP)
 Destroy a ft_spin_rotation_plan.
 
ft_spin_rotation_planft_plan_rotspinsphere (const int n, const int s)
 
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 order m down to order m%2.
 
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 order m%2 up to order m.
 
void ft_execute_sph_hi2lo (const ft_rotation_plan *RP, double *A, double *B, const int M)
 
void ft_execute_sph_lo2hi (const ft_rotation_plan *RP, double *A, double *B, const int M)
 
void ft_execute_sphv_hi2lo (const ft_rotation_plan *RP, double *A, double *B, const int M)
 
void ft_execute_sphv_lo2hi (const ft_rotation_plan *RP, double *A, double *B, const int M)
 
void ft_execute_tri_hi2lo (const ft_rotation_plan *RP, double *A, double *B, const int M)
 
void ft_execute_tri_lo2hi (const ft_rotation_plan *RP, double *A, double *B, const int M)
 
void ft_execute_disk_hi2lo (const ft_rotation_plan *RP, double *A, double *B, const int M)
 
void ft_execute_disk_lo2hi (const ft_rotation_plan *RP, double *A, double *B, const int M)
 
void ft_execute_rectdisk_hi2lo (const ft_rotation_plan *RP, double *A, double *B, const int M)
 
void ft_execute_rectdisk_lo2hi (const ft_rotation_plan *RP, double *A, double *B, const int M)
 
void ft_execute_tet_hi2lo (const ft_rotation_plan *RP1, const ft_rotation_plan *RP2, double *A, const int L, const int M)
 
void ft_execute_tet_lo2hi (const ft_rotation_plan *RP1, const ft_rotation_plan *RP2, double *A, const int L, const int M)
 
void ft_execute_spinsph_hi2lo (const ft_spin_rotation_plan *SRP, ft_complex *A, ft_complex *B, const int M)
 
void ft_execute_spinsph_lo2hi (const ft_spin_rotation_plan *SRP, ft_complex *A, ft_complex *B, const int M)
 
void ft_destroy_harmonic_plan (ft_harmonic_plan *P)
 Destroy a ft_harmonic_plan.
 
ft_harmonic_planft_plan_sph2fourier (const int n)
 Plan a spherical harmonic transform.
 
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.
 
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.
 
void ft_execute_sphv2fourier (const char TRANS, const ft_harmonic_plan *P, double *A, const int N, const int M)
 
void ft_execute_fourier2sphv (const char TRANS, const ft_harmonic_plan *P, double *A, const int N, const int M)
 
ft_harmonic_planft_plan_tri2cheb (const int n, const double alpha, const double beta, const double gamma)
 Plan a triangular harmonic transform.
 
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_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_harmonic_planft_plan_disk2cxf (const int n, const double alpha, const double beta)
 Plan a disk harmonic transform.
 
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.
 
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.
 
ft_harmonic_planft_plan_ann2cxf (const int n, const double alpha, const double beta, const double gamma, const double rho)
 Plan an annulus harmonic transform.
 
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.
 
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.
 
ft_harmonic_planft_plan_rectdisk2cheb (const int n, const double beta)
 Plan a rectangularized disk harmonic transform.
 
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.
 
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_harmonic_planft_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_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.
 
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_spin_harmonic_plan (ft_spin_harmonic_plan *P)
 Destroy a ft_spin_harmonic_plan.
 
ft_spin_harmonic_planft_plan_spinsph2fourier (const int n, const int s)
 Plan a spin-weighted spherical harmonic transform.
 
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.
 
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.
 
int ft_fftw_init_threads (void)
 
void ft_fftw_plan_with_nthreads (const int n)
 
void ft_destroy_sphere_fftw_plan (ft_sphere_fftw_plan *P)
 Destroy a ft_sphere_fftw_plan.
 
ft_sphere_fftw_planft_plan_sph_with_kind (const int N, const int M, const fftw_r2r_kind kind[3][1], const unsigned flags)
 
ft_sphere_fftw_planft_plan_sph_synthesis (const int N, const int M, const unsigned flags)
 Plan FFTW synthesis on the sphere.
 
ft_sphere_fftw_planft_plan_sph_analysis (const int N, const int M, const unsigned flags)
 Plan FFTW analysis on the sphere.
 
ft_sphere_fftw_planft_plan_sphv_synthesis (const int N, const int M, const unsigned flags)
 
ft_sphere_fftw_planft_plan_sphv_analysis (const int N, const int M, const unsigned flags)
 
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_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.
 
void ft_execute_sphv_synthesis (const char TRANS, const ft_sphere_fftw_plan *P, double *X, const int N, const int M)
 
void ft_execute_sphv_analysis (const char TRANS, const ft_sphere_fftw_plan *P, double *X, const int N, const int M)
 
void ft_destroy_triangle_fftw_plan (ft_triangle_fftw_plan *P)
 Destroy a ft_triangle_fftw_plan.
 
ft_triangle_fftw_planft_plan_tri_with_kind (const int N, const int M, const fftw_r2r_kind kind0, const fftw_r2r_kind kind1, const unsigned flags)
 
ft_triangle_fftw_planft_plan_tri_synthesis (const int N, const int M, const unsigned flags)
 Plan FFTW synthesis on the triangle.
 
ft_triangle_fftw_planft_plan_tri_analysis (const int N, const int M, const unsigned flags)
 Plan FFTW analysis on the triangle.
 
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_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_destroy_tetrahedron_fftw_plan (ft_tetrahedron_fftw_plan *P)
 
ft_tetrahedron_fftw_planft_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)
 
ft_tetrahedron_fftw_planft_plan_tet_synthesis (const int N, const int L, const int M, const unsigned flags)
 
ft_tetrahedron_fftw_planft_plan_tet_analysis (const int N, const int L, const int M, const unsigned flags)
 
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)
 
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)
 
void ft_destroy_disk_fftw_plan (ft_disk_fftw_plan *P)
 Destroy a ft_disk_fftw_plan.
 
ft_disk_fftw_planft_plan_disk_with_kind (const int N, const int M, const fftw_r2r_kind kind[3][1], const unsigned flags)
 
ft_disk_fftw_planft_plan_disk_synthesis (const int N, const int M, const unsigned flags)
 Plan FFTW synthesis on the disk.
 
ft_disk_fftw_planft_plan_disk_analysis (const int N, const int M, const unsigned flags)
 Plan FFTW analysis on the disk.
 
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_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_destroy_annulus_fftw_plan (ft_annulus_fftw_plan *P)
 Destroy a ft_annulus_fftw_plan.
 
double ft_get_rho_annulus_fftw_plan (ft_annulus_fftw_plan *P)
 
ft_annulus_fftw_planft_plan_annulus_with_kind (const int N, const int M, const double rho, const fftw_r2r_kind kind[3][1], const unsigned flags)
 
ft_annulus_fftw_planft_plan_annulus_synthesis (const int N, const int M, const double rho, const unsigned flags)
 Plan FFTW synthesis on the annulus.
 
ft_annulus_fftw_planft_plan_annulus_analysis (const int N, const int M, const double rho, const unsigned flags)
 Plan FFTW analysis on the annulus.
 
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_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_destroy_rectdisk_fftw_plan (ft_rectdisk_fftw_plan *P)
 Destroy a ft_rectdisk_fftw_plan.
 
ft_rectdisk_fftw_planft_plan_rectdisk_with_kind (const int N, const int M, const fftw_r2r_kind kind[3][1], const unsigned flags)
 
ft_rectdisk_fftw_planft_plan_rectdisk_synthesis (const int N, const int M, const unsigned flags)
 Plan FFTW synthesis on the rectangularized disk.
 
ft_rectdisk_fftw_planft_plan_rectdisk_analysis (const int N, const int M, const unsigned flags)
 Plan FFTW analysis on the rectangularized disk.
 
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.
 
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_destroy_spinsphere_fftw_plan (ft_spinsphere_fftw_plan *P)
 Destroy a ft_spinsphere_fftw_plan.
 
int ft_get_spin_spinsphere_fftw_plan (const ft_spinsphere_fftw_plan *P)
 
ft_spinsphere_fftw_planft_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)
 
ft_spinsphere_fftw_planft_plan_spinsph_synthesis (const int N, const int M, const int S, const unsigned flags)
 Plan FFTW synthesis on the sphere with spin.
 
ft_spinsphere_fftw_planft_plan_spinsph_analysis (const int N, const int M, const int S, const unsigned flags)
 Plan FFTW analysis on the sphere with spin.
 
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.
 
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.
 
void ft_destroy_gradient_plan (ft_gradient_plan *P)
 
ft_gradient_planft_plan_sph_gradient (const int n)
 
void ft_execute_sph_gradient (ft_gradient_plan *P, double *U, double *Ut, double *Up, const int N, const int M)
 
void ft_execute_sph_curl (ft_gradient_plan *P, double *U, double *Ut, double *Up, const int N, const int M)
 
void ft_destroy_helmholtzhodge_plan (ft_helmholtzhodge_plan *P)
 
ft_helmholtzhodge_planft_plan_sph_helmholtzhodge (const int n)
 
void ft_execute_sph_helmholtzhodge (ft_helmholtzhodge_plan *P, double *U1, double *U2, double *V1, double *V2, const int N, const int M)
 
ft_ZYZR ft_create_ZYZR (ft_orthogonal_transformation Q)
 
void ft_apply_ZYZR (ft_ZYZR Q, ft_orthogonal_transformation *U)
 
void ft_apply_reflection (ft_reflection Q, ft_orthogonal_transformation *U)
 
void ft_execute_sph_polar_rotation (double *A, const int N, const int M, double s, double c)
 
void ft_execute_sph_polar_reflection (double *A, const int N, const int M)
 
void ft_destroy_partial_sph_isometry_plan (ft_partial_sph_isometry_plan *F)
 
void ft_destroy_sph_isometry_plan (ft_sph_isometry_plan *F)
 
ft_partial_sph_isometry_planft_plan_partial_sph_isometry (const int l)
 
ft_sph_isometry_planft_plan_sph_isometry (const int n)
 
void ft_execute_sph_yz_axis_exchange (ft_sph_isometry_plan *J, double *A, const int N, const int M)
 
void ft_execute_sph_isometry (ft_sph_isometry_plan *J, ft_ZYZR Q, double *A, const int N, const int M)
 
void 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)
 
void ft_execute_sph_reflection (ft_sph_isometry_plan *J, ft_reflection W, double *A, const int N, const int M)
 
void ft_execute_sph_orthogonal_transformation (ft_sph_isometry_plan *J, ft_orthogonal_transformation Q, double *A, const int N, const int M)
 

Detailed Description

fasttransforms.h is the main header file for FastTransforms.

Function Documentation

◆ ft_plan_chebyshev_to_jacobi()

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 in double precision so that ft_bfmv converts between expansions:

\[ \sum_{\ell=0}^{n-1} c_\ell^{\mathrm{Chebyshev}} T_\ell(x) = \sum_{\ell=0}^{n-1} c_\ell^{\mathrm{Jacobi}} P_\ell^{(\alpha,\beta)}(x). \]

normcheb and normjac govern the normalizations, either standard ( == 0) or orthonormalized ( == 1).
See also ft_plan_chebyshev_to_jacobif, ft_plan_chebyshev_to_jacobil, and ft_mpfr_plan_chebyshev_to_jacobi.

◆ ft_plan_chebyshev_to_legendre()

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 in double precision so that ft_bfmv converts between expansions:

\[ \sum_{\ell=0}^{n-1} c_\ell^{\mathrm{Chebyshev}} T_\ell(x) = \sum_{\ell=0}^{n-1} c_\ell^{\mathrm{Legendre}} P_\ell(x). \]

normcheb and normleg govern the normalizations, either standard ( == 0) or orthonormalized ( == 1).
See also ft_plan_chebyshev_to_legendref, ft_plan_chebyshev_to_legendrel, and ft_mpfr_plan_chebyshev_to_legendre.

◆ ft_plan_chebyshev_to_ultraspherical()

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 polynomials in double precision so that ft_bfmv converts between expansions:

\[ \sum_{\ell=0}^{n-1} c_\ell^{\mathrm{Chebyshev}} T_\ell(x) = \sum_{\ell=0}^{n-1} c_\ell^{\mathrm{ultraspherical}} C_\ell^{(\lambda)}(x). \]

normcheb and normultra govern the normalizations, either standard ( == 0) or orthonormalized ( == 1).
See also ft_plan_chebyshev_to_ultrasphericalf, ft_plan_chebyshev_to_ultrasphericall, and ft_mpfr_plan_chebyshev_to_ultraspherical.

◆ ft_plan_jacobi_to_chebyshev()

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 in double precision so that ft_bfmv converts between expansions:

\[ \sum_{\ell=0}^{n-1} c_\ell^{\mathrm{Jacobi}} P_\ell^{(\alpha,\beta)}(x) = \sum_{\ell=0}^{n-1} c_\ell^{\mathrm{Chebyshev}} T_\ell(x). \]

normjac and normcheb govern the normalizations, either standard ( == 0) or orthonormalized ( == 1).
See also ft_plan_jacobi_to_chebyshevf, ft_plan_jacobi_to_chebyshevl, and ft_mpfr_plan_jacobi_to_chebyshev.

◆ ft_plan_jacobi_to_jacobi()

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 precision so that ft_bfmv converts between expansions:

\[ \sum_{\ell=0}^{n-1} c_\ell^{(1)} P_\ell^{(\alpha,\beta)}(x) = \sum_{\ell=0}^{n-1} c_\ell^{(2)} P_\ell^{(\gamma,\delta)}(x). \]

norm1 and norm2 govern the normalizations, either standard ( == 0) or orthonormalized ( == 1).
See also ft_plan_jacobi_to_jacobif, ft_plan_jacobi_to_jacobil, and ft_mpfr_plan_jacobi_to_jacobi.

Examples
nonlocaldiffusion.c.

◆ ft_plan_jacobi_to_ultraspherical()

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 polynomials in double precision so that ft_bfmv converts between expansions:

\[ \sum_{\ell=0}^{n-1} c_\ell^{\mathrm{Jacobi}} P_\ell^{(\alpha,\beta)}(x) = \sum_{\ell=0}^{n-1} c_\ell^{\mathrm{ultraspherical}} C_\ell^{(\lambda)}(x). \]

normjac and normultra govern the normalizations, either standard ( == 0) or orthonormalized ( == 1).
See also ft_plan_jacobi_to_ultrasphericalf, ft_plan_jacobi_to_ultrasphericall, and ft_mpfr_plan_jacobi_to_ultraspherical.

◆ ft_plan_laguerre_to_laguerre()

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 precision so that ft_bfmv converts between expansions:

\[ \sum_{\ell=0}^{n-1} c_\ell^{(1)} L_\ell^{(\alpha)}(x) = \sum_{\ell=0}^{n-1} c_\ell^{(2)} L_\ell^{(\beta)}(x). \]

norm1 and norm2 govern the normalizations, either standard ( == 0) or orthonormalized ( == 1).
See also ft_plan_laguerre_to_laguerref, ft_plan_laguerre_to_laguerrel, and ft_mpfr_plan_laguerre_to_laguerre.

Examples
subspaceangles.c.

◆ ft_plan_legendre_to_chebyshev()

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 in double precision so that ft_bfmv converts between expansions:

\[ \sum_{\ell=0}^{n-1} c_\ell^{\mathrm{Legendre}} P_\ell(x) = \sum_{\ell=0}^{n-1} c_\ell^{\mathrm{Chebyshev}} T_\ell(x). \]

normleg and normcheb govern the normalizations, either standard ( == 0) or orthonormalized ( == 1).
See also ft_plan_legendre_to_chebyshevf, ft_plan_legendre_to_chebyshevl, and ft_mpfr_plan_legendre_to_chebyshev.

◆ ft_plan_modified_hermite_to_hermite()

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, orthonormal with respect to a rationally modified weight and orthonormal Hermite polynomials. The rational function is expressed as a ratio of orthonormal Hermite polynomial expansions:

\[ r(x) = \frac{u(x)}{v(x)} = \frac{\displaystyle \sum_{k=0}^{n_u-1} u_k \tilde{H}_k(x)}{\displaystyle \sum_{k=0}^{n_v-1} v_k \tilde{H}_k(x)}. \]

If \(n_v \le 1\), then \(r(x)\) is polynomial and the algorithm is direct and faster. See also ft_plan_modified_hermite_to_hermitef and ft_plan_modified_hermite_to_hermitel.

◆ ft_plan_modified_jacobi_to_jacobi()

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, orthonormal with respect to a rationally modified weight and orthonormal Jacobi polynomials. The rational function is expressed as a ratio of orthonormal Jacobi polynomial expansions:

\[ r(x) = \frac{u(x)}{v(x)} = \frac{\displaystyle \sum_{k=0}^{n_u-1} u_k \tilde{P}_k^{(\alpha,\beta)}(x)}{\displaystyle \sum_{k=0}^{n_v-1} v_k \tilde{P}_k^{(\alpha,\beta)}(x)}. \]

If \(n_v \le 1\), then \(r(x)\) is polynomial and the algorithm is direct and faster. See also ft_plan_modified_jacobi_to_jacobif and ft_plan_modified_jacobi_to_jacobil.

◆ ft_plan_modified_laguerre_to_laguerre()

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, orthonormal with respect to a rationally modified weight and orthonormal Laguerre polynomials. The rational function is expressed as a ratio of orthonormal Laguerre polynomial expansions:

\[ r(x) = \frac{u(x)}{v(x)} = \frac{\displaystyle \sum_{k=0}^{n_u-1} u_k \tilde{L}_k^{(\alpha)}(x)}{\displaystyle \sum_{k=0}^{n_v-1} v_k \tilde{L}_k^{(\alpha)}(x)}. \]

If \(n_v \le 1\), then \(r(x)\) is polynomial and the algorithm is direct and faster. See also ft_plan_modified_laguerre_to_laguerref and ft_plan_modified_laguerre_to_laguerrel.

◆ ft_plan_ultraspherical_to_chebyshev()

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 polynomials in double precision so that ft_bfmv converts between expansions:

\[ \sum_{\ell=0}^{n-1} c_\ell^{\mathrm{ultraspherical}} C_\ell^{(\lambda)}(x) = \sum_{\ell=0}^{n-1} c_\ell^{\mathrm{Chebyshev}} T_\ell(x). \]

normultra and normcheb govern the normalizations, either standard ( == 0) or orthonormalized ( == 1).
See also ft_plan_ultraspherical_to_chebyshevf, ft_plan_ultraspherical_to_chebyshevl, and ft_mpfr_plan_ultraspherical_to_chebyshev.

◆ ft_plan_ultraspherical_to_jacobi()

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 polynomials in double precision so that ft_bfmv converts between expansions:

\[ \sum_{\ell=0}^{n-1} c_\ell^{\mathrm{ultraspherical}} C_\ell^{(\lambda)}(x) = \sum_{\ell=0}^{n-1} c_\ell^{\mathrm{Jacobi}} P_\ell^{(\alpha,\beta)}(x). \]

normultra and normjac govern the normalizations, either standard ( == 0) or orthonormalized ( == 1).
See also ft_plan_ultraspherical_to_jacobif, ft_plan_ultraspherical_to_jacobil, and ft_mpfr_plan_ultraspherical_to_jacobi.

◆ ft_plan_ultraspherical_to_ultraspherical()

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 double precision so that ft_bfmv converts between expansions:

\[ \sum_{\ell=0}^{n-1} c_\ell^{(1)} C_\ell^{(\lambda)}(x) = \sum_{\ell=0}^{n-1} c_\ell^{(2)} C_\ell^{(\mu)}(x). \]

norm1 and norm2 govern the normalizations, either standard ( == 0) or orthonormalized ( == 1).
See also ft_plan_ultraspherical_to_ultrasphericalf, ft_plan_ultraspherical_to_ultrasphericall, and ft_mpfr_plan_ultraspherical_to_ultraspherical.