S2kit  1.1
Toolkit for working with functions defined on the sphere
FST_semi_fly.c File Reference

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...
 

Detailed Description

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.

Function Documentation

◆ ConvOn2SphereSemiFly()

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.

Note
See an example of use in test_conv_semi_fly.c

Memory requirements for Conv2Sphere:

4*bw^2 + 2*bw + 10*bw^2 + 24*bw = 14*bw^2 + 26*bw
Parameters
rdataarray of length 4*bw*bw of the real part of sampled function
idataarray of length 4*bw*bw of the imaginary part of sampled function
rfilterarray of length 4*bw*bw of the real part of sampled filter function
ifilterarray of length 4*bw*bw of the imaginary part of sampled filter function
rresarray of length 4*bw*bw which will contain the real part of result function
iresarray of length 4*bw*bw which will contain the imaginary part of result function
bwbandwidth of problem
workspacespace for computations of size (14*bw*bw)+(26*bw)
Note
We assume that data is real.
This function will do seminaive algorithm for all orders.
If you want to do multiple convolutions, you can freely reuse once allocated workspace. Same for FSTSemiFly(), FZTSemiFly() and InvFSTSemiFly() functions.

Definition at line 454 of file FST_semi_fly.c.

◆ FSTSemiFly()

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:

f(0,0) f(0,1) f(0,2) ... f(0,bw-1)
f(1,1) f(1,2) ... f(1,bw-1)
etc.
f(bw-2,bw-2), f(bw-2,bw-1)
f(bw-1,bw-1)
f(-(bw-1),bw-1)
f(-(bw-2),bw-2), f(-(bw-2),bw-1)
etc.
f(-2,2) ... f(-2,bw-1)
f(-1,1) f(-1,2) ... f(-1,bw-1)

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.

Note
See an example of use in test_s2_semi_fly.c
Parameters
rdataarray of length 4*bw*bw of the real part of the function samples
idataarray of length 4*bw*bw of the imaginary part of the function samples
rcoeffsarray of length bw*bw which will contain the real part of harmonic coefficients in a linearized form
icoeffsarray of length bw*bw which will contain the imaginary part of harmonic coefficients in a linearized form
bwbandwidth of problem
workspacespace for computations of size (10*bw*bw)+(21*bw)
data_formatdetermines the format of the computed data
cutoffdetermines order to switch from seminaive to naive algorithm
DCT_planplan for DLT which is used as argument for call DLTSemi()
FFT_planplan for FFT
weightsarray which is used as argument for call DLTSemi() and DLTNaive()
Note
See more info about DCT_plan and weights in DLTSemi().

Definition at line 67 of file FST_semi_fly.c.

◆ FZTSemiFly()

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.

Parameters
rdataarray of length 4*bw*bw of the real part of the data samples
idataarray of length 4*bw*bw of the imaginary part of the data samples
rresarray of length bw which will contain the real part of harmonic coefficients
iresarray of length bw which will contain the imaginary part of harmonic coefficients
bwbandwidth of problem
workspacespace for computations of size 12*bw
data_formatdetermines the format of the computed data
DCT_planplan for DLT which is used as argument for call DLTSemi()
weightsarray which is used as argument for call DLTSemi()
Note
See more info about DCT_plan and weights in DLTSemi().

Definition at line 390 of file FST_semi_fly.c.

◆ InvFSTSemiFly()

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.

Note
See an example of use in test_s2_semi_fly.c
Parameters
rcoeffsarray of length bw*bw of the real part of harmonic coefficients
icoeffsarray of length bw*bw of the imaginary part of harmonic coefficients
rdataarray of length 4*bw*bw which will contain the real part of transformed result
idataarray of length 4*bw*bw which will contain the imaginary part of transformed result
bwbandwidth of problem
workspacespace for computations of size (10*bw*bw)+(24*bw)
data_formatdetermines the format of the computed data
cutoffdetermines order to switch from seminaive to naive algorithm
inv_DCT_planplan for inverse DLT which is used as argument for call InvDLTSemi()
inv_FFT_planplan for inverse FFT
Note
The real and imaginary part of harmonic coefficients should be stored in rcoeffs and icoeffs before passing to this function in the same order as it was mentioned in FSTSemiFly().
See more info about inv_DCT_plan in InvDLTSemi().

Definition at line 235 of file FST_semi_fly.c.