45 #define SYMBOL_ABEL_POISSON(k,h) (pow(h,k))
48 #define SYMBOL_SINGULARITY(k,h) ((2.0/(2*k+1))*pow(h,k))
53 #define KT_ABEL_POISSON (0)
55 #define KT_SINGULARITY (1)
57 #define KT_LOC_SUPP (2)
59 #define KT_GAUSSIAN (3)
62 enum pvalue {NO = 0, YES = 1, BOTH = 2};
78 static inline double innerProduct(
const double phi1,
const double theta1,
79 const double phi2,
const double theta2)
81 double pi2theta1 = PI2*theta1, pi2theta2 = PI2*theta2;
82 return (cos(pi2theta1)*cos(pi2theta2)
83 + sin(pi2theta1)*sin(pi2theta2)*cos(PI2*(phi1-phi2)));
99 return (1.0/(PI4))*((1.0-h)*(1.0+h))/pow(sqrt(1.0-2.0*h*x+h*h),3.0);
115 return (1.0/(PI2))/sqrt(1.0-2.0*h*x+h*h);
134 return (x<=h)?(0.0):(pow((x-h),lambda));
151 return exp(2.0*sigma*(x-1.0));
164 int main (
int argc,
char **argv)
213 fftw_complex *prec = NULL;
228 fscanf(stdin,
"testcases=%d\n",&tc_max);
229 fprintf(stdout,
"%d\n",tc_max);
232 for (tc = 0; tc < tc_max; tc++)
235 fscanf(stdin,
"nfsft=%d\n",&use_nfsft);
236 fprintf(stdout,
"%d\n",use_nfsft);
240 fscanf(stdin,
"nfft=%d\n",&use_nfft);
241 fprintf(stdout,
"%d\n",use_nfft);
245 fscanf(stdin,
"cutoff=%d\n",&cutoff);
246 fprintf(stdout,
"%d\n",cutoff);
255 fscanf(stdin,
"fpt=%d\n",&use_fpt);
256 fprintf(stdout,
"%d\n",use_fpt);
258 fscanf(stdin,
"threshold=%lf\n",&threshold);
259 fprintf(stdout,
"%lf\n",threshold);
266 threshold = 1000000000000.0;
282 fscanf(stdin,
"kernel=%d\n",&kt);
283 fprintf(stdout,
"%d\n",kt);
286 fscanf(stdin,
"parameter_sets=%d\n",&ip_max);
287 fprintf(stdout,
"%d\n",ip_max);
290 p = (
double**)
nfft_malloc(ip_max*
sizeof(
double*));
295 fscanf(stdin,
"parameters=%d\n",&ipp_max);
296 fprintf(stdout,
"%d\n",ipp_max);
298 for (ip = 0; ip < ip_max; ip++)
301 p[ip] = (
double*)
nfft_malloc(ipp_max*
sizeof(
double));
304 for (ipp = 0; ipp < ipp_max; ipp++)
307 fscanf(stdin,
"%lf\n",&p[ip][ipp]);
308 fprintf(stdout,
"%lf\n",p[ip][ipp]);
313 fscanf(stdin,
"bandwidths=%d\n",&im_max);
314 fprintf(stdout,
"%d\n",im_max);
318 for (im = 0; im < im_max; im++)
321 fscanf(stdin,
"%d\n",&m[im]);
322 fprintf(stdout,
"%d\n",m[im]);
327 fscanf(stdin,
"node_sets=%d\n",&ild_max);
328 fprintf(stdout,
"%d\n",ild_max);
332 for (ild = 0; ild < ild_max; ild++)
338 fscanf(stdin,
"L=%d ",&ld[ild][0]);
339 fprintf(stdout,
"%d\n",ld[ild][0]);
343 fscanf(stdin,
"D=%d ",&ld[ild][1]);
344 fprintf(stdout,
"%d\n",ld[ild][1]);
348 fscanf(stdin,
"compare=%d ",&ld[ild][2]);
349 fprintf(stdout,
"%d\n",ld[ild][2]);
352 if (ld[ild][2] == YES)
355 fscanf(stdin,
"precomputed=%d\n",&ld[ild][3]);
356 fprintf(stdout,
"%d\n",ld[ild][3]);
360 fscanf(stdin,
"repetitions=%d\n",&ld[ild][4]);
361 fprintf(stdout,
"%d\n",ld[ild][4]);
364 if (ld[ild][3] == YES)
367 ld_max_prec =
NFFT_MAX(ld_max_prec,ld[ild][0]*ld[ild][1]);
369 l_max_prec =
NFFT_MAX(l_max_prec,ld[ild][0]);
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));
391 if (precompute == YES)
393 prec = (fftw_complex*)
nfft_malloc(ld_max_prec*
sizeof(fftw_complex));
397 for (l = 0; l < l_max; l++)
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);
405 for (d = 0; d < d_max; d++)
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);
412 nfsft_precompute(m_max,threshold,
413 ((use_nfsft==NO)?(NFSFT_NO_FAST_ALGORITHM):(0U)), 0U);
416 for (ip = 0; ip < ip_max; ip++)
423 for (k = 0; k <= m_max; k++)
424 a[k] = SYMBOL_ABEL_POISSON(k,p[ip][0]);
430 for (k = 0; k <= m_max; k++)
431 a[k] = SYMBOL_SINGULARITY(k,p[ip][0]);
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]);
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];
457 for (k = 0; k <= m_max; k++)
458 a[k] *= (2*k+1)/(PI4);
461 for (ild = 0; ild < ild_max; ild++)
464 if (ld[ild][2] != NO)
468 if (ld[ild][3] != NO)
473 rinc = l_max_prec-ld[ild][0];
476 for (d = 0; d < ld[ild][1]; d++)
479 for (l = 0; l < ld[ild][0]; l++)
483 temp =
innerProduct(eta[2*l],eta[2*l+1],xi[2*d],xi[2*d+1]);
523 for (i = 0; i < ld[ild][4]; i++)
529 rinc = l_max_prec-ld[ild][0];
537 constant = ((p[ip][1]+1)/(PI2*pow(1-p[ip][0],p[ip][1]+1)));
540 for (d = 0; d < ld[ild][1]; d++)
546 for (l = 0; l < ld[ild][0]; l++)
547 f[d] += b[l]*(*ptr++);
559 for (d = 0; d < ld[ild][1]; d++)
565 for (l = 0; l < ld[ild][0]; l++)
566 f[d] += b[l]*(*ptr++);
576 t_dp = nfft_elapsed_seconds(t1,t0);
579 t_dp = t_dp/((double)ld[ild][4]);
594 for (i = 0; i < ld[ild][4]; i++)
602 for (d = 0; d < ld[ild][1]; d++)
608 for (l = 0; l < ld[ild][0]; l++)
612 temp =
innerProduct(eta[2*l],eta[2*l+1],xi[2*d],xi[2*d+1]);
623 for (d = 0; d < ld[ild][1]; d++)
629 for (l = 0; l < ld[ild][0]; l++)
633 temp =
innerProduct(eta[2*l],eta[2*l+1],xi[2*d],xi[2*d+1]);
644 constant = ((p[ip][1]+1)/(PI2*pow(1-p[ip][0],p[ip][1]+1)));
647 for (d = 0; d < ld[ild][1]; d++)
653 for (l = 0; l < ld[ild][0]; l++)
657 temp =
innerProduct(eta[2*l],eta[2*l+1],xi[2*d],xi[2*d+1]);
671 for (d = 0; d < ld[ild][1]; d++)
677 for (l = 0; l < ld[ild][0]; l++)
681 temp =
innerProduct(eta[2*l],eta[2*l+1],xi[2*d],xi[2*d+1]);
693 t_d = nfft_elapsed_seconds(t1,t0);
695 t_d = t_d/((double)ld[ild][4]);
712 for (im = 0; im < im_max; im++)
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 |
726 plan_adjoint.
f_hat = f_hat;
727 plan_adjoint.
x = eta;
732 nfsft_precompute_x(&plan_adjoint);
733 nfsft_precompute_x(&plan);
736 if (use_nfsft == BOTH)
745 for (i = 0; i < ld[ild][4]; i++)
749 nfsft_adjoint_direct(&plan_adjoint);
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];
757 nfsft_trafo_direct(&plan);
763 t_fd = nfft_elapsed_seconds(t1,t0);
766 t_fd = t_fd/((double)ld[ild][4]);
769 if (ld[ild][2] != NO)
772 err_fd = X(error_l_infty_1_complex)(f, f_m, ld[ild][1], b,
794 for (i = 0; i < ld[ild][4]; i++)
800 nfsft_adjoint(&plan_adjoint);
805 nfsft_adjoint_direct(&plan_adjoint);
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];
822 nfsft_trafo_direct(&plan);
830 t_f = nfft_elapsed_seconds(t1,t0);
832 t_fd = nfft_elapsed_seconds(t1,t0);
838 t_f = t_f/((double)ld[ild][4]);
843 t_fd = t_fd/((double)ld[ild][4]);
847 if (ld[ild][2] != NO)
853 err_f = X(error_l_infty_1_complex)(f, f_m, ld[ild][1], b,
859 err_fd = X(error_l_infty_1_complex)(f, f_m, ld[ild][1], b,
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,
869 nfsft_finalize(&plan_adjoint);
870 nfsft_finalize(&plan);
881 if (precompute == YES)
896 for (ild = 0; ild < ild_max; ild++)
904 for (ip = 0; ip < ip_max; ip++)
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 ...
double * x
the nodes for ,
static double gaussianKernel(const double x, const double sigma)
Evaluates the spherical Gaussian kernel at a node .
#define KT_SINGULARITY
Singularity kernel.
static double poissonKernel(const double x, const double h)
Evaluates the Poisson kernel at a node .
fftw_complex * f_hat
Vector of Fourier coefficients, size is N_total * sizeof( fftw_complex )
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...
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 .
static double singularityKernel(const double x, const double h)
Evaluates the singularity kernel at a node .
#define PI
Formerly known to be an irrational number.
#define KT_ABEL_POISSON
Abel-Poisson kernel.
void * nfft_malloc(size_t n)
fftw_complex * f
Vector of samples, size is M_total * sizeof( fftw_complex )
pvalue
Enumeration type for yes/no/both-type parameters.
#define KT_LOC_SUPP
Locally supported kernel.
#define NFFT_MAX(a, b)
Maximum of its two arguments.
#define KT_GAUSSIAN
Gaussian kernel.