NFFT Logo 3.0.2 API Reference
Main Page | Modules | Data Structures | Directories | File List | Data Fields | Globals

iterS2.c

00001 /* $Id$
00002  *
00003  * iterS2 - Iterative reconstruction on the sphere S2
00004  *
00005  * Copyright (C) 2006 Jens Keiner
00006  */
00007 
00014 /* Include standard C headers. */
00015 #include <stdio.h>
00016 #include <stdlib.h>
00017 #include <math.h>
00018 
00019 /* Include NFFT3 library header. */
00020 #include "nfft3.h"
00021 
00022 /* Include NFFT 3 utilities headers. */
00023 #include "util.h"
00024 
00025 #include "legendre.h"
00026 
00028 enum boolean {NO = 0, YES = 1};
00029 
00040 int main (int argc, char **argv)
00041 {
00042   int T;
00043   int N;
00044   int M;
00045   int M2;
00046 
00047   int t;                       /* Index variable for testcases                */
00048   double time;                 /* Time for fast algorithm in seconds          */
00049   double err_f;                /* Error E_infty for fast algorithm            */
00050   nfsft_plan plan;             /* NFSFT plan                                  */
00051   nfsft_plan plan2;            /* NFSFT plan                                  */
00052   infsft_plan iplan;           /* NFSFT plan                                  */
00053   int j;                       /*                                             */
00054   int k;                       /*                                             */
00055   int m;                       /*                                             */
00056   int n;                       /*                                             */
00057   int use_nfsft;               /*                                             */
00058   int use_nfft;                /*                                             */
00059   int use_fpt;                 /*                                             */
00060   int cutoff;                  
00061   double threshold;            
00062   double re;
00063   double im;
00064   double a,b;
00065   double *scratch;
00066   double xs;
00067   double *ys;
00068   double *temp;
00069   double complex *temp2;
00070   int qlength;
00071   double *qweights;
00072   fftw_plan fplan;
00073   fpt_set set;
00074   int npt;
00075   int npt_exp;
00076   double *alpha, *beta, *gamma;
00077 
00078   /* Read the number of testcases. */
00079   fscanf(stdin,"testcases=%d\n",&T);
00080   fprintf(stderr,"%d\n",T);
00081 
00082   /* Process each testcase. */
00083   for (t = 0; t < T; t++)
00084   {
00085     /* Check if the fast transform shall be used. */
00086     fscanf(stdin,"nfsft=%d\n",&use_nfsft);
00087     fprintf(stderr,"%d\n",use_nfsft);
00088     if (use_nfsft != NO)
00089     {
00090       /* Check if the NFFT shall be used. */
00091       fscanf(stdin,"nfft=%d\n",&use_nfft);
00092       fprintf(stderr,"%d\n",use_nfsft);
00093       if (use_nfft != NO)
00094       {
00095         /* Read the cut-off parameter. */
00096         fscanf(stdin,"cutoff=%d\n",&cutoff);
00097         fprintf(stderr,"%d\n",cutoff);
00098       }
00099       else
00100       {
00101         /* TODO remove this */
00102         /* Initialize unused variable with dummy value. */
00103         cutoff = 1;
00104       }
00105       /* Check if the fast polynomial transform shall be used. */
00106       fscanf(stdin,"fpt=%d\n",&use_fpt);
00107       fprintf(stderr,"%d\n",use_fpt);
00108       if (use_fpt != NO)
00109       {
00110         /* Read the NFSFT threshold parameter. */
00111         fscanf(stdin,"threshold=%lf\n",&threshold);
00112         fprintf(stderr,"%lf\n",threshold);
00113       }
00114       else
00115       {
00116         /* TODO remove this */
00117         /* Initialize unused variable with dummy value. */
00118         threshold = 1000.0;
00119       }
00120     }
00121     else
00122     {
00123       /* TODO remove this */
00124       /* Set dummy values. */
00125       use_nfft = NO;
00126       use_fpt = NO;
00127       cutoff = 3;
00128       threshold = 1000.0;
00129     }
00130 
00131     /* Read the bandwidth. */
00132     fscanf(stdin,"bandwidth=%d\n",&N);
00133     fprintf(stderr,"%d\n",N);
00134 
00135     /* Do precomputation. */
00136     nfsft_precompute(N,threshold,
00137       ((use_nfsft==NO)?(NFSFT_NO_FAST_ALGORITHM):(0U/*NFSFT_NO_DIRECT_ALGORITHM*/)), 0U);
00138 
00139     /* Read the number of nodes. */
00140     fscanf(stdin,"nodes=%d\n",&M);
00141     fprintf(stderr,"%d\n",M);
00142 
00143     /* */
00144     if ((N+1)*(N+1) > M)
00145     {
00146       nfft_next_power_of_2_exp(N, &npt, &npt_exp);
00147       fprintf(stderr, "npt = %d, npt_exp = %d\n", npt, npt_exp);
00148       fprintf(stderr,"Optimal interpolation!\n");
00149       scratch = (double*) malloc(4*sizeof(double));
00150       ys = (double*) malloc((N+1)*sizeof(double));
00151       temp = (double*) malloc((2*N+1)*sizeof(double));
00152       temp2 = (double complex*) malloc((N+1)*sizeof(double complex));
00153 
00154       a = 0.0;
00155       for (j = 0; j <= N; j++)
00156       {
00157         xs = 2.0 + (2.0*j)/(N+1);
00158         ys[j] = (2.0-((j == 0)?(1.0):(0.0)))*4.0*nfft_bspline(4,xs,scratch);
00159         //fprintf(stdout,"%3d: g(%le) = %le\n",j,xs,ys[j]);
00160         a += ys[j];
00161       }
00162       //fprintf(stdout,"a = %le\n",a);
00163       for (j = 0; j <= N; j++)
00164       {
00165         ys[j] *= 1.0/a;
00166       }
00167 
00168       qlength = 2*N+1;
00169       qweights = (double*) malloc(qlength*sizeof(double));
00170 
00171       fplan = fftw_plan_r2r_1d(N+1, qweights, qweights, FFTW_REDFT00, 0U);
00172       for (j = 0; j < N+1; j++)
00173       {
00174         qweights[j] = -2.0/(4*j*j-1);
00175       }
00176       fftw_execute(fplan);
00177       qweights[0] *= 0.5;
00178 
00179       for (j = 0; j < N+1; j++)
00180       {
00181         qweights[j] *= 1.0/(2.0*N+1.0);
00182         qweights[2*N+1-1-j] = qweights[j];
00183       }
00184 
00185       fplan = fftw_plan_r2r_1d(2*N+1, temp, temp, FFTW_REDFT00, 0U);
00186       for (j = 0; j <= N; j++)
00187       {
00188         temp[j] = ((j==0 || j == 2*N)?(1.0):(0.5))*ys[j];
00189       }
00190       for (j = N+1; j < 2*N+1; j++)
00191       {
00192         temp[j] = 0.0;
00193       }
00194       fftw_execute(fplan);
00195 
00196       for (j = 0; j < 2*N+1; j++)
00197       {
00198         temp[j] *= qweights[j];
00199       }
00200 
00201       fftw_execute(fplan);
00202 
00203       for (j = 0; j < 2*N+1; j++)
00204       {
00205         temp[j] *= ((j==0 || j == 2*N)?(1.0):(0.5));
00206         if (j <= N)
00207         {
00208           temp2[j] = temp[j];
00209         }
00210       }
00211 
00212       set = fpt_init(0, npt_exp, 0U);
00213 
00214       alpha = (double*) malloc((N+2)*sizeof(double));
00215       beta = (double*) malloc((N+2)*sizeof(double));
00216       gamma = (double*) malloc((N+2)*sizeof(double));
00217 
00218       alpha_al_row(alpha, N, 0);
00219       beta_al_row(beta, N, 0);
00220       gamma_al_row(gamma, N, 0);
00221 
00222       fpt_precompute(set, 0, alpha, beta, gamma, 0, 1000.0);
00223 
00224       fpt_transposed(set,0, temp2, temp2, N, 0U);
00225 
00226       for (j = 0; j <= N; j++)
00227       {
00228         //temp2[j] *= (2*j+1.0)/2.0;
00229         //fprintf(stdout,"temp2[%d] = %le\n",j,cabs(temp2[j]));
00230       }
00231 
00232       fpt_finalize(set);
00233 
00234       free(alpha);
00235       free(beta);
00236       free(gamma);
00237 
00238       fftw_destroy_plan(fplan);
00239 
00240       free(scratch);
00241       free(qweights);
00242       free(ys);
00243       free(temp);
00244       //free(temp2);
00245       //return EXIT_SUCCESS;
00246     }
00247 
00248     /* Init transform plans. */
00249     nfsft_init_guru(&plan, N, M,
00250       ((use_nfft!=0)?(0U):(NFSFT_USE_NDFT)) |
00251       ((use_fpt!=0)?(0U):(NFSFT_USE_DPT)) | NFSFT_MALLOC_F | NFSFT_MALLOC_X |
00252       NFSFT_MALLOC_F_HAT | NFSFT_NORMALIZED | NFSFT_ZERO_F_HAT,
00253       PRE_PHI_HUT | PRE_PSI | FFTW_INIT |
00254       FFT_OUT_OF_PLACE,
00255       cutoff);
00256 
00257     if ((N+1)*(N+1) > M)
00258     {
00259       infsft_init_advanced (&iplan, &plan, CGNE | PRECOMPUTE_DAMP);
00260     }
00261     else
00262     {
00263       infsft_init_advanced (&iplan, &plan, CGNR | PRECOMPUTE_WEIGHT | PRECOMPUTE_DAMP);
00264     }
00265 
00266     /* Read the nodes and function values. */
00267     for (j = 0; j < M; j++)
00268     {
00269       fscanf(stdin,"%le %le %le %le\n",&plan.x[2*j+1],&plan.x[2*j],&re,&im);
00270       plan.x[2*j+1] = plan.x[2*j+1]/(2.0*PI);
00271       plan.x[2*j] = plan.x[2*j]/(2.0*PI);
00272       if (plan.x[2*j] >= 0.5)
00273       {
00274         plan.x[2*j] = plan.x[2*j] - 1;
00275       }
00276       iplan.y[j] = re + I * im;
00277       fprintf(stderr,"%le %le %le %le\n",plan.x[2*j+1],plan.x[2*j],
00278         creal(iplan.y[j]),cimag(iplan.y[j]));
00279     }
00280 
00281     /* Read the number of nodes. */
00282     fscanf(stdin,"nodes_eval=%d\n",&M2);
00283     fprintf(stderr,"%d\n",M2);
00284 
00285     /* Init transform plans. */
00286     nfsft_init_guru(&plan2, N, M2,
00287       ((use_nfft!=0)?(0U):(NFSFT_USE_NDFT)) |
00288       ((use_fpt!=0)?(0U):(NFSFT_USE_DPT)) | NFSFT_MALLOC_F | NFSFT_MALLOC_X |
00289       NFSFT_MALLOC_F_HAT | NFSFT_NORMALIZED | NFSFT_ZERO_F_HAT,
00290       PRE_PHI_HUT | PRE_PSI | FFTW_INIT |
00291       FFT_OUT_OF_PLACE,
00292       cutoff);
00293 
00294     /* Read the nodes and function values. */
00295     for (j = 0; j < M2; j++)
00296     {
00297       fscanf(stdin,"%le %le\n",&plan2.x[2*j+1],&plan2.x[2*j]);
00298       plan2.x[2*j+1] = plan2.x[2*j+1]/(2.0*PI);
00299       plan2.x[2*j] = plan2.x[2*j]/(2.0*PI);
00300       if (plan2.x[2*j] >= 0.5)
00301       {
00302         plan2.x[2*j] = plan2.x[2*j] - 1;
00303       }
00304       fprintf(stderr,"%le %le\n",plan2.x[2*j+1],plan2.x[2*j]);
00305     }
00306 
00307     nfsft_precompute_x(&plan);
00308 
00309     nfsft_precompute_x(&plan2);
00310 
00311     /* Frequency weights. */
00312     if ((N+1)*(N+1) > M)
00313     {
00314       /* Compute Voronoi weights. */
00315       //nfft_voronoi_weights_S2(iplan.w, plan.x, M);
00316 
00317       /* Print out Voronoi weights. */
00318       /*a = 0.0;
00319       for (j = 0; j < plan.M_total; j++)
00320       {
00321         fprintf(stderr,"%le\n",iplan.w[j]);
00322         a += iplan.w[j];
00323       }
00324       fprintf(stderr,"sum = %le\n",a);*/
00325 
00326       for (j = 0; j < plan.N_total; j++)
00327       {
00328         iplan.w_hat[j] = 0.0;
00329       }
00330 
00331       for (k = 0; k <= N; k++)
00332       {
00333         for (j = -k; j <= k; j++)
00334         {
00335           iplan.w_hat[NFSFT_INDEX(k,j,&plan)] = 1.0/(pow(k+1.0,2.0)); /*temp2[j]*/;
00336         }
00337       }
00338     }
00339     else
00340     {
00341       for (j = 0; j < plan.N_total; j++)
00342       {
00343         iplan.w_hat[j] = 0.0;
00344       }
00345 
00346       for (k = 0; k <= N; k++)
00347       {
00348         for (j = -k; j <= k; j++)
00349         {
00350           iplan.w_hat[NFSFT_INDEX(k,j,&plan)] = 1/(pow(k+1.0,2.5));
00351         }
00352       }
00353 
00354       /* Compute Voronoi weights. */
00355       nfft_voronoi_weights_S2(iplan.w, plan.x, M);
00356 
00357       /* Print out Voronoi weights. */
00358       a = 0.0;
00359       for (j = 0; j < plan.M_total; j++)
00360       {
00361         fprintf(stderr,"%le\n",iplan.w[j]);
00362         a += iplan.w[j];
00363       }
00364       fprintf(stderr,"sum = %le\n",a);
00365     }
00366 
00367     fprintf(stderr, "N_total = %d\n", plan.N_total);
00368     fprintf(stderr, "M_total = %d\n", plan.M_total);
00369 
00370     /* init some guess */
00371     for (k = 0; k < plan.N_total; k++)
00372     {
00373       iplan.f_hat_iter[k] = 0.0;
00374     }
00375 
00376     /* inverse trafo */
00377     infsft_before_loop(&iplan);
00378 
00379     /*for (k = 0; k < plan.M_total; k++)
00380     {
00381       printf("%le %le\n",creal(iplan.r_iter[k]),cimag(iplan.r_iter[k]));
00382     }*/
00383 
00384     for (m = 0; m < 29; m++)
00385     {
00386       fprintf(stderr,"Residual ||r||=%e,\n",sqrt(iplan.dot_r_iter));
00387       infsft_loop_one_step(&iplan);
00388     }
00389 
00390     /*NFFT_SWAP_complex(iplan.f_hat_iter, plan.f_hat);
00391     nfsft_trafo(&plan);
00392     NFFT_SWAP_complex(iplan.f_hat_iter, plan.f_hat);
00393 
00394     a = 0.0;
00395     b = 0.0;
00396     for (k = 0; k < plan.M_total; k++)
00397     {
00398       printf("%le %le %le\n",cabs(iplan.y[k]),cabs(plan.f[k]),
00399         cabs(iplan.y[k]-plan.f[k]));
00400       a += cabs(iplan.y[k]-plan.f[k])*cabs(iplan.y[k]-plan.f[k]);
00401       b += cabs(iplan.y[k])*cabs(iplan.y[k]);
00402     }
00403 
00404     fprintf(stderr,"relative error in 2-norm: %le\n",a/b);*/
00405 
00406     NFFT_SWAP_complex(iplan.f_hat_iter, plan2.f_hat);
00407     nfsft_trafo(&plan2);
00408     NFFT_SWAP_complex(iplan.f_hat_iter, plan2.f_hat);
00409     for (k = 0; k < plan2.M_total; k++)
00410     {
00411       fprintf(stdout,"%le\n",cabs(plan2.f[k]));
00412     }
00413 
00414     infsft_finalize (&iplan);
00415 
00416     nfsft_finalize(&plan);
00417 
00418     nfsft_finalize(&plan2);
00419 
00420     /* Delete precomputed data. */
00421     nfsft_forget();
00422 
00423     if ((N+1)*(N+1) > M)
00424     {
00425       free(temp2);
00426     }
00427 
00428   } /* Process each testcase. */
00429 
00430   /* Return exit code for successful run. */
00431   return EXIT_SUCCESS;
00432 }
00433 /* \} */

Generated on 22 Jan 2007 by Doxygen 1.4.1