/**
 * reconstruct makes a 2d-adjoint-nfft
 */
static void reconstruct(char* filename, int N, int M, int weight)
{
  int j;                   /* some variables  */
  double weights;          /* store one weight temporary */
  double real,imag;        /* to read the real and imag part of a complex number */
  nfft_plan my_plan;       /* plan for the two dimensional nfft  */
  FILE* fin;               /* input file  */
  FILE* fweight;           /* input file for the weights */
  FILE *fout_real;         /* output file  */
  FILE *fout_imag;         /* output file  */
  int my_N[2],my_n[2];
  int flags = PRE_PHI_HUT| PRE_PSI |MALLOC_X| MALLOC_F_HAT|
                      MALLOC_F| FFTW_INIT| FFT_OUT_OF_PLACE|
                      FFTW_MEASURE| FFTW_DESTROY_INPUT;

  /* initialise nfft */
  my_N[0]=N; my_n[0]=ceil(N*1.2);
  my_N[1]=N; my_n[1]=ceil(N*1.2);
  nfft_init_guru(&my_plan, 2, my_N, M, my_n, 6,flags,
                      FFTW_MEASURE| FFTW_DESTROY_INPUT);

  fin=fopen(filename,"r");

  fweight=fopen("weights.dat","r");
  for(j=0;j<my_plan.M_total;j++)
  {
    fscanf(fweight,"%le ",&weights);
    fscanf(fin,"%le %le %le %le",&my_plan.x[2*j+0],&my_plan.x[2*j+1],&real,&imag);
    my_plan.f[j] = real + _Complex_I*imag;
    if (weight)
      my_plan.f[j] = my_plan.f[j] * weights;
  }
  fclose(fweight);

  /* precompute psi */
  if(my_plan.nfft_flags & PRE_PSI)
    nfft_precompute_psi(&my_plan);

  /* precompute full psi */
  if(my_plan.nfft_flags & PRE_FULL_PSI)
    nfft_precompute_full_psi(&my_plan);


  /* compute the adjoint nfft */
  nfft_adjoint(&my_plan);

  fout_real=fopen("output_real.dat","w");
  fout_imag=fopen("output_imag.dat","w");

  for (j=0;j<N*N;j++) {
    fprintf(fout_real,"%le ",creal(my_plan.f_hat[j]));
    fprintf(fout_imag,"%le ",cimag(my_plan.f_hat[j]));
  }

  fclose(fin);
  fclose(fout_real);
  fclose(fout_imag);

  nfft_finalize(&my_plan);
}
/**
 * reconstruct makes an 2d-adjoint-nfft for every slice
 */
static void reconstruct(char* filename,int N,int M,int Z, int weight ,fftw_complex *mem)
{
  int j,k,z;               /* some variables  */
  double weights;          /* store one weight temporary */
  double tmp;              /* tmp to read the obsolent z from the input file */
  double real,imag;        /* to read the real and imag part of a complex number */
  nfft_plan my_plan;       /* plan for the two dimensional nfft  */
  int my_N[2],my_n[2];     /* to init the nfft */
  FILE* fin;               /* input file  */
  FILE* fweight;           /* input file for the weights */

  /* initialise my_plan */
  my_N[0]=N; my_n[0]=ceil(N*1.2);
  my_N[1]=N; my_n[1]=ceil(N*1.2);
  nfft_init_guru(&my_plan, 2, my_N, M/Z, my_n, 6, PRE_PHI_HUT| PRE_PSI|
                        MALLOC_X| MALLOC_F_HAT| MALLOC_F|
                        FFTW_INIT| FFT_OUT_OF_PLACE,
                        FFTW_MEASURE| FFTW_DESTROY_INPUT);

  /* precompute lin psi if set */
  if(my_plan.nfft_flags & PRE_LIN_PSI)
    nfft_precompute_lin_psi(&my_plan);

  fin=fopen(filename,"r");

  for(z=0;z<Z;z++) {
    fweight=fopen("weights.dat","r");
    for(j=0;j<my_plan.M_total;j++)
    {
      fscanf(fweight,"%le ",&weights);
      fscanf(fin,"%le %le %le %le %le",
             &my_plan.x[2*j+0],&my_plan.x[2*j+1],&tmp,&real,&imag);
      my_plan.f[j] = real + _Complex_I*imag;
      if(weight)
        my_plan.f[j] = my_plan.f[j] * weights;
    }
    fclose(fweight);

    /* precompute psi if set just one time because the knots equal each slice */
    if(z==0 && my_plan.nfft_flags & PRE_PSI)
      nfft_precompute_psi(&my_plan);

    /* precompute full psi if set just one time because the knots equal each slice */
    if(z==0 && my_plan.nfft_flags & PRE_FULL_PSI)
      nfft_precompute_full_psi(&my_plan);

    /* compute the adjoint nfft */
    nfft_adjoint(&my_plan);

    for(k=0;k<my_plan.N_total;k++) {
      /* write every slice in the memory.
      here we make an fftshift direct */
      mem[(Z*N*N/2+z*N*N+ k)%(Z*N*N)] = my_plan.f_hat[k];
    }
  }
  fclose(fin);

  nfft_finalize(&my_plan);
}
static void reconstruct(char* filename,int N,int M,int iteration , int weight)
{
  int j,k,l;
  double t0, t1;
  double time,min_time,max_time,min_inh,max_inh;
  double t,real,imag;
  double w,epsilon=0.0000003;     /* epsilon is a the break criterium for
                                   the iteration */;
  mri_inh_3d_plan my_plan;
  solver_plan_complex my_iplan;
  FILE* fp,*fw,*fout_real,*fout_imag,*finh,*ftime;
  int my_N[3],my_n[3];
  int flags = PRE_PHI_HUT| PRE_PSI |MALLOC_X| MALLOC_F_HAT|
                      MALLOC_F| FFTW_INIT| FFT_OUT_OF_PLACE;
  unsigned infft_flags = CGNR | PRECOMPUTE_DAMP;

  double Ts;
  double W;
  int N3;
  int m=2;
  double sigma = 1.25;

  ftime=fopen("readout_time.dat","r");
  finh=fopen("inh.dat","r");

  min_time=INT_MAX; max_time=INT_MIN;
  for(j=0;j<M;j++)
  {
    fscanf(ftime,"%le ",&time);
    if(time<min_time)
      min_time = time;
    if(time>max_time)
      max_time = time;
  }

  fclose(ftime);

  Ts=(min_time+max_time)/2.0;


  min_inh=INT_MAX; max_inh=INT_MIN;
  for(j=0;j<N*N;j++)
  {
    fscanf(finh,"%le ",&w);
    if(w<min_inh)
      min_inh = w;
    if(w>max_inh)
      max_inh = w;
  }
  fclose(finh);

  N3=ceil((MAX(fabs(min_inh),fabs(max_inh))*(max_time-min_time)/2.0+m/(2*sigma))*4*sigma);
	/* N3 has to be even */
	if(N3%2!=0)
	  N3++;

  W= MAX(fabs(min_inh),fabs(max_inh))/(0.5-((double) m)/N3);

  my_N[0]=N;my_n[0]=ceil(N*sigma);
  my_N[1]=N; my_n[1]=ceil(N*sigma);
  my_N[2]=N3; my_n[2]=ceil(N3*sigma);

  /* initialise nfft */
  mri_inh_3d_init_guru(&my_plan, my_N, M, my_n, m, sigma, flags,
                      FFTW_MEASURE| FFTW_DESTROY_INPUT);

  if (weight)
    infft_flags = infft_flags | PRECOMPUTE_WEIGHT;

  /* initialise my_iplan, advanced */
  solver_init_advanced_complex(&my_iplan,(nfft_mv_plan_complex*)(&my_plan), infft_flags );

  /* get the weights */
  if(my_iplan.flags & PRECOMPUTE_WEIGHT)
  {
    fw=fopen("weights.dat","r");
    for(j=0;j<my_plan.M_total;j++)
    {
        fscanf(fw,"%le ",&my_iplan.w[j]);
    }
    fclose(fw);
  }

  /* get the damping factors */
  if(my_iplan.flags & PRECOMPUTE_DAMP)
  {
    for(j=0;j<N;j++){
      for(k=0;k<N;k++) {
        int j2= j-N/2;
        int k2= k-N/2;
        double r=sqrt(j2*j2+k2*k2);
        if(r>(double) N/2)
          my_iplan.w_hat[j*N+k]=0.0;
        else
          my_iplan.w_hat[j*N+k]=1.0;
      }
    }
  }

  fp=fopen(filename,"r");
  ftime=fopen("readout_time.dat","r");

  for(j=0;j<my_plan.M_total;j++)
  {
    fscanf(fp,"%le %le %le %le",&my_plan.plan.x[3*j+0],&my_plan.plan.x[3*j+1],&real,&imag);
    my_iplan.y[j]=real+ _Complex_I*imag;
    fscanf(ftime,"%le ",&my_plan.plan.x[3*j+2]);

    my_plan.plan.x[3*j+2] = (my_plan.plan.x[3*j+2]-Ts)*W/N3;
  }
  fclose(fp);
  fclose(ftime);


  finh=fopen("inh.dat","r");
  for(j=0;j<N*N;j++)
  {
    fscanf(finh,"%le ",&my_plan.w[j]);
    my_plan.w[j]/=W;
  }
  fclose(finh);


  if(my_plan.plan.flags & PRE_PSI) {
    nfft_precompute_psi(&my_plan.plan);
  }
  if(my_plan.plan.flags & PRE_FULL_PSI) {
      nfft_precompute_full_psi(&my_plan.plan);
  }

  /* init some guess */
  for(j=0;j<my_plan.N_total;j++)
  {
    my_iplan.f_hat_iter[j]=0.0;
  }

  t0 = nfft_clock_gettime_seconds();

  /* inverse trafo */
  solver_before_loop_complex(&my_iplan);
  for(l=0;l<iteration;l++)
  {
    /* break if dot_r_iter is smaller than epsilon*/
    if(my_iplan.dot_r_iter<epsilon)
    break;
    fprintf(stderr,"%e,  %i of %i\n",sqrt(my_iplan.dot_r_iter),
    l+1,iteration);
    solver_loop_one_step_complex(&my_iplan);
  }

  t1 = nfft_clock_gettime_seconds();
  t = t1-t0;

  fout_real=fopen("output_real.dat","w");
  fout_imag=fopen("output_imag.dat","w");

  for (j=0;j<N*N;j++) {
    /* Verschiebung wieder herausrechnen */
    my_iplan.f_hat_iter[j]*=cexp(-2.0*_Complex_I*M_PI*Ts*my_plan.w[j]*W);

    fprintf(fout_real,"%le ",creal(my_iplan.f_hat_iter[j]));
    fprintf(fout_imag,"%le ",cimag(my_iplan.f_hat_iter[j]));
  }

  fclose(fout_real);
  fclose(fout_imag);
  solver_finalize_complex(&my_iplan);
  mri_inh_3d_finalize(&my_plan);
}
Exemple #4
0
/** precomputation for fastsum */
void fastsum_precompute(fastsum_plan *ths)
{
  int j,k,t;
  int n_total;
  ticks t0, t1;

  ths->MEASURE_TIME_t[0] = 0.0;
  ths->MEASURE_TIME_t[1] = 0.0;
  ths->MEASURE_TIME_t[2] = 0.0;
  ths->MEASURE_TIME_t[3] = 0.0;

#ifdef MEASURE_TIME
  t0 = getticks();
#endif


  if (ths->flags & NEARFIELD_BOXES)
  {
    BuildBox(ths);
  }
  else
  {
    /** sort source knots */
    BuildTree(ths->d,0,ths->x,ths->alpha,ths->N_total);
  }

#ifdef MEASURE_TIME
  t1 = getticks();
  ths->MEASURE_TIME_t[3] += nfft_elapsed_seconds(t1,t0);
#endif


#ifdef MEASURE_TIME
  t0 = getticks();
#endif
  /** precompute spline values for near field*/
  if (!(ths->flags & EXACT_NEARFIELD))
  {
    if (ths->d==1)
      #pragma omp parallel for default(shared) private(k)
      for (k=-ths->Ad/2-2; k <= ths->Ad/2+2; k++)
        ths->Add[k+ths->Ad/2+2] = regkern1(ths->k, ths->eps_I*(double)k/ths->Ad*2, ths->p, ths->kernel_param, ths->eps_I, ths->eps_B);
    else
      #pragma omp parallel for default(shared) private(k)
      for (k=0; k <= ths->Ad+2; k++)
        ths->Add[k] = regkern3(ths->k, ths->eps_I*(double)k/ths->Ad, ths->p, ths->kernel_param, ths->eps_I, ths->eps_B);
  }
#ifdef MEASURE_TIME
  t1 = getticks();
  ths->MEASURE_TIME_t[0] += nfft_elapsed_seconds(t1,t0);
#endif


#ifdef MEASURE_TIME
  t0 = getticks();
#endif
  /** init NFFT plan for transposed transform in first step*/
  for (k=0; k<ths->mv1.M_total; k++)
    for (t=0; t<ths->mv1.d; t++)
      ths->mv1.x[ths->mv1.d*k+t] = - ths->x[ths->mv1.d*k+t];  /* note the factor -1 for transposed transform instead of adjoint*/

  /** precompute psi, the entries of the matrix B */
  if(ths->mv1.nfft_flags & PRE_LIN_PSI)
    nfft_precompute_lin_psi(&(ths->mv1));

  if(ths->mv1.nfft_flags & PRE_PSI)
    nfft_precompute_psi(&(ths->mv1));

  if(ths->mv1.nfft_flags & PRE_FULL_PSI)
    nfft_precompute_full_psi(&(ths->mv1));
#ifdef MEASURE_TIME
  t1 = getticks();
  ths->MEASURE_TIME_t[1] += nfft_elapsed_seconds(t1,t0);
#endif

  /** init Fourier coefficients */
  for(k=0; k<ths->mv1.M_total;k++)
    ths->mv1.f[k] = ths->alpha[k];

#ifdef MEASURE_TIME
  t0 = getticks();
#endif
  /** init NFFT plan for transform in third step*/
  for (j=0; j<ths->mv2.M_total; j++)
    for (t=0; t<ths->mv2.d; t++)
      ths->mv2.x[ths->mv2.d*j+t] = - ths->y[ths->mv2.d*j+t];  /* note the factor -1 for conjugated transform instead of standard*/

  /** precompute psi, the entries of the matrix B */
  if(ths->mv2.nfft_flags & PRE_LIN_PSI)
    nfft_precompute_lin_psi(&(ths->mv2));

  if(ths->mv2.nfft_flags & PRE_PSI)
    nfft_precompute_psi(&(ths->mv2));

  if(ths->mv2.nfft_flags & PRE_FULL_PSI)
    nfft_precompute_full_psi(&(ths->mv2));
#ifdef MEASURE_TIME
  t1 = getticks();
  ths->MEASURE_TIME_t[2] += nfft_elapsed_seconds(t1,t0);
#endif


#ifdef MEASURE_TIME
  t0 = getticks();
#endif
  /** precompute Fourier coefficients of regularised kernel*/
  n_total = 1;
  for (t=0; t<ths->d; t++)
    n_total *= ths->n;

  #pragma omp parallel for default(shared) private(j,k,t)
  for (j=0; j<n_total; j++)
  {
    if (ths->d==1)
      ths->b[j] = regkern1(ths->k, (double)j / (ths->n) - 0.5, ths->p, ths->kernel_param, ths->eps_I, ths->eps_B)/n_total;
    else
    {
      k=j;
      ths->b[j]=0.0;
      for (t=0; t<ths->d; t++)
      {
        ths->b[j] += ((double)(k % (ths->n)) / (ths->n) - 0.5) * ((double)(k % (ths->n)) / (ths->n) - 0.5);
        k = k / (ths->n);
      }
      ths->b[j] = regkern3(ths->k, sqrt(ths->b[j]), ths->p, ths->kernel_param, ths->eps_I, ths->eps_B)/n_total;
    }
  }

  nfft_fftshift_complex(ths->b, ths->mv1.d, ths->mv1.N);
  fftw_execute(ths->fft_plan);
  nfft_fftshift_complex(ths->b, ths->mv1.d, ths->mv1.N);
#ifdef MEASURE_TIME
  t1 = getticks();
  ths->MEASURE_TIME_t[0] += nfft_elapsed_seconds(t1,t0);
#endif
}
/**
 * reconstruct makes an inverse 3d-nfft
 */
static void reconstruct(char* filename,int N,int M,int Z,int iteration, int weight)
{
  int j,k,z,l;                  /* some variables  */
  double real,imag;             /* to read the real and imag part of a complex number */
  nfft_plan my_plan;            /* plan for the two dimensional nfft  */
  solver_plan_complex my_iplan;          /* plan for the two dimensional infft */
  FILE* fin;                    /* input file                         */
  FILE* fout_real;              /* output file (real part) */
  FILE* fout_imag;              /* output file (imag part) */
  int my_N[3],my_n[3];          /* to init the nfft */
  double epsilon=0.0000003;     /* tmp to read the obsolent z from 700.acs
                                   epsilon is a the break criterion for
                                   the iteration */
  unsigned infft_flags = CGNR | PRECOMPUTE_DAMP;  /* flags for the infft */

  /* initialise my_plan, specific.
     we don't precompute psi */
  my_N[0]=Z; my_n[0]=ceil(Z*1.2);
  my_N[1]=N; my_n[1]=ceil(N*1.2);
  my_N[2]=N; my_n[2]=ceil(N*1.2);
  nfft_init_guru(&my_plan, 3, my_N, M, my_n, 6,
                      PRE_PHI_HUT| PRE_PSI |MALLOC_X| MALLOC_F_HAT|
                      MALLOC_F| FFTW_INIT| FFT_OUT_OF_PLACE,
                      FFTW_MEASURE| FFTW_DESTROY_INPUT);

  /* precompute lin psi */
  if(my_plan.nfft_flags & PRE_LIN_PSI)
    nfft_precompute_lin_psi(&my_plan);

  if (weight)
    infft_flags = infft_flags | PRECOMPUTE_WEIGHT;

  /* initialise my_iplan, advanced */
  solver_init_advanced_complex(&my_iplan,(nfft_mv_plan_complex*)(&my_plan), infft_flags );

  /* get the weights */
  if(my_iplan.flags & PRECOMPUTE_WEIGHT)
  {
    fin=fopen("weights.dat","r");
    for(j=0;j<M;j++)
    {
      fscanf(fin,"%le ",&my_iplan.w[j]);
    }
    fclose(fin);
  }

  /* get the damping factors */
  if(my_iplan.flags & PRECOMPUTE_DAMP)
  {
    for(j=0;j<N;j++){
      for(k=0;k<N;k++) {
        for(z=0;z<N;z++) {
        int j2= j-N/2;
        int k2= k-N/2;
        int z2= z-N/2;
        double r=sqrt(j2*j2+k2*k2+z2*z2);
        if(r>(double) N/2)
          my_iplan.w_hat[z*N*N+j*N+k]=0.0;
        else
          my_iplan.w_hat[z*N*N+j*N+k]=1.0;
        }
      }
    }
  }

  /* open the input file */
  fin=fopen(filename,"r");

  /* open the output files */
  fout_real=fopen("output_real.dat","w");
  fout_imag=fopen("output_imag.dat","w");

  /* read x,y,freal and fimag from the knots */
  for(j=0;j<M;j++)
  {
    fscanf(fin,"%le %le %le %le %le ",&my_plan.x[3*j+1],&my_plan.x[3*j+2], &my_plan.x[3*j+0],
    &real,&imag);
    my_iplan.y[j] = real + _Complex_I*imag;
  }

  /* precompute psi */
  if(my_plan.nfft_flags & PRE_PSI)
    nfft_precompute_psi(&my_plan);

  /* precompute full psi */
  if(my_plan.nfft_flags & PRE_FULL_PSI)
    nfft_precompute_full_psi(&my_plan);

  /* init some guess */
  for(k=0;k<my_plan.N_total;k++)
    my_iplan.f_hat_iter[k]=0.0;

  /* inverse trafo */
  solver_before_loop_complex(&my_iplan);
  for(l=0;l<iteration;l++)
  {
    /* break if dot_r_iter is smaller than epsilon*/
    if(my_iplan.dot_r_iter<epsilon)
      break;
    fprintf(stderr,"%e,  %i of %i\n",sqrt(my_iplan.dot_r_iter),
    l+1,iteration);
    solver_loop_one_step_complex(&my_iplan);
  }

  for(l=0;l<Z;l++)
  {
    for(k=0;k<N*N;k++)
    {
      /* write every Layer in the files */
      fprintf(fout_real,"%le ",creal(my_iplan.f_hat_iter[ k+N*N*l ]));
      fprintf(fout_imag,"%le ",cimag(my_iplan.f_hat_iter[ k+N*N*l ]));
    }
    fprintf(fout_real,"\n");
    fprintf(fout_imag,"\n");
  }

  fclose(fout_real);
  fclose(fout_imag);

  solver_finalize_complex(&my_iplan);
  nfft_finalize(&my_plan);
}
Exemple #6
0
/** inverse NFFT-based mpolar FFT */
static int inverse_mpolar_fft(fftw_complex *f, int T, int R, fftw_complex *f_hat, int NN, int max_i, int m)
{
  ticks t0, t1;
  int j,k;                              /**< index for nodes and freqencies   */
  nfft_plan my_nfft_plan;               /**< plan for the nfft-2D             */
  solver_plan_complex my_infft_plan;             /**< plan for the inverse nfft        */

  double *x, *w;                        /**< knots and associated weights     */
  int l;                                /**< index for iterations             */

  int N[2],n[2];
  int M;                                /**< number of knots                  */

  N[0]=NN; n[0]=2*N[0];                 /**< oversampling factor sigma=2      */
  N[1]=NN; n[1]=2*N[1];                 /**< oversampling factor sigma=2      */

  x = (double *)nfft_malloc(5*T*R/2*(sizeof(double)));
  if (x==NULL)
    return -1;

  w = (double *)nfft_malloc(5*T*R/4*(sizeof(double)));
  if (w==NULL)
    return -1;

  /** init two dimensional NFFT plan */
  M=mpolar_grid(T,R,x,w);
  nfft_init_guru(&my_nfft_plan, 2, N, M, n, m,
                  PRE_PHI_HUT| PRE_PSI| MALLOC_X | MALLOC_F_HAT| MALLOC_F| FFTW_INIT | FFT_OUT_OF_PLACE,
                  FFTW_MEASURE| FFTW_DESTROY_INPUT);

  /** init two dimensional infft plan */
    solver_init_advanced_complex(&my_infft_plan,(nfft_mv_plan_complex*)(&my_nfft_plan), CGNR | PRECOMPUTE_WEIGHT );

  /** init nodes, given samples and weights */
  for(j=0;j<my_nfft_plan.M_total;j++)
  {
    my_nfft_plan.x[2*j+0] = x[2*j+0];
    my_nfft_plan.x[2*j+1] = x[2*j+1];
    my_infft_plan.y[j]    = f[j];
    my_infft_plan.w[j]    = w[j];
  }

  /** precompute psi, the entries of the matrix B */
  if(my_nfft_plan.nfft_flags & PRE_LIN_PSI)
    nfft_precompute_lin_psi(&my_nfft_plan);

  if(my_nfft_plan.nfft_flags & PRE_PSI)
    nfft_precompute_psi(&my_nfft_plan);

  if(my_nfft_plan.nfft_flags & PRE_FULL_PSI)
    nfft_precompute_full_psi(&my_nfft_plan);


 /** initialise damping factors */
 if(my_infft_plan.flags & PRECOMPUTE_DAMP)
   for(j=0;j<my_nfft_plan.N[0];j++)
     for(k=0;k<my_nfft_plan.N[1];k++)
     {
        my_infft_plan.w_hat[j*my_nfft_plan.N[1]+k]=
          (sqrt(pow(j-my_nfft_plan.N[0]/2,2)+pow(k-my_nfft_plan.N[1]/2,2))>(my_nfft_plan.N[0]/2)?0:1);
     }

  /** initialise some guess f_hat_0 */
  for(k=0;k<my_nfft_plan.N_total;k++)
      my_infft_plan.f_hat_iter[k] = 0.0 + _Complex_I*0.0;

  t0 = getticks();

  /** solve the system */
  solver_before_loop_complex(&my_infft_plan);

  if (max_i<1)
  {
    l=1;
    for(k=0;k<my_nfft_plan.N_total;k++)
      my_infft_plan.f_hat_iter[k] = my_infft_plan.p_hat_iter[k];
  }
  else
  {
    for(l=1;l<=max_i;l++)
    {
      solver_loop_one_step_complex(&my_infft_plan);
    }
  }

  t1 = getticks();
  GLOBAL_elapsed_time = nfft_elapsed_seconds(t1,t0);

  /** copy result */
  for(k=0;k<my_nfft_plan.N_total;k++)
    f_hat[k] = my_infft_plan.f_hat_iter[k];

  /** finalise the plans and free the variables */
  solver_finalize_complex(&my_infft_plan);
  nfft_finalize(&my_nfft_plan);
  nfft_free(x);
  nfft_free(w);

  return EXIT_SUCCESS;
}
Exemple #7
0
/** NFFT-based mpolar FFT */
static int mpolar_fft(fftw_complex *f_hat, int NN, fftw_complex *f, int T, int R, int m)
{
  ticks t0, t1;
  int j,k;                              /**< index for nodes and freqencies   */
  nfft_plan my_nfft_plan;               /**< plan for the nfft-2D             */

  double *x, *w;                        /**< knots and associated weights     */

  int N[2],n[2];
  int M;                                /**< number of knots                  */

  N[0]=NN; n[0]=2*N[0];                 /**< oversampling factor sigma=2      */
  N[1]=NN; n[1]=2*N[1];                 /**< oversampling factor sigma=2      */

  x = (double *)nfft_malloc(5*T*R/2*(sizeof(double)));
  if (x==NULL)
    return -1;

  w = (double *)nfft_malloc(5*T*R/4*(sizeof(double)));
  if (w==NULL)
    return -1;

  /** init two dimensional NFFT plan */
  M=mpolar_grid(T,R,x,w);
  nfft_init_guru(&my_nfft_plan, 2, N, M, n, m,
                  PRE_PHI_HUT| PRE_PSI| MALLOC_X | MALLOC_F_HAT| MALLOC_F| FFTW_INIT | FFT_OUT_OF_PLACE,
                  FFTW_MEASURE| FFTW_DESTROY_INPUT);

  /** init nodes from mpolar grid*/

  for(j=0;j<my_nfft_plan.M_total;j++)
  {
    my_nfft_plan.x[2*j+0] = x[2*j+0];
    my_nfft_plan.x[2*j+1] = x[2*j+1];
  }

  /** precompute psi, the entries of the matrix B */
  if(my_nfft_plan.nfft_flags & PRE_LIN_PSI)
    nfft_precompute_lin_psi(&my_nfft_plan);

  if(my_nfft_plan.nfft_flags & PRE_PSI)
    nfft_precompute_psi(&my_nfft_plan);

  if(my_nfft_plan.nfft_flags & PRE_FULL_PSI)
    nfft_precompute_full_psi(&my_nfft_plan);

  /** init Fourier coefficients from given image */
  for(k=0;k<my_nfft_plan.N_total;k++)
    my_nfft_plan.f_hat[k] = f_hat[k];

  t0 = getticks();

  /** NFFT-2D */
  nfft_trafo(&my_nfft_plan);

  t1 = getticks();
  GLOBAL_elapsed_time = nfft_elapsed_seconds(t1,t0);

  /** copy result */
  for(j=0;j<my_nfft_plan.M_total;j++)
    f[j] = my_nfft_plan.f[j];

  /** finalise the plans and free the variables */
  nfft_finalize(&my_nfft_plan);
  nfft_free(x);
  nfft_free(w);

  return EXIT_SUCCESS;
}
Exemple #8
0
/** computes the inverse discrete Radon transform of Rf
 *  on the grid given by gridfcn() with T angles and R offsets
 *  by a NFFT-based CG-type algorithm
 */
int Inverse_Radon_trafo(int (*gridfcn)(), int T, int S, double *Rf, int NN, double *f, int max_i)
{
  int j,k;                              /**< index for nodes and freqencies   */
  nfft_plan my_nfft_plan;               /**< plan for the nfft-2D             */
  solver_plan_complex my_infft_plan;             /**< plan for the inverse nfft        */

  fftw_complex *fft;                    /**< variable for the fftw-1Ds        */
  fftw_plan my_fftw_plan;               /**< plan for the fftw-1Ds            */

  int t,r;                              /**< index for directions and offsets */
  double *x, *w;                        /**< knots and associated weights     */
  int l;                                /**< index for iterations             */

  int N[2],n[2];
  int M=T*S;

  N[0]=NN; n[0]=2*N[0];
  N[1]=NN; n[1]=2*N[1];

  fft = (fftw_complex *)nfft_malloc(S*sizeof(fftw_complex));
  my_fftw_plan = fftw_plan_dft_1d(S,fft,fft,FFTW_FORWARD,FFTW_MEASURE);

  x = (double *)nfft_malloc(2*T*S*(sizeof(double)));
  if (x==NULL)
    return -1;

  w = (double *)nfft_malloc(T*S*(sizeof(double)));
  if (w==NULL)
    return -1;

  /** init two dimensional NFFT plan */
  nfft_init_guru(&my_nfft_plan, 2, N, M, n, 4,
                  PRE_PHI_HUT| PRE_PSI| MALLOC_X | MALLOC_F_HAT| MALLOC_F| FFTW_INIT | FFT_OUT_OF_PLACE,
                  FFTW_MEASURE| FFTW_DESTROY_INPUT);

  /** init two dimensional infft plan */
  solver_init_advanced_complex(&my_infft_plan,(nfft_mv_plan_complex*)(&my_nfft_plan), CGNR | PRECOMPUTE_WEIGHT);

  /** init nodes and weights of grid*/
  gridfcn(T,S,x,w);
  for(j=0;j<my_nfft_plan.M_total;j++)
  {
    my_nfft_plan.x[2*j+0] = x[2*j+0];
    my_nfft_plan.x[2*j+1] = x[2*j+1];
    if (j%S)
      my_infft_plan.w[j]    = w[j];
    else
      my_infft_plan.w[j]    = 0.0;
  }

  /** precompute psi, the entries of the matrix B */
  if(my_nfft_plan.nfft_flags & PRE_LIN_PSI)
    nfft_precompute_lin_psi(&my_nfft_plan);

  if(my_nfft_plan.nfft_flags & PRE_PSI)
    nfft_precompute_psi(&my_nfft_plan);

  if(my_nfft_plan.nfft_flags & PRE_FULL_PSI)
    nfft_precompute_full_psi(&my_nfft_plan);

  /** compute 1D-ffts and init given samples and weights */
  for(t=0; t<T; t++)
  {
/*    for(r=0; r<R/2; r++)
       fft[r] = cexp(I*KPI*r)*Rf[t*R+(r+R/2)];
      for(r=0; r<R/2; r++)
       fft[r+R/2] = cexp(I*KPI*r)*Rf[t*R+r];
 */

    for(r=0; r<S; r++)
      fft[r] = Rf[t*S+r] + _Complex_I*0.0;

    nfft_fftshift_complex(fft, 1, &S);
    fftw_execute(my_fftw_plan);
    nfft_fftshift_complex(fft, 1, &S);

    my_infft_plan.y[t*S] = 0.0;
    for(r=-S/2+1; r<S/2; r++)
      my_infft_plan.y[t*S+(r+S/2)] = fft[r+S/2]/KERNEL(r);
  }

  /** initialise some guess f_hat_0 */
  for(k=0;k<my_nfft_plan.N_total;k++)
    my_infft_plan.f_hat_iter[k] = 0.0 + _Complex_I*0.0;

  /** solve the system */
  solver_before_loop_complex(&my_infft_plan);

  if (max_i<1)
  {
    l=1;
    for(k=0;k<my_nfft_plan.N_total;k++)
      my_infft_plan.f_hat_iter[k] = my_infft_plan.p_hat_iter[k];
  }
  else
  {
    for(l=1;l<=max_i;l++)
    {
      solver_loop_one_step_complex(&my_infft_plan);
      /*if (sqrt(my_infft_plan.dot_r_iter)<=1e-12) break;*/
    }
  }
  /*printf("after %d iteration(s): weighted 2-norm of original residual vector = %g\n",l-1,sqrt(my_infft_plan.dot_r_iter));*/

  /** copy result */
  for(k=0;k<my_nfft_plan.N_total;k++)
    f[k] = creal(my_infft_plan.f_hat_iter[k]);

  /** finalise the plans and free the variables */
  fftw_destroy_plan(my_fftw_plan);
  nfft_free(fft);
  solver_finalize_complex(&my_infft_plan);
  nfft_finalize(&my_nfft_plan);
  nfft_free(x);
  nfft_free(w);
  return 0;
}