NFFT  3.4.1
Macros | Typedefs | Enumerations | Functions
Util - Auxiliary functions

This module implements frequently used utility functions. More...

Macros

#define CONCAT(prefix, name)   prefix ## name
 
#define Y(name)   CONCAT(nfft_,name)
 
#define FFTW(name)   CONCAT(fftw_,name)
 
#define NFFT(name)   CONCAT(nfft_,name)
 
#define NFCT(name)   CONCAT(nfct_,name)
 
#define NFST(name)   CONCAT(nfst_,name)
 
#define NFSFT(name)   CONCAT(nfsft_,name)
 
#define SOLVER(name)   CONCAT(solver_,name)
 
#define X(name)   Y(name)
 
#define STRINGIZEx(x)   #x
 
#define STRINGIZE(x)   STRINGIZEx(x)
 
#define K(x)   ((R) x)
 
#define DK(name, value)   const R name = K(value)
 
#define KPI   K(3.1415926535897932384626433832795028841971693993751)
 
#define K2PI   K(6.2831853071795864769252867665590057683943387987502)
 
#define K4PI   K(12.5663706143591729538505735331180115367886775975004)
 
#define KE   K(2.7182818284590452353602874713526624977572470937000)
 
#define IF(x, a, b)   ((x)?(a):(b))
 
#define MIN(a, b)   (((a)<(b))?(a):(b))
 
#define MAX(a, b)   (((a)>(b))?(a):(b))
 
#define ABS(x)   (((x)>K(0.0))?(x):(-(x)))
 
#define SIGN(a)   (((a)>=0)?1:-1)
 
#define SIGN(a)   (((a)>=0)?1:-1)
 
#define SIGNF(a)   IF((a)<K(0.0),K(-1.0),K(1.0))
 
#define SIZE(x)   sizeof(x)/sizeof(x[0])
 
#define CSWAP(x, y)
 Swap two vectors. More...
 
#define RSWAP(x, y)
 Swap two vectors. More...
 
#define PHI_HUT(n, k, d)   (Y(bessel_i0)((R)(ths->m) * SQRT(ths->b[d] * ths->b[d] - (K(2.0) * KPI * (R)(k) / (R)(n)) * (K(2.0) * KPI * (R)(k) / (R)(n)))))
 
#define PHI(n, x, d)
 
#define WINDOW_HELP_INIT
 
#define WINDOW_HELP_FINALIZE   {Y(free)(ths->b);}
 
#define WINDOW_HELP_ESTIMATE_m   8
 
#define COPYSIGN   copysign
 
#define NEXTAFTER   nextafter
 
#define MKNAN   nan
 
#define CEIL   ceil
 
#define FLOOR   floor
 
#define NEARBYINT   nearbyint
 
#define RINT   rint
 
#define ROUND   round
 
#define LRINT   lrint
 
#define LROUND   lround
 
#define LLRINT   llrint
 
#define LLROUND   llround
 
#define TRUNC   trunc
 
#define FMOD   fmod
 
#define REMAINDER   remainder
 
#define REMQUO   remquo
 
#define FDIM   fdim
 
#define FMAX   fmax
 
#define FMIN   fmin
 
#define FFMA   fma
 
#define FABS   fabs
 
#define SQRT   sqrt
 
#define CBRT   cbrt
 
#define HYPOT   hypot
 
#define EXP   exp
 
#define EXP2   exp2
 
#define EXPM1   expm1
 
#define LOG   log
 
#define LOG2   log2
 
#define LOG10   log10
 
#define LOG1P   log1p
 
#define LOGB   logb
 
#define ILOGB   ilogb
 
#define MODF   modf
 
#define FREXP   frexp
 
#define LDEXP   ldexp
 
#define SCALBN   scalbn
 
#define SCALBLN   scalbln
 
#define POW   pow
 
#define COS   cos
 
#define SIN   sin
 
#define TAN   tan
 
#define COSH   cosh
 
#define SINH   sinh
 
#define TANH   tanh
 
#define ACOS   acos
 
#define ASIN   asin
 
#define ATAN   atan
 
#define ATAN2   atan2
 
#define ACOSH   acosh
 
#define ASINH   asinh
 
#define ATANH   atanh
 
#define TGAMMA   tgamma
 
#define LGAMMA   lgamma
 
#define J0   j0
 
#define J1   j1
 
#define JN   jn
 
#define Y0   y0
 
#define Y1   y1
 
#define YN   yn
 
#define ERF   erf
 
#define ERFC   erfc
 
#define CREAL   creal
 
#define CIMAG   cimag
 
#define CABS   cabs
 
#define CARG   carg
 
#define CONJ   conj
 
#define CPROJ   cproj
 
#define CSQRT   csqrt
 
#define CEXP   cexp
 
#define CLOG   clog
 
#define CPOW   cpow
 
#define CSIN   csin
 
#define CCOS   ccos
 
#define CTAN   ctan
 
#define CASIN   casin
 
#define CACOS   cacos
 
#define CATAN   catan
 
#define CSINH   csinh
 
#define CCOSH   ccosh
 
#define CTANH   ctanh
 
#define CASINH   casinh
 
#define CACOSH   cacosh
 
#define CATANH   catanh
 
#define MANT_DIG   DBL_MANT_DIG
 
#define MIN_EXP   DBL_MIN_EXP
 
#define MAX_EXP   DBL_MAX_EXP
 
#define FLTROUND   0.0
 
#define R_RADIX   FLT_RADIX
 
#define II   _Complex_I
 
#define __FGS__   "lg"
 
#define __FES__   "lE"
 
#define __FE__   "% 20.16lE"
 
#define __FI__   "%lf"
 
#define __FIS__   "lf"
 
#define __FR__   "%le"
 
#define TRUE   1
 
#define FALSE   0
 
#define __D__   "%td"
 
#define UNUSED(x)   (void)x
 Dummy use of unused parameters to silence compiler warnings.
 
#define UNUSED(x)   (void)x
 Dummy use of unused parameters to silence compiler warnings.
 
#define STACK_MALLOC(T, p, x)   p = (T)alloca(x)
 
#define STACK_FREE(x)   /* Nothing. Cleanup done automatically. */
 
#define TIC(a)
 Timing, method works since the inaccurate timer is updated mostly in the measured function. More...
 
#define TOC(a)
 
#define TIC_FFTW(a)
 
#define TOC_FFTW(a)
 
#define CK(ex)   (void)((ex) || (Y(assertion_failed)(#ex, __LINE__, __FILE__), 0))
 
#define A(ex)   /* nothing */
 

Typedefs

typedef double R
 
typedef double _Complex C
 
typedef ptrdiff_t INT
 

Enumerations

enum  float_property {
  NFFT_SAFE__MIN = 1, NFFT_BASE = 2, NFFT_PRECISION = 3, NFFT_MANT_DIG = 4,
  NFFT_FLTROUND = 5, NFFT_E_MIN = 6, NFFT_R_MIN = 7, NFFT_E_MAX = 8,
  NFFT_R_MAX = 9
}
 

Functions

INT nfft_m2K (const INT m)
 
nfft_elapsed_seconds (ticks t1, ticks t0)
 Return number of elapsed seconds between two time points. More...
 
nfft_sinc (R x)
 
nfft_lambda (R z, R eps)
 
nfft_lambda2 (R mu, R nu)
 
nfft_bessel_i0 (R x)
 
nfft_bsplines (const INT, const R x)
 
nfft_float_property (float_property)
 
nfft_prod_real (R *vec, INT d)
 
INT nfft_log2i (const INT m)
 
void nfft_next_power_of_2_exp (const INT N, INT *N2, INT *t)
 
void nfft_next_power_of_2_exp_int (const int N, int *N2, int *t)
 
nfft_error_l_infty_double (const R *x, const R *y, const INT n)
 
nfft_error_l_infty_1_double (const R *x, const R *y, const INT n, const R *z, const INT m)
 
nfft_error_l_2_complex (const C *x, const C *y, const INT n)
 
nfft_error_l_2_double (const R *x, const R *y, const INT n)
 
void nfft_sort_node_indices_radix_msdf (INT n, INT *keys0, INT *keys1, INT rhigh)
 
void nfft_sort_node_indices_radix_lsdf (INT n, INT *keys0, INT *keys1, INT rhigh)
 
void nfft_assertion_failed (const char *s, int line, const char *file)
 
nfft_dot_double (R *x, INT n)
 Computes the inner/dot product $x^H x$. More...
 
nfft_dot_w_complex (C *x, R *w, INT n)
 Computes the weighted inner/dot product $x^H (w \odot x)$. More...
 
nfft_dot_w_double (R *x, R *w, INT n)
 Computes the weighted inner/dot product $x^H (w \odot x)$. More...
 
nfft_dot_w_w2_complex (C *x, R *w, R *w2, INT n)
 Computes the weighted inner/dot product $x^H (w\odot w2\odot w2 \odot x)$. More...
 
nfft_dot_w2_complex (C *x, R *w2, INT n)
 Computes the weighted inner/dot product $x^H (w2\odot w2 \odot x)$. More...
 
void nfft_cp_complex (C *x, C *y, INT n)
 Copies $x \leftarrow y$. More...
 
void nfft_cp_double (R *x, R *y, INT n)
 Copies $x \leftarrow y$. More...
 
void nfft_cp_a_complex (C *x, R a, C *y, INT n)
 Copies $x \leftarrow a y$. More...
 
void nfft_cp_a_double (R *x, R a, R *y, INT n)
 Copies $x \leftarrow a y$. More...
 
void nfft_cp_w_complex (C *x, R *w, C *y, INT n)
 Copies $x \leftarrow w\odot y$. More...
 
void nfft_cp_w_double (R *x, R *w, R *y, INT n)
 Copies $x \leftarrow w\odot y$. More...
 
void nfft_upd_axpy_double (R *x, R a, R *y, INT n)
 Updates $x \leftarrow a x + y$. More...
 
void nfft_upd_xpay_complex (C *x, R a, C *y, INT n)
 Updates $x \leftarrow x + a y$. More...
 
void nfft_upd_xpay_double (R *x, R a, R *y, INT n)
 Updates $x \leftarrow x + a y$. More...
 
void nfft_upd_axpby_complex (C *x, R a, C *y, R b, INT n)
 Updates $x \leftarrow a x + b y$. More...
 
void nfft_upd_axpby_double (R *x, R a, R *y, R b, INT n)
 Updates $x \leftarrow a x + b y$. More...
 
void nfft_upd_xpawy_complex (C *x, R a, R *w, C *y, INT n)
 Updates $x \leftarrow x + a w\odot y$. More...
 
void nfft_upd_xpawy_double (R *x, R a, R *w, R *y, INT n)
 Updates $x \leftarrow x + a w\odot y$. More...
 
void nfft_upd_axpwy_complex (C *x, R a, R *w, C *y, INT n)
 Updates $x \leftarrow a x + w\odot y$. More...
 
void nfft_upd_axpwy_double (R *x, R a, R *w, R *y, INT n)
 Updates $x \leftarrow a x + w\odot y$. More...
 
void nfft_voronoi_weights_1d (R *w, R *x, const INT M)
 
nfft_modified_fejer (const INT N, const INT kk)
 Compute damping factor for modified Fejer kernel: /f${2}{N}(1-{|2k+1|}{N})/f$.
 
nfft_modified_jackson2 (const INT N, const INT kk)
 Compute damping factor for modified Jackson kernel. More...
 
nfft_modified_jackson4 (const INT N, const INT kk)
 Compute damping factor for modified generalised Jackson kernel. More...
 
nfft_modified_sobolev (const R mu, const INT kk)
 Compute damping factor for modified Sobolev kernel. More...
 
nfft_modified_multiquadric (const R mu, const R c, const INT kk)
 Comput damping factor for modified multiquadric kernel. More...
 

Detailed Description

This module implements frequently used utility functions.

In particular, this includes simple measurement of resources, evaluation of window functions, vector routines for basic linear algebra tasks, and computation of weights for the inverse transforms.

Macro Definition Documentation

#define CSWAP (   x,
 
)
Value:
{C* NFFT_SWAP_temp__; \
NFFT_SWAP_temp__=(x); (x)=(y); (y)=NFFT_SWAP_temp__;}

Swap two vectors.

Definition at line 139 of file infft.h.

Referenced by main(), solver_loop_one_step_cgnr_complex(), solver_loop_one_step_landweber_complex(), solver_loop_one_step_steepest_descent_complex(), and taylor_time_accuracy().

#define RSWAP (   x,
 
)
Value:
{R* NFFT_SWAP_temp__; NFFT_SWAP_temp__=(x); \
(x)=(y); (y)=NFFT_SWAP_temp__;}

Swap two vectors.

Definition at line 143 of file infft.h.

Referenced by solver_loop_one_step_cgnr_double(), solver_loop_one_step_landweber_double(), and solver_loop_one_step_steepest_descent_double().

#define PHI (   n,
  x,
 
)
Value:
( (((R)(ths->m) * (R)(ths->m) - (x) * (R)(n) * (x) * (R)(n)) > K(0.0)) \
? SINH(ths->b[d] * SQRT((R)(ths->m) * (R)(ths->m) - (x) * (R)(n) * (x) * (R)(n))) \
/ (KPI * SQRT((R)(ths->m) * (R)(ths->m) - (x) * (R)(n) * (x) * (R)(n))) \
: ((((R)(ths->m) * (R)(ths->m) - (x) * (R)(n) * (x) * (R)(n)) < K(0.0)) \
? SIN(ths->b[d] * SQRT((x) * (R)(n) * (x) * (R)(n) - (R)(ths->m) * (R)(ths->m))) \
/ (KPI * SQRT((x) * (R)(n) * (x) * (R)(n) - (R)(ths->m) * (R)(ths->m))) \
: ths->b[d] / KPI))

Definition at line 209 of file infft.h.

#define WINDOW_HELP_INIT
Value:
{ \
int WINDOW_idx; \
ths->b = (R*) Y(malloc)((size_t)(ths->d) * sizeof(R)); \
for (WINDOW_idx = 0; WINDOW_idx < ths->d; WINDOW_idx++) \
ths->b[WINDOW_idx] = (KPI * (K(2.0) - K(1.0) / ths->sigma[WINDOW_idx])); \
}

Definition at line 216 of file infft.h.

#define TIC (   a)

Timing, method works since the inaccurate timer is updated mostly in the measured function.

For small times not every call of the measured function will also produce a 'unit' time step. Measuring the fftw might cause a wrong output vector due to the repeated ffts.

Definition at line 1397 of file infft.h.

Function Documentation

R nfft_elapsed_seconds ( ticks  t1,
ticks  t0 
)

Return number of elapsed seconds between two time points.

Referenced by fastsum_precompute_source_nodes(), fastsum_trafo(), main(), nfsft_adjoint(), and nfsft_trafo().

R nfft_dot_double ( R *  x,
INT  n 
)

Computes the inner/dot product $x^H x$.

R nfft_dot_w_complex ( C *  x,
R *  w,
INT  n 
)

Computes the weighted inner/dot product $x^H (w \odot x)$.

R nfft_dot_w_double ( R *  x,
R *  w,
INT  n 
)

Computes the weighted inner/dot product $x^H (w \odot x)$.

R nfft_dot_w_w2_complex ( C *  x,
R *  w,
R *  w2,
INT  n 
)

Computes the weighted inner/dot product $x^H (w\odot w2\odot w2 \odot x)$.

R nfft_dot_w2_complex ( C *  x,
R *  w2,
INT  n 
)

Computes the weighted inner/dot product $x^H (w2\odot w2 \odot x)$.

void nfft_cp_complex ( C *  x,
C *  y,
INT  n 
)

Copies $x \leftarrow y$.

void nfft_cp_double ( R *  x,
R *  y,
INT  n 
)

Copies $x \leftarrow y$.

void nfft_cp_a_complex ( C *  x,
a,
C *  y,
INT  n 
)

Copies $x \leftarrow a y$.

void nfft_cp_a_double ( R *  x,
a,
R *  y,
INT  n 
)

Copies $x \leftarrow a y$.

void nfft_cp_w_complex ( C *  x,
R *  w,
C *  y,
INT  n 
)

Copies $x \leftarrow w\odot y$.

void nfft_cp_w_double ( R *  x,
R *  w,
R *  y,
INT  n 
)

Copies $x \leftarrow w\odot y$.

void nfft_upd_axpy_double ( R *  x,
a,
R *  y,
INT  n 
)

Updates $x \leftarrow a x + y$.

void nfft_upd_xpay_complex ( C *  x,
a,
C *  y,
INT  n 
)

Updates $x \leftarrow x + a y$.

void nfft_upd_xpay_double ( R *  x,
a,
R *  y,
INT  n 
)

Updates $x \leftarrow x + a y$.

void nfft_upd_axpby_complex ( C *  x,
a,
C *  y,
b,
INT  n 
)

Updates $x \leftarrow a x + b y$.

void nfft_upd_axpby_double ( R *  x,
a,
R *  y,
b,
INT  n 
)

Updates $x \leftarrow a x + b y$.

void nfft_upd_xpawy_complex ( C *  x,
a,
R *  w,
C *  y,
INT  n 
)

Updates $x \leftarrow x + a w\odot y$.

void nfft_upd_xpawy_double ( R *  x,
a,
R *  w,
R *  y,
INT  n 
)

Updates $x \leftarrow x + a w\odot y$.

void nfft_upd_axpwy_complex ( C *  x,
a,
R *  w,
C *  y,
INT  n 
)

Updates $x \leftarrow a x + w\odot y$.

void nfft_upd_axpwy_double ( R *  x,
a,
R *  w,
R *  y,
INT  n 
)

Updates $x \leftarrow a x + w\odot y$.

R nfft_modified_jackson2 ( const INT  N,
const INT  kk 
)

Compute damping factor for modified Jackson kernel.

R nfft_modified_jackson4 ( const INT  N,
const INT  kk 
)

Compute damping factor for modified generalised Jackson kernel.

R nfft_modified_sobolev ( const R  mu,
const INT  kk 
)

Compute damping factor for modified Sobolev kernel.

R nfft_modified_multiquadric ( const R  mu,
const R  c,
const INT  kk 
)

Comput damping factor for modified multiquadric kernel.