S2kit  1.1
Toolkit for working with functions defined on the sphere
test_s2_semi_memo_inv.c
Go to the documentation of this file.
1 
58 #include <math.h>
59 #include <stdio.h>
60 #include <stdlib.h>
61 #include <string.h>
62 #include <time.h>
63 
64 #include <fftw3.h>
65 
66 #include "s2kit/FST_semi_memo.h"
67 #include "s2kit/cospml.h"
68 #include "s2kit/weights.h"
69 
70 #include "util/csecond.h"
71 
72 int main(int argc, char** argv) {
73  if (argc < 4) {
74  fprintf(stdout, "Usage: test_s2_semi_memo_inv coeffs_file output_file bw\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
82  DataFormat data_format = COMPLEX;
83 
84  double* workspace = (double*)malloc(sizeof(double) * ((8 * (bw * bw)) + (10 * bw)));
85  double* seminaive_naive_tablespace = (double*)malloc(
86  sizeof(double) * (Reduced_Naive_TableSize(bw, cutoff) + Reduced_SpharmonicTableSize(bw, cutoff)));
87  double* trans_seminaive_naive_tablespace = (double*)malloc(
88  sizeof(double) * (Reduced_Naive_TableSize(bw, cutoff) + Reduced_SpharmonicTableSize(bw, cutoff)));
89 
90  // precompute the Legendres (that's what memo suffix for)
91  fprintf(stdout, "Generating seminaive_naive tables...\n");
92  double** seminaive_naive_table = SemiNaive_Naive_Pml_Table(bw, cutoff, seminaive_naive_tablespace, workspace);
93 
94  fprintf(stdout, "Generating trans_seminaive_naive tables...\n");
95  double** trans_seminaive_naive_table = Transpose_SemiNaive_Naive_Pml_Table(
96  seminaive_naive_table, bw, cutoff, trans_seminaive_naive_tablespace, workspace);
97 
98  double* weights = (double*)malloc(sizeof(double) * 4 * bw);
99  double* rdata = (double*)malloc(sizeof(double) * (size * size));
100  double* idata = (double*)malloc(sizeof(double) * (size * size));
101 
102  // Make inverse DCT plan. Note that I will be using the GURU interface to execute these plans within the routines
103  fftw_plan inv_DCT_plan = fftw_plan_r2r_1d(2 * bw, weights, rdata, FFTW_REDFT01, FFTW_ESTIMATE);
104 
105  // fftw "preamble"
106  // Note that FFT plans assumes that I'm working with a transposed array, e.g. the inputs for a length 2*bw transform
107  // are placed every 2*bw apart, the output will be consecutive entries in the array
108 
109  int rank = 1;
110  fftw_iodim dims[rank];
111  dims[0].n = 2 * bw;
112  dims[0].is = 2 * bw;
113  dims[0].os = 1;
114 
115  int howmany_rank = 1;
116  fftw_iodim howmany_dims[howmany_rank];
117  howmany_dims[0].n = 2 * bw;
118  howmany_dims[0].is = 1;
119  howmany_dims[0].os = 2 * bw;
120 
121  fftw_plan inv_FFT_plan = fftw_plan_guru_split_dft(rank, dims, howmany_rank, howmany_dims, rdata, idata, workspace,
122  workspace + (4 * bw * bw), FFTW_ESTIMATE);
123 
124  GenerateWeightsForDLT(bw, weights);
125 
126  double* rcoeffs = (double*)malloc(sizeof(double) * (bw * bw));
127  double* icoeffs = (double*)malloc(sizeof(double) * (bw * bw));
128 
129  // read coefficients
130  FILE* fp = fopen(argv[1], "r");
131  for (int i = 0; i < bw * bw; ++i) {
132  fscanf(fp, "%lf", rcoeffs + i); // the real part of the coefficient
133  fscanf(fp, "%lf", icoeffs + i); // the imaginary part
134  }
135  fclose(fp);
136 
137  double time_start = csecond();
138  // inverse spherical transform
139  InvFSTSemiMemo(rcoeffs, icoeffs, rdata, idata, bw, trans_seminaive_naive_table, workspace, data_format, cutoff,
140  &inv_DCT_plan, &inv_FFT_plan);
141 
142  fprintf(stderr, "inv time \t = %.4e\n", csecond() - time_start);
143  fprintf(stdout, "about to write out samples\n");
144 
145  fp = fopen(argv[2], "w");
146  for (int i = 0; i < size * size; ++i)
147  fprintf(fp, "%.15f\n%.15f\n", rdata[i], idata[i]);
148  fclose(fp);
149 
150  fprintf(stdout, "finished writing samples\n");
151 
152  fftw_destroy_plan(inv_FFT_plan);
153  fftw_destroy_plan(inv_DCT_plan);
154 
155  free(trans_seminaive_naive_table);
156  free(seminaive_naive_table);
157  free(trans_seminaive_naive_tablespace);
158  free(seminaive_naive_tablespace);
159  free(workspace);
160  free(weights);
161  free(idata);
162  free(rdata);
163  free(icoeffs);
164  free(rcoeffs);
165 
166  return 0;
167 }
Reduced_Naive_TableSize
int Reduced_Naive_TableSize(const int, const int)
Definition: cospml.c:434
weights.h
Transpose_SemiNaive_Naive_Pml_Table
double ** Transpose_SemiNaive_Naive_Pml_Table(double **, const int, const int, double *, double *)
Definition: cospml.c:490
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
main
int main(int argc, char **argv)
Definition: test_s2_semi_memo_inv.c:72
cospml.h
csecond.h
csecond
double csecond()
Definition: csecond.c:32
InvFSTSemiMemo
void InvFSTSemiMemo(double *, double *, double *, double *, const int, double **, double *, DataFormat, const int, fftw_plan *, fftw_plan *)
Computes inverse spherical harmonic transform.
Definition: FST_semi_memo.c:228
FST_semi_memo.h