NFFT Logo 3.0 API Reference
Main Page | Modules | Data Structures | Directories | File List | Data Fields | Globals

fpt.c File Reference

Implementation file for the FPT module. More...

#include <math.h>
#include <string.h>
#include <stdlib.h>
#include <stdbool.h>
#include "nfft3.h"
#include "util.h"

Go to the source code of this file.

Data Structures

struct  fpt_step_
 Holds data for a single multiplication step in the cascade summation. More...
struct  fpt_data_
 Holds data for a single cascade summation. More...
struct  fpt_set_s_
 Holds data for a set of cascade summations. More...

Defines

#define K_START_TILDE(x, y)   (NFFT_MAX(NFFT_MIN(x,y-2),0))
 Computes the minimum degree at the top of a cascade.
#define K_END_TILDE(x, y)   NFFT_MIN(x,y-1)
 Computes the maximum degree at the top of a cascade.
#define FIRST_L(x, y)   ((int)floor((x)/(double)y))
 Computes the index of the first block of four functions in a cascade level.
#define LAST_L(x, y)   ((int)ceil(((x)+1)/(double)y)-1)
 Computes the index of the last block of four functions in a cascade level.
#define N_TILDE(y)   (y-1)
#define IS_SYMMETRIC(x, y, z)   (x >= ((y-1.0)/z))
#define ABUVXPWY_SYMMETRIC(NAME, S1, S2)
#define ABUVXPWY_SYMMETRIC_1(NAME, S1)
#define ABUVXPWY_SYMMETRIC_2(NAME, S1)
#define AUVXPWY_SYMMETRIC(NAME, S1, S2)
#define FPT_DO_STEP(NAME, M1_FUNCTION, M2_FUNCTION)
#define FPT_DO_STEP_TRANSPOSED(NAME, M1_FUNCTION, M2_FUNCTION)

Typedefs

typedef fpt_step_ fpt_step
 Holds data for a single multiplication step in the cascade summation.
typedef fpt_data_ fpt_data
 Holds data for a single cascade summation.
typedef fpt_set_s_ fpt_set_s
 Holds data for a set of cascade summations.

Functions

void abuvxpwy (double a, double b, double complex *u, double complex *x, double *v, double complex *y, double *w, int n)
void abuvxpwy_symmetric1 (double a, double b, double complex *u, double complex *x, double *v, double complex *y, double *w, int n)
void abuvxpwy_symmetric2 (double a, double b, double complex *u, double complex *x, double *v, double complex *y, double *w, int n)
void abuvxpwy_symmetric1_1 (double a, double b, double complex *u, double complex *x, double *v, double complex *y, double *w, int n)
void abuvxpwy_symmetric1_2 (double a, double b, double complex *u, double complex *x, double *v, double complex *y, double *w, int n)
void abuvxpwy_symmetric2_1 (double a, double b, double complex *u, double complex *x, double *v, double complex *y, double *w, int n)
void abuvxpwy_symmetric2_2 (double a, double b, double complex *u, double complex *x, double *v, double complex *y, double *w, int n)
void auvxpwy (double a, double complex *u, double complex *x, double *v, double complex *y, double *w, int n)
void auvxpwy_symmetric (double a, double complex *u, double complex *x, double *v, double complex *y, double *w, int n)
void fpt_do_step (double complex *a, double complex *b, double *a11, double *a12, double *a21, double *a22, double gamma, int tau, fpt_set set)
void fpt_do_step_symmetric (double complex *a, double complex *b, double *a11, double *a12, double *a21, double *a22, double gamma, int tau, fpt_set set)
void fpt_do_step_symmetric_u (double complex *a, double complex *b, double *a11, double *a12, double *a21, double *a22, double gamma, int tau, fpt_set set)
void fpt_do_step_symmetric_l (double complex *a, double complex *b, double *a11, double *a12, double *a21, double *a22, double gamma, int tau, fpt_set set)
void fpt_do_step_transposed (double complex *a, double complex *b, double *a11, double *a12, double *a21, double *a22, double gamma, int tau, fpt_set set)
void fpt_do_step_transposed_symmetric (double complex *a, double complex *b, double *a11, double *a12, double *a21, double *a22, double gamma, int tau, fpt_set set)
void fpt_do_step_transposed_symmetric_u (double complex *a, double complex *b, double *a11, double *a12, double *a21, double *a22, double gamma, int tau, fpt_set set)
void fpt_do_step_transposed_symmetric_l (double complex *a, double complex *b, double *a11, double *a12, double *a21, double *a22, double gamma, int tau, fpt_set set)
void eval_clenshaw (const double *x, double *y, int size, int k, const double *alpha, const double *beta, const double *gamma)
int eval_clenshaw_thresh (const double *x, double *y, int size, int k, const double *alpha, const double *beta, const double *gamma, const double threshold)
void eval_sum_clenshaw (int N, int M, double complex *a, double *x, double complex *y, double complex *temp, double *alpha, double *beta, double *gamma, double lambda)
 Clenshaw algorithm.
void eval_sum_clenshaw_transposed (int N, int M, double complex *a, double *x, double complex *y, double complex *temp, double *alpha, double *beta, double *gamma, double lambda)
 Clenshaw algorithm.
fpt_set fpt_init (const int M, const int t, const unsigned int flags)
 Initializes a set of precomputed data for DPT transforms of equal length.
void fpt_precompute (fpt_set set, const int m, const double *alpha, const double *beta, const double *gamma, int k_start, const double threshold)
 Computes the data required for a single DPT transform.
void dpt_trafo (fpt_set set, const int m, const double complex *x, double complex *y, const int k_end, const unsigned int flags)
 Computes a single DPT transform.
void fpt_trafo (fpt_set set, const int m, const double complex *x, double complex *y, const int k_end, const unsigned int flags)
 Computes a single DPT transform.
void dpt_transposed (fpt_set set, const int m, double complex *x, double complex *y, const int k_end, const unsigned int flags)
 Computes a single DPT transform.
void fpt_transposed (fpt_set set, const int m, double complex *x, const double complex *y, const int k_end, const unsigned int flags)
 Computes a single DPT transform.
void fpt_finalize (fpt_set set)


Detailed Description

Implementation file for the FPT module.

Author:
Jens Keiner

Definition in file fpt.c.


Define Documentation

#define ABUVXPWY_SYMMETRIC NAME,
S1,
S2   ) 
 

Value:

inline void NAME(double a, double b, double complex* u, double complex* x, double* v, \
  double complex* y, double* w, int n) \
{ \
  int l; \
  double complex *u_ptr, *x_ptr, *y_ptr; \
  double *v_ptr, *w_ptr; \
  \
  u_ptr = u; \
  x_ptr = x; \
  v_ptr = v; \
  y_ptr = y; \
  w_ptr = w; \
  \
  for (l = 0; l < n/2; l++) \
  { \
    *u_ptr++ = a * (b * (*v_ptr++) * (*x_ptr++) + (*w_ptr++) * (*y_ptr++)); \
  } \
  v_ptr--; \
  w_ptr--; \
  for (l = 0; l < n/2; l++) \
  { \
    *u_ptr++ = a * (b * S1 * (*v_ptr--) * (*x_ptr++) + S2 * (*w_ptr--) * (*y_ptr++)); \
  } \
}

Definition at line 133 of file fpt.c.

#define ABUVXPWY_SYMMETRIC_1 NAME,
S1   ) 
 

Value:

inline void NAME(double a, double b, double complex* u, double complex* x, double* v, \
  double complex* y, double* w, int n) \
{ \
  int l; \
  double complex *u_ptr, *x_ptr, *y_ptr; \
  double *v_ptr, *w_ptr; \
  \
  u_ptr = u; \
  x_ptr = x; \
  v_ptr = v; \
  y_ptr = y; \
  w_ptr = w; \
  \
  for (l = 0; l < n/2; l++) \
  { \
    *u_ptr++ = a * (b * (*v_ptr++) * (*x_ptr++) + (*w_ptr++) * (*y_ptr++)); \
  } \
  v_ptr--; \
  /*w_ptr--;*/ \
  for (l = 0; l < n/2; l++) \
  { \
    *u_ptr++ = a * (b * S1 * (*v_ptr--) * (*x_ptr++) + (*w_ptr++) * (*y_ptr++)); \
  } \
}

Definition at line 162 of file fpt.c.

#define ABUVXPWY_SYMMETRIC_2 NAME,
S1   ) 
 

Value:

inline void NAME(double a, double b, double complex* u, double complex* x, double* v, \
  double complex* y, double* w, int n) \
{ \
  int l; \
  double complex *u_ptr, *x_ptr, *y_ptr; \
  double *v_ptr, *w_ptr; \
  \
  u_ptr = u; \
  x_ptr = x; \
  v_ptr = v; \
  y_ptr = y; \
  w_ptr = w; \
  \
  for (l = 0; l < n/2; l++) \
  { \
    *u_ptr++ = a * (b * (*v_ptr++) * (*x_ptr++) + (*w_ptr++) * (*y_ptr++)); \
  } \
  /*v_ptr--;*/ \
  w_ptr--; \
  for (l = 0; l < n/2; l++) \
  { \
    *u_ptr++ = a * (b * (*v_ptr++) * (*x_ptr++) + S1 * (*w_ptr--) * (*y_ptr++)); \
  } \
}

Definition at line 191 of file fpt.c.

#define FPT_DO_STEP_TRANSPOSED NAME,
M1_FUNCTION,
M2_FUNCTION   ) 
 

Value:

inline void NAME(double complex  *a, double complex *b, double *a11, \
  double *a12, double *a21, double *a22, double gamma, int tau, fpt_set set) \
{ \ \
  int length = 1<<(tau+1); \ \
  double norm = 1.0/(length<<1); \
  \
  /* Compute function values from Chebyshev-coefficients using a DCT-III. */ \
  fftw_execute_r2r(set->plans_dct3[tau-1],(double*)a,(double*)a); \
  fftw_execute_r2r(set->plans_dct3[tau-1],(double*)b,(double*)b); \
  \
  /* Perform matrix multiplication. */ \
  M1_FUNCTION(norm,gamma,set->z,a,a11,b,a21,length); \
  M2_FUNCTION(norm,gamma,b,a,a12,b,a22,length); \
  memcpy(a,set->z,length*sizeof(double complex)); \
  \
  /* Compute Chebyshev-coefficients using a DCT-II. */ \
  fftw_execute_r2r(set->plans_dct2[tau-1],(double*)a,(double*)a); \
  fftw_execute_r2r(set->plans_dct2[tau-1],(double*)b,(double*)b); \
}

Definition at line 333 of file fpt.c.


Function Documentation

void fpt_do_step double complex *  a,
double complex *  b,
double *  a11,
double *  a12,
double *  a21,
double *  a22,
double  gamma,
int  tau,
fpt_set  set
[inline]
 

The length of the coefficient arrays.

Twice the length of the coefficient arrays.

Definition at line 328 of file fpt.c.

Referenced by fpt_trafo().

void fpt_do_step_symmetric double complex *  a,
double complex *  b,
double *  a11,
double *  a12,
double *  a21,
double *  a22,
double  gamma,
int  tau,
fpt_set  set
[inline]
 

The length of the coefficient arrays.

Twice the length of the coefficient arrays.

Definition at line 329 of file fpt.c.

Referenced by fpt_trafo().

void fpt_do_step_symmetric_u double complex *  a,
double complex *  b,
double *  a11,
double *  a12,
double *  a21,
double *  a22,
double  gamma,
int  tau,
fpt_set  set
[inline]
 

The length of the coefficient arrays.

Twice the length of the coefficient arrays.

Definition at line 330 of file fpt.c.

Referenced by fpt_trafo().

void fpt_do_step_symmetric_l double complex *  a,
double complex *  b,
double *  a11,
double *  a12,
double *  a21,
double *  a22,
double  gamma,
int  tau,
fpt_set  set
[inline]
 

The length of the coefficient arrays.

Twice the length of the coefficient arrays.

Definition at line 331 of file fpt.c.

Referenced by fpt_trafo().

void fpt_do_step_transposed double complex *  a,
double complex *  b,
double *  a11,
double *  a12,
double *  a21,
double *  a22,
double  gamma,
int  tau,
fpt_set  set
[inline]
 

The length of the coefficient arrays.

Twice the length of the coefficient arrays.

Definition at line 356 of file fpt.c.

Referenced by fpt_transposed().

void fpt_do_step_transposed_symmetric double complex *  a,
double complex *  b,
double *  a11,
double *  a12,
double *  a21,
double *  a22,
double  gamma,
int  tau,
fpt_set  set
[inline]
 

The length of the coefficient arrays.

Twice the length of the coefficient arrays.

Definition at line 357 of file fpt.c.

Referenced by fpt_transposed().

void fpt_do_step_transposed_symmetric_u double complex *  a,
double complex *  b,
double *  a11,
double *  a12,
double *  a21,
double *  a22,
double  gamma,
int  tau,
fpt_set  set
[inline]
 

The length of the coefficient arrays.

Twice the length of the coefficient arrays.

Definition at line 358 of file fpt.c.

Referenced by fpt_transposed().

void fpt_do_step_transposed_symmetric_l double complex *  a,
double complex *  b,
double *  a11,
double *  a12,
double *  a21,
double *  a22,
double  gamma,
int  tau,
fpt_set  set
[inline]
 

The length of the coefficient arrays.

Twice the length of the coefficient arrays.

Definition at line 359 of file fpt.c.

Referenced by fpt_transposed().

void eval_sum_clenshaw int  N,
int  M,
double complex *  a,
double *  x,
double complex *  y,
double complex *  temp,
double *  alpha,
double *  beta,
double *  gamma,
double  lambda
 

Clenshaw algorithm.

Evaluates a sum of real-valued functions $P_k : \mathbb{R} \rightarrow \mathbb{R}$

\[ f(x) = \sum_{k=0}^N a_k P_k(x) \quad (N \in \mathbb{N}_0) \]

obeying a three-term recurrence relation

\[ P_{k+1}(x) = (alpha_k * x + beta_k)*P_{k}(x) + gamma_k P_{k-1}(x) \quad (alpha_k, beta_k, gamma_k \in \mathbb{R},\; k \ge 0) \]

with initial conditions $P_{-1}(x) := 0$ , $P_0(x) := \lambda$ for given double complex coefficients $\left(a_k\right)_{k=0}^N \in \mathbb{C}^{N+1}$ at given nodes $\left(x_j\right)_{j=0}^M \in \mathbb{R}^{M+1}$ , $M \in \mathbb{N}_0$ .

Definition at line 477 of file fpt.c.

Referenced by dpt_trafo().

void eval_sum_clenshaw_transposed int  N,
int  M,
double complex *  a,
double *  x,
double complex *  y,
double complex *  temp,
double *  alpha,
double *  beta,
double *  gamma,
double  lambda
 

Clenshaw algorithm.

Evaluates a sum of real-valued functions $P_k : \mathbb{R} \rightarrow \mathbb{R}$

\[ f(x) = \sum_{k=0}^N a_k P_k(x) \quad (N \in \mathbb{N}_0) \]

obeying a three-term recurrence relation

\[ P_{k+1}(x) = (alpha_k * x + beta_k)*P_{k}(x) + gamma_k P_{k-1}(x) \quad (alpha_k, beta_k, gamma_k \in \mathbb{R},\; k \ge 0) \]

with initial conditions $P_{-1}(x) := 0$ , $P_0(x) := \lambda$ for given double complex coefficients $\left(a_k\right)_{k=0}^N \in \mathbb{C}^{N+1}$ at given nodes $\left(x_j\right)_{j=0}^M \in \mathbb{R}^{M+1}$ , $M \in \mathbb{N}_0$ .

Definition at line 543 of file fpt.c.

Referenced by dpt_transposed().


Generated on 1 Nov 2006 by Doxygen 1.4.4