My Project
so3_test_csv.c
Go to the documentation of this file.
1 // S03 package to perform Wigner transform on the rotation group SO(3)
2 // Copyright (C) 2013 Martin Büttner and Jason McEwen
3 // See LICENSE.txt for license details
4 
39 #include <stdio.h>
40 #include <stdlib.h>
41 #include <complex.h>
42 #include <time.h>
43 
44 #include <so3.h>
45 
46 #define NREPEAT 10
47 #define MIN(a,b) ((a < b) ? (a) : (b))
48 #define MAX(a,b) ((a > b) ? (a) : (b))
49 
50 double get_max_error(complex double *expected, complex double *actual, int n);
51 double ran2_dp(int idum);
52 void so3_test_gen_flmn_complex(complex double *flmn, const so3_parameters_t *parameters, int seed);
53 void so3_test_gen_flmn_real(complex double *flmn, const so3_parameters_t *parameters, int seed);
54 
55 int main(int argc, char **argv)
56 {
57  so3_parameters_t parameters = {};
58  int Lmax, L, useLasN, N, L0;
59  complex double *flmn_orig, *flmn_syn;
60  complex double *f;
61  double *f_real;
62  int seed;
63  clock_t time_start, time_end;
64  int i, real;
65  int flmn_size;
66 
67  const char *n_order_str[SO3_N_ORDER_SIZE];
68  const char *storage_mode_str[SO3_STORAGE_SIZE];
69  const char *n_mode_str[SO3_N_MODE_SIZE];
70  const char *reality_str[2];
71  const char *sampling_str[2];
72 
73  double min_duration_inverse;
74  double min_duration_forward;
75  double avg_duration_inverse;
76  double avg_duration_forward;
77  double avg_error;
78 
79  n_order_str[SO3_N_ORDER_ZERO_FIRST] = "n = 0 first";
80  n_order_str[SO3_N_ORDER_NEGATIVE_FIRST] = "n = -N+1 first";
81 
82  storage_mode_str[SO3_STORAGE_PADDED] = "padded storage";
83  storage_mode_str[SO3_STORAGE_COMPACT] = "compact storage";
84 
85  n_mode_str[SO3_N_MODE_ALL] = "all n";
86  n_mode_str[SO3_N_MODE_EVEN] = "only even n";
87  n_mode_str[SO3_N_MODE_ODD] = "only odd n";
88  n_mode_str[SO3_N_MODE_MAXIMUM] = "only |n| = N-1";
89 
90  reality_str[0] = "complex";
91  reality_str[1] = "real";
92 
93  sampling_str[0] = "MW";
94  sampling_str[1] = "MWSS";
95 
96  // Parse command line arguments
97  N = Lmax = 64;
98  useLasN = 1; // true
99  L0 = 0;
100  seed = 1;
101  if (argc > 1)
102  Lmax = atoi(argv[1]);
103  if (argc > 2)
104  L0 = atoi(argv[2]);
105  if (argc > 3)
106  {
107  useLasN = 0; // false
108  N = atoi(argv[3]);
109  parameters.N = N;
110  }
111 
112  parameters.L0 = L0;
113  parameters.verbosity = 0;
114 
115  // (2*N-1)*L*L is the largest number of flmn ever needed. For more
116  // compact storage modes, only part of the memory will be used.
117  flmn_orig = malloc((2*N-1)*Lmax*Lmax * sizeof *flmn_orig);
118  SO3_ERROR_MEM_ALLOC_CHECK(flmn_orig);
119  flmn_syn = malloc((2*N-1)*Lmax*Lmax * sizeof *flmn_syn);
120  SO3_ERROR_MEM_ALLOC_CHECK(flmn_syn);
121 
122  // We only need (2*L) * (L+1) * (2*N-1) samples for MW symmetric sampling.
123  // For the usual MW sampling, only part of the memory will be used.
124  f = malloc((2*Lmax)*(Lmax+1)*(2*N-1) * sizeof *f);
125  SO3_ERROR_MEM_ALLOC_CHECK(f);
126  f_real = malloc((2*Lmax)*(Lmax+1)*(2*N-1) * sizeof *f_real);
127  SO3_ERROR_MEM_ALLOC_CHECK(f_real);
128 
129  // Output header row
130  printf("reality;L;N;avg_duration_inverse;avg_duration_forward;min_duration_inverse;min_duration_forward;avg_error\n");
131 
132  parameters.sampling_scheme = SO3_SAMPLING_MW;
133  parameters.n_order = SO3_N_ORDER_ZERO_FIRST;
134  parameters.storage = SO3_STORAGE_PADDED;
135  parameters.n_mode = SO3_N_MODE_ALL;
136 
137  // real == 0 --> complex signal
138  // real == 1 --> real signal
139  for (real = 0; real < 2; ++real)
140  {
141  parameters.reality = real;
142 
143  L = 1;
144  while(L <= Lmax)
145  {
146  if (L <= L0 || (!useLasN && L < N))
147  {
148  L *= 2;
149  continue;
150  }
151 
152  parameters.L = L;
153  if (useLasN)
154  {
155  N = L;
156  parameters.N = L;
157  }
158 
159  min_duration_inverse = 0.0;
160  min_duration_forward = 0.0;
161  avg_duration_inverse = 0.0;
162  avg_duration_forward = 0.0;
163  avg_error = 0.0;
164 
165  for (i = 0; i < NREPEAT; ++i)
166  {
167  int j;
168  double duration;
169 
170  // Reset output array
171  for (j = 0; j < (2*N-1)*L*L; ++j)
172  flmn_syn[j] = 0.0;
173 
174  if (real) so3_test_gen_flmn_real(flmn_orig, &parameters, seed);
175  else so3_test_gen_flmn_complex(flmn_orig, &parameters, seed);
176 
177  time_start = clock();
178  if (real) so3_core_inverse_via_ssht_real(f_real, flmn_orig, &parameters);
179  else so3_core_inverse_via_ssht(f, flmn_orig, &parameters);
180  time_end = clock();
181 
182  duration = (time_end - time_start) / (double)CLOCKS_PER_SEC;
183  avg_duration_inverse = avg_duration_inverse + duration / NREPEAT;
184  if (!i || duration < min_duration_inverse)
185  min_duration_inverse = duration;
186 
187  time_start = clock();
188  if (real) so3_core_forward_via_ssht_real(flmn_syn, f_real, &parameters);
189  else so3_core_forward_via_ssht(flmn_syn, f, &parameters);
190  time_end = clock();
191 
192  duration = (time_end - time_start) / (double)CLOCKS_PER_SEC;
193  avg_duration_forward = avg_duration_forward + duration / NREPEAT;
194  if (!i || duration < min_duration_forward)
195  min_duration_forward = duration;
196 
197  flmn_size = so3_sampling_flmn_size(&parameters);
198  avg_error += get_max_error(flmn_orig, flmn_syn, flmn_size)/NREPEAT;
199  }
200 
201  printf("%d;%d;%d;%f;%f;%f;%f;%e\n",
202  real,
203  L,
204  N,
205  avg_duration_inverse,
206  avg_duration_forward,
207  min_duration_inverse,
208  min_duration_forward,
209  avg_error);
210 
211  L *= 2;
212  }
213  }
214 
215  free(flmn_orig);
216  free(flmn_syn);
217  free(f);
218  free(f_real);
219 
220  return 0;
221 }
222 
223 double get_max_error(complex double *expected, complex double *actual, int n)
224 {
225  int i;
226  double error, maxError = 0;
227 
228  for (i = 0; i < n; ++i)
229  {
230  error = cabs(expected[i] - actual[i]);
231  maxError = MAX(error, maxError);
232  }
233 
234  return maxError;
235 }
236 
254  complex double *flmn,
255  const so3_parameters_t *parameters,
256  int seed)
257 {
258  int L0, L, N;
259  int i, el, m, n, n_start, n_stop, n_inc, ind;
260 
261  L0 = parameters->L0;
262  L = parameters->L;
263  N = parameters->N;
264 
265  for (i = 0; i < (2*N-1)*L*L; ++i)
266  flmn[i] = 0.0;
267 
268  switch (parameters->n_mode)
269  {
270  case SO3_N_MODE_ALL:
271  n_start = -N+1;
272  n_stop = N-1;
273  n_inc = 1;
274  break;
275  case SO3_N_MODE_EVEN:
276  n_start = ((N-1) % 2 == 0) ? -N+1 : -N+2;
277  n_stop = ((N-1) % 2 == 0) ? N-1 : N-2;
278  n_inc = 2;
279  break;
280  case SO3_N_MODE_ODD:
281  n_start = ((N-1) % 2 != 0) ? -N+1 : -N+2;
282  n_stop = ((N-1) % 2 != 0) ? N-1 : N-2;
283  n_inc = 2;
284  break;
285  case SO3_N_MODE_MAXIMUM:
286  n_start = -N+1;
287  n_stop = N-1;
288  n_inc = 2*N - 2;
289  break;
290  default:
291  SO3_ERROR_GENERIC("Invalid n-mode.");
292  }
293 
294  for (n = n_start; n <= n_stop; n += n_inc)
295  {
296  for (el = MAX(L0, abs(n)); el < L; ++el)
297  {
298  for (m = -el; m <= el; ++m)
299  {
300  so3_sampling_elmn2ind(&ind, el, m, n, parameters);
301  flmn[ind] = (2.0*ran2_dp(seed) - 1.0) + I * (2.0*ran2_dp(seed) - 1.0);
302  }
303  }
304  }
305 }
306 
327  complex double *flmn,
328  const so3_parameters_t *parameters,
329  int seed)
330 {
331  int L0, L, N;
332  int i, el, m, n, n_start, n_stop, n_inc, ind;
333  double real, imag;
334 
335  L0 = parameters->L0;
336  L = parameters->L;
337  N = parameters->N;
338 
339  for (i = 0; i < (2*N-1)*L*L; ++i)
340  flmn[i] = 0.0;
341 
342  switch (parameters->n_mode)
343  {
344  case SO3_N_MODE_ALL:
345  n_start = 0;
346  n_stop = N-1;
347  n_inc = 1;
348  break;
349  case SO3_N_MODE_EVEN:
350  n_start = 0;
351  n_stop = ((N-1) % 2 == 0) ? N-1 : N-2;
352  n_inc = 2;
353  break;
354  case SO3_N_MODE_ODD:
355  n_start = 1;
356  n_stop = ((N-1) % 2 != 0) ? N-1 : N-2;
357  n_inc = 2;
358  break;
359  case SO3_N_MODE_MAXIMUM:
360  n_start = N-1;
361  n_stop = N-1;
362  n_inc = 1;
363  break;
364  default:
365  SO3_ERROR_GENERIC("Invalid n-mode.");
366  }
367 
368  for (n = n_start; n <= n_stop; n += n_inc)
369  {
370  if (n == 0)
371  {
372  // Fill fl00 with random real values
373  for (el = L0; el < L; ++el)
374  {
375  so3_sampling_elmn2ind_real(&ind, el, 0, 0, parameters);
376  flmn[ind] = (2.0*ran2_dp(seed) - 1.0);
377  }
378 
379  // Fill fl+-m0 with conjugated random values
380 
381  for (el = L0; el < L; ++el)
382  {
383  for (m = 1; m <= el; ++m)
384  {
385  real = (2.0*ran2_dp(seed) - 1.0);
386  imag = (2.0*ran2_dp(seed) - 1.0);
387  so3_sampling_elmn2ind_real(&ind, el, m, 0, parameters);
388  flmn[ind] = real + imag * I;
389  so3_sampling_elmn2ind_real(&ind, el, -m, 0, parameters);
390  flmn[ind] = real - imag * I;
391  if (m % 2)
392  flmn[ind] = - real + imag * I;
393  else
394  flmn[ind] = real - imag * I;
395  }
396  }
397  }
398  else
399  {
400  for (el = MAX(L0, n); el < L; ++el)
401  {
402  for (m = -el; m <= el; ++m)
403  {
404 
405  so3_sampling_elmn2ind_real(&ind, el, m, n, parameters);
406  flmn[ind] = (2.0*ran2_dp(seed) - 1.0) + I * (2.0*ran2_dp(seed) - 1.0);
407  }
408  }
409  }
410  }
411 }
412 
426 double ran2_dp(int idum)
427 {
428  int IM1=2147483563,IM2=2147483399,IMM1=IM1-1,
429  IA1=40014,IA2=40692,IQ1=53668,IQ2=52774,IR1=12211,IR2=3791,
430  NTAB=32,NDIV=1+IMM1/NTAB;
431 
432  double AM=1./IM1,EPS=1.2e-7,RNMX=1.-EPS;
433  int j,k;
434  static int iv[32],iy,idum2 = 123456789;
435  // N.B. in C static variables are initialised to 0 by default.
436 
437  if (idum <= 0) {
438  idum= (-idum>1 ? -idum : 1); // max(-idum,1);
439  idum2=idum;
440  for(j=NTAB+8;j>=1;j--) {
441  k=idum/IQ1;
442  idum=IA1*(idum-k*IQ1)-k*IR1;
443  if (idum < 0) idum=idum+IM1;
444  if (j < NTAB) iv[j-1]=idum;
445  }
446  iy=iv[0];
447  }
448  k=idum/IQ1;
449  idum=IA1*(idum-k*IQ1)-k*IR1;
450  if (idum < 0) idum=idum+IM1;
451  k=idum2/IQ2;
452  idum2=IA2*(idum2-k*IQ2)-k*IR2;
453  if (idum2 < 0) idum2=idum2+IM2;
454  j=1+iy/NDIV;
455  iy=iv[j-1]-idum2;
456  iv[j-1]=idum;
457  if(iy < 1)iy=iy+IMM1;
458  return (AM*iy < RNMX ? AM*iy : RNMX); // min(AM*iy,RNMX);
459 }
so3_test_gen_flmn_real
void so3_test_gen_flmn_real(complex double *flmn, const so3_parameters_t *parameters, int seed)
Definition: so3_test_csv.c:326
so3_sampling_elmn2ind_real
void so3_sampling_elmn2ind_real(int *ind, int el, int m, int n, const so3_parameters_t *parameters)
Definition: so3_sampling.c:569
so3_core_inverse_via_ssht
void so3_core_inverse_via_ssht(complex double *f, const complex double *flmn, const so3_parameters_t *parameters)
Definition: so3_core.c:48
MAX
#define MAX(a, b)
Definition: so3_test_csv.c:48
so3_test_gen_flmn_complex
void so3_test_gen_flmn_complex(complex double *flmn, const so3_parameters_t *parameters, int seed)
Definition: so3_test_csv.c:253
NREPEAT
#define NREPEAT
Definition: so3_test_csv.c:46
so3_core_forward_via_ssht_real
void so3_core_forward_via_ssht_real(complex double *flmn, const double *f, const so3_parameters_t *parameters)
Definition: so3_core.c:700
so3_core_forward_via_ssht
void so3_core_forward_via_ssht(complex double *flmn, const complex double *f, const so3_parameters_t *parameters)
Definition: so3_core.c:253
so3_sampling_flmn_size
int so3_sampling_flmn_size(const so3_parameters_t *parameters)
Definition: so3_sampling.c:314
so3_core_inverse_via_ssht_real
void so3_core_inverse_via_ssht_real(double *f, const complex double *flmn, const so3_parameters_t *parameters)
Definition: so3_core.c:475
get_max_error
double get_max_error(complex double *expected, complex double *actual, int n)
Definition: so3_test_csv.c:223
ran2_dp
double ran2_dp(int idum)
Definition: so3_test_csv.c:426
main
int main(int argc, char **argv)
Definition: so3_test_csv.c:55
so3_sampling_elmn2ind
void so3_sampling_elmn2ind(int *ind, int el, int m, int n, const so3_parameters_t *parameters)
Definition: so3_sampling.c:367