FNFT
Files | Classes | Functions
Fast nonlinear Fourier transforms
Collaboration diagram for Fast nonlinear Fourier transforms:

Files

file  fnft_kdv_discretization_t.h
 Lists discretizations for the Korteweg-de Vries equation.
 
file  fnft_kdvv.h
 Fast nonlinear Fourier transform for the vanishing Korteweg-de Vries equation.
 
file  fnft_nse_discretization_t.h
 Lists discretizations for the nonlinear Schroedinger equation.
 
file  fnft_nsep.h
 Fast nonlinear Fourier transform for the periodic nonlinear Schroedinger equation.
 
file  fnft_nsev.h
 Fast nonlinear Fourier transform for the vanishing nonlinear Schroedinger equation.
 
file  fnft_version.h
 Provides a function to query the version of the FNFT library.
 

Classes

struct  fnft_kdvv_opts_t
 Stores additional options for the routine fnft_kdvv. More...
 
struct  fnft_nsep_opts_t
 Stores additional options for the routine fnft_nsep. More...
 
struct  fnft_nsev_opts_t
 Stores additional options for the routine fnft_nsev. More...
 

Functions

fnft_kdvv_opts_t fnft_kdvv_default_opts ()
 Creates a new options variable for fnft_kdvv with default settings. More...
 
FNFT_INT fnft_kdvv (const FNFT_UINT D, FNFT_COMPLEX *const q, FNFT_REAL const *const T, const FNFT_UINT M, FNFT_COMPLEX *const contspec, FNFT_REAL const *const XI, FNFT_UINT *const K_ptr, FNFT_COMPLEX *const bound_states, FNFT_COMPLEX *const normconsts_or_residues, fnft_kdvv_opts_t *opts_ptr)
 Fast nonlinear Fourier transform for the Korteweg-de Vries equation with vanishing boundary conditions. More...
 
fnft_nsep_opts_t fnft_nsep_default_opts ()
 Creates a new options variable for fnft_nsep with default settings. More...
 
FNFT_INT fnft_nsep (const FNFT_UINT D, FNFT_COMPLEX const *const q, FNFT_REAL const *const T, FNFT_UINT *const K_ptr, FNFT_COMPLEX *const main_spec, FNFT_UINT *const M_ptr, FNFT_COMPLEX *const aux_spec, FNFT_REAL *const sheet_indices, const FNFT_INT kappa, fnft_nsep_opts_t *opts)
 Fast nonlinear Fourier transform for the nonlinear Schroedinger equation with periodic boundary conditions. More...
 
fnft_nsev_opts_t fnft_nsev_default_opts ()
 Creates a new options variable for fnft_nsev with default settings. More...
 
FNFT_UINT fnft_nsev_max_K (const FNFT_UINT D, fnft_nsev_opts_t const *const opts)
 Returns the maximum number of bound states that can be detected by fnft_nsev. More...
 
FNFT_INT fnft_nsev (const FNFT_UINT D, FNFT_COMPLEX *const q, FNFT_REAL const *const T, const FNFT_UINT M, FNFT_COMPLEX *const contspec, FNFT_REAL const *const XI, FNFT_UINT *const K_ptr, FNFT_COMPLEX *const bound_states, FNFT_COMPLEX *const normconsts_or_residues, const FNFT_INT kappa, fnft_nsev_opts_t *opts)
 Fast nonlinear Fourier transform for the nonlinear Schroedinger equation with vanishing boundary conditions. More...
 
FNFT_INT fnft_version (FNFT_UINT *const major, FNFT_UINT *const minor, FNFT_UINT *const patch, char suffix[FNFT_VERSION_SUFFIX_MAXLEN+1])
 

Detailed Description

Function Documentation

◆ fnft_kdvv()

FNFT_INT fnft_kdvv ( const FNFT_UINT  D,
FNFT_COMPLEX *const  q,
FNFT_REAL const *const  T,
const FNFT_UINT  M,
FNFT_COMPLEX *const  contspec,
FNFT_REAL const *const  XI,
FNFT_UINT *const  K_ptr,
FNFT_COMPLEX *const  bound_states,
FNFT_COMPLEX *const  normconsts_or_residues,
fnft_kdvv_opts_t opts_ptr 
)

Fast nonlinear Fourier transform for the Korteweg-de Vries equation with vanishing boundary conditions.

This routine computes the nonlinear Fourier transform for the Korteweg-de Vries equation

\[ q_x + 6qq_{t} + q_{ttt}=0, \quad q=q(x,t), \]

of Gardner et al. (Phys. Rev. Lett., 1967) for initial conditions with vanishing boundaries

\[ \lim_{t\to \pm \infty }q(x_0,t) = 0 \text{ sufficiently rapidly.} \]

Parameters
[in]DNumber of samples.
[in]qArray of length D, contains samples \( q(t_n)=q(x_0, t_n) \), where \( t_n = T[0] + n(T[1]-T[0])/(D-1) \) and \(n=0,1,\dots,D-1\), of the to-be-transformed signal in ascending order (i.e., \( q(t_0), q(t_1), \dots, q(t_{D-1}) \)).
[in]TArray of length 2, contains the position in time of the first and of the last sample. It should be T[0]<T[1].
[in]MNumber of points at which the continuous spectrum (aka reflection coefficient) should be computed.
[out]contspecArray of length M in which the routine will store the desired samples \( R(\xi_m) \) of the continuous spectrum (aka reflection coefficient) in ascending order, where \( \xi_m = XI[0]+m(XI[1]-XI[0])/(M-1) \) and \(m=0,1,\dots,M-1\). Has to be preallocated by the user.
[in]XIArray of length 2, contains the position of the first and the last sample of the continuous spectrum. It should be XI[0]<XI[1]. Can also be NULL if contspec==NULL.
[in,out]K_ptrNot yet implemented. Pass NULL.
[out]bound_statesNot yet implemented. Pass NULL.
[out]normconsts_or_residuesNot yet implemented. Pass NULL.
[in]opts_ptrPointer to a fnft_kdvv_opts_t object. The object can be used to modify the behavior of the routine. Use the routine fnft_kdvv_default_opts to generate such an object and modify as desired. It is also possible to pass NULL, in which case the routine will use the default options. The user is reponsible to freeing the object after the routine has returned.
Returns
FNFT_SUCCESS or one of the FNFT_EC_... error codes defined in fnft_errwarn.h.

◆ fnft_kdvv_default_opts()

fnft_kdvv_opts_t fnft_kdvv_default_opts ( )

Creates a new options variable for fnft_kdvv with default settings.

Returns
A fnft_kdvv_opts_t object with default options.

◆ fnft_nsep()

FNFT_INT fnft_nsep ( const FNFT_UINT  D,
FNFT_COMPLEX const *const  q,
FNFT_REAL const *const  T,
FNFT_UINT *const  K_ptr,
FNFT_COMPLEX *const  main_spec,
FNFT_UINT *const  M_ptr,
FNFT_COMPLEX *const  aux_spec,
FNFT_REAL *const  sheet_indices,
const FNFT_INT  kappa,
fnft_nsep_opts_t opts 
)

Fast nonlinear Fourier transform for the nonlinear Schroedinger equation with periodic boundary conditions.



This routine computes the nonlinear Fourier transform for the nonlinear Schroedinger equation

\[ iq_x + q_{tt} \pm 2q|q|^2=0, \quad q=q(x,t), \]

of Kotlyarov and Its (Problems of Mathemetical Physics and Functional Analysis, 1976) for initial conditions with periodic boundary conditions,

\[ q(x_0,t + L) = q(x_0,t). \]


The main references for the numerical algorithm are:

  • Wahls and Poor, "Fast numerical nonlinear Fourier transforms," IEEE Trans. Inform. Theor. 61(12), 2015.
  • Prins and Wahls, "Higher order exponential splittings for the fast non-linear Fourier transform of the KdV equation," to appear in Proc. ICASSP 2018.
Parameters
[in]DNumber of samples. Has to be a power of two.
[in]qArray of length D, contains samples \( q(t_n)=q(x_0, t_n) \), where \( t_n = T[0] + n*L/D \), where L=T[2]-T[1] is the period and \(n=0,1,\dots,D-1\), of the to-be-transformed signal in ascending order (i.e., \( q(t_0), q(t_1), \dots, q(t_{D-1}) \))
[in]TArray of length 2. T[0] is the position in time of the first sample. T[2] is the beginning of the next period. (The location of the last sample is thus t_{D-1}=T[2]-L/D.) It should be T[0]<T[1].
[in,out]K_ptrUpon entry, *K_ptr should contain the length of the array main_spec. Upon return, *K_ptr contains the number of actually detected points in the main spectrum. If the length of the array main_spec was not sufficient to store all of the detected points in the main spectrum, a warning is printed and as many points as possible are returned instead. Note that in order to skip the computation of the main spectrum completely, it is not sufficient to pass *K_ptr==NULL. Instead, pass main_spec==NULL.
[out]main_specArray. Upon return, the routine has stored the detected main specrum points (i.e., the points for which the trace of the monodromy matrix is either +2 or -2) in the first *K_ptr entries of this array. If NULL is passed instead, the main spectrum will not be computed. Has to be preallocated by the user. The user can choose an arbitrary length. Typically, D is a good choice.
[in,out]M_ptrUpon entry, *M_ptr should contain the length of the array aux_spec. Upon return, *M_ptr contains the number of actually detected points in the auxiliary spectrum. If the length of the array aux_spec was not sufficient to store all of the detected points in the auxiliary spectrum, a warning is printed and as many points as possible are returned instead. Note that in order to skip the computation of the auxiliary spectrum completely, it is not sufficient to pass *M_ptr==NULL. Instead, pass aux_spec==NULL.
[out]aux_specArray. Upon return, the routine has stored the detected auxiliary specrum points (i.e., the points for which the upper right element of the monodromy matrix is zero) in the first *M_ptr entries of this array. If NULL is passed instead, the auxiliary spectrum will not be computed. Has to be preallocated by the user. The user can choose an arbitrary length. Typically, D is a good choice.
[in]sheet_indicesNot yet implemented. Pass NULL.
[in]kappa=+1 for the focusing nonlinear Schroedinger equation, =-1 for the defocusing one
[in]optsPointer to a fnft_nsep_opts_t object. The object can be used to modify the behavior of the routine. Use the routine fnft_nsep_default_opts to generate such an object and modify as desired. It is also possible to pass NULL, in which case the routine will use the default options. The user is reponsible to freeing the object after the routine has returned.
Returns
FNFT_SUCCESS or one of the FNFT_EC_... error codes defined in fnft_errwarn.h.

◆ fnft_nsep_default_opts()

fnft_nsep_opts_t fnft_nsep_default_opts ( )

Creates a new options variable for fnft_nsep with default settings.

Returns
A fnft_nsep_opts_t object with the following options.
localization = fnft_nsep_loc_MIXED
filtering = fnft_nsep_filt_AUTO
max_evals = 20
bounding_box[0] = -FNFT_INF
bounding_box[1] = FNFT_INF
bounding_box[2] = -FNFT_INF
bounding_box[3] = FNFT_INF
normalization_flag = 1
discretization = fnft_nse_discretization_2SPLIT2A
floquet_range = {-1, 1}
floquet_nvals = 2
Dsub = 0

◆ fnft_nsev()

FNFT_INT fnft_nsev ( const FNFT_UINT  D,
FNFT_COMPLEX *const  q,
FNFT_REAL const *const  T,
const FNFT_UINT  M,
FNFT_COMPLEX *const  contspec,
FNFT_REAL const *const  XI,
FNFT_UINT *const  K_ptr,
FNFT_COMPLEX *const  bound_states,
FNFT_COMPLEX *const  normconsts_or_residues,
const FNFT_INT  kappa,
fnft_nsev_opts_t opts 
)

Fast nonlinear Fourier transform for the nonlinear Schroedinger equation with vanishing boundary conditions.

This routine computes the nonlinear Fourier transform for the nonlinear Schroedinger equation

\[ iq_x + q_{tt} \pm 2q|q|^2=0, \quad q=q(x,t), \]

of Zakharov and Shabat (Soviet. Phys. JTEP 31(1), 1972) for initial conditions with vanishing boundaries

\[ \lim_{t\to \pm \infty }q(x_0,t) = 0 \text{ sufficiently rapidly.} \]


The main references are:

The routine also utilizes ideas from the following papers:

Parameters
[in]DNumber of samples
[in]qArray of length D, contains samples \( q(t_n)=q(x_0, t_n) \), where \( t_n = T[0] + n(T[1]-T[0])/(D-1) \) and \(n=0,1,\dots,D-1\), of the to-be-transformed signal in ascending order (i.e., \( q(t_0), q(t_1), \dots, q(t_{D-1}) \))
[in]TArray of length 2, contains the position in time of the first and of the last sample. It should be T[0]<T[1].
[in]MNumber of points at which the continuous spectrum (aka reflection coefficient) should be computed.
[out]contspecArray of length M in which the routine will store the desired samples \( r(\xi_m) \) of the continuous spectrum (aka reflection coefficient) in ascending order, where \( \xi_m = XI[0]+(XI[1]-XI[0])/(M-1) \) and \(m=0,1,\dots,M-1\). Has to be preallocated by the user. If NULL is passed instead, the continuous spectrum will not be computed. By changing the options, it is also possible to compute the values of \( a(\xi) \) and \( b(\xi) \) instead. In that case, twice the amount of memory has to be allocated.
[in]XIArray of length 2, contains the position of the first and the last sample of the continuous spectrum. It should be XI[0]<XI[1]. Can also be NULL if contspec==NULL.
[in,out]K_ptrUpon entry, *K_ptr should contain the length of the array bound_states. Upon return, *K_ptr contains the number of actually detected bound states. If the length of the array bound_states was not sufficient to store all of the detected bound states, a warning is printed and as many bound states as possible are returned instead. Note that in order to skip the computation of the bound states completely, it is not sufficient to pass *K_ptr==0. Instead, one needs to pass bound_states==NULL.
[out]bound_statesArray. Upon return, the routine has stored the detected bound states (aka eigenvalues) in the first *K_ptr entries of this array. If NULL is passed instead, the discrete spectrum will not be computed. Has to be preallocated by the user. The user can choose an arbitrary length. Typically, D is a good choice.
[out]normconsts_or_residuesArray of the same length as bound_states. Upon return, the routine has stored the residues (aka spectral amplitudes) \(\rho_k = b_k\big/ \frac{da(\lambda_k)}{d\lambda}\) in the first *K_ptr entries of this array. By passing a proper opts, it is also possible to store the norming constants (the ' \( b_k \)') or both. Has to be pre-allocated by the user. If NULL is passed instead, the residues will not be computed.
[in]kappa=+1 for the focusing nonlinear Schroedinger equation, =-1 for the defocusing one
[in]optsPointer to a fnft_nsev_opts_t object. The object can be used to modify the behavior of the routine. Use the routine fnft_nsev_default_opts to generate such an object and modify as desired. It is also possible to pass NULL, in which case the routine will use the default options. The user is reponsible to freeing the object after the routine has returned.
Returns
FNFT_SUCCESS or one of the FNFT_EC_... error codes defined in fnft_errwarn.h.

◆ fnft_nsev_default_opts()

fnft_nsev_opts_t fnft_nsev_default_opts ( )

Creates a new options variable for fnft_nsev with default settings.

Returns
A fnft_nsev_opts_t object with the following options.
bound_state_filtering = fnft_nsev_bsfilt_FULL
bound_state_localization = fnft_nsev_bsloc_SUBSAMPLE_AND_REFINE
niter = 10
discspec_type = fnft_nsev_dstype_NORMING_CONSTANTS
contspec_type = fnft_nsev_cstype_REFLECTION_COEFFICIENT
normalization_flag = 1
discretization = fnft_nse_discretization_2SPLIT4B

◆ fnft_nsev_max_K()

FNFT_UINT fnft_nsev_max_K ( const FNFT_UINT  D,
fnft_nsev_opts_t const *const  opts 
)

Returns the maximum number of bound states that can be detected by fnft_nsev.

Parameters
[in]DNumber of samples that will be passed to fnft_nsev. Should be larger than zero.
[in]optsOptions that will be passed to fnft_nsev. If NULL is passed, the default options will be used.
Returns
Returns the maximum number of bound states or zero on error.

◆ fnft_version()

FNFT_INT fnft_version ( FNFT_UINT *const  major,
FNFT_UINT *const  minor,
FNFT_UINT *const  patch,
char  suffix[FNFT_VERSION_SUFFIX_MAXLEN+1] 
)

Provides the version of the FNFT library.

The version is of the form $major.$minor.$patch$suffix, e.g., "0.2.1-al". The function simply returns the corresponding values defined fnft_config.h.

Parameters
[out]major*major is overwritten with the value of FNFT_VERSION_MAJOR
[out]minor*minor is overwritten with the value of FNFT_VERSION_MINOR
[out]patch*patch is overwritten with the value of FNFT_VERSION_PATCH
[out]suffixsuffix is overwritten with the value of FNFT_VERSION_SUFFIX
Returns
FNFT_SUCCESS or one of the FNFT_EC_... error codes defined in fnft_errwarn.h.