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

Source code for seminaive discrete Legendre transform. More...

#include "s2kit/seminaive.h"
#include <math.h>
#include <stdio.h>
#include <string.h>
#include <fftw3.h>
#include "s2kit/cospml.h"

Go to the source code of this file.

Functions

void InvDLTSemi (double *coeffs, const int bw, const int m, double *result, double *trans_cos_pml_table, double *sin_values, double *workspace, fftw_plan *plan)
 Computes the inverse Legendre transform using the transposed seminaive algorithm. More...
 
void DLTSemi (double *data, const int bw, const int m, double *result, double *workspace, double *cos_pml_table, double *weights, fftw_plan *plan)
 Computes the Legendre transform of data using the seminaive algorithm. More...
 

Detailed Description

Source code for seminaive discrete Legendre transform.

Source code for computing the Legendre transform where projections are carried out in cosine space, i.e., the seminaive algorithm. For a description, see the related paper or Sean's thesis. // TODO link

Definition in file seminaive.c.

Function Documentation

◆ DLTSemi()

void DLTSemi ( double *  data,
const int  bw,
const int  m,
double *  result,
double *  workspace,
double *  cos_pml_table,
double *  weights,
fftw_plan *  plan 
)

Computes the Legendre transform of data using the seminaive algorithm.

Note
See an example of use in test_DLT_semi.c

The forward transform looks like:

l = PCWf

where f is the data vector, W is a quadrature matrix, C is a cosine transform matrix, P is a matrix full of coefficients of the cosine series representation of each Pml function P(m,m) P(m,m+1) ... P(m,bw-1), and l is the (associated) Legendre series representation of f.

Parameters
dataarray of size 2*bw containing function to be transformed, assumes sampling at Chebyshev nodes
bwbandwidth of the problem
morder of the problem
resultarray of length bw for returning computed Legendre coefficients, contains bw-m coeffs, with the <f,P(m,m)> coefficient located in result[0]
workspacespace for computations of size 4*bw
cos_pml_tablearray containing the cosine series coefficients of the Pmls (or Gmls) for this problem
weightsarray of size 4*bw which holds the weights for both even (starting at weights[0]) and odd (weights[2*bw]) transforms
planpointer to fftw_plan with input array being weighted_data and output being cos_data
Note
This function was written to be a part of the full spherical harmonic transform, so a lot of precomputation has been assumed.
cos_pml_table should be generated with GenerateCosPmlTable(). The offset for a particular Pml can be found with TableOffset(), the size of the table - with the TableSize(). Since the cosine series are always zero-striped, the zeroes have been removed.
weights should be generated by GenerateWeightsForDLT().
Warning
We use the guru interface to execute inverse DFT so we can apply the same plan to different arrays. Thus, the plan should be:
fftw_plan_r2r_1d(2*bw, weighted_data, cos_data, FFTW_REDFT10, FFTW_ESTIMATE)

Definition at line 153 of file seminaive.c.

◆ InvDLTSemi()

void InvDLTSemi ( double *  coeffs,
const int  bw,
const int  m,
double *  result,
double *  trans_cos_pml_table,
double *  sin_values,
double *  workspace,
fftw_plan *  plan 
)

Computes the inverse Legendre transform using the transposed seminaive algorithm.

Note
See an example of use in test_DLT_semi.c

Because the Legendre transform is orthogonal, the inverse can be computed by transposing the matrix formulation of the problem.

(see the description of DLTSemi() for forward Legendre transform)

Then to do the inverse:

f = trans(C) trans(P) l

so we need to transpose the matrix P from the forward transform and then do a cosine series evaluation. No quadrature matrix is necessary. If order m is odd, then there is also a sin factor that needs to be taken into account.

Parameters
coeffsarray of length bw-m containing associated Legendre series coefficients, assumed that first entry contains the P(m,m) coefficient
bwproblem bandwidth
morder of the associated Legendre functions
resultarray of 2*bw samples representing the evaluation of the Legendre series at 2*bw Chebyshev nodes
trans_cos_pml_tablearray of the linearized form of trans(P) from the description of this function
sin_valuesarray is used when m is odd in generation of the values in trans(P)
workspacespace for computations of size 2*bw
planpointer to fftw_plan with input array being fcos and output being result
Note
This function was written to be a part of the full spherical harmonic transform, so a lot of precomputation has been assumed.
trans_cos_pml_table should be generated by TransposeCosPmlTable().
Warning
We use the guru interface to execute inverse DFT so we can apply the same plan to different arrays. Thus, the plan should be:
fftw_plan_r2r_1d(2*bw, fcos, result, FFTW_REDFT01, FFTW_ESTIMATE)

Definition at line 56 of file seminaive.c.