NFFT Logo 3.1.0 API Reference

Fast summation
[Applications]


Modules

 fastsum_matlab
 fastsum_test

Data Structures

struct  fastsum_plan
 plan for fast summation algorithm More...

Defines

#define EXACT_NEARFIELD   (1U<< 0)
 Include header for C99 complex datatype.

Functions

double fak (int n)
 factorial
double binom (int n, int m)
 binomial coefficient
double BasisPoly (int m, int r, double xx)
 basis polynomial for regularized kernel
double regkern (double _Complex(*kernel)(double, int, const double *), double xx, int p, const double *param, double a, double b)
 regularized kernel with K_I arbitrary and K_B smooth to zero
double regkern1 (double _Complex(*kernel)(double, int, const double *), double xx, int p, const double *param, double a, double b)
 regularized kernel with K_I arbitrary and K_B periodized (used in 1D)
double regkern2 (double _Complex(*kernel)(double, int, const double *), double xx, int p, const double *param, double a, double b)
 regularized kernel for even kernels with K_I even and K_B mirrored
double regkern3 (double _Complex(*kernel)(double, int, const double *), double xx, int p, const double *param, double a, double b)
 regularized kernel for even kernels with K_I even and K_B mirrored smooth to K(1/2) (used in dD, d>1)
double kubintkern (double x, double *Add, int Ad, double a)
 cubic spline interpolation in near field with even kernels
double kubintkern1 (double x, double *Add, int Ad, double a)
 cubic spline interpolation in near field with arbitrary kernels
void quicksort (int d, int t, double *x, double _Complex *alpha, int N)
 quicksort algorithm for source knots and associated coefficients
void BuildTree (int d, int t, double *x, double _Complex *alpha, int N)
 recursive sort of source knots dimension by dimension to get tree structure
double _Complex SearchTree (int d, int t, double *x, double _Complex *alpha, double *xmin, double *xmax, int N, double _Complex(*kernel)(double, int, const double *), const double *param, int Ad, double *Add, int p, unsigned flags)
 fast search in tree of source knots for near field computation
void fastsum_init_guru (fastsum_plan *ths, int d, int N_total, int M_total, double _Complex(*kernel)(), double *param, unsigned flags, int nn, int m, int p, double eps_I, double eps_B)
 initialization of fastsum plan
void fastsum_finalize (fastsum_plan *ths)
 finalization of fastsum plan
void fastsum_exact (fastsum_plan *ths)
 direct computation of sums
void fastsum_precompute (fastsum_plan *ths)
 precomputation for fastsum
void fastsum_trafo (fastsum_plan *ths)
 fast NFFT-based summation
double _Complex gaussian (double x, int der, const double *param)
double _Complex multiquadric (double x, int der, const double *param)
double _Complex inverse_multiquadric (double x, int der, const double *param)
double _Complex logarithm (double x, int der, const double *param)
double _Complex thinplate_spline (double x, int der, const double *param)
double _Complex one_over_square (double x, int der, const double *param)
double _Complex one_over_modulus (double x, int der, const double *param)
double _Complex one_over_x (double x, int der, const double *param)
double _Complex inverse_multiquadric3 (double x, int der, const double *param)
double _Complex sinc_kernel (double x, int der, const double *param)
double _Complex cosc (double x, int der, const double *param)
double _Complex kcot (double x, int der, const double *param)
double _Complex one_over_cube (double x, int der, const double *param)

Define Documentation

#define EXACT_NEARFIELD   (1U<< 0)

Include header for C99 complex datatype.

Include header for utils from NFFT3 library. Include header for NFFT3 library. Constant symbols

Definition at line 51 of file fastsum.h.

Referenced by fastsum_finalize(), fastsum_init_guru(), fastsum_precompute(), FieldsSearchTree(), and SearchTree().


Function Documentation

void fastsum_init_guru ( fastsum_plan ths,
int  d,
int  N_total,
int  M_total,
double _Complex(*)()  kernel,
double *  param,
unsigned  flags,
int  nn,
int  m,
int  p,
double  eps_I,
double  eps_B 
)

initialization of fastsum plan

initialize fast summation plan

Parameters:
ths The pointer to a fastsum plan.
d The dimension of the problem.
N_total The number of source knots x.
M_total The number of target knots y.
kernel The kernel function.
param The parameters for the kernel function.
flags Fastsum flags.
nn The expansion degree.
m The cut-off parameter for the NFFT.
p The degree of smoothness.
eps_I The inner boundary.
eps_B the outer boundary.

Definition at line 390 of file fastsum.c.

References fastsum_plan::Ad, fastsum_plan::Add, fastsum_plan::alpha, fastsum_plan::b, fastsum_plan::d, fastsum_plan::eps_B, fastsum_plan::eps_I, EXACT_NEARFIELD, fastsum_plan::f, FFT_OUT_OF_PLACE, fastsum_plan::fft_plan, FFTW_INIT, fastsum_plan::flags, fastsum_plan::kernel, fastsum_plan::kernel_param, fastsum_plan::M_total, MALLOC_F, MALLOC_F_HAT, MALLOC_X, fastsum_plan::mv1, fastsum_plan::mv2, fastsum_plan::n, fastsum_plan::N_total, nfft_init_guru(), nfft_malloc(), fastsum_plan::p, PRE_PHI_HUT, PRE_PSI, fastsum_plan::x, and fastsum_plan::y.

Referenced by fields_fastsum_init(), and main().

void fastsum_finalize ( fastsum_plan ths  ) 

finalization of fastsum plan

finalize plan

Parameters:
ths The pointer to a fastsum plan.

Definition at line 456 of file fastsum.c.

References fastsum_plan::Add, fastsum_plan::alpha, fastsum_plan::b, EXACT_NEARFIELD, fastsum_plan::f, fastsum_plan::fft_plan, fastsum_plan::flags, fastsum_plan::mv1, fastsum_plan::mv2, nfft_finalize(), nfft_free(), fastsum_plan::x, and fastsum_plan::y.

Referenced by fields_fastsum_finalize(), and main().

void fastsum_exact ( fastsum_plan ths  ) 

direct computation of sums

direct summation

Parameters:
ths The pointer to a fastsum plan.

Definition at line 474 of file fastsum.c.

References fastsum_plan::alpha, fastsum_plan::d, fastsum_plan::f, fastsum_plan::kernel, fastsum_plan::kernel_param, fastsum_plan::M_total, fastsum_plan::N_total, fastsum_plan::x, and fastsum_plan::y.

Referenced by main().

void fastsum_precompute ( fastsum_plan ths  ) 

precomputation for fastsum

sort source nodes, precompute Fourier coefficients, etc.

Parameters:
ths The pointer to a fastsum plan.

Definition at line 500 of file fastsum.c.

References fastsum_plan::Ad, fastsum_plan::Add, fastsum_plan::alpha, fastsum_plan::b, BuildTree(), nfft_plan::d, fastsum_plan::d, fastsum_plan::eps_B, fastsum_plan::eps_I, EXACT_NEARFIELD, nfft_plan::f, fastsum_plan::fft_plan, fastsum_plan::flags, fastsum_plan::kernel, fastsum_plan::kernel_param, nfft_plan::M_total, fastsum_plan::mv1, fastsum_plan::mv2, nfft_plan::N, fastsum_plan::n, fastsum_plan::N_total, nfft_fftshift_complex(), nfft_plan::nfft_flags, nfft_precompute_full_psi(), nfft_precompute_lin_psi(), nfft_precompute_psi(), fastsum_plan::p, PRE_FULL_PSI, PRE_LIN_PSI, PRE_PSI, regkern1(), regkern3(), nfft_plan::x, fastsum_plan::x, and fastsum_plan::y.

Referenced by fields_fastsum_precompute(), fields_fastsum_simple(), and main().

void fastsum_trafo ( fastsum_plan ths  ) 

fast NFFT-based summation

fast NFFT-based summation algorithm

Parameters:
ths The pointer to a fastsum plan.

Definition at line 583 of file fastsum.c.

References fastsum_plan::Ad, fastsum_plan::Add, fastsum_plan::alpha, fastsum_plan::b, fastsum_plan::d, fastsum_plan::eps_I, nfft_plan::f, fastsum_plan::f, nfft_plan::f_hat, fastsum_plan::flags, fastsum_plan::kernel, fastsum_plan::kernel_param, fastsum_plan::M_total, fastsum_plan::mv1, fastsum_plan::mv2, fastsum_plan::N_total, nfft_plan::N_total, nfft_adjoint(), nfft_malloc(), nfft_trafo(), fastsum_plan::p, SearchTree(), fastsum_plan::x, and fastsum_plan::y.

Referenced by fields_fastsum_simple(), and main().


Generated on 19 Mar 2009 by Doxygen 1.5.3