NFFT Logo 3.2.4
fastsumS2.c
1 /*
2  * Copyright (c) 2002, 2012 Jens Keiner, Stefan Kunis, Daniel Potts
3  *
4  * This program is free software; you can redistribute it and/or modify it under
5  * the terms of the GNU General Public License as published by the Free Software
6  * Foundation; either version 2 of the License, or (at your option) any later
7  * version.
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
12  * details.
13  *
14  * You should have received a copy of the GNU General Public License along with
15  * this program; if not, write to the Free Software Foundation, Inc., 51
16  * Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
17  */
18 
19 /* $Id: fastsumS2.c 3896 2012-10-10 12:19:26Z tovo $ */
20 
26 #include "config.h"
27 
28 /* standard headers */
29 #include <stdio.h>
30 #include <stdlib.h>
31 #include <math.h>
32 #include <float.h>
33 #ifdef HAVE_COMPLEX_H
34 #include <complex.h>
35 #endif
36 
37 /* NFFT3 header */
38 #include "nfft3.h"
39 
40 /* NFFT3 utilities */
41 #include "nfft3util.h"
42 #include "infft.h"
43 
44 /* Fourier-Legendre coefficients for Abel-Poisson kernel */
45 #define SYMBOL_ABEL_POISSON(k,h) (pow(h,k))
46 
47 /* Fourier-Legendre coefficients for singularity kernel */
48 #define SYMBOL_SINGULARITY(k,h) ((2.0/(2*k+1))*pow(h,k))
49 
50 /* Flags for the different kernel functions */
51 
53 #define KT_ABEL_POISSON (0)
54 
55 #define KT_SINGULARITY (1)
56 
57 #define KT_LOC_SUPP (2)
58 
59 #define KT_GAUSSIAN (3)
60 
62 enum pvalue {NO = 0, YES = 1, BOTH = 2};
63 
78 static inline double innerProduct(const double phi1, const double theta1,
79  const double phi2, const double theta2)
80 {
81  double pi2theta1 = PI2*theta1, pi2theta2 = PI2*theta2;
82  return (cos(pi2theta1)*cos(pi2theta2)
83  + sin(pi2theta1)*sin(pi2theta2)*cos(PI2*(phi1-phi2)));
84 }
85 
97 static inline double poissonKernel(const double x, const double h)
98 {
99  return (1.0/(PI4))*((1.0-h)*(1.0+h))/pow(sqrt(1.0-2.0*h*x+h*h),3.0);
100 }
101 
113 static inline double singularityKernel(const double x, const double h)
114 {
115  return (1.0/(PI2))/sqrt(1.0-2.0*h*x+h*h);
116 }
117 
131 static inline double locallySupportedKernel(const double x, const double h,
132  const double lambda)
133 {
134  return (x<=h)?(0.0):(pow((x-h),lambda));
135 }
136 
149 static inline double gaussianKernel(const double x, const double sigma)
150 {
151  return exp(2.0*sigma*(x-1.0));
152 }
153 
164 int main (int argc, char **argv)
165 {
166  double **p; /* The array containing the parameter sets *
167  * for the kernel functions */
168  int *m; /* The array containing the cut-off degrees M */
169  int **ld; /* The array containing the numbers of source *
170  * and target nodes, L and D */
171  int ip; /* Index variable for p */
172  int im; /* Index variable for m */
173  int ild; /* Index variable for l */
174  int ipp; /* Index for kernel parameters */
175  int ip_max; /* The maximum index for p */
176  int im_max; /* The maximum index for m */
177  int ild_max; /* The maximum index for l */
178  int ipp_max; /* The maximum index for ip */
179  int tc_max; /* The number of testcases */
180  int m_max; /* The maximum cut-off degree M for the *
181  * current dataset */
182  int l_max; /* The maximum number of source nodes L for *
183  * the current dataset */
184  int d_max; /* The maximum number of target nodes D for *
185  * the current dataset */
186  long ld_max_prec; /* The maximum number of source and target *
187  * nodes for precomputation multiplied */
188  long l_max_prec; /* The maximum number of source nodes for *
189  * precomputation */
190  int tc; /* Index variable for testcases */
191  int kt; /* The kernel function */
192  int cutoff; /* The current NFFT cut-off parameter */
193  double threshold; /* The current NFSFT threshold parameter */
194  double t_d; /* Time for direct algorithm in seconds */
195  double t_dp; /* Time for direct algorithm with *
196  precomputation in seconds */
197  double t_fd; /* Time for fast direct algorithm in seconds */
198  double t_f; /* Time for fast algorithm in seconds */
199  double temp; /* */
200  double err_f; /* Error E_infty for fast algorithm */
201  double err_fd; /* Error E_\infty for fast direct algorithm */
202  ticks t0, t1; /* */
203  int precompute = NO; /* */
204  fftw_complex *ptr; /* */
205  double* steed; /* */
206  fftw_complex *b; /* The weights (b_l)_{l=0}^{L-1} */
207  fftw_complex *f_hat; /* The spherical Fourier coefficients */
208  fftw_complex *a; /* The Fourier-Legendre coefficients */
209  double *xi; /* Target nodes */
210  double *eta; /* Source nodes */
211  fftw_complex *f_m; /* Approximate function values */
212  fftw_complex *f; /* Exact function values */
213  fftw_complex *prec = NULL; /* */
214  nfsft_plan plan; /* NFSFT plan */
215  nfsft_plan plan_adjoint; /* adjoint NFSFT plan */
216  int i; /* */
217  int k; /* */
218  int n; /* */
219  int d; /* */
220  int l; /* */
221  int use_nfsft; /* */
222  int use_nfft; /* */
223  int use_fpt; /* */
224  int rinc; /* */
225  double constant; /* */
226 
227  /* Read the number of testcases. */
228  fscanf(stdin,"testcases=%d\n",&tc_max);
229  fprintf(stdout,"%d\n",tc_max);
230 
231  /* Process each testcase. */
232  for (tc = 0; tc < tc_max; tc++)
233  {
234  /* Check if the fast transform shall be used. */
235  fscanf(stdin,"nfsft=%d\n",&use_nfsft);
236  fprintf(stdout,"%d\n",use_nfsft);
237  if (use_nfsft != NO)
238  {
239  /* Check if the NFFT shall be used. */
240  fscanf(stdin,"nfft=%d\n",&use_nfft);
241  fprintf(stdout,"%d\n",use_nfft);
242  if (use_nfft != NO)
243  {
244  /* Read the cut-off parameter. */
245  fscanf(stdin,"cutoff=%d\n",&cutoff);
246  fprintf(stdout,"%d\n",cutoff);
247  }
248  else
249  {
250  /* TODO remove this */
251  /* Initialize unused variable with dummy value. */
252  cutoff = 1;
253  }
254  /* Check if the fast polynomial transform shall be used. */
255  fscanf(stdin,"fpt=%d\n",&use_fpt);
256  fprintf(stdout,"%d\n",use_fpt);
257  /* Read the NFSFT threshold parameter. */
258  fscanf(stdin,"threshold=%lf\n",&threshold);
259  fprintf(stdout,"%lf\n",threshold);
260  }
261  else
262  {
263  /* TODO remove this */
264  /* Set dummy values. */
265  cutoff = 3;
266  threshold = 1000000000000.0;
267  }
268 
269  /* Initialize bandwidth bound. */
270  m_max = 0;
271  /* Initialize source nodes bound. */
272  l_max = 0;
273  /* Initialize target nodes bound. */
274  d_max = 0;
275  /* Initialize source nodes bound for precomputation. */
276  l_max_prec = 0;
277  /* Initialize source and target nodes bound for precomputation. */
278  ld_max_prec = 0;
279 
280  /* Read the kernel type. This is one of KT_ABEL_POISSON, KT_SINGULARITY,
281  * KT_LOC_SUPP and KT_GAUSSIAN. */
282  fscanf(stdin,"kernel=%d\n",&kt);
283  fprintf(stdout,"%d\n",kt);
284 
285  /* Read the number of parameter sets. */
286  fscanf(stdin,"parameter_sets=%d\n",&ip_max);
287  fprintf(stdout,"%d\n",ip_max);
288 
289  /* Allocate memory for pointers to parameter sets. */
290  p = (double**) nfft_malloc(ip_max*sizeof(double*));
291 
292  /* We now read in the parameter sets. */
293 
294  /* Read number of parameters. */
295  fscanf(stdin,"parameters=%d\n",&ipp_max);
296  fprintf(stdout,"%d\n",ipp_max);
297 
298  for (ip = 0; ip < ip_max; ip++)
299  {
300  /* Allocate memory for the parameters. */
301  p[ip] = (double*) nfft_malloc(ipp_max*sizeof(double));
302 
303  /* Read the parameters. */
304  for (ipp = 0; ipp < ipp_max; ipp++)
305  {
306  /* Read the next parameter. */
307  fscanf(stdin,"%lf\n",&p[ip][ipp]);
308  fprintf(stdout,"%lf\n",p[ip][ipp]);
309  }
310  }
311 
312  /* Read the number of cut-off degrees. */
313  fscanf(stdin,"bandwidths=%d\n",&im_max);
314  fprintf(stdout,"%d\n",im_max);
315  m = (int*) nfft_malloc(im_max*sizeof(int));
316 
317  /* Read the cut-off degrees. */
318  for (im = 0; im < im_max; im++)
319  {
320  /* Read cut-off degree. */
321  fscanf(stdin,"%d\n",&m[im]);
322  fprintf(stdout,"%d\n",m[im]);
323  m_max = NFFT_MAX(m_max,m[im]);
324  }
325 
326  /* Read number of node specifications. */
327  fscanf(stdin,"node_sets=%d\n",&ild_max);
328  fprintf(stdout,"%d\n",ild_max);
329  ld = (int**) nfft_malloc(ild_max*sizeof(int*));
330 
331  /* Read the run specification. */
332  for (ild = 0; ild < ild_max; ild++)
333  {
334  /* Allocate memory for the run parameters. */
335  ld[ild] = (int*) nfft_malloc(5*sizeof(int));
336 
337  /* Read number of source nodes. */
338  fscanf(stdin,"L=%d ",&ld[ild][0]);
339  fprintf(stdout,"%d\n",ld[ild][0]);
340  l_max = NFFT_MAX(l_max,ld[ild][0]);
341 
342  /* Read number of target nodes. */
343  fscanf(stdin,"D=%d ",&ld[ild][1]);
344  fprintf(stdout,"%d\n",ld[ild][1]);
345  d_max = NFFT_MAX(d_max,ld[ild][1]);
346 
347  /* Determine whether direct and fast algorithm shall be compared. */
348  fscanf(stdin,"compare=%d ",&ld[ild][2]);
349  fprintf(stdout,"%d\n",ld[ild][2]);
350 
351  /* Check if precomputation for the direct algorithm is used. */
352  if (ld[ild][2] == YES)
353  {
354  /* Read whether the precomputed version shall also be used. */
355  fscanf(stdin,"precomputed=%d\n",&ld[ild][3]);
356  fprintf(stdout,"%d\n",ld[ild][3]);
357 
358  /* Read the number of repetitions over which measurements are
359  * averaged. */
360  fscanf(stdin,"repetitions=%d\n",&ld[ild][4]);
361  fprintf(stdout,"%d\n",ld[ild][4]);
362 
363  /* Update ld_max_prec and l_max_prec. */
364  if (ld[ild][3] == YES)
365  {
366  /* Update ld_max_prec. */
367  ld_max_prec = NFFT_MAX(ld_max_prec,ld[ild][0]*ld[ild][1]);
368  /* Update l_max_prec. */
369  l_max_prec = NFFT_MAX(l_max_prec,ld[ild][0]);
370  /* Turn on the precomputation for the direct algorithm. */
371  precompute = YES;
372  }
373  }
374  else
375  {
376  /* Set default value for the number of repetitions. */
377  ld[ild][4] = 1;
378  }
379  }
380 
381  /* Allocate memory for data structures. */
382  b = (fftw_complex*) nfft_malloc(l_max*sizeof(fftw_complex));
383  eta = (double*) nfft_malloc(2*l_max*sizeof(double));
384  f_hat = (fftw_complex*) nfft_malloc(NFSFT_F_HAT_SIZE(m_max)*sizeof(fftw_complex));
385  a = (fftw_complex*) nfft_malloc((m_max+1)*sizeof(fftw_complex));
386  xi = (double*) nfft_malloc(2*d_max*sizeof(double));
387  f_m = (fftw_complex*) nfft_malloc(d_max*sizeof(fftw_complex));
388  f = (fftw_complex*) nfft_malloc(d_max*sizeof(fftw_complex));
389 
390  /* Allocate memory for precomputed data. */
391  if (precompute == YES)
392  {
393  prec = (fftw_complex*) nfft_malloc(ld_max_prec*sizeof(fftw_complex));
394  }
395 
396  /* Generate random source nodes and weights. */
397  for (l = 0; l < l_max; l++)
398  {
399  b[l] = (((double)rand())/RAND_MAX) - 0.5;
400  eta[2*l] = (((double)rand())/RAND_MAX) - 0.5;
401  eta[2*l+1] = acos(2.0*(((double)rand())/RAND_MAX) - 1.0)/(PI2);
402  }
403 
404  /* Generate random target nodes. */
405  for (d = 0; d < d_max; d++)
406  {
407  xi[2*d] = (((double)rand())/RAND_MAX) - 0.5;
408  xi[2*d+1] = acos(2.0*(((double)rand())/RAND_MAX) - 1.0)/(PI2);
409  }
410 
411  /* Do precomputation. */
412  nfsft_precompute(m_max,threshold,
413  ((use_nfsft==NO)?(NFSFT_NO_FAST_ALGORITHM):(0U/*NFSFT_NO_DIRECT_ALGORITHM*/)), 0U);
414 
415  /* Process all parameter sets. */
416  for (ip = 0; ip < ip_max; ip++)
417  {
418  /* Compute kernel coeffcients up to the maximum cut-off degree m_max. */
419  switch (kt)
420  {
421  case KT_ABEL_POISSON:
422  /* Compute Fourier-Legendre coefficients for the Poisson kernel. */
423  for (k = 0; k <= m_max; k++)
424  a[k] = SYMBOL_ABEL_POISSON(k,p[ip][0]);
425  break;
426 
427  case KT_SINGULARITY:
428  /* Compute Fourier-Legendre coefficients for the singularity
429  * kernel. */
430  for (k = 0; k <= m_max; k++)
431  a[k] = SYMBOL_SINGULARITY(k,p[ip][0]);
432  break;
433 
434  case KT_LOC_SUPP:
435  /* Compute Fourier-Legendre coefficients for the locally supported
436  * kernel. */
437  a[0] = 1.0;
438  if (1 <= m_max)
439  a[1] = ((p[ip][1]+1+p[ip][0])/(p[ip][1]+2.0))*a[0];
440  for (k = 2; k <= m_max; k++)
441  a[k] = (1.0/(k+p[ip][1]+1))*((2*k-1)*p[ip][0]*a[k-1] -
442  (k-p[ip][1]-2)*a[k-2]);
443  break;
444 
445  case KT_GAUSSIAN:
446  /* Fourier-Legendre coefficients */
447  steed = (double*) nfft_malloc((m_max+1)*sizeof(double));
448  nfft_smbi(2.0*p[ip][0],0.5,m_max+1,2,steed);
449  for (k = 0; k <= m_max; k++)
450  a[k] = PI2*(sqrt(PI/p[ip][0]))*steed[k];
451 
452  nfft_free(steed);
453  break;
454  }
455 
456  /* Normalize Fourier-Legendre coefficients. */
457  for (k = 0; k <= m_max; k++)
458  a[k] *= (2*k+1)/(PI4);
459 
460  /* Process all node sets. */
461  for (ild = 0; ild < ild_max; ild++)
462  {
463  /* Check if the fast algorithm shall be used. */
464  if (ld[ild][2] != NO)
465  {
466  /* Check if the direct algorithm with precomputation should be
467  * tested. */
468  if (ld[ild][3] != NO)
469  {
470  /* Get pointer to start of data. */
471  ptr = prec;
472  /* Calculate increment from one row to the next. */
473  rinc = l_max_prec-ld[ild][0];
474 
475  /* Process al target nodes. */
476  for (d = 0; d < ld[ild][1]; d++)
477  {
478  /* Process all source nodes. */
479  for (l = 0; l < ld[ild][0]; l++)
480  {
481  /* Compute inner product between current source and target
482  * node. */
483  temp = innerProduct(eta[2*l],eta[2*l+1],xi[2*d],xi[2*d+1]);
484 
485  /* Switch by the kernel type. */
486  switch (kt)
487  {
488  case KT_ABEL_POISSON:
489  /* Evaluate the Poisson kernel for the current value. */
490  *ptr++ = poissonKernel(temp,p[ip][0]);
491  break;
492 
493  case KT_SINGULARITY:
494  /* Evaluate the singularity kernel for the current
495  * value. */
496  *ptr++ = singularityKernel(temp,p[ip][0]);
497  break;
498 
499  case KT_LOC_SUPP:
500  /* Evaluate the localized kernel for the current
501  * value. */
502  *ptr++ = locallySupportedKernel(temp,p[ip][0],p[ip][1]);
503  break;
504 
505  case KT_GAUSSIAN:
506  /* Evaluate the spherical Gaussian kernel for the current
507  * value. */
508  *ptr++ = gaussianKernel(temp,p[ip][0]);
509  break;
510  }
511  }
512  /* Increment pointer for next row. */
513  ptr += rinc;
514  }
515 
516  /* Initialize cumulative time variable. */
517  t_dp = 0.0;
518 
519  /* Initialize time measurement. */
520  t0 = getticks();
521 
522  /* Cycle through all runs. */
523  for (i = 0; i < ld[ild][4]; i++)
524  {
525 
526  /* Reset pointer to start of precomputed data. */
527  ptr = prec;
528  /* Calculate increment from one row to the next. */
529  rinc = l_max_prec-ld[ild][0];
530 
531  /* Check if the localized kernel is used. */
532  if (kt == KT_LOC_SUPP)
533  {
534  /* Perform final summation */
535 
536  /* Calculate the multiplicative constant. */
537  constant = ((p[ip][1]+1)/(PI2*pow(1-p[ip][0],p[ip][1]+1)));
538 
539  /* Process all target nodes. */
540  for (d = 0; d < ld[ild][1]; d++)
541  {
542  /* Initialize function value. */
543  f[d] = 0.0;
544 
545  /* Process all source nodes. */
546  for (l = 0; l < ld[ild][0]; l++)
547  f[d] += b[l]*(*ptr++);
548 
549  /* Multiply with the constant. */
550  f[d] *= constant;
551 
552  /* Proceed to next row. */
553  ptr += rinc;
554  }
555  }
556  else
557  {
558  /* Process all target nodes. */
559  for (d = 0; d < ld[ild][1]; d++)
560  {
561  /* Initialize function value. */
562  f[d] = 0.0;
563 
564  /* Process all source nodes. */
565  for (l = 0; l < ld[ild][0]; l++)
566  f[d] += b[l]*(*ptr++);
567 
568  /* Proceed to next row. */
569  ptr += rinc;
570  }
571  }
572  }
573 
574  /* Calculate the time needed. */
575  t1 = getticks();
576  t_dp = nfft_elapsed_seconds(t1,t0);
577 
578  /* Calculate average time needed. */
579  t_dp = t_dp/((double)ld[ild][4]);
580  }
581  else
582  {
583  /* Initialize cumulative time variable with dummy value. */
584  t_dp = -1.0;
585  }
586 
587  /* Initialize cumulative time variable. */
588  t_d = 0.0;
589 
590  /* Initialize time measurement. */
591  t0 = getticks();
592 
593  /* Cycle through all runs. */
594  for (i = 0; i < ld[ild][4]; i++)
595  {
596  /* Switch by the kernel type. */
597  switch (kt)
598  {
599  case KT_ABEL_POISSON:
600 
601  /* Process all target nodes. */
602  for (d = 0; d < ld[ild][1]; d++)
603  {
604  /* Initialize function value. */
605  f[d] = 0.0;
606 
607  /* Process all source nodes. */
608  for (l = 0; l < ld[ild][0]; l++)
609  {
610  /* Compute the inner product for the current source and
611  * target nodes. */
612  temp = innerProduct(eta[2*l],eta[2*l+1],xi[2*d],xi[2*d+1]);
613 
614  /* Evaluate the Poisson kernel for the current value and add
615  * to the result. */
616  f[d] += b[l]*poissonKernel(temp,p[ip][0]);
617  }
618  }
619  break;
620 
621  case KT_SINGULARITY:
622  /* Process all target nodes. */
623  for (d = 0; d < ld[ild][1]; d++)
624  {
625  /* Initialize function value. */
626  f[d] = 0.0;
627 
628  /* Process all source nodes. */
629  for (l = 0; l < ld[ild][0]; l++)
630  {
631  /* Compute the inner product for the current source and
632  * target nodes. */
633  temp = innerProduct(eta[2*l],eta[2*l+1],xi[2*d],xi[2*d+1]);
634 
635  /* Evaluate the Poisson kernel for the current value and add
636  * to the result. */
637  f[d] += b[l]*singularityKernel(temp,p[ip][0]);
638  }
639  }
640  break;
641 
642  case KT_LOC_SUPP:
643  /* Calculate the multiplicative constant. */
644  constant = ((p[ip][1]+1)/(PI2*pow(1-p[ip][0],p[ip][1]+1)));
645 
646  /* Process all target nodes. */
647  for (d = 0; d < ld[ild][1]; d++)
648  {
649  /* Initialize function value. */
650  f[d] = 0.0;
651 
652  /* Process all source nodes. */
653  for (l = 0; l < ld[ild][0]; l++)
654  {
655  /* Compute the inner product for the current source and
656  * target nodes. */
657  temp = innerProduct(eta[2*l],eta[2*l+1],xi[2*d],xi[2*d+1]);
658 
659  /* Evaluate the Poisson kernel for the current value and add
660  * to the result. */
661  f[d] += b[l]*locallySupportedKernel(temp,p[ip][0],p[ip][1]);
662  }
663 
664  /* Multiply result with constant. */
665  f[d] *= constant;
666  }
667  break;
668 
669  case KT_GAUSSIAN:
670  /* Process all target nodes. */
671  for (d = 0; d < ld[ild][1]; d++)
672  {
673  /* Initialize function value. */
674  f[d] = 0.0;
675 
676  /* Process all source nodes. */
677  for (l = 0; l < ld[ild][0]; l++)
678  {
679  /* Compute the inner product for the current source and
680  * target nodes. */
681  temp = innerProduct(eta[2*l],eta[2*l+1],xi[2*d],xi[2*d+1]);
682  /* Evaluate the Poisson kernel for the current value and add
683  * to the result. */
684  f[d] += b[l]*gaussianKernel(temp,p[ip][0]);
685  }
686  }
687  break;
688  }
689  }
690 
691  /* Calculate and add the time needed. */
692  t1 = getticks();
693  t_d = nfft_elapsed_seconds(t1,t0);
694  /* Calculate average time needed. */
695  t_d = t_d/((double)ld[ild][4]);
696  }
697  else
698  {
699  /* Initialize cumulative time variable with dummy value. */
700  t_d = -1.0;
701  t_dp = -1.0;
702  }
703 
704  /* Initialize error and cumulative time variables for the fast
705  * algorithm. */
706  err_fd = -1.0;
707  err_f = -1.0;
708  t_fd = -1.0;
709  t_f = -1.0;
710 
711  /* Process all cut-off bandwidths. */
712  for (im = 0; im < im_max; im++)
713  {
714  /* Init transform plans. */
715  nfsft_init_guru(&plan_adjoint, m[im],ld[ild][0],
716  ((use_nfft!=0)?(0U):(NFSFT_USE_NDFT)) |
717  ((use_fpt!=0)?(0U):(NFSFT_USE_DPT)),
718  PRE_PHI_HUT | PRE_PSI | FFTW_INIT |
719  FFT_OUT_OF_PLACE, cutoff);
720  nfsft_init_guru(&plan,m[im],ld[ild][1],
721  ((use_nfft!=0)?(0U):(NFSFT_USE_NDFT)) |
722  ((use_fpt!=0)?(0U):(NFSFT_USE_DPT)),
723  PRE_PHI_HUT | PRE_PSI | FFTW_INIT |
724  FFT_OUT_OF_PLACE,
725  cutoff);
726  plan_adjoint.f_hat = f_hat;
727  plan_adjoint.x = eta;
728  plan_adjoint.f = b;
729  plan.f_hat = f_hat;
730  plan.x = xi;
731  plan.f = f_m;
732  nfsft_precompute_x(&plan_adjoint);
733  nfsft_precompute_x(&plan);
734 
735  /* Check if direct algorithm shall also be tested. */
736  if (use_nfsft == BOTH)
737  {
738  /* Initialize cumulative time variable. */
739  t_fd = 0.0;
740 
741  /* Initialize time measurement. */
742  t0 = getticks();
743 
744  /* Cycle through all runs. */
745  for (i = 0; i < ld[ild][4]; i++)
746  {
747 
748  /* Execute adjoint direct NDSFT transformation. */
749  nfsft_adjoint_direct(&plan_adjoint);
750 
751  /* Multiplication with the Fourier-Legendre coefficients. */
752  for (k = 0; k <= m[im]; k++)
753  for (n = -k; n <= k; n++)
754  f_hat[NFSFT_INDEX(k,n,&plan_adjoint)] *= a[k];
755 
756  /* Execute direct NDSFT transformation. */
757  nfsft_trafo_direct(&plan);
758 
759  }
760 
761  /* Calculate and add the time needed. */
762  t1 = getticks();
763  t_fd = nfft_elapsed_seconds(t1,t0);
764 
765  /* Calculate average time needed. */
766  t_fd = t_fd/((double)ld[ild][4]);
767 
768  /* Check if error E_infty should be computed. */
769  if (ld[ild][2] != NO)
770  {
771  /* Compute the error E_infinity. */
772  err_fd = X(error_l_infty_1_complex)(f, f_m, ld[ild][1], b,
773  ld[ild][0]);
774  }
775  }
776 
777  /* Check if the fast NFSFT algorithm shall also be tested. */
778  if (use_nfsft != NO)
779  {
780  /* Initialize cumulative time variable for the NFSFT algorithm. */
781  t_f = 0.0;
782  }
783  else
784  {
785  /* Initialize cumulative time variable for the direct NDSFT
786  * algorithm. */
787  t_fd = 0.0;
788  }
789 
790  /* Initialize time measurement. */
791  t0 = getticks();
792 
793  /* Cycle through all runs. */
794  for (i = 0; i < ld[ild][4]; i++)
795  {
796  /* Check if the fast NFSFT algorithm shall also be tested. */
797  if (use_nfsft != NO)
798  {
799  /* Execute the adjoint NFSFT transformation. */
800  nfsft_adjoint(&plan_adjoint);
801  }
802  else
803  {
804  /* Execute the adjoint direct NDSFT transformation. */
805  nfsft_adjoint_direct(&plan_adjoint);
806  }
807 
808  /* Multiplication with the Fourier-Legendre coefficients. */
809  for (k = 0; k <= m[im]; k++)
810  for (n = -k; n <= k; n++)
811  f_hat[NFSFT_INDEX(k,n,&plan_adjoint)] *= a[k];
812 
813  /* Check if the fast NFSFT algorithm shall also be tested. */
814  if (use_nfsft != NO)
815  {
816  /* Execute the NFSFT transformation. */
817  nfsft_trafo(&plan);
818  }
819  else
820  {
821  /* Execute the NDSFT transformation. */
822  nfsft_trafo_direct(&plan);
823  }
824  }
825 
826  /* Check if the fast NFSFT algorithm has been used. */
827  t1 = getticks();
828 
829  if (use_nfsft != NO)
830  t_f = nfft_elapsed_seconds(t1,t0);
831  else
832  t_fd = nfft_elapsed_seconds(t1,t0);
833 
834  /* Check if the fast NFSFT algorithm has been used. */
835  if (use_nfsft != NO)
836  {
837  /* Calculate average time needed. */
838  t_f = t_f/((double)ld[ild][4]);
839  }
840  else
841  {
842  /* Calculate average time needed. */
843  t_fd = t_fd/((double)ld[ild][4]);
844  }
845 
846  /* Check if error E_infty should be computed. */
847  if (ld[ild][2] != NO)
848  {
849  /* Check if the fast NFSFT algorithm has been used. */
850  if (use_nfsft != NO)
851  {
852  /* Compute the error E_infinity. */
853  err_f = X(error_l_infty_1_complex)(f, f_m, ld[ild][1], b,
854  ld[ild][0]);
855  }
856  else
857  {
858  /* Compute the error E_infinity. */
859  err_fd = X(error_l_infty_1_complex)(f, f_m, ld[ild][1], b,
860  ld[ild][0]);
861  }
862  }
863 
864  /* Print out the error measurements. */
865  fprintf(stdout,"%e\n%e\n%e\n%e\n%e\n%e\n\n",t_d,t_dp,t_fd,t_f,err_fd,
866  err_f);
867 
868  /* Finalize the NFSFT plans */
869  nfsft_finalize(&plan_adjoint);
870  nfsft_finalize(&plan);
871  } /* for (im = 0; im < im_max; im++) - Process all cut-off
872  * bandwidths.*/
873  } /* for (ild = 0; ild < ild_max; ild++) - Process all node sets. */
874  } /* for (ip = 0; ip < ip_max; ip++) - Process all parameter sets. */
875 
876  /* Delete precomputed data. */
877  nfsft_forget();
878 
879  /* Check if memory for precomputed data of the matrix K has been
880  * allocated. */
881  if (precompute == YES)
882  {
883  /* Free memory for precomputed matrix K. */
884  nfft_free(prec);
885  }
886  /* Free data arrays. */
887  nfft_free(f);
888  nfft_free(f_m);
889  nfft_free(xi);
890  nfft_free(eta);
891  nfft_free(a);
892  nfft_free(f_hat);
893  nfft_free(b);
894 
895  /* Free memory for node sets. */
896  for (ild = 0; ild < ild_max; ild++)
897  nfft_free(ld[ild]);
898  nfft_free(ld);
899 
900  /* Free memory for cut-off bandwidths. */
901  nfft_free(m);
902 
903  /* Free memory for parameter sets. */
904  for (ip = 0; ip < ip_max; ip++)
905  nfft_free(p[ip]);
906  nfft_free(p);
907  } /* for (tc = 0; tc < tc_max; tc++) - Process each testcase. */
908 
909  /* Return exit code for successful run. */
910  return EXIT_SUCCESS;
911 }
912 /* \} */
int nfft_smbi(const double x, const double alpha, const int nb, const int ize, double *b)
Calculates the modified bessel function , possibly scaled by , for real non-negative with ...
Definition: util.c:950
double * x
the nodes for ,
Definition: nfft3.h:570
static double gaussianKernel(const double x, const double sigma)
Evaluates the spherical Gaussian kernel at a node .
Definition: fastsumS2.c:149
#define KT_SINGULARITY
Singularity kernel.
Definition: fastsumS2.c:55
static double poissonKernel(const double x, const double h)
Evaluates the Poisson kernel at a node .
Definition: fastsumS2.c:97
fftw_complex * f_hat
Vector of Fourier coefficients, size is N_total * sizeof( fftw_complex )
Definition: nfft3.h:570
static double innerProduct(const double phi1, const double theta1, const double phi2, const double theta2)
Computes the standard inner product between two vectors on the unit sphere given in spherical coord...
Definition: fastsumS2.c:78
void nfft_free(void *p)
Header file for utility functions used by the nfft3 library.
static double locallySupportedKernel(const double x, const double h, const double lambda)
Evaluates the locally supported kernel at a node .
Definition: fastsumS2.c:131
static double singularityKernel(const double x, const double h)
Evaluates the singularity kernel at a node .
Definition: fastsumS2.c:113
#define PI
Formerly known to be an irrational number.
Definition: nfft3util.h:62
#define KT_ABEL_POISSON
Abel-Poisson kernel.
Definition: fastsumS2.c:53
void * nfft_malloc(size_t n)
fftw_complex * f
Vector of samples, size is M_total * sizeof( fftw_complex )
Definition: nfft3.h:570
pvalue
Enumeration type for yes/no/both-type parameters.
Definition: fastsumS2.c:62
#define KT_LOC_SUPP
Locally supported kernel.
Definition: fastsumS2.c:57
#define NFFT_MAX(a, b)
Maximum of its two arguments.
Definition: nfft3util.h:68
#define KT_GAUSSIAN
Gaussian kernel.
Definition: fastsumS2.c:59

Generated on Thu May 7 2015 by Doxygen 1.8.5