NFFT Logo 3.1.0 API Reference

nfsoft.c

00001 /*
00002  * Copyright (c) 2002, 2009 Jens Keiner, Stefan Kunis, Daniel Potts
00003  *
00004  * This program is free software; you can redistribute it and/or modify it under
00005  * the terms of the GNU General Public License as published by the Free Software
00006  * Foundation; either version 2 of the License, or (at your option) any later
00007  * version.
00008  *
00009  * This program is distributed in the hope that it will be useful, but WITHOUT
00010  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
00011  * FOR A PARTICULAR PURPOSE.  See the GNU General Public License for more
00012  * details.
00013  *
00014  * You should have received a copy of the GNU General Public License along with
00015  * this program; if not, write to the Free Software Foundation, Inc., 51
00016  * Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
00017  */
00018 
00019 /* $Id: nfsoft.c 3119 2009-03-16 21:30:59Z vollrath $ */
00020 
00021 #include <stdio.h>
00022 #include <math.h>
00023 #include <string.h>
00024 #include <stdlib.h>
00025 #include <complex.h>
00026 #include "nfft3.h"
00027 #include "util.h"
00028 #include "infft.h"
00029 #include "wigner.h"
00030 
00031 #define DEFAULT_NFFT_CUTOFF    6
00032 #define FPT_THRESHOLD          1000
00033 
00034 void nfsoft_init(nfsoft_plan *plan, int N, int M)
00035 {
00036   nfsoft_init_advanced(plan, N, M, NFSOFT_MALLOC_X | NFSOFT_MALLOC_F
00037       | NFSOFT_MALLOC_F_HAT);
00038 }
00039 
00040 void nfsoft_init_advanced(nfsoft_plan *plan, int N, int M,
00041     unsigned int nfsoft_flags)
00042 {
00043   nfsoft_init_guru(plan, N, M, nfsoft_flags, PRE_PHI_HUT | PRE_PSI | MALLOC_X
00044       | MALLOC_F_HAT | MALLOC_F | FFTW_INIT | FFT_OUT_OF_PLACE,
00045       DEFAULT_NFFT_CUTOFF, FPT_THRESHOLD);
00046 }
00047 
00048 void nfsoft_init_guru(nfsoft_plan *plan, int B, int M,
00049     unsigned int nfsoft_flags, unsigned int nfft_flags, int nfft_cutoff,
00050     int fpt_kappa)
00051 {
00052   int N[3];
00053   int n[3];
00054 
00055   N[0] = 2* B + 2;
00056   N[1] = 2* B + 2;
00057   N[2] = 2* B + 2;
00058 
00059   n[0] = 8* B ;
00060   n[1] = 8* B ;
00061   n[2] = 8* B ;
00062 
00063   nfft_init_guru(&plan->nfft_plan, 3, N, M, n, nfft_cutoff, nfft_flags,
00064       FFTW_ESTIMATE | FFTW_DESTROY_INPUT);
00065 
00066   if ((plan->nfft_plan).nfft_flags & PRE_LIN_PSI)
00067   {
00068     nfft_precompute_lin_psi(&(plan->nfft_plan));
00069   }
00070 
00071   plan->N_total = B;
00072   plan->M_total = M;
00073   plan->fpt_kappa = fpt_kappa;
00074   plan->flags = nfsoft_flags;
00075 
00076   if (plan->flags & NFSOFT_MALLOC_F_HAT)
00077   {
00078     plan->f_hat = (C*) nfft_malloc((B + 1) * (4* (B +1)*(B+1)-1)/3*sizeof(C));
00079   }
00080 
00081   if (plan->f_hat == NULL ) printf("Allocation failed!\n");
00082 
00083   if (plan->flags & NFSOFT_MALLOC_X)
00084   {
00085     plan->x = (R*) nfft_malloc(plan->M_total*3*sizeof(R));
00086   }
00087   if (plan->flags & NFSOFT_MALLOC_F)
00088   {
00089     plan->f = (C*) nfft_malloc(plan->M_total*sizeof(C));
00090   }
00091 
00092   if (plan->x == NULL ) printf("Allocation failed!\n");
00093   if (plan->f == NULL ) printf("Allocation failed!\n");
00094 
00095   plan->wig_coeffs = (C*) nfft_malloc((nfft_next_power_of_2(B)+1)*sizeof(C));
00096 
00097   plan->cheby = (C*) nfft_malloc((2*B+2)*sizeof(C));
00098   plan->aux = (C*) nfft_malloc((2*B+4)*sizeof(C));
00099 
00100   if (plan->wig_coeffs == NULL ) printf("Allocation failed!\n");
00101   if (plan->cheby == NULL ) printf("Allocation failed!\n");
00102   if (plan->aux == NULL ) printf("Allocation failed!\n");
00103 
00104   plan->mv_trafo = (void (*) (void* ))nfsoft_trafo;
00105   plan->mv_adjoint = (void (*) (void* ))nfsoft_adjoint;
00106 }
00107 
00108 static void c2e(nfsoft_plan *my_plan, int even)
00109 {
00110   int j, N;
00111 
00113   N = 2* (my_plan ->N_total+1);
00114 
00116   my_plan->cheby[my_plan->N_total+1] = my_plan->wig_coeffs[0];
00117   my_plan->cheby[0]=0.0;
00118 
00119   for (j=1;j<my_plan->N_total+1;j++)
00120   {
00121     my_plan->cheby[my_plan->N_total+1+j]=0.5* my_plan->wig_coeffs[j];
00122     my_plan->cheby[my_plan->N_total+1-j]=0.5* my_plan->wig_coeffs[j];
00123   }
00124 
00125   C *aux= (C*) nfft_malloc((N+2)*sizeof(C));
00126 
00127   for(j=1;j<N;j++)
00128   aux[j]=my_plan->cheby[j];
00129 
00130   aux[0]=0.;
00131   aux[N]=0.;
00132 
00133   if (even>0)
00134   {
00135     my_plan->cheby[0]=(C) (-1.)/(2.*_Complex_I) * aux[1];
00136     for (j=1;j<N;j++)
00137     {
00138       my_plan->cheby[j]=(1./(2.*_Complex_I)*(aux[j+1]-aux[j-1]));
00139     }
00140 
00141   }
00142   free(aux);
00143   aux = NULL;
00144 }
00145 
00146 static fpt_set SO3_fpt_init(int l, fpt_set set, unsigned int flags, int kappa)
00147 {
00148   int N, t, k_start, k_end, k, m;
00149   int glo = 0;
00150   R *alpha, *beta, *gamma;
00151 
00153   if (flags & NFSOFT_USE_DPT)
00154   {
00155     if (l < 2)
00156       N = 2;
00157     else
00158       N = l;
00159 
00160     t = (int) log2(nfft_next_power_of_2(N));
00161 
00162   }
00163   else
00164   {
00166     if (l < 2)
00167       N = 2;
00168     else
00169       N = nfft_next_power_of_2(l);
00170 
00171     t = (int) log2(N);
00172   }
00173 
00175   alpha = (R*) nfft_malloc((N + 2) * sizeof(R));
00176   beta = (R*) nfft_malloc((N + 2) * sizeof(R));
00177   gamma = (R*) nfft_malloc((N + 2) * sizeof(R));
00178 
00180   if (flags & NFSOFT_NO_STABILIZATION)
00181   {
00182     set = fpt_init((2* N + 1) * (2* N + 1), t, 0U | FPT_NO_STABILIZATION);
00183   }
00184   else
00185   {
00186     set = fpt_init((2* N + 1) * (2* N + 1), t, 0U);
00187   }
00188 
00189   for (k = -N; k <= N; k++)
00190     for (m = -N; m <= N; m++)
00191     {
00193       k_start = (ABS(k) >= ABS(m)) ? ABS(k) : ABS(m);
00194       k_end = N;
00195 
00196       SO3_alpha_row(alpha, N, k, m);
00197       SO3_beta_row(beta, N, k, m);
00198       SO3_gamma_row(gamma, N, k, m);
00199 
00200       fpt_precompute(set, glo, alpha, beta, gamma, k_start, kappa);
00201       glo++;
00202     }
00203 
00204   free(alpha);
00205   free(beta);
00206   free(gamma);
00207   alpha = NULL;
00208   beta = NULL;
00209   gamma = NULL;
00210 
00211   return set;
00212 }
00213 
00214 void SO3_fpt(C *coeffs, fpt_set set, int l, int k, int m, unsigned int flags)
00215 {
00216   int N;
00218   C* x;
00220   C* y;
00221 
00222   int trafo_nr; 
00223   int k_start, k_end, j;
00224   int function_values = 0;
00225 
00227   if (flags & NFSOFT_USE_DPT)
00228   {
00229     N = l;
00230     if (l < 2)
00231       N = 2;
00232   }
00233   else
00234   {
00235     if (l < 2)
00236       N = 2;
00237     else
00238       N = nfft_next_power_of_2(l);
00239   }
00240 
00242   k_start = (ABS(k) >= ABS(m)) ? ABS(k) : ABS(m);
00243   k_end = N;
00244   trafo_nr = (N + k) * (2* N + 1) + (m + N);
00245 
00247   x = (C*) nfft_malloc((k_end + 1) * sizeof(C));
00248 
00249   //for (j = 0; j <= k_end-k_start; j++)
00250   //{
00251   // x[j+k_start] = coeffs[j];
00252   //}
00253 
00254   for (j = 0; j <= l - k_start; j++)
00255   {
00256     x[j + k_start] = coeffs[j];
00257   }
00258   for (j = l - k_start + 1; j <= k_end - k_start; j++)
00259   {
00260     x[j + k_start] = K(0.0);
00261   }
00262 
00264   y = (C*) nfft_malloc((k_end + 1) * sizeof(C));
00265 
00266   if (flags & NFSOFT_USE_DPT)
00267   { 
00268     dpt_trafo(set, trafo_nr, &x[k_start], y, k_end, 0U
00269         | (function_values ? FPT_FUNCTION_VALUES : 0U));
00270   }
00271   else
00272   { 
00273     fpt_trafo(set, trafo_nr, &x[k_start], y, k_end, 0U
00274         | (function_values ? FPT_FUNCTION_VALUES : 0U));
00275   }
00276 
00278   for (j = 0; j <= l; j++)
00279   {
00280     coeffs[j] = y[j];
00281   }
00282 
00285   free(x);
00286   free(y);
00287   x = NULL;
00288   y = NULL;
00289 }
00290 
00291 void SO3_fpt_transposed(C *coeffs, fpt_set set, int l, int k, int m,
00292     unsigned int flags)
00293 {
00294   int N, k_start, k_end, j;
00295   int trafo_nr; 
00296   int function_values = 0;
00298   C* x;
00300   C* y;
00301 
00304   if (flags & NFSOFT_USE_DPT)
00305   {
00306     N = l;
00307     if (l < 2)
00308       N = 2;
00309   }
00310   else
00311   {
00312     if (l < 2)
00313       N = 2;
00314     else
00315       N = nfft_next_power_of_2(l);
00316   }
00317 
00319   k_start = (ABS(k) >= ABS(m)) ? ABS(k) : ABS(m);
00320   k_end = N;
00321   trafo_nr = (N + k) * (2* N + 1) + (m + N);
00322 
00324   y = (C*) nfft_malloc((k_end + 1) * sizeof(C));
00326   x = (C*) nfft_malloc((k_end + 1) * sizeof(C));
00327 
00328   for (j = 0; j <= l; j++)
00329   {
00330     y[j] = coeffs[j];
00331   }
00332   for (j = l + 1; j <= k_end; j++)
00333   {
00334     y[j] = K(0.0);
00335   }
00336 
00337   if (flags & NFSOFT_USE_DPT)
00338   {
00339     dpt_transposed(set, trafo_nr, &x[k_start], y, k_end, 0U
00340         | (function_values ? FPT_FUNCTION_VALUES : 0U));
00341   }
00342   else
00343   {
00344     fpt_transposed(set, trafo_nr, &x[k_start], y, k_end, 0U
00345         | (function_values ? FPT_FUNCTION_VALUES : 0U));
00346   }
00347 
00348   for (j = 0; j <= l; j++)
00349   {
00350     coeffs[j] = x[j];
00351   }
00352 
00354   free(x);
00355   free(y);
00356   x = NULL;
00357   y = NULL;
00358 }
00359 
00360 void nfsoft_precompute(nfsoft_plan *plan3D)
00361 {
00362   int j;
00363   int N = plan3D->N_total;
00364   int M = plan3D->M_total;
00365 
00368   for (j = 0; j < M; j++)
00369   {
00370     plan3D->nfft_plan.x[3* j ] = plan3D->x[3* j + 2];
00371     plan3D->nfft_plan.x[3* j + 1] = plan3D->x[3* j ];
00372     plan3D->nfft_plan.x[3* j + 2] = plan3D->x[3* j + 1];
00373   }
00374 
00375   for (j = 0; j < 3* plan3D ->nfft_plan.M_total; j++)
00376   {
00377     plan3D->nfft_plan.x[j] = plan3D->nfft_plan.x[j] * (1 / (2* PI ));
00378   }
00379 
00380   if ((plan3D->nfft_plan).nfft_flags & FG_PSI)
00381   {
00382     nfft_precompute_one_psi(&(plan3D->nfft_plan));
00383   }
00384   if ((plan3D->nfft_plan).nfft_flags & PRE_PSI)
00385   {
00386     nfft_precompute_one_psi(&(plan3D->nfft_plan));
00387   }
00388 
00390   plan3D->fpt_set = SO3_fpt_init(N, plan3D->fpt_set, plan3D->flags,
00391       plan3D->fpt_kappa);
00392 
00393   for (j = 0; j < plan3D->nfft_plan.N_total; j++)
00394     plan3D->nfft_plan.f_hat[j] = 0.0;
00395 
00396 }
00397 
00398 void nfsoft_trafo(nfsoft_plan *plan3D)
00399 {
00400   int i, j, m, k, max, glo1, glo2;
00401 
00402   i = 0;
00403   glo1 = 0;
00404   glo2 = 0;
00405 
00406   int N = plan3D->N_total;
00407   int M = plan3D->M_total;
00408 
00410   if (N == 0)
00411   {
00412     for (j = 0; j < M; j++)
00413       plan3D->f[j] = plan3D->f_hat[0];
00414     return;
00415   }
00416 
00417   for (j = 0; j < plan3D->nfft_plan.N_total; j++)
00418     plan3D->nfft_plan.f_hat[j] = 0.0;
00419 
00420   for (k = -N; k <= N; k++)
00421   {
00422     for (m = -N; m <= N; m++)
00423     {
00424 
00425       max = (ABS(m) > ABS(k) ? ABS(m) : ABS(k));
00426 
00427       for (j = 0; j <= N - max; j++)
00428       {
00429         plan3D->wig_coeffs[j] = plan3D->f_hat[glo1];
00430 
00431         if ((plan3D->flags & NFSOFT_NORMALIZED))
00432         {
00433           plan3D->wig_coeffs[j] = plan3D->wig_coeffs[j] * (1. / (2. * PI))
00434               * SQRT(0.5 * (2. * (max + j) + 1.));
00435         }
00436 
00437         if ((plan3D->flags & NFSOFT_REPRESENT))
00438         {
00439           if ((k < 0) && (k % 2))
00440           {
00441             plan3D->wig_coeffs[j] = plan3D->wig_coeffs[j] * (-1);
00442           }
00443           if ((m < 0) && (m % 2))
00444             plan3D->wig_coeffs[j] = plan3D->wig_coeffs[j] * (-1);
00445         }
00446 
00447         glo1++;
00448       }
00449 
00450       for (j = N - max + 1; j < nfft_next_power_of_2(N) + 1; j++)
00451         plan3D->wig_coeffs[j] = 0.0;
00452       //fprintf(stdout,"\n k= %d, m= %d \n",k,m);
00453       SO3_fpt(plan3D->wig_coeffs, plan3D->fpt_set, N, k, m, plan3D->flags);
00454 
00455       c2e(plan3D, ABS((k + m) % 2));
00456 
00457       for (i = 1; i <= 2* plan3D ->N_total + 2; i++)
00458       {
00459         plan3D->nfft_plan.f_hat[NFSOFT_INDEX(k, m, i - N - 1, N) - 1]
00460             = plan3D->cheby[i - 1];
00461         //fprintf(stdout,"%f \t", plan3D->nfft_plan.f_hat[NFSOFT_INDEX(k,m,i-N-1,N)-1]);
00462         //fprintf(stdout,"another index: %d for k=%d,m=%d,l=%d,N=%d \n", NFSOFT_INDEX(k,m,i-N-1,N)-1,k,m,i-N-1,N);
00463       }
00464 
00465     }
00466   }
00467 
00468   if (plan3D->flags & NFSOFT_USE_NDFT)
00469   {
00470     ndft_trafo(&(plan3D->nfft_plan));
00471   }
00472   else
00473   {
00474     nfft_trafo(&(plan3D->nfft_plan));
00475   }
00476 
00477   for (j = 0; j < plan3D->M_total; j++)
00478     plan3D->f[j] = plan3D->nfft_plan.f[j];
00479 
00480 }
00481 
00482 static void e2c(nfsoft_plan *my_plan, int even)
00483 {
00484   int N;
00485   int j;
00486 
00488   N = 2* (my_plan ->N_total+1);
00489   //nfft_vpr_complex(my_plan->cheby,N+1,"chebychev");
00490 
00491 
00492       if (even>0)
00493       {
00494         //my_plan->aux[N-1]= -1/(2*I)* my_plan->cheby[N-2];
00495         my_plan->aux[0]= 1/(2*_Complex_I)*my_plan->cheby[1];
00496 
00497         for(j=1;j<N-1;j++)
00498         {
00499           my_plan->aux[j]=1/(2*_Complex_I)*(my_plan->cheby[j+1]-my_plan->cheby[j-1]);
00500 }
00501 my_plan->aux[N-1]=1/(2*_Complex_I)*(-my_plan->cheby[j-1]);
00502 
00503 
00504 for(j=0;j<N;j++)
00505 {
00506 my_plan->cheby[j]= my_plan->aux[j];
00507 }
00508 }
00509 
00510 my_plan->wig_coeffs[0]=my_plan->cheby[my_plan->N_total+1];
00511 
00512 for(j=1;j<=my_plan->N_total;j++)
00513 {
00514 my_plan->wig_coeffs[j]=0.5*(my_plan->cheby[my_plan->N_total+j+1]+my_plan->cheby[my_plan->N_total+1-j]);
00515 }
00516 
00517 
00518 
00519 //nfft_vpr_complex(my_plan->wig_coeffs,my_plan->N_total,"chebychev ");
00520 
00521 }
00522 
00523 void nfsoft_adjoint(nfsoft_plan *plan3D)
00524 {
00525   int i, j, m, k, max, glo1, glo2;
00526 
00527   i = 0;
00528   glo1 = 0;
00529   glo2 = 0;
00530 
00531   int N = plan3D->N_total;
00532   int M = plan3D->M_total;
00533 
00534   //nothing much to be done for polynomial degree 0
00535   if (N == 0)
00536   {
00537     plan3D->f_hat[0]=0;
00538     for (j = 0; j < M; j++)
00539       plan3D->f_hat[0] += plan3D->f[j];
00540     return;
00541   }
00542 
00543   for (j = 0; j < M; j++)
00544   {
00545     plan3D->nfft_plan.f[j] = plan3D->f[j];
00546   }
00547 
00548   if (plan3D->flags & NFSOFT_USE_NDFT)
00549   {
00550     ndft_adjoint(&(plan3D->nfft_plan));
00551   }
00552   else
00553   {
00554     nfft_adjoint(&(plan3D->nfft_plan));
00555   }
00556 
00557   //nfft_vpr_complex(plan3D->nfft_plan.f_hat,plan3D->nfft_plan.N_total,"all results");
00558 
00559   glo1 = 0;
00560 
00561   for (k = -N; k <= N; k++)
00562   {
00563     for (m = -N; m <= N; m++)
00564     {
00565 
00566       max = (ABS(m) > ABS(k) ? ABS(m) : ABS(k));
00567 
00568       for (i = 1; i < 2* plan3D ->N_total + 3; i++)
00569       {
00570         plan3D->cheby[i - 1] = plan3D->nfft_plan.f_hat[NFSOFT_INDEX(k, m, i - N
00571             - 1, N) - 1];
00572       }
00573 
00574       //fprintf(stdout,"k=%d,m=%d \n",k,m);
00575       //nfft_vpr_complex(plan3D->cheby,2*plan3D->N_total+2,"euler");
00576       e2c(plan3D, ABS((k + m) % 2));
00577 
00578       //nfft_vpr_complex(plan3D->wig_coeffs,plan3D->N_total+1,"chebys");
00579       SO3_fpt_transposed(plan3D->wig_coeffs, plan3D->fpt_set, N, k, m,
00580           plan3D->flags);
00581       //nfft_vpr_complex(plan3D->wig_coeffs,plan3D->N_total+1,"wigners");
00582       //  SO3_fpt_transposed(plan3D->wig_coeffs,N,k,m,plan3D->flags,plan3D->fpt_kappa);
00583 
00584 
00585       for (j = max; j <= N; j++)
00586       {
00587         if ((plan3D->flags & NFSOFT_REPRESENT))
00588         {
00589           if ((k < 0) && (k % 2))
00590           {
00591             plan3D->wig_coeffs[j] = -plan3D->wig_coeffs[j];
00592           }
00593           if ((m < 0) && (m % 2))
00594             plan3D->wig_coeffs[j] = -plan3D->wig_coeffs[j];
00595         }
00596 
00597         plan3D->f_hat[glo1] = plan3D->wig_coeffs[j];
00598 
00599         if ((plan3D->flags & NFSOFT_NORMALIZED))
00600         {
00601           plan3D->f_hat[glo1] = plan3D->f_hat[glo1] * (1 / (2. * PI)) * SQRT(
00602               0.5 * (2. * (j) + 1.));
00603         }
00604 
00605         glo1++;
00606       }
00607 
00608     }
00609   }
00610 }
00611 
00612 void nfsoft_finalize(nfsoft_plan *plan)
00613 {
00614   /* Finalise the nfft plan. */
00615   nfft_finalize(&plan->nfft_plan);
00616   free(plan->wig_coeffs);
00617   free(plan->cheby);
00618   free(plan->aux);
00619 
00620   fpt_finalize(plan->fpt_set);
00621   plan->fpt_set = NULL;
00622 
00623   if (plan->flags & NFSOFT_MALLOC_F_HAT)
00624   {
00625     //fprintf(stderr,"deallocating f_hat\n");
00626     free(plan->f_hat);
00627   }
00628 
00629   /* De-allocate memory for samples, if neccesary. */
00630   if (plan->flags & NFSOFT_MALLOC_F)
00631   {
00632     //fprintf(stderr,"deallocating f\n");
00633     free(plan->f);
00634   }
00635 
00636   /* De-allocate memory for nodes, if neccesary. */
00637   if (plan->flags & NFSOFT_MALLOC_X)
00638   {
00639     //fprintf(stderr,"deallocating x\n");
00640     free(plan->x);
00641   }
00642 }
00643 
00644 int posN(int n, int m, int B)
00645 {
00646   int pos;
00647 
00648   if (n > -B)
00649     pos = posN(n - 1, m, B) + B + 1 - MAX(ABS(m), ABS(n - 1));
00650   else
00651     pos = 0;
00652   //(n > -B? pos=posN(n-1,m,B)+B+1-MAX(ABS(m),ABS(n-1)): pos= 0)
00653   return pos;
00654 }
00655 

Generated on 19 Mar 2009 by Doxygen 1.5.3