Zephyr Scientific Library (zscilib)
|
Go to the documentation of this file.
22 #ifndef ZEPHYR_INCLUDE_ZSL_MATRICES_H_
23 #define ZEPHYR_INCLUDE_ZSL_MATRICES_H_
41 #define EEIGENSIZE (100)
43 #define ECOMPLEXVAL (101)
61 #define ZSL_MATRIX_DEF(name, m, n); \
62 zsl_real_t name ## _mtx[m * n]; \
63 struct zsl_mtx name = { \
66 .data = name ## _mtx \
157 struct zsl_mtx *mc,
size_t i,
size_t j);
594 #ifndef CONFIG_ZSL_SINGLE_PRECISION
694 struct zsl_mtx *mi,
size_t i,
size_t j);
873 #ifndef CONFIG_ZSL_SINGLE_PRECISION
890 #ifndef CONFIG_ZSL_SINGLE_PRECISION
913 #ifndef CONFIG_ZSL_SINGLE_PRECISION
940 #ifndef CONFIG_ZSL_SINGLE_PRECISION
956 struct zsl_mtx *v,
size_t iter);
959 #ifndef CONFIG_ZSL_SINGLE_PRECISION
int zsl_mtx_augm_diag(struct zsl_mtx *m, struct zsl_mtx *maug)
Augments the input square matrix with additional rows and columns, based on the size 'diff' between m...
int zsl_mtx_mult_d(struct zsl_mtx *ma, struct zsl_mtx *mb)
Multiplies matrix 'ma' by 'mb', assigning the output to 'ma'. Matrices 'ma' and 'mb' must be compatib...
int zsl_mtx_svd(struct zsl_mtx *m, struct zsl_mtx *u, struct zsl_mtx *e, struct zsl_mtx *v, size_t iter)
Performs singular value decomposition, converting input matrix 'm' into matrices 'u',...
int zsl_mtx_sum_rows_scaled_d(struct zsl_mtx *m, size_t i, size_t j, zsl_real_t s)
This function takes the coefficients of row 'j' and multiplies them by scalar 's',...
int zsl_mtx_scalar_mult_row_d(struct zsl_mtx *m, size_t i, zsl_real_t s)
Multiplies the elements of row 'i' in matrix 'm' by scalar 's'.
int zsl_mtx_set(struct zsl_mtx *m, size_t i, size_t j, zsl_real_t x)
Sets a single value at the specified row (i) and column (j).
int zsl_mtx_get_row(struct zsl_mtx *m, size_t i, zsl_real_t *v)
Gets the contents of row 'i' from matrix 'm', assigning the array of data to 'v'.
int zsl_mtx_norm_elem_d(struct zsl_mtx *m, struct zsl_mtx *mi, size_t i, size_t j)
Normalises elements in matrix m such that the element at position (i, j) is equal to 1....
int zsl_mtx_cholesky(struct zsl_mtx *m, struct zsl_mtx *l)
Calculates the Cholesky decomposition of a symmetric square matrix using the Cholesky–Crout algorithm...
int zsl_mtx_get(struct zsl_mtx *m, size_t i, size_t j, zsl_real_t *x)
Gets a single value from the specified row (i) and column (j).
int(* zsl_mtx_unary_fn_t)(struct zsl_mtx *m, size_t i, size_t j)
Function prototype called when applying a set of component-wise unary operations to a matrix via zsl_...
Represents a m x n matrix, with data stored in row-major order.
int zsl_mtx_qrd(struct zsl_mtx *m, struct zsl_mtx *q, struct zsl_mtx *r, bool hessenberg)
If 'hessenberg' is set to false, this function performs the QR decomposition, which is a factorisatio...
int zsl_mtx_max(struct zsl_mtx *m, zsl_real_t *x)
Traverses the matrix elements to find the maximum element value.
int(* zsl_mtx_init_entry_fn_t)(struct zsl_mtx *m, size_t i, size_t j)
Function prototype called when populating a matrix via zsl_mtx_init.
int zsl_mtx_reduce_iter(struct zsl_mtx *m, struct zsl_mtx *mred)
Reduces the number of rows/columns in the input square matrix 'm' to match the shape of 'mred',...
@ ZSL_MTX_BINARY_OP_GREAT
int zsl_mtx_sub_d(struct zsl_mtx *ma, struct zsl_mtx *mb)
Subtracts matrix 'mb' from 'ma', assigning the output to 'ma'. Matrices 'ma', and 'mb' must be identi...
int zsl_mtx_min(struct zsl_mtx *m, zsl_real_t *x)
Traverses the matrix elements to find the minimum element value.
int zsl_mtx_eigenvectors(struct zsl_mtx *m, struct zsl_mtx *mev, size_t iter, bool orthonormal)
Calcualtes the set of eigenvectors for input matrix 'm', using the specified number of iterations to ...
enum zsl_mtx_binary_op zsl_mtx_binary_op_t
Component-wise binary operations.
int zsl_mtx_mult(struct zsl_mtx *ma, struct zsl_mtx *mb, struct zsl_mtx *mc)
Multiplies matrix 'ma' by 'mb', assigning the output to 'mc'. Matrices 'ma' and 'mb' must be compatib...
@ ZSL_MTX_BINARY_OP_EQUAL
int zsl_mtx_get_col(struct zsl_mtx *m, size_t j, zsl_real_t *v)
Gets the contents of column 'j' from matrix 'm', assigning the array of data to 'v'.
int zsl_mtx_deter_3x3(struct zsl_mtx *m, zsl_real_t *d)
Calculates the determinant of the input 3x3 matrix 'm'.
int zsl_mtx_min_idx(struct zsl_mtx *m, size_t *i, size_t *j)
Traverses the matrix elements to find the (i,j) index of the minimum element value....
int zsl_mtx_gauss_reduc(struct zsl_mtx *m, struct zsl_mtx *mi, struct zsl_mtx *mg)
Given matrix 'm', puts the matrix into echelon form using Gauss-Jordan reduction.
int zsl_mtx_qrd_iter(struct zsl_mtx *m, struct zsl_mtx *mout, size_t iter)
Computes recursively the QR decompisition method to put the input square matrix into upper triangular...
enum zsl_mtx_unary_op zsl_mtx_unary_op_t
Component-wise unary operations.
int zsl_mtx_balance(struct zsl_mtx *m, struct zsl_mtx *mout)
Balances the square matrix 'm', a process in which the eigenvalues of the output matrix are the same ...
int zsl_mtx_deter(struct zsl_mtx *m, zsl_real_t *d)
Calculates the determinant of the input square matrix 'm'.
int zsl_mtx_add_d(struct zsl_mtx *ma, struct zsl_mtx *mb)
Adds matrices 'ma' and 'mb', assigning the output to 'ma'. Matrices 'ma', and 'mb' must be identicall...
int zsl_mtx_eigenvalues(struct zsl_mtx *m, struct zsl_vec *v, size_t iter)
Calculates the eigenvalues for input matrix 'm' using QR decomposition recursively....
int zsl_mtx_gauss_elim(struct zsl_mtx *m, struct zsl_mtx *mg, struct zsl_mtx *mi, size_t i, size_t j)
Given the element (i,j) in matrix 'm', this function performs gaussian elimination by adding row 'i' ...
int zsl_mtx_entry_fn_empty(struct zsl_mtx *m, size_t i, size_t j)
Assigns a zero-value to all entries in the matrix.
int zsl_mtx_vec_wedge(struct zsl_mtx *m, struct zsl_vec *v)
Calculates the wedge product of n-1 vectors of size n, which are the rows of the matrix 'm'....
int zsl_mtx_sum_rows_d(struct zsl_mtx *m, size_t i, size_t j)
Adds the values of row 'j' to row 'i' in matrix 'm'. This operation is destructive for row 'i'.
int zsl_mtx_init(struct zsl_mtx *m, zsl_mtx_init_entry_fn_t entry_fn)
Initialises matrix 'm' using the specified entry function to assign values.
int zsl_mtx_pinv(struct zsl_mtx *m, struct zsl_mtx *pinv, size_t iter)
Performs the pseudo-inverse (aka pinv or Moore-Penrose inverse) on input matrix 'm'.
int zsl_mtx_add(struct zsl_mtx *ma, struct zsl_mtx *mb, struct zsl_mtx *mc)
Adds matrices 'ma' and 'mb', assigning the output to 'mc'. Matrices 'ma', 'mb' and 'mc' must all be i...
int zsl_mtx_max_idx(struct zsl_mtx *m, size_t *i, size_t *j)
Traverses the matrix elements to find the (i,j) index of the maximum element value....
@ ZSL_MTX_BINARY_OP_EXPON
int zsl_mtx_cols_norm(struct zsl_mtx *m, struct zsl_mtx *mnorm)
Updates the values of every column vector in input matrix 'm' to have unitary values.
int zsl_mtx_trans(struct zsl_mtx *ma, struct zsl_mtx *mb)
Transposes the matrix 'ma' into matrix 'mb'. Note that output matrix 'mb' must have 'ma->sz_rows' col...
int zsl_mtx_entry_fn_identity(struct zsl_mtx *m, size_t i, size_t j)
Sets the value to '1.0' if the entry is on the diagonal (row=col), otherwise '0.0'.
int zsl_mtx_entry_fn_random(struct zsl_mtx *m, size_t i, size_t j)
Sets the value to a random number between -1.0 and 1.0.
int(* zsl_mtx_binary_fn_t)(struct zsl_mtx *ma, struct zsl_mtx *mb, struct zsl_mtx *mc, size_t i, size_t j)
Function prototype called when applying a set of component-wise binary operations using a pair of sym...
int zsl_mtx_gauss_elim_d(struct zsl_mtx *m, struct zsl_mtx *mi, size_t i, size_t j)
Given the element (i,j) in matrix 'm', this function performs gaussian elimination by adding row 'i' ...
int zsl_mtx_binary_func(struct zsl_mtx *ma, struct zsl_mtx *mb, struct zsl_mtx *mc, zsl_mtx_binary_fn_t fn)
Applies a component-wise binary operztion on every coefficient in symmetrical matrices 'ma' and 'mb',...
int zsl_mtx_copy(struct zsl_mtx *mdest, struct zsl_mtx *msrc)
Copies the contents of matrix 'msrc' into matrix 'mdest'.
API header file for vectors in zscilib.
int zsl_mtx_norm_elem(struct zsl_mtx *m, struct zsl_mtx *mn, struct zsl_mtx *mi, size_t i, size_t j)
Normalises elements in matrix m such that the element at position (i, j) is equal to 1....
bool zsl_mtx_is_sym(struct zsl_mtx *m)
Checks if the square input matrix is symmetric.
int zsl_mtx_gram_schmidt(struct zsl_mtx *m, struct zsl_mtx *mort)
Performs the Gram-Schmidt algorithm on the set of column vectors in matrix 'm'. This algorithm calcul...
int zsl_mtx_set_row(struct zsl_mtx *m, size_t i, zsl_real_t *v)
Sets the contents of row 'i' in matrix 'm', assigning the values found in array 'v'.
int zsl_mtx_inv_3x3(struct zsl_mtx *m, struct zsl_mtx *mi)
Calculates the inverse of 3x3 matrix 'm'. If the determinant of 'm' is zero, an identity matrix will ...
int zsl_mtx_print(struct zsl_mtx *m)
Printf the supplied matrix using printf in a human-readable manner.
bool zsl_mtx_is_equal(struct zsl_mtx *ma, struct zsl_mtx *mb)
Checks if two matrices are identical in shape and content.
bool zsl_mtx_is_notneg(struct zsl_mtx *m)
Checks if all elements in matrix m are >= zero.
API header file for zscilib.
int zsl_mtx_scalar_mult_d(struct zsl_mtx *m, zsl_real_t s)
Multiplies all elements in matrix 'm' by scalar value 's'.
int zsl_mtx_unary_func(struct zsl_mtx *m, zsl_mtx_unary_fn_t fn)
Applies a unary function on every coefficient in matrix 'm', using the specified 'zsl_mtx_apply_unary...
int zsl_mtx_householder(struct zsl_mtx *m, struct zsl_mtx *h, bool hessenberg)
Calculates the householder reflection of 'm'. Used as part of QR decomposition. When 'hessenberg' is ...
int zsl_mtx_from_arr(struct zsl_mtx *m, zsl_real_t *a)
Converts an array of values into a matrix.
int zsl_mtx_set_col(struct zsl_mtx *m, size_t j, zsl_real_t *v)
Sets the contents of column 'j' in matrix 'm', assigning the values found in array 'v'.
@ ZSL_MTX_UNARY_OP_NEGATIVE
int zsl_mtx_adjoint(struct zsl_mtx *m, struct zsl_mtx *ma)
Calculates the ajoint matrix, based on the input square matrix 'm'.
int zsl_mtx_inv(struct zsl_mtx *m, struct zsl_mtx *mi)
Calculates the inverse of square matrix 'm'.
@ ZSL_MTX_BINARY_OP_NEQUAL
@ ZSL_MTX_UNARY_OP_INCREMENT
int zsl_mtx_adjoint_3x3(struct zsl_mtx *m, struct zsl_mtx *ma)
Calculates the ajoint matrix, based on the input 3x3 matrix 'm'.
int zsl_mtx_sub(struct zsl_mtx *ma, struct zsl_mtx *mb, struct zsl_mtx *mc)
Subtracts matrices 'mb' from 'ma', assigning the output to 'mc'. Matrices 'ma', 'mb' and 'mc' must al...
@ ZSL_MTX_UNARY_OP_DECREMENT
int zsl_mtx_reduce(struct zsl_mtx *m, struct zsl_mtx *mr, size_t i, size_t j)
Removes row 'i' and column 'j' from square matrix 'm', assigning the remaining elements in the matrix...
zsl_mtx_unary_op
Component-wise unary operations.
zsl_mtx_binary_op
Component-wise binary operations.