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

Routines to perform convolutions on the 2-sphere using a combination of seminaive and naive algorithms. More...

#include "s2kit/FST_semi_memo.h"
#include <math.h>
#include <stdlib.h>
#include <string.h>
#include <fftw3.h>
#include "s2kit/cospml.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 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. More...
 
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. More...
 
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. 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. More...
 

Detailed Description

Routines to perform convolutions on the 2-sphere using a combination of seminaive and naive algorithms.

Assumes that all precomputed data is already in memory (e.g. tables are generated).
For descriptions on calling these functions, see the documentation preceding each function.

Definition in file FST_semi_memo.c.

Function Documentation

◆ ConvOn2SphereSemiMemo()

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.

Note
See an example of use in test_conv_semi_memo.c

Conv2Sphere requires memory for spharmonic tables, local workspace and workspace for FSTSemiMemo() (reuses in InvFSTSemiMemo()):

legendre_size = Reduced_Naive_TableSize(bw, cutoff) + Reduced_SpharmonicTableSize(bw,cutoff)
2*legendre_size + (8*(bw*bw)+10*bw) + (4*(bw*bw)+2*bw) = 2*legendre_size + (12*(bw*bw)+12*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 (2*legendre_size)+(12*(bw*bw)+12*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 and generated spharmonic Pml tables. Same for FSTSemiMemo(), FZTSemiMemo() and InvFSTSemiMemo() functions.

Definition at line 439 of file FST_semi_memo.c.

◆ FSTSemiMemo()

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:

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_memo.c and test_s2_semi_memo_fwd.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
seminaive_naive_tablepre-generated spharmonic Pml table
workspacespace for computations of size (8*bw^2)+(7*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
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.
See more info about DCT_plan and weights in DLTSemi().

Definition at line 68 of file FST_semi_memo.c.

◆ FZTSemiMemo()

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.

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
cos_pml_tablepre-generated table of Legendre coefficients of P(0,l) functions
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
cos_pml_table should be generated by GenerateCosPmlTable() with m = 0.
See more info about DCT_plan and weights in DLTSemi().

Definition at line 374 of file FST_semi_memo.c.

◆ InvFSTSemiMemo()

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.

Note
See an example of use in test_s2_semi_memo.c and test_s2_semi_memo_inv.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
transpose_seminaive_naive_tablepre-generated transposed spharmonic Pml table
workspacespace for computations of size (8*bw^2)+(10*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 FSTSemiMemo().
transpose_seminaive_naive_table should be generated by Transpose_Spharmonic_Pml_Table().
See more info about inv_DCT_plan in InvDLTSemi().

Definition at line 228 of file FST_semi_memo.c.

Reduced_Naive_TableSize
int Reduced_Naive_TableSize(const int, const int)
Definition: cospml.c:434
Reduced_SpharmonicTableSize
int Reduced_SpharmonicTableSize(const int, const int)
Definition: cospml.c:107