|
S2kit
1.1
Toolkit for working with functions defined on the sphere
|
Go to the source code of this file.
Functions | |
| void | FSTSemiMemo (double *, double *, double *, double *, const int, double **, double *, DataFormat, const int, fftw_plan *, fftw_plan *, double *) |
| Computes spherical harmonic transform using the seminaive and naive algorithms. More... | |
| void | InvFSTSemiMemo (double *, double *, double *, double *, const int, double **, double *, DataFormat, const int, fftw_plan *, fftw_plan *) |
| Computes inverse spherical harmonic transform. More... | |
| void | FZTSemiMemo (double *, double *, double *, double *, const int, double *, double *, DataFormat, fftw_plan *, double *) |
| Computes zonal harmonic transform using seminaive algorithm. More... | |
| void | ConvOn2SphereSemiMemo (double *, double *, double *, double *, double *, double *, const int, double *) |
| Convolves two functions defined on the 2-sphere. More... | |
| void ConvOn2SphereSemiMemo | ( | double * | rdata, |
| double * | idata, | ||
| double * | rfilter, | ||
| double * | ifilter, | ||
| double * | rres, | ||
| double * | ires, | ||
| const int | bw, | ||
| double * | workspace | ||
| ) |
Convolves two functions defined on the 2-sphere.
Uses seminaive algorithms for spherical harmonic transforms.
Conv2Sphere requires memory for spharmonic tables, local workspace and workspace for FSTSemiMemo() (reuses in InvFSTSemiMemo()):
| rdata | array of length 4*bw*bw of the real part of sampled function |
| idata | array of length 4*bw*bw of the imaginary part of sampled function |
| rfilter | array of length 4*bw*bw of the real part of sampled filter function |
| ifilter | array of length 4*bw*bw of the imaginary part of sampled filter function |
| rres | array of length 4*bw*bw which will contain the real part of result function |
| ires | array of length 4*bw*bw which will contain the imaginary part of result function |
| bw | bandwidth of problem |
| workspace | space for computations of size (2*legendre_size)+(12*(bw*bw)+12*bw) |
Definition at line 439 of file FST_semi_memo.c.
| void FSTSemiMemo | ( | double * | rdata, |
| double * | idata, | ||
| double * | rcoeffs, | ||
| double * | icoeffs, | ||
| const int | bw, | ||
| double ** | seminaive_naive_table, | ||
| double * | workspace, | ||
| DataFormat | data_format, | ||
| const int | cutoff, | ||
| fftw_plan * | DCT_plan, | ||
| fftw_plan * | FFT_plan, | ||
| double * | weights | ||
| ) |
Computes spherical harmonic transform using the seminaive and naive algorithms.
Output ordering of coeffs (in rcoeffs and icoeffs) f(m,l) is:
This only requires an array of size bw*bw. If zero-padding is used to make the indexing nice, then you need a an (2bw-1)*bw array, but that is not done here. Because of the amount of space necessary for doing large transforms, it is important not to use any more than necessary.
| rdata | array of length 4*bw*bw of the real part of the function samples |
| idata | array of length 4*bw*bw of the imaginary part of the function samples |
| rcoeffs | array of length bw*bw which will contain the real part of harmonic coefficients in a linearized form |
| icoeffs | array of length bw*bw which will contain the imaginary part of harmonic coefficients in a linearized form |
| bw | bandwidth of problem |
| seminaive_naive_table | pre-generated spharmonic Pml table |
| workspace | space for computations of size (8*bw^2)+(7*bw) |
| data_format | determines the format of the computed data |
| cutoff | determines order to switch from seminaive to naive algorithm |
| DCT_plan | plan for DLT which is used as argument for call DLTSemi() |
| FFT_plan | plan for FFT |
| weights | array which is used as argument for call DLTSemi() and DLTNaive() |
seminaive_naive_table should be generated by Spharmonic_Pml_Table(). This table can be re-used in the inverse transform and, for example, in series of convolutions. DCT_plan and weights in DLTSemi(). Definition at line 68 of file FST_semi_memo.c.
| void FZTSemiMemo | ( | double * | rdata, |
| double * | idata, | ||
| double * | rres, | ||
| double * | ires, | ||
| const int | bw, | ||
| double * | cos_pml_table, | ||
| double * | workspace, | ||
| const DataFormat | data_format, | ||
| fftw_plan * | DCT_plan, | ||
| double * | weights | ||
| ) |
Computes zonal harmonic transform using seminaive algorithm.
This transform is used in convolutions. Only computes spherical harmonics for m=0.
| rdata | array of length 4*bw*bw of the real part of the data samples |
| idata | array of length 4*bw*bw of the imaginary part of the data samples |
| rres | array of length bw which will contain the real part of harmonic coefficients |
| ires | array of length bw which will contain the imaginary part of harmonic coefficients |
| bw | bandwidth of problem |
| cos_pml_table | pre-generated table of Legendre coefficients of P(0,l) functions |
| workspace | space for computations of size 12*bw |
| data_format | determines the format of the computed data |
| DCT_plan | plan for DLT which is used as argument for call DLTSemi() |
| weights | array which is used as argument for call DLTSemi() |
cos_pml_table should be generated by GenerateCosPmlTable() with m = 0. DCT_plan and weights in DLTSemi(). Definition at line 374 of file FST_semi_memo.c.
| void InvFSTSemiMemo | ( | double * | rcoeffs, |
| double * | icoeffs, | ||
| double * | rdata, | ||
| double * | idata, | ||
| const int | bw, | ||
| double ** | transpose_seminaive_naive_table, | ||
| double * | workspace, | ||
| DataFormat | data_format, | ||
| const int | cutoff, | ||
| fftw_plan * | inv_DCT_plan, | ||
| fftw_plan * | inv_FFT_plan | ||
| ) |
Computes inverse spherical harmonic transform.
| rcoeffs | array of length bw*bw of the real part of harmonic coefficients |
| icoeffs | array of length bw*bw of the imaginary part of harmonic coefficients |
| rdata | array of length 4*bw*bw which will contain the real part of transformed result |
| idata | array of length 4*bw*bw which will contain the imaginary part of transformed result |
| bw | bandwidth of problem |
| transpose_seminaive_naive_table | pre-generated transposed spharmonic Pml table |
| workspace | space for computations of size (8*bw^2)+(10*bw) |
| data_format | determines the format of the computed data |
| cutoff | determines order to switch from seminaive to naive algorithm |
| inv_DCT_plan | plan for inverse DLT which is used as argument for call InvDLTSemi() |
| inv_FFT_plan | plan for inverse FFT |
rcoeffs and icoeffs before passing to this function in the same order as it was mentioned in FSTSemiMemo(). transpose_seminaive_naive_table should be generated by Transpose_Spharmonic_Pml_Table(). inv_DCT_plan in InvDLTSemi(). Definition at line 228 of file FST_semi_memo.c.