S2kit
1.1
Toolkit for working with functions defined on the sphere
|
Routines to perform spherical harmonic transforms and convolutions on the 2-sphere using a combination of semi-naive and naive algorithms. More...
#include "s2kit/FST_semi_fly.h"
#include <math.h>
#include <stdlib.h>
#include <string.h>
#include <fftw3.h>
#include "s2kit/cospml.h"
#include "s2kit/pml.h"
#include "s2kit/naive.h"
#include "s2kit/seminaive.h"
#include "s2kit/weights.h"
#include "s2kit/chebyshev_nodes.h"
#include "s2kit/util.h"
Go to the source code of this file.
Functions | |
void | FSTSemiFly (double *rdata, double *idata, double *rcoeffs, double *icoeffs, const int bw, 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. More... | |
void | InvFSTSemiFly (double *rcoeffs, double *icoeffs, double *rdata, double *idata, const int bw, double *workspace, DataFormat data_format, const int cutoff, fftw_plan *inv_DCT_plan, fftw_plan *inv_FFT_plan) |
Computes inverse spherical harmonic transform. More... | |
void | FZTSemiFly (double *rdata, double *idata, double *rres, double *ires, const int bw, double *workspace, DataFormat data_format, fftw_plan *DCT_plan, double *weights) |
Computes zonal harmonic transform using seminaive algorithm. More... | |
void | ConvOn2SphereSemiFly (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. More... | |
Routines to perform spherical harmonic transforms and convolutions on the 2-sphere using a combination of semi-naive and naive algorithms.
Just like FST_semi_memo.c, except that these routines compute associated Legendre functions on the fly.
For descriptions on calling these functions, see the documentation preceding each function.
Definition in file FST_semi_fly.c.
void ConvOn2SphereSemiFly | ( | 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.
Memory requirements for Conv2Sphere:
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 (14*bw*bw)+(26*bw) |
Definition at line 454 of file FST_semi_fly.c.
void FSTSemiFly | ( | double * | rdata, |
double * | idata, | ||
double * | rcoeffs, | ||
double * | icoeffs, | ||
const int | bw, | ||
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 |
workspace | space for computations of size (10*bw*bw)+(21*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() |
DCT_plan
and weights
in DLTSemi(). Definition at line 67 of file FST_semi_fly.c.
void FZTSemiFly | ( | double * | rdata, |
double * | idata, | ||
double * | rres, | ||
double * | ires, | ||
const int | bw, | ||
double * | workspace, | ||
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 |
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() |
DCT_plan
and weights
in DLTSemi(). Definition at line 390 of file FST_semi_fly.c.
void InvFSTSemiFly | ( | double * | rcoeffs, |
double * | icoeffs, | ||
double * | rdata, | ||
double * | idata, | ||
const int | bw, | ||
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 |
workspace | space for computations of size (10*bw*bw)+(24*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 FSTSemiFly(). inv_DCT_plan
in InvDLTSemi(). Definition at line 235 of file FST_semi_fly.c.