S2kit  1.1
Toolkit for working with functions defined on the sphere
test_s2_semi_memo_fwd.c
Go to the documentation of this file.
1 
56 #include <math.h>
57 #include <stdio.h>
58 #include <stdlib.h>
59 #include <string.h>
60 #include <time.h>
61 
62 #include <fftw3.h>
63 
64 #include "s2kit/FST_semi_memo.h"
65 #include "s2kit/cospml.h"
66 #include "s2kit/weights.h"
67 #include "s2kit/util.h"
68 
69 #include "util/csecond.h"
70 
71 int main(int argc, char** argv) {
72  if (argc < 4) {
73  fprintf(stdout, "Usage: test_s2_semi_memo_fwd sample_file output_file bw "
74  "[output_format]\n");
75  exit(0);
76  }
77 
78  int bw = atoi(argv[3]);
79 
80  int size = 2 * bw;
81  int cutoff = bw; // seminaive all orders // TODO?
82  DataFormat data_format = COMPLEX;
83 
84  double* workspace = (double*)malloc(sizeof(double) * ((8 * (bw * bw)) + (7 * bw)));
85  double* seminaive_naive_tablespace = (double*)malloc(
86  sizeof(double) * (Reduced_Naive_TableSize(bw, cutoff) + Reduced_SpharmonicTableSize(bw, cutoff)));
87 
88  // precompute the Legendres (that's what memo suffix for)
89  fprintf(stdout, "Generating seminaive_naive tables...\n");
90  double** seminaive_naive_table = SemiNaive_Naive_Pml_Table(bw, cutoff, seminaive_naive_tablespace, workspace);
91 
92  double* weights = (double*)malloc(sizeof(double) * 4 * bw);
93  double* rdata = (double*)malloc(sizeof(double) * (size * size));
94  double* idata = (double*)malloc(sizeof(double) * (size * size));
95 
96  // Make DCT plan. Note that I will be using the GURU interface to execute these plans within the routines
97  fftw_plan DCT_plan = fftw_plan_r2r_1d(2 * bw, weights, rdata, FFTW_REDFT10, FFTW_ESTIMATE);
98 
99  // fftw "preamble"
100  // Note that FFT plan places the output in a transposed array
101  int rank = 1;
102  fftw_iodim dims[rank];
103  dims[0].n = 2 * bw;
104  dims[0].is = 1;
105  dims[0].os = 2 * bw;
106 
107  int howmany_rank = 1;
108  fftw_iodim howmany_dims[howmany_rank];
109  howmany_dims[0].n = 2 * bw;
110  howmany_dims[0].is = 2 * bw;
111  howmany_dims[0].os = 1;
112 
113  fftw_plan FFT_plan = fftw_plan_guru_split_dft(rank, dims, howmany_rank, howmany_dims, rdata, idata, workspace,
114  workspace + (4 * bw * bw), FFTW_ESTIMATE);
115 
116  GenerateWeightsForDLT(bw, weights);
117 
118  // read samples
119  FILE* fp = fopen(argv[1], "r");
120  for (int i = 0; i < size * size; ++i) {
121  fscanf(fp, "%lf", rdata + i); // the real part of the sample
122  fscanf(fp, "%lf", idata + i); // the imaginary part
123  }
124  fclose(fp);
125 
126  double* rcoeffs = (double*)malloc(sizeof(double) * (bw * bw));
127  double* icoeffs = (double*)malloc(sizeof(double) * (bw * bw));
128 
129  double time_start = csecond();
130  // forward spherical transform
131  FSTSemiMemo(rdata, idata, rcoeffs, icoeffs, bw, seminaive_naive_table, workspace, data_format, cutoff, &DCT_plan,
132  &FFT_plan, weights);
133 
134  fprintf(stdout, "forward time \t = %.4e\n", csecond() - time_start);
135  fprintf(stdout, "about to write out coefficients\n");
136 
137  // choose format for writing coefficients
138  int order = 0; // code format
139  if (argc == 5)
140  order = atoi(argv[4]); // human-readable format
141 
142  fp = fopen(argv[2], "w");
143  if (order == 0)
144  for (int i = 0; i < bw * bw; ++i)
145  fprintf(fp, "%.15f\n%.15f\n", rcoeffs[i], icoeffs[i]);
146  else
147  for (int l = 0; l < bw; ++l)
148  for (int m = -l; m < l + 1; ++m) {
149  int index = IndexOfHarmonicCoeff(m, l, bw);
150  fprintf(fp, "l = %d\t m = %d\t %.15f + %.15f I\n", l, m, rcoeffs[index], icoeffs[index]);
151  }
152  fclose(fp);
153 
154  fprintf(stdout, "finished writing coefficients\n");
155 
156  fftw_destroy_plan(FFT_plan);
157  fftw_destroy_plan(DCT_plan);
158 
159  free(workspace);
160  free(seminaive_naive_table);
161  free(seminaive_naive_tablespace);
162  free(weights);
163  free(icoeffs);
164  free(rcoeffs);
165  free(idata);
166  free(rdata);
167 
168  return 0;
169 }
FSTSemiMemo
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.
Definition: FST_semi_memo.c:68
Reduced_Naive_TableSize
int Reduced_Naive_TableSize(const int, const int)
Definition: cospml.c:434
weights.h
IndexOfHarmonicCoeff
int IndexOfHarmonicCoeff(const int, const int, const int)
Gives the position of the coefficient f-hat(m,l) in the one-row array.
Definition: util.c:42
COMPLEX
Definition: util.h:11
GenerateWeightsForDLT
void GenerateWeightsForDLT(const int, double *)
Generates weights for both even and odd order Legendre transforms for a given bandwidth.
Definition: weights.c:32
SemiNaive_Naive_Pml_Table
double ** SemiNaive_Naive_Pml_Table(const int, const int, double *, double *)
Definition: cospml.c:451
Reduced_SpharmonicTableSize
int Reduced_SpharmonicTableSize(const int, const int)
Definition: cospml.c:107
DataFormat
DataFormat
Result data format for FST functions.
Definition: util.h:10
cospml.h
csecond.h
csecond
double csecond()
Definition: csecond.c:32
main
int main(int argc, char **argv)
Definition: test_s2_semi_memo_fwd.c:71
util.h
FST_semi_memo.h