35 #define DEFAULT_NFFT_CUTOFF 6
36 #define FPT_THRESHOLD 1000
38 static fpt_set SO3_fpt_init(
int l,
unsigned int flags,
int kappa);
42 nfsoft_init_advanced(plan, N, M, NFSOFT_MALLOC_X | NFSOFT_MALLOC_F
43 | NFSOFT_MALLOC_F_HAT);
46 void nfsoft_init_advanced(
nfsoft_plan *plan,
int N,
int M,
47 unsigned int nfsoft_flags)
49 nfsoft_init_guru(plan, N, M, nfsoft_flags, PRE_PHI_HUT | PRE_PSI | MALLOC_X
50 | MALLOC_F_HAT | MALLOC_F | FFTW_INIT | FFT_OUT_OF_PLACE,
51 DEFAULT_NFFT_CUTOFF, FPT_THRESHOLD);
54 void nfsoft_init_guru(
nfsoft_plan *plan,
int B,
int M,
55 unsigned int nfsoft_flags,
unsigned int nfft_flags,
int nfft_cutoff,
69 nfft_init_guru(&plan->
p_nfft, 3, N, M, n, nfft_cutoff, nfft_flags,
70 FFTW_ESTIMATE | FFTW_DESTROY_INPUT);
72 if ((plan->
p_nfft).nfft_flags & PRE_LIN_PSI)
79 plan->fpt_kappa = fpt_kappa;
80 plan->
flags = nfsoft_flags;
82 if (plan->
flags & NFSOFT_MALLOC_F_HAT)
85 if (plan->
f_hat == NULL ) printf(
"Allocation failed!\n");
88 if (plan->
flags & NFSOFT_MALLOC_X)
91 if (plan->
x == NULL ) printf(
"Allocation failed!\n");
93 if (plan->
flags & NFSOFT_MALLOC_F)
96 if (plan->
f == NULL ) printf(
"Allocation failed!\n");
103 if (plan->
wig_coeffs == NULL ) printf(
"Allocation failed!\n");
104 if (plan->
cheby == NULL ) printf(
"Allocation failed!\n");
105 if (plan->
aux == NULL ) printf(
"Allocation failed!\n");
107 plan->
mv_trafo = (void (*) (
void* ))nfsoft_trafo;
108 plan->
mv_adjoint = (void (*) (
void* ))nfsoft_adjoint;
123 my_plan->
cheby[0]=0.0;
125 for (j=1;j<my_plan->
N_total+1;j++)
134 aux[j]=my_plan->
cheby[j];
141 my_plan->
cheby[0]=(C) (-1.)/(2.*_Complex_I) * aux[1];
144 my_plan->
cheby[j]=(1./(2.*_Complex_I)*(aux[j+1]-aux[j-1]));
153 static fpt_set SO3_fpt_init(
int l,
unsigned int flags,
int kappa)
156 int N, t, k_start, k_end, k, m;
161 if (flags & NFSOFT_USE_DPT)
168 t = (int) log2(X(next_power_of_2)(N));
177 N = X(next_power_of_2)(l);
188 if (flags & NFSOFT_NO_STABILIZATION)
190 set = fpt_init((2* N + 1) * (2* N + 1), t, 0U | FPT_NO_STABILIZATION);
194 set = fpt_init((2* N + 1) * (2* N + 1), t, 0U);
197 for (k = -N; k <= N; k++)
198 for (m = -N; m <= N; m++)
201 k_start = (ABS(k) >= ABS(m)) ? ABS(k) : ABS(m);
208 fpt_precompute(set, glo, alpha, beta, gamma, k_start, kappa);
222 static fpt_set SO3_single_fpt_init(
int l,
int k,
int m,
unsigned int flags,
int kappa)
224 int N, t, k_start, k_end;
229 if (flags & NFSOFT_USE_DPT)
236 t = (int) log2(X(next_power_of_2)(N));
245 N = X(next_power_of_2)(l);
257 unsigned int fptflags = 0U
258 | IF(flags & NFSOFT_USE_DPT,FPT_NO_FAST_ALGORITHM,IF(t > 1,FPT_NO_DIRECT_ALGORITHM,0U))
259 | IF(flags & NFSOFT_NO_STABILIZATION,FPT_NO_STABILIZATION,0U);
260 set = fpt_init(1, t, fptflags);
264 k_start = (ABS(k) >= ABS(m)) ? ABS(k) : ABS(m);
277 fpt_precompute(set, 0, alpha, beta, gamma, k_start, kappa);
289 void SO3_fpt(C *coeffs,
fpt_set set,
int l,
int k,
int m,
unsigned int flags)
298 int k_start, k_end, j;
299 int function_values = 0;
302 if (flags & NFSOFT_USE_DPT)
313 N = X(next_power_of_2)(l);
317 k_start = (ABS(k) >= ABS(m)) ? ABS(k) : ABS(m);
319 trafo_nr = (N + k) * (2* N + 1) + (m + N);
324 for (j = 0; j <= k_end; j++)
328 for (j = 0; j <= l - k_start; j++)
330 x[j + k_start] = coeffs[j];
332 for (j = l - k_start + 1; j <= k_end - k_start; j++)
334 x[j + k_start] = K(0.0);
340 if (flags & NFSOFT_USE_DPT)
342 fpt_trafo_direct(set, trafo_nr, &x[k_start], y, k_end, 0U
343 | (function_values ? FPT_FUNCTION_VALUES : 0U));
347 fpt_trafo(set, trafo_nr, &x[k_start], y, k_end, 0U
348 | (function_values ? FPT_FUNCTION_VALUES : 0U));
352 for (j = 0; j <= l; j++)
365 void SO3_fpt_transposed(C *coeffs,
fpt_set set,
int l,
int k,
int m,
368 int N, k_start, k_end, j;
370 int function_values = 0;
378 if (flags & NFSOFT_USE_DPT)
389 N = X(next_power_of_2)(l);
393 k_start = (ABS(k) >= ABS(m)) ? ABS(k) : ABS(m);
395 trafo_nr = (N + k) * (2* N + 1) + (m + N);
402 for (j = 0; j <= l; j++)
406 for (j = l + 1; j <= k_end; j++)
411 if (flags & NFSOFT_USE_DPT)
413 fpt_transposed_direct(set, trafo_nr, &x[k_start], y, k_end, 0U
414 | (function_values ? FPT_FUNCTION_VALUES : 0U));
418 fpt_transposed(set, trafo_nr, &x[k_start], y, k_end, 0U
419 | (function_values ? FPT_FUNCTION_VALUES : 0U));
422 for (j = 0; j <= l; j++)
442 for (j = 0; j < M; j++)
444 plan3D->
p_nfft.
x[3* j ] = plan3D->
x[3* j + 2];
445 plan3D->
p_nfft.
x[3* j + 1] = plan3D->
x[3* j ];
446 plan3D->
p_nfft.
x[3* j + 2] = plan3D->
x[3* j + 1];
454 if ((plan3D->
p_nfft).nfft_flags & FG_PSI)
456 nfft_precompute_one_psi(&(plan3D->
p_nfft));
458 if ((plan3D->
p_nfft).nfft_flags & PRE_PSI)
460 nfft_precompute_one_psi(&(plan3D->
p_nfft));
467 int i, j, m, k, max, glo1, glo2;
479 for (j = 0; j < M; j++)
480 plan3D->
f[j] = plan3D->
f_hat[0];
487 for (k = -N; k <= N; k++)
489 for (m = -N; m <= N; m++)
492 max = (ABS(m) > ABS(k) ? ABS(m) : ABS(k));
494 for (j = 0; j <= N - max; j++)
498 if ((plan3D->
flags & NFSOFT_NORMALIZED))
501 * SQRT(0.5 * (2. * (max + j) + 1.));
504 if ((plan3D->
flags & NFSOFT_REPRESENT))
506 if ((k < 0) && (k % 2))
510 if ((m < 0) && (m % 2))
521 for (j = N - max + 1; j < X(next_power_of_2)(N) + 1; j++)
526 c2e(plan3D, ABS((k + m) % 2));
528 for (i = 1; i <= 2* plan3D ->
N_total + 2; i++)
530 plan3D->
p_nfft.
f_hat[NFSOFT_INDEX(k, m, i - N - 1, N) - 1]
531 = plan3D->
cheby[i - 1];
539 if (plan3D->
flags & NFSOFT_USE_NDFT)
541 nfft_trafo_direct(&(plan3D->
p_nfft));
548 for (j = 0; j < plan3D->
M_total; j++)
566 my_plan->
aux[0]= 1/(2*_Complex_I)*my_plan->
cheby[1];
570 my_plan->
aux[j]=1/(2*_Complex_I)*(my_plan->
cheby[j+1]-my_plan->
cheby[j-1]);
572 my_plan->
aux[N-1]=1/(2*_Complex_I)*(-my_plan->
cheby[j-1]);
583 for(j=1;j<=my_plan->
N_total;j++)
596 int i, j, m, k, max, glo1, glo2;
609 for (j = 0; j < M; j++)
610 plan3D->
f_hat[0] += plan3D->
f[j];
614 for (j = 0; j < M; j++)
619 if (plan3D->
flags & NFSOFT_USE_NDFT)
621 nfft_adjoint_direct(&(plan3D->
p_nfft));
625 nfft_adjoint(&(plan3D->
p_nfft));
632 for (k = -N; k <= N; k++)
634 for (m = -N; m <= N; m++)
637 max = (ABS(m) > ABS(k) ? ABS(m) : ABS(k));
639 for (i = 1; i < 2* plan3D ->
N_total + 3; i++)
647 e2c(plan3D, ABS((k + m) % 2));
656 for (j = max; j <= N; j++)
658 if ((plan3D->
flags & NFSOFT_REPRESENT))
660 if ((k < 0) && (k % 2))
664 if ((m < 0) && (m % 2))
674 if ((plan3D->
flags & NFSOFT_NORMALIZED))
676 plan3D->
f_hat[glo1] = plan3D->
f_hat[glo1] * (1 / (2. *
PI)) * SQRT(
677 0.5 * (2. * (j) + 1.));
690 nfft_finalize(&plan->
p_nfft);
698 if (plan->
flags & NFSOFT_MALLOC_F_HAT)
705 if (plan->
flags & NFSOFT_MALLOC_F)
712 if (plan->
flags & NFSOFT_MALLOC_X)
719 int posN(
int n,
int m,
int B)
724 pos = posN(n - 1, m, B) + B + 1 - MAX(ABS(m), ABS(n - 1));
static void c2e(nfsft_plan *plan)
Converts coefficients with , from a linear combination of Chebyshev polynomials to coefficients m...
double * gamma
Precomputed recursion coefficients /f$^n/f$ for /f$k = 0,/ldots, N_{{max}}; n=-k,/ldots,k/f$ of associated Legendre-functions /f$P_k^n/f$.
int M_total
Total number of samples.
fftw_complex * f_hat
Vector of Fourier coefficients, size is N_total * sizeof( fftw_complex )
fftw_complex * aux
used when converting Chebychev to Fourier coeffcients
fftw_complex * f
Vector of samples, size is M_total * sizeof( fftw_complex )
fftw_complex * f_hat
Vector of Fourier coefficients, size is N_total * sizeof( fftw_complex )
int N_total
Total number of Fourier coefficients.
void SO3_gamma_row(double *gamma, int N, int m, int n)
Compute three-term-recurrence coefficients of Wigner-d functions for all degrees ...
fpt_set set
Structure for discrete polynomial transform (DPT)
fftw_complex * wig_coeffs
contains a set of SO(3) Fourier coefficients for fixed orders m and n
void SO3_beta_row(double *beta, int N, int m, int n)
Compute three-term-recurrence coefficients of Wigner-d functions for all degrees ...
void nfft_precompute_lin_psi(nfft_plan *ths)
create a lookup table, but NOT for each node good idea K=2^xx TODO: estimate K, call from init assume...
Holds data for a set of cascade summations.
fftw_complex * f
Vector of samples, size is M_total * sizeof( fftw_complex )
void(* mv_adjoint)(void *)
Pointer to the own adjoint.
void nfft_trafo(nfft_plan *ths)
user routines
void SO3_alpha_row(double *alpha, int N, int m, int n)
Compute three-term-recurrence coefficients of Wigner-d functions for all degrees ...
int M_total
Total number of samples.
Header file for utility functions used by the nfft3 library.
void(* mv_trafo)(void *)
Pointer to the own transform.
fpt_set internal_fpt_set
the internal FPT plan
fftw_complex * cheby
contains a set of Chebychev coefficients for fixed orders m and n
#define PI
Formerly known to be an irrational number.
void * nfft_malloc(size_t n)
double * alpha
Precomputed recursion coefficients /f$^n/f$ for /f$k = 0,/ldots, N_{{max}}; n=-k,/ldots,k/f$ of associated Legendre-functions /f$P_k^n/f$.
double * beta
Precomputed recursion coefficients /f$^n/f$ for /f$k = 0,/ldots, N_{{max}}; n=-k,/ldots,k/f$ of associated Legendre-functions /f$P_k^n/f$.
Header file for functions related to Wigner-d/D functions.
double * x
Nodes in time/spatial domain, size is doubles.
int N_total
Total number of Fourier coefficients.
nfft_plan p_nfft
the internal NFFT plan
unsigned int flags
the planner flags