My Project
ssht_sampling.c
Go to the documentation of this file.
1 // SSHT package to perform spin spherical harmonic transforms
2 // Copyright (C) 2011 Jason McEwen
3 // See LICENSE.txt for license details
4 
14 #include <complex.h>
15 #include <math.h>
16 #include "ssht_types.h"
17 
18 #define ABS(x) ((x) >= 0 ? (x) : -(x))
19 
20 void gauleg(double x1, double x2, double *x, double *w, int n);
21 
22 
23 //============================================================================
24 // Sampling weights
25 //============================================================================
26 
27 
37 SSHT_COMPLEX(double) ssht_sampling_weight_mw(int p) {
38 
39  if (p == 1) {
40  return I * SSHT_PION2;
41  }
42  else if (p == -1) {
43  return - I * SSHT_PION2;
44  }
45  else if (p % 2 == 0) {
46  return 2.0 / (1.0 - p*p);
47  }
48  else {
49  return 0.0;
50  }
51 
52 }
53 
54 
65 double ssht_sampling_weight_dh(double theta_t, int L) {
66 
67  double w;
68  int k;
69 
70  w = 0.0;
71  for (k=0; k<=L-1; k++)
72  w += sin((2.0*k+1.0)*theta_t) / (double)(2.0*k+1.0);
73 
74  w *= 2.0/((double) L) * sin(theta_t);
75 
76  return w;
77 
78 }
79 
80 
94 void ssht_sampling_gl_thetas_weights(double *thetas, double *weights, int L) {
95 
96  int t;
97  double tmp;
98 
99  gauleg(-1.0, 1.0, thetas, weights, L);
100 
101  for (t=0; t<=L-1; t++)
102  thetas[t] = acos(thetas[t]);
103 
104  // Reorder thetas array.
105  for (t=0; t<=(L-1)/2; t++) {
106  tmp = thetas[t];
107  thetas[t] = thetas[L-1-t];
108  thetas[L-1-t] = tmp;
109  }
110 
111 }
112 
113 
129 void gauleg(double x1, double x2, double *x, double *w, int n) {
130 
131  double EPS = 1e-14;
132  int i,j,m;
133  double p1,p2,p3,pp,xl,xm,z,z1;
134 
135  m=(n+1)/2;
136  xm=0.5*(x2+x1);
137  xl=0.5*(x2-x1);
138 
139  for (i=1; i<=m; i++) {
140  z=cos(3.141592654*(i-.25)/(n+.5));
141  do {
142  p1=1.0;
143  p2=0.0;
144  for (j=1; j<=n; j++) {
145  p3=p2;
146  p2=p1;
147  p1=((2.0*j-1.0)*z*p2-(j-1.0)*p3)/j;
148  }
149  pp=n*(z*p1-p2)/(z*z-1.0);
150  z1=z;
151  z=z1-p1/pp;
152  } while (ABS(z-z1) > EPS);
153  x[i-1]=xm-xl*z;
154  x[n+1-i-1]=xm+xl*z;
155  w[i-1]=2.0*xl/((1.0-z*z)*pp*pp);
156  w[n+1-i-1]=w[i-1];
157  }
158 
159 }
160 
161 
162 //============================================================================
163 // Sampling relations
164 //============================================================================
165 
166 
179 
180  return L;
181 
182 }
183 
184 
197 double ssht_sampling_mw_t2theta(int t, int L) {
198 
199  return (2.0*t + 1.0) * SSHT_PI / (2.0*L - 1.0);
200 
201 }
202 
203 
213 
214  return 2*L-1;
215 
216 }
217 
218 
231 double ssht_sampling_mw_p2phi(int p, int L) {
232 
233  return 2.0 * p * SSHT_PI / (2.0*L - 1.0);
234 
235 }
236 
237 
249 int ssht_sampling_mw_n(int L) {
250 
251  return (L-1)*(2*L-1) + 1;
252 
253 }
254 
255 
269 double ssht_sampling_mw_ss_t2theta(int t, int L) {
270 
271  return 2.0 * t * SSHT_PI / (2.0 * L);
272 
273 }
274 
275 
289 
290  return L+1;
291 
292 }
293 
294 
307 double ssht_sampling_mw_ss_p2phi(int p, int L) {
308 
309  return 2.0 * p * SSHT_PI / (2.0*L);
310 
311 }
312 
313 
324 
325  return 2*L;
326 
327 }
328 
329 
343 
344  return (L-1)*(2*L) + 2;
345 
346 }
347 
348 
361 double ssht_sampling_dh_t2theta(int t, int L) {
362 
363  return (2.0*t + 1.0) * SSHT_PI / (4.0*L);
364 
365 }
366 
367 
377 
378  return 2*L;
379 
380 }
381 
382 
395 double ssht_sampling_dh_p2phi(int p, int L) {
396 
397  return 2.0 * p * SSHT_PI / (2.0*L - 1.0);
398 
399 }
400 
401 
411 
412  return 2*L-1;
413 
414 }
415 
416 
425 int ssht_sampling_dh_n(int L) {
426 
427  return 2*L*(2*L-1);
428 
429 }
430 
431 
441 
442  return L;
443 
444 }
445 
446 
459 double ssht_sampling_gl_p2phi(int p, int L) {
460 
461  return 2.0 * p * SSHT_PI / (2.0*L - 1.0);
462 
463 }
464 
465 
475 
476  return 2*L-1;
477 
478 }
479 
480 
489 int ssht_sampling_gl_n(int L) {
490 
491  return L*(2*L-1);
492 
493 }
494 
495 
496 
gauleg
void gauleg(double x1, double x2, double *x, double *w, int n)
Definition: ssht_sampling.c:129
SSHT_PION2
#define SSHT_PION2
Definition: ssht_types.h:39
ssht_sampling_dh_n
int ssht_sampling_dh_n(int L)
Definition: ssht_sampling.c:425
SSHT_PI
#define SSHT_PI
Definition: ssht_types.h:38
SSHT_COMPLEX
SSHT_COMPLEX(double)
Definition: ssht_sampling.c:37
ssht_sampling_gl_n
int ssht_sampling_gl_n(int L)
Definition: ssht_sampling.c:489
ssht_sampling_mw_t2theta
double ssht_sampling_mw_t2theta(int t, int L)
Definition: ssht_sampling.c:197
ssht_sampling_gl_nphi
int ssht_sampling_gl_nphi(int L)
Definition: ssht_sampling.c:474
ssht_sampling_gl_p2phi
double ssht_sampling_gl_p2phi(int p, int L)
Definition: ssht_sampling.c:459
ssht_sampling_dh_t2theta
double ssht_sampling_dh_t2theta(int t, int L)
Definition: ssht_sampling.c:361
ssht_sampling_mw_ntheta
int ssht_sampling_mw_ntheta(int L)
Definition: ssht_sampling.c:178
ssht_sampling_mw_p2phi
double ssht_sampling_mw_p2phi(int p, int L)
Definition: ssht_sampling.c:231
ssht_sampling_mw_ss_ntheta
int ssht_sampling_mw_ss_ntheta(int L)
Definition: ssht_sampling.c:288
ssht_sampling_mw_n
int ssht_sampling_mw_n(int L)
Definition: ssht_sampling.c:249
ABS
#define ABS(x)
Definition: ssht_sampling.c:18
ssht_sampling_mw_nphi
int ssht_sampling_mw_nphi(int L)
Definition: ssht_sampling.c:212
ssht_sampling_weight_dh
double ssht_sampling_weight_dh(double theta_t, int L)
Definition: ssht_sampling.c:65
ssht_sampling_gl_thetas_weights
void ssht_sampling_gl_thetas_weights(double *thetas, double *weights, int L)
Definition: ssht_sampling.c:94
ssht_sampling_mw_ss_p2phi
double ssht_sampling_mw_ss_p2phi(int p, int L)
Definition: ssht_sampling.c:307
ssht_sampling_gl_ntheta
int ssht_sampling_gl_ntheta(int L)
Definition: ssht_sampling.c:440
ssht_sampling_mw_ss_t2theta
double ssht_sampling_mw_ss_t2theta(int t, int L)
Definition: ssht_sampling.c:269
ssht_sampling_mw_ss_nphi
int ssht_sampling_mw_ss_nphi(int L)
Definition: ssht_sampling.c:323
ssht_sampling_dh_ntheta
int ssht_sampling_dh_ntheta(int L)
Definition: ssht_sampling.c:376
ssht_sampling_mw_ss_n
int ssht_sampling_mw_ss_n(int L)
Definition: ssht_sampling.c:342
ssht_sampling_dh_p2phi
double ssht_sampling_dh_p2phi(int p, int L)
Definition: ssht_sampling.c:395
ssht_sampling_dh_nphi
int ssht_sampling_dh_nphi(int L)
Definition: ssht_sampling.c:410
ssht_types.h