SPRAAK
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
Data Structures | Macros | Typedefs | Enumerations | Functions | Variables
API_LowLvl::lib::math
+ Collaboration diagram for API_LowLvl::lib::math:

Data Structures

union  SprFMIntAsFlt
 
union  SprFMFltAsInt
 
struct  SprVarlist
 

Macros

#define SprFMFltU
 unaligned float vector of optimal vector size More...
 
#define SprFMFltA
 aligned float vector of optimal vector size More...
 
#define SprFMIntU
 unaligned int vector of optimal vector size More...
 
#define SprFMIntA
 aligned int vector of optimal vector size More...
 
#define SprFMUIntU
 unaligned unsigned int vector of optimal vector size More...
 
#define SprFMUIntA
 aligned unsigned int vector of optimal vector size More...
 

Typedefs

typedef float SprFevalReal
 
typedef struct SprRandGen SprRandGen
 

Enumerations

enum  {
  SPR_FEVAL_IS_VEC, SPR_FEVAL_IS_STR, SPR_FEVAL_IS_FUN1, SPR_FEVAL_IS_IFX1,
  SPR_FEVAL_IS_IFX2, SPR_FEVAL_IS_FUN2, SPR_FEVAL_IS_FUN3, SPR_FEVAL_IS_FUN,
  SPR_FEVAL_IS_UNDEF, SPR_FEVAL_IS_SCRATCH, SPR_FEVAL_IS_CONST, SPR_FEVAL_IS_PROTECT
}
 
enum  {
  SPR_FEVAL_VAR_UNDEF, SPR_FEVAL_VAR_VAL, SPR_FEVAL_VAR_VEC, SPR_FEVAL_VAR_STR,
  SPR_FEVAL_VAR_FUN1, SPR_FEVAL_VAR_FUN2, SPR_FEVAL_VAR_FUN3
}
 
enum  { SPR_FEVAL_NO_VEC, SPR_FEVAL_NO_IGN, SPR_FEVAL_NO_END, SPR_FEVAL_COMPILE }
 explanation: see feval More...
 
enum  {
  SPR_FEVAL_SUBST_VAR, SPR_FEVAL_SUBST_ENV, SPR_FEVAL_SUBST_ENV_BR, SPR_FEVAL_SUBST_ENV_EXT,
  SPR_FEVAL_SUBST_ENV_EXPR, SPR_FEVAL_SUBST_FREE, SPR_FEVAL_SUBST_DYNSTR, SPR_FEVAL_SUBST_KEEP,
  SPR_FEVAL_SUBST_KEEP_ERR
}
 

Functions

float spr_math_fast_erfinv (float x)
 
double spr_math_erfinv (double x)
 
double spr_math_normcdf_inv (double p)
 
double spr_math_normcdf_inv2 (double r)
 
SprFMFltA spr_fm_vx_max (SprFMFltA x, SprFMFltA y)
 Element wise maximum. More...
 
SprFMFltA spr_fm_vx_min (SprFMFltA x, SprFMFltA y)
 Element wise minimum. More...
 
SprFMFltA spr_fm_vx_exp (SprFMFltA x)
 
float spr_fm_exp (float x)
 
SprFMFltA spr_fm_vx_log (SprFMFltA p)
 
float spr_fm_log (float p)
 
SprFMFltA spr_fm_vx_pow (SprFMFltA x, SprFMFltA y)
 
float spr_fm_pow (float x, float y)
 
SprFMFltA spr_fm_vx_logadd_ (SprFMFltA x, SprFMFltA y)
 
float spr_fm_logadd_ (float x, float y)
 
SprFMFltA spr_fm_vx_logadd (SprFMFltA x, SprFMFltA y)
 
float spr_fm_logadd (float x, float y)
 
float spr_fm_vx_sum (SprFMFltA x)
 
SprFMFltA spr_fm_vx_getN (const float *x, int N)
 
float * spr_fm_vec_exp (float *restrict vec, int N)
 
float * spr_fm_vec_log (float *restrict vec, int N)
 
float * spr_fm_vec_log_offs (float *restrict vec, int N, float offs)
 
float * spr_fm_vec_log_mul (float *restrict vec, int N, float fac)
 
float * spr_fm_vec_log_offs_mul (float *restrict vec, int N, float offs, float fac)
 
SprVarlistspr_feval_var_init (SprVarlist *var_ptr, const char *name)
 
SprVarlistspr_feval_var_clear (SprVarlist *var_ptr)
 
int spr_feval_result_set (SprVarlist *res_ptr, int type,...)
 
SprVarlistspr_feval_var_new (const char *name, int type,...)
 
SprVarlistspr_feval_var_add (SprVarlist *varlist, SprVarlist *argptr)
 Add an item to the list of known variables. More...
 
void spr_feval_var_free (SprVarlist *var_ptr)
 
SprVarlistspr_feval_var_rm (SprVarlist *varlist, SprVarlist *var_ptr)
 
SprVarlistspr_feval_varlist_free (SprVarlist *varlist, SprVarlist *upto)
 
SprVarlistspr_feval_var_get (SprVarlist *varlist, const char *name, int length)
 
int spr_feval_var_set (SprVarlist *var_ptr, int type,...)
 
SprVarlistspr_feval_var_assign (SprVarlist *varlist, const char *name, int type,...)
 
void spr_feval_whos (SprStream *fd, SprVarlist *varlist, const char *name)
 
void spr_feval_var_print (SprStream *fd, const SprVarlist *var_ptr, const char *rm)
 Print a variable. More...
 
char * spr_feval (const char *eval_str, const char *eval_end, SprVarlist *varlist, SprVarlist *res, int flags)
 
SprFevalReal spr_feval_num (const char *eval_str, const char *eval_end, SprVarlist *varlist, int *err_ret)
 
char * spr_feval_var_substitute (const char *str, SprVarlist *list, int flags)
 
float * spr_f32_fht1 (float *restrict x, int N)
 
float * spr_f32_fht2 (float *restrict x, int N)
 
float * spr_f32_fht4 (float *restrict x, int N)
 
float * spr_f32_fft_c2c (float *restrict x, int N)
 
float * spr_f32_ifft_c2c (float *restrict x, int N)
 
float * spr_f32_fft_r2c (float *restrict x, int N)
 
float * spr_f32_ifft_c2r (float *restrict x, int N)
 
double * spr_f64_fht1 (double *restrict x, int N)
 
double * spr_f64_fht2 (double *restrict x, int N)
 
double * spr_f64_fht4 (double *restrict x, int N)
 
double * spr_f64_fft_c2c (double *restrict x, int N)
 
double * spr_f64_ifft_c2c (double *restrict x, int N)
 
double * spr_f64_fft_r2c (double *restrict x, int N)
 
double * spr_f64_ifft_c2r (double *restrict x, int N)
 
double ** spr_lda_stats_alloc (double **moments, int nclasses, int vlen)
 
double ** spr_lda_stats_diag_alloc (double **moments, int nclasses, int vlen)
 
void spr_lda_stats_update (double *restrict moments, const SprMathReal *data, int vlen, float weight)
 Update the 2nd order statistics for a given vector (and class). More...
 
void spr_lda_stats_diag_update (double *restrict moments, const SprMathReal *data, int vlen, float weight)
 
void spr_lda_calc_class_covar (double **moments, int vlen, int nclasses, double Beta, SprStatsStats *within_class_covar, SprStatsStats *between_class_covar)
 
void spr_lda_calc_transf (SprStatsStats *lda_transf, int vlen, SprStatsStats *within_class_covar, SprStatsStats *between_class_covar)
 
int spr_math_mat_print (SprStream *fd, const void *A, int nrows, int ncols, const SprDT *dt, const char *info)
 
float * spr_mat_f32_A_plus_scaleB_insitu (float *restrict mA, int size, const float *mB, float facB)
 
float * spr_mat_f32_scaleA_plus_B (float *restrict mC, int size, const float *mA, const float *mB, float facA)
 
float * spr_mat_f32_scaleA_plus_scaleB (float *restrict mC, int size, const float *mA, const float *mB, float facA, float facB)
 
float * spr_mat_f32_scaleA_insitu (float *restrict mA, float facA, int size)
 
float * spr_mat_f32_scaleA (float *restrict mC, int size, const float *mA, float facA)
 
float * spr_mat_f32_A_plus_B (float *restrict mC, int size, const float *mA, const float *mB)
 
float * spr_mat_f32_A_plus_B_insitu (float *restrict mA, int size, const float *mB)
 
float * spr_mat_f32_A_minus_B (float *restrict mC, int size, const float *mA, const float *mB)
 
float * spr_mat_f32_A_minus_B_insitu (float *restrict mA, int size, const float *mB)
 
float * spr_mat_f32_A_dotmul_B (float *restrict mC, int size, const float *mA, const float *mB)
 
float * spr_mat_f32_A_dotdiv_B (float *restrict mC, int size, const float *mA, const float *mB)
 
double spr_mat_f32_A_inprod_B (const float *mA, int size, const float *mB)
 
double spr_mat_f32_rowA_inprod_colB (const float *mA, int size, const float *mB, int ncolB)
 
float * spr_mat_f32_A_mul_B (float *restrict mC, const float *mA, int nrowA, int ncolA, const float *mB, int ncolB)
 
float * spr_mat_f32_At_mul_B (float *restrict mC, const float *mA, int nrowA, int ncolA, const float *mB, int ncolB)
 
float * spr_mat_f32_A_mul_Bt (float *restrict mC, const float *mA, int nrowA, int ncolA, const float *mB, int nrowB)
 
float * spr_mat_f32_At_mul_Bt (float *restrict mC, const float *mA, int nrowA, int ncolA, const float *mB, int nrowB)
 
float * spr_mat_f32_colA (float *restrict mC, const float *mA, int nrowA, int ncolA)
 
float * spr_mat_f32_At (float *restrict mC, const float *mA, int nrowA, int ncolA)
 
double spr_mat_f32_detA (float *restrict mA, float *restrict copy, int size)
 
double spr_mat_f32_norm2A (const float *mA, int size)
 
double spr_mat_f32_prodA (const float *mA, int size)
 
double * spr_mat_f64_A_plus_scaleB_insitu (double *restrict mA, int size, const double *mB, double facB)
 
double * spr_mat_f64_scaleA_plus_B (double *restrict mC, int size, const double *mA, const double *mB, double facA)
 
double * spr_mat_f64_scaleA_plus_scaleB (double *restrict mC, int size, const double *mA, const double *mB, double facA, double facB)
 
double * spr_mat_f64_scaleA_insitu (double *restrict mA, double facA, int size)
 
double * spr_mat_f64_scaleA (double *restrict mC, int size, const double *mA, double facA)
 
double * spr_mat_f64_A_plus_B (double *restrict mC, int size, const double *mA, const double *mB)
 
double * spr_mat_f64_A_plus_B_insitu (double *restrict mA, int size, const double *mB)
 
double * spr_mat_f64_A_minus_B (double *restrict mC, int size, const double *mA, const double *mB)
 
double * spr_mat_f64_A_minus_B_insitu (double *restrict mA, int size, const double *mB)
 
double * spr_mat_f64_A_dotmul_B (double *restrict mC, int size, const double *mA, const double *mB)
 
double * spr_mat_f64_A_dotdiv_B (double *restrict mC, int size, const double *mA, const double *mB)
 
double spr_mat_f64_A_inprod_B (const double *mA, int size, const double *mB)
 
double spr_mat_f64_rowA_inprod_colB (const double *mA, int size, const double *mB, int ncolB)
 
double * spr_mat_f64_A_mul_B (double *restrict mC, const double *mA, int nrowA, int ncolA, const double *mB, int ncolB)
 
double * spr_mat_f64_At_mul_B (double *restrict mC, const double *mA, int nrowA, int ncolA, const double *mB, int ncolB)
 
double * spr_mat_f64_A_mul_Bt (double *restrict mC, const double *mA, int nrowA, int ncolA, const double *mB, int nrowB)
 
double * spr_mat_f64_At_mul_Bt (double *restrict mC, const double *mA, int nrowA, int ncolA, const double *mB, int nrowB)
 
double * spr_mat_f64_colA (double *restrict mC, const double *mA, int nrowA, int ncolA)
 
double * spr_mat_f64_At (double *restrict mC, const double *mA, int nrowA, int ncolA)
 
double spr_mat_f64_detA (double *restrict mA, double *restrict copy, int size)
 
double spr_mat_f64_norm2A (const double *mA, int size)
 
double spr_mat_f64_prodA (const double *mA, int size)
 
SprMathReal * spr_math_omat_to_nmat (int nrow, int ncol, const SprMathReal *const *from, SprMathReal *to)
 
SprMathReal ** spr_math_nmat_to_omat (int nrow, int ncol, const SprMathReal *from, SprMathReal **to)
 
void spr_math_fmat_to_fmat_trunc (const SprMathReal *from, int from_ncol, SprMathReal *restrict to, int to_nrow, int to_ncol)
 
void spr_math_utmat_to_fmat (int vlen, const SprMathReal *from, SprMathReal *to)
 
void spr_math_utmat_to_fmat_trunc (int vlen, const SprMathReal *from, int nlen, SprMathReal *restrict to)
 
void spr_math_utmat_to_utmat_trunc (int vlen, const SprMathReal *from, int nlen, SprMathReal *restrict to)
 
void spr_math_fmat_to_utmat (int vlen, const SprMathReal *from, SprMathReal *to)
 
SprMathReal * spr_math_utmat_to_diag (int vlen, const SprMathReal *from, SprMathReal *restrict to)
 
SprMathReal * spr_math_fmat_to_diag (int vlen, const SprMathReal *from, SprMathReal *restrict to)
 
void spr_math_v_lin_comb (SprMathReal *vec1, int vlen, SprMathReal *vec2, double c, double s)
 
SprMathReal * spr_math_hess_fact_t (SprMathReal *A, int vlen, SprMathReal *Q)
 
SprMathReal * spr_math_diag_eig_t (SprMathReal *diag0, int vlen, SprMathReal *diag1, SprMathReal *eigvec)
 
SprMathReal * spr_math_symm_eig_t (SprMathReal *A, int vlen, SprMathReal *eigvec, SprMathReal *eigval)
 
SprMathReal * spr_math_scale_eig_vec_t (int vlen, SprMathReal *eig_vec, SprMathReal *eig_val)
 
int spr_math_io_nmat_write (const char *fname, const void *mat, int Nrows, int Ncols, const SprDT *dt)
 
void * spr_math_io_nmat_read (const char *fname, void *mat, int *rows, int *cols, const SprDT *dt)
 
SprMathReal * spr_math_lu_full (SprMathReal *mA, int nrow, int ncol, int *pivot)
 
double spr_math_lu_full_det (const SprMathReal *LU, int vlen, const int *pivot)
 
SprMathReal * spr_math_lu_l1_inv (SprMathReal *res, SprMathReal *L1, int vlen)
 
SprMathReal * spr_math_lu_l1_complete (SprMathReal *restrict L1, int vlen)
 
SprMathReal * spr_math_lu_l1_solve (SprMathReal *L1, SprMathReal *B, int nrow, int ncol)
 
SprMathReal * spr_math_lu_u_solve (SprMathReal *U, SprMathReal *B, int nrow, int ncol)
 
SprMathReal * spr_math_lu_full_column_pivot (SprMathReal *A, int nrow, int ncol, int *pivot)
 Pivot columns. More...
 
SprMathReal * spr_math_lu_full_column_pivot_r (SprMathReal *A, int nrow, int ncol, int *pivot)
 Pivot columns (reverse order). More...
 
SprMathReal * spr_math_lu_full_row_pivot (SprMathReal *A, int nrow, int ncol, int *pivot)
 Pivot rows. More...
 
SprMathReal * spr_math_lu_full_row_pivot_r (SprMathReal *A, int nrow, int ncol, int *pivot)
 Pivot rows (reverse order). More...
 
SprMathReal * spr_math_lu_full_inv (SprMathReal *res, SprMathReal *LU, int vlen, int *pivot)
 
SprMathReal * spr_math_lu_full_solve (SprMathReal *LU, SprMathReal *B, int nrow, int ncol, int *pivot)
 
SprMathReal * spr_math_lu_symm (SprMathReal *restrict mA, int vlen, int *restrict pvt)
 
double spr_math_lu_symm_det (SprMathReal *LU, int vlen)
 
SprMathReal * spr_math_lu_symm_pivot (SprMathReal *A, int vlen, const int *pivot)
 Pivot rows and columns. More...
 
SprMathReal * spr_math_lu_symm_pivot_r (SprMathReal *A, int vlen, const int *pivot)
 Pivot rows and columns (reverse order). More...
 
SprMathReal * spr_math_lu_symm_u1_inv (SprMathReal *restrict U1, int vlen)
 
void spr_math_lu_symm_u1_scale (SprMathReal *U1, int vlen)
 
SprMathReal * spr_math_lu_symm_inv (SprMathReal *LU, int vlen, const int *pvt)
 
double spr_math_lu_symm_mul_vt_A_v (const SprMathReal *mA, const SprMathReal *vec, int vlen)
 
void * spr_math_matrix_alloc (size_t rows, size_t cols, size_t el_size)
 
void * spr_math_matrix_free (void *mat, int rows)
 
void * spr_math_matrix_read (char *fname, void *mat, int *rows, int *cols, const SprDT *dt)
 
int spr_math_matrix_write (const char *fname, const void *mat, const SprDT *dt, int Nrows, int Ncols)
 
void spr_rand_init (SprRandGen *restrict rg_state, const uint32_t *val, int Nval)
 
uint32_t spr_rand_u32 (SprRandGen *restrict rg_state)
 
uint64_t spr_rand_u64 (SprRandGen *restrict rg_state)
 
float spr_rand_f32 (SprRandGen *restrict rg_state)
 
double spr_rand_f64q (SprRandGen *restrict rg_state)
 Calculate the next 32 (pseudo) random bits and update the state of the random generator rg_state. More...
 
double spr_rand_f64 (SprRandGen *restrict rg_state)
 
double spr_randn_f64 (SprRandGen *restrict rg_state)
 
unsigned int spr_rand_modN (unsigned int n, SprRandGen *restrict rg_state)
 
float * spr_rand_vec_f32 (float x[], int n, float rmax, SprRandGen *restrict rg_state)
 
float * spr_randn_vec_f32 (float x[], int n, SprRandGen *restrict rg_state)
 
int spr_rand_m (int indx[], int n, int m, SprRandGen *restrict rg_state)
 
float spr_math_vector_dist2 (int n, const float *x, const float *y)
 
float spr_math_vec_dist (int n, const float *x, const float *y)
 
void spr_math_vec_pme (int n, float *x, float *y, float *z)
 
void spr_math_vec_init (int n, float *x, float a)
 Initializes the first n elements of x[] with a. More...
 
void spr_math_vec_add (int n, float *x, float *y, float *z)
 
int spr_math_vec_findmax (int n, const float *x)
 out: index of highest value More...
 
int spr_math_vec_findmin (int n, const float *x)
 out: index of lowest value More...
 

Variables

const SprFMIntA spr_fm_msk_abs
 
const SprFMIntA spr_fm_msk_sgn
 
const SprFMFltA spr_fm_vec_ones
 
const SprFMIntA spr_fm_vec_pos
 
const char spr_f32_fft_implementation []
 
const char spr_f64_fft_implementation []
 
SprRandGen spr_rand_default
 

Detailed Description

Macro Definition Documentation

#define SprFMFltU

unaligned float vector of optimal vector size

#define SprFMFltA

aligned float vector of optimal vector size

#define SprFMIntU

unaligned int vector of optimal vector size

#define SprFMIntA

aligned int vector of optimal vector size

#define SprFMUIntU

unaligned unsigned int vector of optimal vector size

#define SprFMUIntA

aligned unsigned int vector of optimal vector size

Typedef Documentation

typedef float SprFevalReal
typedef struct SprRandGen SprRandGen

Enumeration Type Documentation

anonymous enum
Enumerator
SPR_FEVAL_IS_VEC 
SPR_FEVAL_IS_STR 
SPR_FEVAL_IS_FUN1 
SPR_FEVAL_IS_IFX1 
SPR_FEVAL_IS_IFX2 
SPR_FEVAL_IS_FUN2 
SPR_FEVAL_IS_FUN3 
SPR_FEVAL_IS_FUN 
SPR_FEVAL_IS_UNDEF 
SPR_FEVAL_IS_SCRATCH 
SPR_FEVAL_IS_CONST 
SPR_FEVAL_IS_PROTECT 
anonymous enum
Enumerator
SPR_FEVAL_VAR_UNDEF 
SPR_FEVAL_VAR_VAL 
SPR_FEVAL_VAR_VEC 
SPR_FEVAL_VAR_STR 
SPR_FEVAL_VAR_FUN1 
SPR_FEVAL_VAR_FUN2 
SPR_FEVAL_VAR_FUN3 
anonymous enum

explanation: see feval

Enumerator
SPR_FEVAL_NO_VEC 
SPR_FEVAL_NO_IGN 
SPR_FEVAL_NO_END 
SPR_FEVAL_COMPILE 
anonymous enum
Enumerator
SPR_FEVAL_SUBST_VAR 
SPR_FEVAL_SUBST_ENV 
SPR_FEVAL_SUBST_ENV_BR 
SPR_FEVAL_SUBST_ENV_EXT 
SPR_FEVAL_SUBST_ENV_EXPR 
SPR_FEVAL_SUBST_FREE 
SPR_FEVAL_SUBST_DYNSTR 
SPR_FEVAL_SUBST_KEEP 
SPR_FEVAL_SUBST_KEEP_ERR 

Function Documentation

float spr_math_fast_erfinv ( float  x)

Calculate the inverse error function. The fast version is only correct up to 6 digits.

double spr_math_erfinv ( double  x)

Calculate the inverse error function. Uses fast_erfinv and performs a two steps of Newton-Raphson correction to achieve full accuracy.

double spr_math_normcdf_inv ( double  p)

Calculate the inverse of the cumulative density function of the normal distribution, i.e. produces the normal deviate z corresponding to a given lower tail area of p; z is accurate to about 16 digits. Based on: Algorithm AS241 Appl. Statist. (1988) Vol. 37, No. 3

double spr_math_normcdf_inv2 ( double  r)

Calculate normcdf_inv(1-r) for r in the range ]0.0,0.5] (handy to create normal distributed random numbers).

See Also
spr_math_normcdf_inv() for details
Note
this routine does not check if the input argument r falls into the correct range (and the range does not include the value 0.0)
SprFMFltA spr_fm_vx_max ( SprFMFltA  x,
SprFMFltA  y 
)

Element wise maximum.

Returns
max(x,y)
SprFMFltA spr_fm_vx_min ( SprFMFltA  x,
SprFMFltA  y 
)

Element wise minimum.

Returns
min(x,y)
SprFMFltA spr_fm_vx_exp ( SprFMFltA  x)

Element wise exponent. The computations are approximate: MSSE=8.67096e-08, Emax=3.49754e-07 (x=4.51795578).

Returns
exp(x)
float spr_fm_exp ( float  x)

Fast approximate exponent. For detail see spr_fm_vx_exp().

Returns
exp(x)
SprFMFltA spr_fm_vx_log ( SprFMFltA  p)

Element wise natural logarithm. The computations are approximate: MSSE=6.27598e-08, Emax=4.87045e-07 (p=1.15748596).

Returns
log(|p|+eps)
float spr_fm_log ( float  p)

Fast approximate natural logarithm. For detail see spr_fm_vx_log().

Returns
log(|x|+eps)
SprFMFltA spr_fm_vx_pow ( SprFMFltA  x,
SprFMFltA  y 
)

Element wise power function. The computations are approximate: MSSE=7.85975e-07, Emax=3.94553e-06 (x=y=12.4452667).

Returns
|x|^y
float spr_fm_pow ( float  x,
float  y 
)

Fast approximate power function. For detail see spr_fm_vx_pow().

Returns
|x|^y
SprFMFltA spr_fm_vx_logadd_ ( SprFMFltA  x,
SprFMFltA  y 
)

Element wise computation of log(exp(x)+exp(y)) = max(x,y) + log(1+exp(-abs(x-y))). The computations are approximate: MSSE=1.1100e-06, Emax=1.9408e-06 (abs(x-y)=4.0967).

Returns
log(exp(x)+exp(y))
float spr_fm_logadd_ ( float  x,
float  y 
)

Fast approximate computation of log(exp(x)+exp(y)) = max(x,y) + log(1+exp(-abs(x-y))). For detail see spr_fm_vx_logadd().

Returns
log(exp(x)+exp(y))
SprFMFltA spr_fm_vx_logadd ( SprFMFltA  x,
SprFMFltA  y 
)

Element wise computation of log(exp(x)+exp(y)) = max(x,y) + log(1+exp(-abs(x-y))). The computations are approximate: MSSE=3.03756e-07, Emax=1.97258e-06 (abs(x-y)=0.0529937744).

Returns
log(exp(x)+exp(y))
float spr_fm_logadd ( float  x,
float  y 
)

Fast approximate computation of log(exp(x)+exp(y)) = max(x,y) + log(1+exp(-abs(x-y))). For detail see spr_fm_vx_logadd().

Returns
log(exp(x)+exp(y))
float spr_fm_vx_sum ( SprFMFltA  x)
Returns
x[0]+x[1]+...
SprFMFltA spr_fm_vx_getN ( const float *  x,
int  N 
)

Get N consequetive floats (N in [1,SPR_FM_VEC_SZ-1]) and fill up the remaining elements with 0.0.

Returns
[x[0],...,c[N-1],0.0,...]
float* spr_fm_vec_exp ( float *restrict  vec,
int  N 
)

Compute vec[i]=exp(vec[i]), i=0...N-1. The exp is approximate. For detail see spr_fm_vx_exp().

float* spr_fm_vec_log ( float *restrict  vec,
int  N 
)

Compute vec[i]=log(|vec[i]|+eps), i=0...N-1. The log is approximate. For detail see spr_fm_vx_log().

float* spr_fm_vec_log_offs ( float *restrict  vec,
int  N,
float  offs 
)

Compute vec[i]=log(|vec[i]+offs|), i=0...N-1. The log is approximate. For detail see spr_fm_vx_log().

float* spr_fm_vec_log_mul ( float *restrict  vec,
int  N,
float  fac 
)

Compute vec[i]=log(|vec[i]|)*fac, i=0...N-1. The log is approximate. For detail see spr_fm_vx_log().

float* spr_fm_vec_log_offs_mul ( float *restrict  vec,
int  N,
float  offs,
float  fac 
)

Compute vec[i]=log(|vec[i]+offs|)*fac, i=0...N-1. The log is approximate. For detail see spr_fm_vx_log().

SprVarlist* spr_feval_var_init ( SprVarlist var_ptr,
const char *  name 
)

Initialize the variable var_ptr to an undefined (SPR_FEVAL_VAR_UNDEF) value. Both the variable itself (var_ptr) and the specified name name must be a dynamically allocated if the variable is going to be freeed. Otherwise, the the spr_feval_var_clear() routine should be used to free the content of the variable.

Returns
The variable itself.
Parameters
var_ptrthe variable to initialize
namethe name of the variable (1)
SprVarlist* spr_feval_var_clear ( SprVarlist var_ptr)

Set the variable var_ptr to an undefined (SPR_FEVAL_VAR_UNDEF) value.

Returns
The variable itself.
Parameters
var_ptrthe variable to clear
int spr_feval_result_set ( SprVarlist res_ptr,
int  type,
  ... 
)

Used to set the result in 'feval'-functions written in C (e.g. strcmp).

See Also
spr_feval_var_new() for the arguments.
Returns
(-1) on error, 0 otherwise.
Note
The variable must exist (i.e. no NULL pointer) and should be empty (the information is not deallocated).
Parameters
res_ptrthe variable to set
typethe type of the result
SprVarlist* spr_feval_var_new ( const char *  name,
int  type,
  ... 
)

Add an item to the list of known variables.
Allowed types:

SPR_FEVAL_VAR_UNDEF
just defines the variable, no value set
extra arguments: /
SPR_FEVAL_VAR_VAL
a single value
extra arguments: the value (real)
SPR_FEVAL_VAR_VEC
a vector
extra arguments:
the length of the vector (int)
the vector values (real*)
SPR_FEVAL_VAR_STR
a string
extra arguments: the string (char*)
SPR_FEVAL_VAR_FUN1
a simpel function, binary encoded extra arguments: address of the function (double (*fun)(double))
SPR_FEVAL_VAR_FUN2
a complex function, binary encoded
extra arguments: address of the function (int (*fun)(SprVarlist *res,SprVarlist *args))
SPR_FEVAL_VAR_FUN3
a complex function, textual description
Not implemented yet.

Allowed flags:

SPR_FEVAL_IS_CONST
variable is a constant
SPR_FEVAL_IS_PROTECT
don't allow assignements to the variable
Note
The name is allocated. String values are duplicated. Vectors and arrays must be allocated, and are owned by the variable from here one.
Parameters
namethe name of the new variable
typethe type of the result
SprVarlist* spr_feval_var_add ( SprVarlist varlist,
SprVarlist argptr 
)

Add an item to the list of known variables.

Parameters
varlistvariable list
argptrvariable to add to the list
void spr_feval_var_free ( SprVarlist var_ptr)

Free a single variable, the variable may not be part of a list.

Note
The vector values are deallocated also.
Parameters
var_ptrvariable to free
SprVarlist* spr_feval_var_rm ( SprVarlist varlist,
SprVarlist var_ptr 
)

Remove the variable var_ptr from the variable list varlist.

Returns
The update variable list.
Parameters
varlistlist containing all variables
var_ptrvariable to free
SprVarlist* spr_feval_varlist_free ( SprVarlist varlist,
SprVarlist upto 
)

Free variables up to a given previous point. To be able to limit a list to the current length later on, the current list pointer has to be remembered so it can be used in a free_vars() call. If upto is NULL, the entire list will be freed.

Note
Vector values are deallocated.
Parameters
varlistvariable list
uptoremove up to where
SprVarlist* spr_feval_var_get ( SprVarlist varlist,
const char *  name,
int  length 
)

Find an item in the list of known variables. If length is not -1, name is only compared up to length characters. The last added variable (see add_var) that matches is used returned. If no matching variable is not found, a NULL is returned

Parameters
varlistvariable list
namename of the requested variable
lengthlength of the name string
int spr_feval_var_set ( SprVarlist var_ptr,
int  type,
  ... 
)

Set a variable to a given value.

See Also
spr_feval_var_new() and spr_feval_var_get().
Parameters
var_ptrvariable to set
typethe new type of the variable
SprVarlist* spr_feval_var_assign ( SprVarlist varlist,
const char *  name,
int  type,
  ... 
)

Set a variable in a variable list to a given value.

Returns
If the variable is found, its adress is returned, else NULL is returned.
See Also
spr_feval_var_new() and spr_feval_var_get().
Parameters
varlistvariable list
namename of the requested variable
typethe new type of the variable
void spr_feval_whos ( SprStream fd,
SprVarlist varlist,
const char *  name 
)

If name equals NULL, show all variables and the amount of memory they use. If name is not NULL, then only show that variable.

Parameters
fddestination file
varlistvariable list
void spr_feval_var_print ( SprStream fd,
const SprVarlist var_ptr,
const char *  rm 
)

Print a variable.

Parameters
fddestination file
var_ptrvariable to print
rmtext to print as right margin
char* spr_feval ( const char *  eval_str,
const char *  eval_end,
SprVarlist varlist,
SprVarlist res,
int  flags 
)

Evaluates the floating point (vector) expression given in eval_str. If eval_end is not NULL, the evaluation will end at eval_end. Local variables are specified in the variable list varlist. Global variables are specified in the global variable spr_feval_global_var. The result of the expression evaluation is returned in res. If res is a NULL pointer, the result is not returned, so the expression should contain assignements to local variables for the function call to have any effect.

For the flags argument, the following flag can be combined:

SPR_FEVAL_NO_VEC
Vectors are not allowed as result.
SPR_FEVAL_NO_END
Evaluate up to the first non allowable character, but gives no error:
e.g. spr_feval("(expr) then,...") will stop at the 't' of 'then'.
SPR_FEVAL_COMPILE
Compile instead of evaluating (NYI).
Returns
The pointer up to where string was processed. If an error occurs, a (NULL) is returned.
Parameters
eval_strstring containing the expression
eval_endstop evaluating at this position
varlistlist with local variables
reslocation to store result
flagsadditional flags
SprFevalReal spr_feval_num ( const char *  eval_str,
const char *  eval_end,
SprVarlist varlist,
int *  err_ret 
)

Evaluate the scalar expression eval_str. If eval_end is not NULL, the evaluation will end at eval_end. The result is guaranteed to be a scalar. The err_ret pointer may be a NULL pointer.

See Also
spr_feval() for more info.
Parameters
eval_strstring to evaluate
eval_endstop evaluating at this position
varlistlist with local variables
err_retflag to set on error
char* spr_feval_var_substitute ( const char *  str,
SprVarlist list,
int  flags 
)

Substitute variable names and environment variables ($name, or ${name}) found in the string str with their value. If the SPR_FEVAL_SUBST_DYNSTR bit flags is set in flags, both input and output string must/will be of the type SprDynStr. If the SPR_FEVAL_SUBST_VAR is set in flags, the variables listed in list are substituted wherever they occur in str. If the SPR_FEVAL_SUBST_ENV bit is set in flags, environment variables (starting with a '$') are substituted as well. If the SPR_FEVAL_SUBST_ENV_BR bit is set in flags, environment variables must be indicated with a ${name} construct. If the SPR_FEVAL_SUBST_ENV_EXT bit is set in flags, a small number of substitution expressions (see below) is supported. If the SPR_FEVAL_SUBST_ENV_EXPR bit is set in flags, the $(some_expression) construct to evaluate an expression (using constants and variables from list) is supported. Unknown variable names / environment variables are copied to the output string unchanged. The original string is dealocated if the SPR_FEVAL_SUBST_FREE bit flags is set in flags (except when SPR_FEVAL_SUBST_KEEP or SPR_FEVAL_SUBST_KEEP_ERR is set and the original string is returned). If the SPR_FEVAL_SUBST_KEEP bit flags is set in flags, the input string pointer will be returned if no variable substitutions were performed. If the SPR_FEVAL_SUBST_KEEP_ERR bit flags is set in flags, the input string pointer will be returned if an error occurred.

Returns
The string with all known names substituted.
Note
the following env. var. substition constructs are supported (behavior as in bash)
  • $var, ${var}
  • ${#var}, ${#var[@]}, ${#var[*]}
  • ${var[ndx]}
      with ndx a numerical expression using global variables (but no env. vars.)
  • ${var:<pos0>}, ${var:pos0:len}
    &nbsp;&nbsp;with pos0 and len numerical expressions using global variables (but no env. vars.)
  • ${var-alt_txt}, ${var:-alt_txt}
  • ${var::str}, ${var##str}, ${var#*str}, ${var##*str}
  • ${varstr}, ${var%%str}, ${varstr*}, ${var%%str*}
Note
the $() construct
  • uses floating point evaluation
  • knows vector expressions
  • will be printed in list format (a b c ...) by default; use $([]) to print vectors in vector format ([a,b,c,...])
Parameters
strinput string
listvariable list
flagssubstitute env. variables also
float* spr_f32_fht1 ( float *restrict  x,
int  N 
)

In-situ fast Hartely transform of a vector x of size N. The size N must be a power of 2 and must be larger than or equal to 4.

Returns
the vector x (in-situ tranformed) on success or NULL or error (invalid size).
float* spr_f32_fht2 ( float *restrict  x,
int  N 
)

Dual in-situ fast Hartely transform of a N*2 matrix x. The elements in x are stored in an interleaving fashion, i.e. x[0][0], x[0][1], x[1][0], x[1][1], x[2][0], ..., x[N-1][0], x[N-1][1]. The size N must be a power of 2 and must be larger than or equal to 4.

Returns
the matrix x (in-situ tranformed) on success or NULL or error (invalid size).
float* spr_f32_fht4 ( float *restrict  x,
int  N 
)

Quad in-situ fast Hartely transform of a N*4 matrix x. The elements in x are stored in an interleaving fashion, i.e. x[0][0...3], x[1][0...3], ..., x[N-1][0...3]. The size N must be a power of 2 and must be larger than or equal to 4.

Returns
the matrix x (in-situ tranformed) on success or NULL or error (invalid size).
float* spr_f32_fft_c2c ( float *restrict  x,
int  N 
)

In-situ FFT of a complex vector x. N, the number of (complex) datapoints, must be a power of 2 and must be larger than or equal to 4.

Returns
x on success and NULL or error (invalid size).
Note
this routine is not all that optimized since rarely needed.
float* spr_f32_ifft_c2c ( float *restrict  x,
int  N 
)

In-situ inverse FFT of a complex vector x. N, the number of (complex) datapoints, must be a power of 2 and must be larger than or equal to 4.

Returns
x on success and NULL or error (invalid size).
Note
this routine is not all that optimized since rarely needed.
float* spr_f32_fft_r2c ( float *restrict  x,
int  N 
)

In-situ FFT of a real real vector x. The resulting vector contains complex values. N, the number of (complex) datapoints, must be a power of 2 and must be larger than or equal to 8.

Note
Since N real points are converted to N/2+1 complex values, the vector x must have space for 2 extra elements (the values at the Nyquist frequency).
Returns
x on success and NULL or error (invalid size).
float* spr_f32_ifft_c2r ( float *restrict  x,
int  N 
)

In-situ inverse FFT of a complex symmetric vector x. The resulting vector contains real values. N, the resulting number of (real) datapoints, must be a power of 2 and must be larger than or equal to 8.

Note
The input vector x must contain N/2+1 complex values and hence must have space for N+2 elements.
Returns
x on success and NULL or error (invalid size).
double* spr_f64_fht1 ( double *restrict  x,
int  N 
)

In-situ fast Hartely transform of a vector x of size N. The size N must be a power of 2 and must be larger than or equal to 4.

Returns
the vector x (in-situ tranformed) on success or NULL or error (invalid size).
double* spr_f64_fht2 ( double *restrict  x,
int  N 
)

Dual in-situ fast Hartely transform of a N*2 matrix x. The elements in x are stored in an interleaving fashion, i.e. x[0][0], x[0][1], x[1][0], x[1][1], x[2][0], ..., x[N-1][0], x[N-1][1]. The size N must be a power of 2 and must be larger than or equal to 4.

Returns
the matrix x (in-situ tranformed) on success or NULL or error (invalid size).
double* spr_f64_fht4 ( double *restrict  x,
int  N 
)

Quad in-situ fast Hartely transform of a N*4 matrix x. The elements in x are stored in an interleaving fashion, i.e. x[0][0...3], x[1][0...3], ..., x[N-1][0...3]. The size N must be a power of 2 and must be larger than or equal to 4.

Returns
the matrix x (in-situ tranformed) on success or NULL or error (invalid size).
double* spr_f64_fft_c2c ( double *restrict  x,
int  N 
)

In-situ FFT of a complex vector x. N, the number of (complex) datapoints, must be a power of 2 and must be larger than or equal to 4.

Returns
x on success and NULL or error (invalid size).
Note
this routine is not all that optimized since rarely needed.
double* spr_f64_ifft_c2c ( double *restrict  x,
int  N 
)

In-situ inverse FFT of a complex vector x. N, the number of (complex) datapoints, must be a power of 2 and must be larger than or equal to 4.

Returns
x on success and NULL or error (invalid size).
Note
this routine is not all that optimized since rarely needed.
double* spr_f64_fft_r2c ( double *restrict  x,
int  N 
)

In-situ FFT of a real real vector x. The resulting vector contains complex values. N, the number of (complex) datapoints, must be a power of 2 and must be larger than or equal to 8.

Note
Since N real points are converted to N/2+1 complex values, the vector x must have space for 2 extra elements (the values at the Nyquist frequency).
Returns
x on success and NULL or error (invalid size).
double* spr_f64_ifft_c2r ( double *restrict  x,
int  N 
)

In-situ inverse FFT of a complex symmetric vector x. The resulting vector contains real values. N, the resulting number of (real) datapoints, must be a power of 2 and must be larger than or equal to 8.

Note
The input vector x must contain N/2+1 complex values and hence must have space for N+2 elements.
Returns
x on success and NULL or error (invalid size).
double** spr_lda_stats_alloc ( double **  moments,
int  nclasses,
int  vlen 
)

Allocate and/or initialize the array that will contain the within class statistics.

double** spr_lda_stats_diag_alloc ( double **  moments,
int  nclasses,
int  vlen 
)

Allocate and/or initialize the array that will contain the within class statistics.

void spr_lda_stats_update ( double *restrict  moments,
const SprMathReal *  data,
int  vlen,
float  weight 
)

Update the 2nd order statistics for a given vector (and class).

void spr_lda_stats_diag_update ( double *restrict  moments,
const SprMathReal *  data,
int  vlen,
float  weight 
)

Update the 2nd order statistics for a given vector (and class). Only the diagonal elements of the covariance matrix are updated.

void spr_lda_calc_class_covar ( double **  moments,
int  vlen,
int  nclasses,
double  Beta,
SprStatsStats within_class_covar,
SprStatsStats between_class_covar 
)

Calculate the within class covariance and the between class covariance.

void spr_lda_calc_transf ( SprStatsStats lda_transf,
int  vlen,
SprStatsStats within_class_covar,
SprStatsStats between_class_covar 
)

Calculate the LDA-transform by solving the (symmetric) generalized eigenvlue problem. The eigenvalues are sorted from largest to smallest. Both the within_class_covar and between_class_covar matrices are used as scratch matrices. Their values are invalid after calling this routine. "The resulting eigenvalues are stored in" lda_transf->mean "The resulting tranformation matrix is stored in" lda_transf->covar

int spr_math_mat_print ( SprStream fd,
const void *  A,
int  nrows,
int  ncols,
const SprDT dt,
const char *  info 
)

Print matrix A of size nrowsxncols type dt to the stream fd. The info string is printed before the matrix (if not NULL).

Returns
0 on success and -1 on error.
Parameters
fddestination file
Amatrix to print
nrowsnumber of rows
ncolsnumber of columns
dtdata type
infointro string
float* spr_mat_f32_A_plus_scaleB_insitu ( float *restrict  mA,
int  size,
const float *  mB,
float  facB 
)

Add a scaled version of a second matrix/vector to a given matrix/vector.

mA += mB*facB
Returns
the matrix/vector mA.
float* spr_mat_f32_scaleA_plus_B ( float *restrict  mC,
int  size,
const float *  mA,
const float *  mB,
float  facA 
)

Calculates the sum of a vector/matrix with a scaled vector/matrix:

mC = mA*facA + mB

with mA, mB and mC of size size (size = nrow x ncol).

Returns
the destination matrix mC.
float* spr_mat_f32_scaleA_plus_scaleB ( float *restrict  mC,
int  size,
const float *  mA,
const float *  mB,
float  facA,
float  facB 
)

Calculates the sum of two scaled a vectors/matrices:

mC = mA*facA + mB*facB

with mA, mB and mC of size size (size = nrow x ncol).

Returns
the destination matrix mC.
float* spr_mat_f32_scaleA_insitu ( float *restrict  mA,
float  facA,
int  size 
)

Scale a matrix/vector mA of size size (multiply all elements with a constant).

mA *= facA
Returns
the vector mA.
float* spr_mat_f32_scaleA ( float *restrict  mC,
int  size,
const float *  mA,
float  facA 
)

Scale the elements in a vector/matrix:

mC = mA*facA

with mA and mC of size size (size = nrow x ncol).

Returns
the destination matrix mC.
float* spr_mat_f32_A_plus_B ( float *restrict  mC,
int  size,
const float *  mA,
const float *  mB 
)

Calculate the sum of two vectors/matrices:

mC = mA + mB

with mA, mB and mC of size size (size = nrow x ncol).

Returns
the destination matrix mC.
float* spr_mat_f32_A_plus_B_insitu ( float *restrict  mA,
int  size,
const float *  mB 
)

Calculate the sum of two vectors/matrices:

mA += mB

with mA and mB of size size (size = nrow x ncol).

Returns
the destination matrix mA.
float* spr_mat_f32_A_minus_B ( float *restrict  mC,
int  size,
const float *  mA,
const float *  mB 
)

Calculate the difference of two vectors/matrices:

mC = mA - mB

with mA, mB and mC of size size (size = nrow x ncol).

Returns
the destination matrix mC.
float* spr_mat_f32_A_minus_B_insitu ( float *restrict  mA,
int  size,
const float *  mB 
)

Calculate the difference of two vectors/matrices:

mA -= mB

with mA and mB of size size (size = nrow x ncol).

Returns
the destination matrix mA.
float* spr_mat_f32_A_dotmul_B ( float *restrict  mC,
int  size,
const float *  mA,
const float *  mB 
)

Calculate the point-wise multiplication of two vectors/matrices:

mC = mA .* mB

with mA, mB and mC of size size (size = nrow x ncol).

Returns
the destination matrix mC.
float* spr_mat_f32_A_dotdiv_B ( float *restrict  mC,
int  size,
const float *  mA,
const float *  mB 
)

Calculate the point-wise division of two vectors/matrices:

mC = mA ./ mB

with mA, mB and mC of size size (size = nrow x ncol).

Returns
the destination matrix mC.
double spr_mat_f32_A_inprod_B ( const float *  mA,
int  size,
const float *  mB 
)

Calculate the in-product of two vectors/matrices:

sum(mA.*mB)
Returns
the in-product.
double spr_mat_f32_rowA_inprod_colB ( const float *  mA,
int  size,
const float *  mB,
int  ncolB 
)

Calculate the in-product of a vector (or first row of a matrix) mA and a column of a second matrix mB:

sum(mA.*mB[:,0])

Another column in the matrix mB or row in matrix mA can be selected by adding the index of the column/row to the matrix pointer(s):

spr_mat_f32_rowA_inprod_colB(mA+irow*size,size,mB+icol,ncolB);
Returns
the in-product.
float* spr_mat_f32_A_mul_B ( float *restrict  mC,
const float *  mA,
int  nrowA,
int  ncolA,
const float *  mB,
int  ncolB 
)

Calculate the product of two matrices:

mC = mA * mB

with mA of size nrowA*ncolA, mB of size ncolA*ncolB and the resulting matrix mC of size nrowA*ncolB.

Returns
the destination matrix mC.
float* spr_mat_f32_At_mul_B ( float *restrict  mC,
const float *  mA,
int  nrowA,
int  ncolA,
const float *  mB,
int  ncolB 
)

Calculate the product of two matrices, the first matrix being transposed:

mC = mA' * mB

with mA of size nrowA*ncolA, mB of size nrowA*ncolB and the resulting matrix mC of size ncolA*ncolB.

Returns
the destination matrix mC.
float* spr_mat_f32_A_mul_Bt ( float *restrict  mC,
const float *  mA,
int  nrowA,
int  ncolA,
const float *  mB,
int  nrowB 
)

Calculate the product of two matrices, the second matrix being transposed:

mC = mA * mB'

with mA of size nrowA*ncolA, mB of size nrowB*ncolA and the resulting matrix mC of size nrowA*nrowB.

Returns
the destination matrix mC.
float* spr_mat_f32_At_mul_Bt ( float *restrict  mC,
const float *  mA,
int  nrowA,
int  ncolA,
const float *  mB,
int  nrowB 
)

Calculate the product of two matrices, both matrices being transposed:

mC = mA' * mB'

with mA of size nrowA*ncolA, mB of size nrowB*nrowA and the resulting matrix mC of size ncolA*nrowB.

Returns
the destination matrix mC.
float* spr_mat_f32_colA ( float *restrict  mC,
const float *  mA,
int  nrowA,
int  ncolA 
)

Extract the first column from a matrix and copy it to a vector.

mC = mA[:,0]

If another column has to be extracted, the routine could be called as follows:

spr_mat_f32_colA(mC,mA+icol,nrowA,ncolA);
Returns
the destination matrix mC.
float* spr_mat_f32_At ( float *restrict  mC,
const float *  mA,
int  nrowA,
int  ncolA 
)

Transpose a matrix:

mC = mA'
Note
This routine does not work in situ!
Returns
the destination matrix mC.
double spr_mat_f32_detA ( float *restrict  mA,
float *restrict  copy,
int  size 
)

Calculate the determinant of a square matrix mA of size size*size via pivoting and column elimination.

Note
The routine destroys the contents of the matrix mA. If the copy argument is not NULL, a temporary copy is made and this copy is then used, hence preserving the content of mA.
Returns
the determinant.
double spr_mat_f32_norm2A ( const float *  mA,
int  size 
)

Calculate the Frobenius norm of a vector/matrix mA of size size:

sqrt(sum(mA.^2))
Returns
the Frobenius norm.
double spr_mat_f32_prodA ( const float *  mA,
int  size 
)

Calculate the product of the elements in a vector/matrix mA of size size:

prod(mA)
Returns
the product of the elements.
double* spr_mat_f64_A_plus_scaleB_insitu ( double *restrict  mA,
int  size,
const double *  mB,
double  facB 
)

Add a scaled version of a second matrix/vector to a given matrix/vector.

mA += mB*facB
Returns
the matrix/vector mA.
double* spr_mat_f64_scaleA_plus_B ( double *restrict  mC,
int  size,
const double *  mA,
const double *  mB,
double  facA 
)

Calculates the sum of a vector/matrix with a scaled vector/matrix:

mC = mA*facA + mB

with mA, mB and mC of size size (size = nrow x ncol).

Returns
the destination matrix mC.
double* spr_mat_f64_scaleA_plus_scaleB ( double *restrict  mC,
int  size,
const double *  mA,
const double *  mB,
double  facA,
double  facB 
)

Calculates the sum of two scaled a vectors/matrices:

mC = mA*facA + mB*facB

with mA, mB and mC of size size (size = nrow x ncol).

Returns
the destination matrix mC.
double* spr_mat_f64_scaleA_insitu ( double *restrict  mA,
double  facA,
int  size 
)

Scale a matrix/vector mA of size size (multiply all elements with a constant).

mA *= facA
Returns
the vector mA.
double* spr_mat_f64_scaleA ( double *restrict  mC,
int  size,
const double *  mA,
double  facA 
)

Scale the elements in a vector/matrix:

mC = mA*facA

with mA and mC of size size (size = nrow x ncol).

Returns
the destination matrix mC.
double* spr_mat_f64_A_plus_B ( double *restrict  mC,
int  size,
const double *  mA,
const double *  mB 
)

Calculate the sum of two vectors/matrices:

mC = mA + mB

with mA, mB and mC of size size (size = nrow x ncol).

Returns
the destination matrix mC.
double* spr_mat_f64_A_plus_B_insitu ( double *restrict  mA,
int  size,
const double *  mB 
)

Calculate the sum of two vectors/matrices:

mA += mB

with mA and mB of size size (size = nrow x ncol).

Returns
the destination matrix mA.
double* spr_mat_f64_A_minus_B ( double *restrict  mC,
int  size,
const double *  mA,
const double *  mB 
)

Calculate the difference of two vectors/matrices:

mC = mA - mB

with mA, mB and mC of size size (size = nrow x ncol).

Returns
the destination matrix mC.
double* spr_mat_f64_A_minus_B_insitu ( double *restrict  mA,
int  size,
const double *  mB 
)

Calculate the difference of two vectors/matrices:

mA -= mB

with mA and mB of size size (size = nrow x ncol).

Returns
the destination matrix mA.
double* spr_mat_f64_A_dotmul_B ( double *restrict  mC,
int  size,
const double *  mA,
const double *  mB 
)

Calculate the point-wise multiplication of two vectors/matrices:

mC = mA .* mB

with mA, mB and mC of size size (size = nrow x ncol).

Returns
the destination matrix mC.
double* spr_mat_f64_A_dotdiv_B ( double *restrict  mC,
int  size,
const double *  mA,
const double *  mB 
)

Calculate the point-wise division of two vectors/matrices:

mC = mA ./ mB

with mA, mB and mC of size size (size = nrow x ncol).

Returns
the destination matrix mC.
double spr_mat_f64_A_inprod_B ( const double *  mA,
int  size,
const double *  mB 
)

Calculate the in-product of two vectors/matrices:

sum(mA.*mB)
Returns
the in-product.
double spr_mat_f64_rowA_inprod_colB ( const double *  mA,
int  size,
const double *  mB,
int  ncolB 
)

Calculate the in-product of a vector (or first row of a matrix) mA and a column of a second matrix mB:

sum(mA.*mB[:,0])

Another column in the matrix mB or row in matrix mA can be selected by adding the index of the column/row to the matrix pointer(s):

spr_mat_f64_rowA_inprod_colB(mA+irow*size,size,mB+icol,ncolB);
Returns
the in-product.
double* spr_mat_f64_A_mul_B ( double *restrict  mC,
const double *  mA,
int  nrowA,
int  ncolA,
const double *  mB,
int  ncolB 
)

Calculate the product of two matrices:

mC = mA * mB

with mA of size nrowA*ncolA, mB of size ncolA*ncolB and the resulting matrix mC of size nrowA*ncolB.

Returns
the destination matrix mC.
double* spr_mat_f64_At_mul_B ( double *restrict  mC,
const double *  mA,
int  nrowA,
int  ncolA,
const double *  mB,
int  ncolB 
)

Calculate the product of two matrices, the first matrix being transposed:

mC = mA' * mB

with mA of size nrowA*ncolA, mB of size nrowA*ncolB and the resulting matrix mC of size ncolA*ncolB.

Returns
the destination matrix mC.
double* spr_mat_f64_A_mul_Bt ( double *restrict  mC,
const double *  mA,
int  nrowA,
int  ncolA,
const double *  mB,
int  nrowB 
)

Calculate the product of two matrices, the second matrix being transposed:

mC = mA * mB'

with mA of size nrowA*ncolA, mB of size nrowB*ncolA and the resulting matrix mC of size nrowA*nrowB.

Returns
the destination matrix mC.
double* spr_mat_f64_At_mul_Bt ( double *restrict  mC,
const double *  mA,
int  nrowA,
int  ncolA,
const double *  mB,
int  nrowB 
)

Calculate the product of two matrices, both matrices being transposed:

mC = mA' * mB'

with mA of size nrowA*ncolA, mB of size nrowB*nrowA and the resulting matrix mC of size ncolA*nrowB.

Returns
the destination matrix mC.
double* spr_mat_f64_colA ( double *restrict  mC,
const double *  mA,
int  nrowA,
int  ncolA 
)

Extract the first column from a matrix and copy it to a vector.

mC = mA[:,0]

If another column has to be extracted, the routine could be called as follows:

spr_mat_f64_colA(mC,mA+icol,nrowA,ncolA);
Returns
the destination matrix mC.
double* spr_mat_f64_At ( double *restrict  mC,
const double *  mA,
int  nrowA,
int  ncolA 
)

Transpose a matrix:

mC = mA'
Note
This routine does not work in situ!
Returns
the destination matrix mC.
double spr_mat_f64_detA ( double *restrict  mA,
double *restrict  copy,
int  size 
)

Calculate the determinant of a square matrix mA of size size*size via pivoting and column elimination.

Note
The routine destroys the contents of the matrix mA. If the copy argument is not NULL, a temporary copy is made and this copy is then used, hence preserving the content of mA.
Returns
the determinant.
double spr_mat_f64_norm2A ( const double *  mA,
int  size 
)

Calculate the Frobenius norm of a vector/matrix mA of size size:

sqrt(sum(mA.^2))
Returns
the Frobenius norm.
double spr_mat_f64_prodA ( const double *  mA,
int  size 
)

Calculate the product of the elements in a vector/matrix mA of size size:

prod(mA)
Returns
the product of the elements.
SprMathReal* spr_math_omat_to_nmat ( int  nrow,
int  ncol,
const SprMathReal *const *  from,
SprMathReal *  to 
)

Convert from the 'old' matrix format (array of pointers) to the new matrix format (one continuous block, row0, row1, ... order). If to is a NULL pointer, a new one will be allocated.

SprMathReal** spr_math_nmat_to_omat ( int  nrow,
int  ncol,
const SprMathReal *  from,
SprMathReal **  to 
)

Convert from the 'new' row0, row1, ... order matrix format to the old matrix format (array of pointers). If to is a NULL pointer, a new one will be allocated.

void spr_math_fmat_to_fmat_trunc ( const SprMathReal *  from,
int  from_ncol,
SprMathReal *restrict  to,
int  to_nrow,
int  to_ncol 
)

Select the left uppercorner sub-matrix from a large matrix. A different sub-matrix can be selected by adding an offset to the from matrix, for example: fmat_to_fmat_trunc(from+irow*from_ncol+icol,...). This operation cannot be done in situ.

void spr_math_utmat_to_fmat ( int  vlen,
const SprMathReal *  from,
SprMathReal *  to 
)

Convert from an upper triangular matrix format to a full matrix format (row0, row1, ... order). This operation can be done in situ.

void spr_math_utmat_to_fmat_trunc ( int  vlen,
const SprMathReal *  from,
int  nlen,
SprMathReal *restrict  to 
)

Convert from an upper triangular matrix format to a full matrix format (row0, row1, ... order). Only the left upper corner nlen elements are filled in. This operation cannot be done in situ.

void spr_math_utmat_to_utmat_trunc ( int  vlen,
const SprMathReal *  from,
int  nlen,
SprMathReal *restrict  to 
)

Convert an upper triangular matrix to a new (smaller) upper triangular matrix. Only the left upper corner nlen elements are copied in. This operation cannot be done in situ.

void spr_math_fmat_to_utmat ( int  vlen,
const SprMathReal *  from,
SprMathReal *  to 
)

Convert from a full matrix format (row0, row1, ... order) to an upper triangular matrix format. This operation can be done in situ.

SprMathReal* spr_math_utmat_to_diag ( int  vlen,
const SprMathReal *  from,
SprMathReal *restrict  to 
)

Extract the diagonal elements from an upper triangular matrix in compact format. This operation can be done in situ.

SprMathReal* spr_math_fmat_to_diag ( int  vlen,
const SprMathReal *  from,
SprMathReal *restrict  to 
)

Extract the diagonal elements from a full matrix format. This operation can be done in situ.

void spr_math_v_lin_comb ( SprMathReal *  vec1,
int  vlen,
SprMathReal *  vec2,
double  c,
double  s 
)

Make an orthogonal decomposition of two vectors:

vec1_out = c*vec1+s*vec2;
vec2_out = c*vec2-s*vec1;
SprMathReal* spr_math_hess_fact_t ( SprMathReal *  A,
int  vlen,
SprMathReal *  Q 
)

Compute Hessenberg factorisation in compact form and, if requested, construct the Hessenberg orthogonalising matrix Q', i.e. Hess(A') = Q'.A'.Q Note: the factorisation is performed in situ.

Returns
The factorized matrix on success or (NULL) on faillure.
SprMathReal* spr_math_diag_eig_t ( SprMathReal *  diag0,
int  vlen,
SprMathReal *  diag1,
SprMathReal *  eigvec 
)

Finds eigenvalues of symmetric tridiagonal matrices. The matrix is represented by its diag entries and sub- & super-diag entries. The eigenvalues are returned in diag0. The eigenvectors are returned in eigvec, each ROW contains one eigenvector.

SprMathReal* spr_math_symm_eig_t ( SprMathReal *  A,
int  vlen,
SprMathReal *  eigvec,
SprMathReal *  eigval 
)

Computes the eigenvalues and eigenvectors of a dense symmetric matrix. The eigenvalues are stored in eigval. The eigenvector are stored in eigvec, each ROW contains one eigenvector.

Returns
The eigenvalues on success, and a NULL on faillure
Note
A MUST be symmetric on entry.
SprMathReal* spr_math_scale_eig_vec_t ( int  vlen,
SprMathReal *  eig_vec,
SprMathReal *  eig_val 
)
int spr_math_io_nmat_write ( const char *  fname,
const void *  mat,
int  Nrows,
int  Ncols,
const SprDT dt 
)

Write the matrix mat of size NrowsxNcols and type dt to stream fname.

Returns
0 on success and (-1) on failure.
void* spr_math_io_nmat_read ( const char *  fname,
void *  mat,
int *  rows,
int *  cols,
const SprDT dt 
)

Read a matrix of size (*rows)x(*cols) and type dt from the file fname into the matrix mat. If *rows and/or *cols not (-1), they will be compared with the values specified in the file, a value of (-1) avoids the test and fills in the correct values. If mat is a NULL pointer, a matrix of the correct size will be allocated.

Returns
the matrix read from file or NULL on failure
SprMathReal* spr_math_lu_full ( SprMathReal *  mA,
int  nrow,
int  ncol,
int *  pivot 
)

Calculates the in-situ LU decomposition of a matrix. The routine assumes that nrow <= ncol. If the pivot argument is not NULL, optimal row pivoting is done and the permutation is stored in the array addressed by pivot. A special value 'LUPivot' can be used to request optimal pivoting without storing the permutation.
The routine returns the matrix mA. The upper triangular part (including the diagonal) of mA contains the matrix U, the lower triangular part contains the matrix L (the diagonal elements of L equal to 1.0).

double spr_math_lu_full_det ( const SprMathReal *  LU,
int  vlen,
const int *  pivot 
)

Calculated the determinant of a matrix. Note that the LU-decomposition should be provided, and not the original matrix. If the pivot argument is not given (LUPivot), then the absolute value of the determinant is determined.

SprMathReal* spr_math_lu_l1_inv ( SprMathReal *  res,
SprMathReal *  L1,
int  vlen 
)

Calculates the inverse of a lower triagular matrix with ones on the diagonal. The ones on the diagonal and the zeros in the upper triangulart part of the result matrix are not set. This operation can be done in-situ.

SprMathReal* spr_math_lu_l1_complete ( SprMathReal *restrict  L1,
int  vlen 
)

Set the diagonal elements to 1.0 and the upper diagonal elements to 0.0. Returns the matrix L1.

SprMathReal* spr_math_lu_l1_solve ( SprMathReal *  L1,
SprMathReal *  B,
int  nrow,
int  ncol 
)

Calculate the solution of 'L1*X = B' in-situ (in B) with L1 a lower triangular matrix with ones on the diagonal. Returns the matrix B, which now contains the solution.

SprMathReal* spr_math_lu_u_solve ( SprMathReal *  U,
SprMathReal *  B,
int  nrow,
int  ncol 
)

Calculate the solution of 'U*X = B' in-situ (in B), with U an upper triangular matrix. Returns the matrix B, which now contains the solution.

SprMathReal* spr_math_lu_full_column_pivot ( SprMathReal *  A,
int  nrow,
int  ncol,
int *  pivot 
)

Pivot columns.

SprMathReal* spr_math_lu_full_column_pivot_r ( SprMathReal *  A,
int  nrow,
int  ncol,
int *  pivot 
)

Pivot columns (reverse order).

SprMathReal* spr_math_lu_full_row_pivot ( SprMathReal *  A,
int  nrow,
int  ncol,
int *  pivot 
)

Pivot rows.

SprMathReal* spr_math_lu_full_row_pivot_r ( SprMathReal *  A,
int  nrow,
int  ncol,
int *  pivot 
)

Pivot rows (reverse order).

SprMathReal* spr_math_lu_full_inv ( SprMathReal *  res,
SprMathReal *  LU,
int  vlen,
int *  pivot 
)

Calculated the inverse of a matrix. Note that the LU-decomposition should be provided, and not the original matrix. If the pivot argument is not given (LUPivot), then the inverse will not be determined on a column permutation. This routine cannot be done in-situ. Returns the inverse.

SprMathReal* spr_math_lu_full_solve ( SprMathReal *  LU,
SprMathReal *  B,
int  nrow,
int  ncol,
int *  pivot 
)

Calculate the solution of 'L*U*X = B' in-situ (in B) with L a lower triangular matrix of size <nrow x nrow> with ones on the diagonal, U an upper triangular matrix of size <nrow x nrow>, and B a full matrix of size <nrow x ncol>. Returns the matrix B, which now contains the solution.

SprMathReal* spr_math_lu_symm ( SprMathReal *restrict  mA,
int  vlen,
int *restrict  pvt 
)

In-situ LU-decomposition of a symmetrix matrix stored in an upper triangular format. The LU-decompostion of a symmetrix matrix < A = L*D*L' > results in a lower diagonal matrix L with ones on the diagonal, and a diagonal matrix D. Both matrices are combined in one upper triangular matrix (transpose of L) and stored on the same location as the original matrix. If the pvt argument is not NULL, a good permutation of rows and columns is performed and stored in the pvt array. A special value 'LUPivot' is also allowed to indicate that permutation is allowed, but should not be stored. Returns the pointer to the LU-decomposition.

double spr_math_lu_symm_det ( SprMathReal *  LU,
int  vlen 
)

Calculate the determinant of a symmetric matrix. Note that the symmetric LU-decomposition should be provided, and not the original matrix.

SprMathReal* spr_math_lu_symm_pivot ( SprMathReal *  A,
int  vlen,
const int *  pivot 
)

Pivot rows and columns.

SprMathReal* spr_math_lu_symm_pivot_r ( SprMathReal *  A,
int  vlen,
const int *  pivot 
)

Pivot rows and columns (reverse order).

SprMathReal* spr_math_lu_symm_u1_inv ( SprMathReal *restrict  U1,
int  vlen 
)

In-situ calculation of the inverse of an upper triangular matrix in compact format and with ones on the diagonal. The diagonal elements are left unchanged.

void spr_math_lu_symm_u1_scale ( SprMathReal *  U1,
int  vlen 
)

Calculate < U1 * |D|^(-1/2) >, with U1 an upper triangular matrix in compact format and with ones on the diagonal, and D a diagonal matrix (stored on the diagonal of the input matrix U1).

SprMathReal* spr_math_lu_symm_inv ( SprMathReal *  LU,
int  vlen,
const int *  pvt 
)

In-situ calculation of the inverse of a positive deterministic symmetric matrix. Note that the symmetric LU-decomposition should be provided, and not the original matrix. If permutations were allowed during the LU-decomposition, the resulting permutation indexes should be provided also.

double spr_math_lu_symm_mul_vt_A_v ( const SprMathReal *  mA,
const SprMathReal *  vec,
int  vlen 
)

Calculate and return vec'*mA*vec with mA a symmetric matrix.

Returns
vec'*mA*vec
void* spr_math_matrix_alloc ( size_t  rows,
size_t  cols,
size_t  el_size 
)

Allocate a matrix with rows rows and cols columns of elements with size el_size.

Returns
The newly allocated matrix or NULL on error.
void* spr_math_matrix_free ( void *  mat,
int  rows 
)

Free a matrix with rows rows. If mat is NULL nothing is done.

Returns
NULL.
void* spr_math_matrix_read ( char *  fname,
void *  mat,
int *  rows,
int *  cols,
const SprDT dt 
)

Reads rows rows of cols columns from the stream fname into the matrix mat. If *rows and/or *cols differe from (-1), their values will be checked agains the the values of DIM1 and DIM2 in the file header. When their values equal -1, the correct value (from the header) will be filled in. If mat is a NULL pointer, the matrix will be allocated. If the specified data type dt doesn't agree with the one in the file this function fails.

Returns
The pointer to the filled in matrix on success and NULL on failure.
int spr_math_matrix_write ( const char *  fname,
const void *  mat,
const SprDT dt,
int  Nrows,
int  Ncols 
)

Write the given NrowsxNcols matrix mat to the file fname. Note that the type of the elements dt must be one of the standard types (I32, F32, ...).

Returns
0 on success and -1 on failure.
void spr_rand_init ( SprRandGen *restrict  rg_state,
const uint32_t *  val,
int  Nval 
)

Initialize a random generator. Up till 32 state values val[] are used (more can be specified, but the excess elements are simply ignored). If Nval is negative, the current time and data is used to initialize the random generator.

uint32_t spr_rand_u32 ( SprRandGen *restrict  rg_state)

Get the next 32 (pseudo) random bits and update the state of the random generator rg_state.

Returns
a (pseudo) random 32 bits unsigned integer value.
uint64_t spr_rand_u64 ( SprRandGen *restrict  rg_state)

Get the next 2x32 (pseudo) random bits and update the state of the random generator rg_state.

Returns
a (pseudo) random 64 bits unsigned integer value.
float spr_rand_f32 ( SprRandGen *restrict  rg_state)

Get the next 32 (pseudo) random bits and update the state of the random generator rg_state.

Returns
a (pseudo) random float value in the range [0,1[.
double spr_rand_f64q ( SprRandGen *restrict  rg_state)

Calculate the next 32 (pseudo) random bits and update the state of the random generator rg_state.

Returns
a (pseudo) random float value in the range [0,1[.
Note
only the 32 most significant bits in the mantise are non-zero

see spr_rand_f64() for more precision.

double spr_rand_f64 ( SprRandGen *restrict  rg_state)

Get the next 2x32 (pseudo) random bits and update the state of the random generator rg_state.

Returns
a (pseudo) random float value in the range [0,1[.
double spr_randn_f64 ( SprRandGen *restrict  rg_state)

Get the next 2x32 (pseudo) random bits and update the state of the random generator rg_state.

Returns
a (pseudo) random float value which is normal distributed (but clipped to the region [-9.1553,9.1553] due to working with 64bits).
unsigned int spr_rand_modN ( unsigned int  n,
SprRandGen *restrict  rg_state 
)

Generate a uniformly distributed random unsigned integer on the interval [0,n[ with no bias.

Returns
A random unsigned integer.
float* spr_rand_vec_f32 ( float  x[],
int  n,
float  rmax,
SprRandGen *restrict  rg_state 
)

Create a vector of uniformly distributed float numbers in the range [0,rmax[ by drawing random bits from the random generator rg_state.

Returns
pointer to the vector x.
float* spr_randn_vec_f32 ( float  x[],
int  n,
SprRandGen *restrict  rg_state 
)

Create a vector of normal distributed float numbers (mean 0, variance 1.0) by drawing random bits from the random generator rg_state.

Returns
pointer to the vector x.
Note
the values are clipped to the region [-9.1553,9.1553] due to working with 64bits.
int spr_rand_m ( int  indx[],
int  n,
int  m,
SprRandGen *restrict  rg_state 
)

Creates an ordered random index vector of length n in the range [0 .. m-1] by drawing random bits from the random generator rg_state. All numbers are unique. If n > m, only m numbers will be returned!

Returns
the number of random numbers, i.e. n or m when (n > m) on success and (-1) on failure (out of memory).
float spr_math_vector_dist2 ( int  n,
const float *  x,
const float *  y 
)

Calculate the square of the sum of the differences of the first n components of the arrays x[] and y[].

float spr_math_vec_dist ( int  n,
const float *  x,
const float *  y 
)

Calculate the square root of the sum of the squares of the differences of the first n components of the arrays x[] and y[].

void spr_math_vec_pme ( int  n,
float *  x,
float *  y,
float *  z 
)

Calculates the pairwise products of the first n components of x[] and y[] and put the result in z[]; z[] may be the same as x[] or y[].

void spr_math_vec_init ( int  n,
float *  x,
float  a 
)

Initializes the first n elements of x[] with a.

void spr_math_vec_add ( int  n,
float *  x,
float *  y,
float *  z 
)

Calculate the pairwise sums of the first n components of x and y and puts the result in z; z may be the same as x or y.

Parameters
nin: nr of elements
xin: vector x
yin: vector y
zout: sum of x and y
int spr_math_vec_findmax ( int  n,
const float *  x 
)

out: index of highest value

Returns
the index of the maximal parameter of x[0...n-1].
Parameters
nin: nr of elements of array
xin: vector
int spr_math_vec_findmin ( int  n,
const float *  x 
)

out: index of lowest value

Returns
the index of the minimal parameter of x[0...n-1].
Parameters
nin: nr of elements of array
xin: vector

Variable Documentation

const SprFMIntA spr_fm_msk_abs
const SprFMIntA spr_fm_msk_sgn
const SprFMFltA spr_fm_vec_ones
const SprFMIntA spr_fm_vec_pos
const char spr_f32_fft_implementation[]
const char spr_f64_fft_implementation[]
SprRandGen spr_rand_default